Functions¶
Solution of least squares by the normal equations
1 2 3 4 5 6 7 8 9 10 11 12 13 14""" lsnormal(A, b) Solve a linear least-squares problem by the normal equations. Returns the minimizer of ||b-Ax||. """ function lsnormal(A, b) N = A' * A z = A' * b R = cholesky(N).U w = forwardsub(R', z) # solve R'z=c x = backsub(R, w) # solve Rx=z return x end
About the code
The syntax on line 9 is a field reference to extract the matrix we want from the structure returned by cholesky.
Solution of least squares by QR factorization
1 2 3 4 5 6 7 8 9 10 11 12""" lsqrfact(A, b) Solve a linear least-squares problem by QR factorization. Returns the minimizer of ||b-Ax||. """ function lsqrfact(A, b) Q, R = qr(A) c = Q' * b x = backsub(R, c) return x end
QR factorization by Householder reflections
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27""" qrfact(A) QR factorization by Householder reflections. Returns Q and R. """ function qrfact(A) m, n = size(A) Qt = diagm(ones(m)) R = float(copy(A)) for k in 1:n z = R[k:m, k] w = [-sign(z[1]) * norm(z) - z[1]; -z[2:end]] nrmw = norm(w) if nrmw < eps() continue # already in place; skip this iteration end v = w / nrmw # Apply the reflection to each relevant column of R and Q for j in k:n R[k:m, j] -= v * (2 * (v' * R[k:m, j])) end for j in 1:m Qt[k:m, j] -= v * (2 * (v' * Qt[k:m, j])) end end return Qt', triu(R) end
Examples¶
3.1 Fitting functions to data¶
Example 3.1.1
Here are 5-year averages of the worldwide temperature anomaly as compared to the 1951–1980 average (source: NASA).
year = 1955:5:2000
temp = [ -0.0480, -0.0180, -0.0360, -0.0120, -0.0040,
0.1180, 0.2100, 0.3320, 0.3340, 0.4560 ]
scatter(year, temp, label="data",
xlabel="year", ylabel="anomaly (degrees C)",
legend=:bottomright)A polynomial interpolant can be used to fit the data. Here we build one using a Vandermonde matrix. First, though, we express time as decades since 1950, as it improves the condition number of the matrix.
t = @. (year - 1950) / 10
n = length(t)
V = [ t[i]^j for i in 1:n, j in 0:n-1 ]
c = V \ temp10-element Vector{Float64}:
-14.114000001832462
76.36173810552113
-165.45597224550528
191.96056669514388
-133.27347224319684
58.015577787494486
-15.962888891734785
2.6948063497166928
-0.2546666667177082
0.010311111113288083The coefficients in vector c are used to create a polynomial. Then we create a function that evaluates the polynomial after changing the time variable as we did for the Vandermonde matrix.
Tip
If you plot a function, then the points are chosen automatically to make a smooth curve.
using Polynomials, Plots
p = Polynomial(c)
f = yr -> p((yr - 1950) / 10)
plot!(f, 1955, 2000, label="interpolant")As you can see, the interpolant does represent the data, in a sense. However it’s a crazy-looking curve for the application. Trying too hard to reproduce all the data exactly is known as overfitting.
Example 3.1.2
Here are the 5-year temperature averages again.
year = 1955:5:2000
t = @. (year - 1950) / 10
temp = [ -0.0480, -0.0180, -0.0360, -0.0120, -0.0040,
0.1180, 0.2100, 0.3320, 0.3340, 0.4560 ]10-element Vector{Float64}:
-0.048
-0.018
-0.036
-0.012
-0.004
0.118
0.21
0.332
0.334
0.456The standard best-fit line results from using a linear polynomial that meets the least-squares criterion.
Tip
Backslash solves overdetermined linear systems in a least-squares sense.
V = [ t.^0 t ] # Vandermonde-ish matrix
@show size(V)
c = V \ temp
p = Polynomial(c)size(V) = (10, 2)
f = yr -> p((yr - 1955) / 10)
scatter(year, temp, label="data",
xlabel="year", ylabel="anomaly (degrees C)", leg=:bottomright)
plot!(f, 1955, 2000, label="linear fit")If we use a global cubic polynomial, the points are fit more closely.
V = [ t[i]^j for i in 1:length(t), j in 0:3 ]
@show size(V);size(V) = (10, 4)
Now we solve the new least-squares problem to redefine the fitting polynomial.
Tip
The definition of f above is in terms of p. When p is changed, then f calls the new version.
p = Polynomial( V \ temp )
plot!(f, 1955, 2000, label="cubic fit")If we were to continue increasing the degree of the polynomial, the residual at the data points would get smaller, but overfitting would increase.
Example 3.1.3
a = [1/k^2 for k=1:100]
s = cumsum(a) # cumulative summation
p = @. sqrt(6*s)
scatter(1:100, p;
title="Sequence convergence",
xlabel=L"k", ylabel=L"p_k")This graph suggests that maybe , but it’s far from clear how close the sequence gets. It’s more informative to plot the sequence of errors, . By plotting the error sequence on a log-log scale, we can see a nearly linear relationship.
ϵ = @. abs(π - p) # error sequence
scatter(1:100, ϵ;
title="Convergence of errors",
xaxis=(:log10,L"k"), yaxis=(:log10,"error"))The straight line on the log-log scale suggests a power-law relationship where , or .
k = 1:100
V = [ k.^0 log.(k) ] # fitting matrix
c = V \ log.(ϵ) # coefficients of linear fit2-element Vector{Float64}:
-0.1823752497282998
-0.967410323312793In terms of the parameters and used above, we have
a, b = exp(c[1]), c[2];
@show b;b = -0.967410323312793
It’s tempting to conjecture that the slope asymptotically. Here is how the numerical fit compares to the original convergence curve.
plot!(k, a * k.^b, l=:dash, label="power-law fit")3.2 The normal equations¶
Example 3.2.1
Because the functions , , and 1 are linearly dependent, we should find that the following matrix is somewhat ill-conditioned.
Tip
The local variable scoping rule for loops applies to comprehensions as well.
t = range(0, 3, 400)
f = [ x -> sin(x)^2, x -> cos((1 + 1e-7) * x)^2, x -> 1. ]
A = [ f(t) for t in t, f in f ]
@show κ = cond(A);κ = cond(A) = 1.8253225428206295e7
Now we set up an artificial linear least-squares problem with a known exact solution that actually makes the residual zero.
x = [1., 2, 1]
b = A * x;Using backslash to find the least-squares solution, we get a relative error that is well below times machine epsilon.
x_BS = A \ b
@show observed_error = norm(x_BS - x) / norm(x);
@show error_bound = κ * eps();observed_error = norm(x_BS - x) / norm(x) = 1.7978987436098622e-10
error_bound = κ * eps() = 4.053030228813602e-9
If we formulate and solve via the normal equations, we get a much larger relative error. With , we may not be left with more than about 2 accurate digits.
N = A' * A
x_NE = N \ (A'*b)
@show observed_err = norm(x_NE - x) / norm(x);
@show digits = -log10(observed_err);observed_err = norm(x_NE - x) / norm(x) = 0.036205099396071284
digits = -(log10(observed_err)) = 1.441230255886605
3.3 The QR factorization¶
Example 3.3.1
Julia provides access to both the thin and full forms of the QR factorization.
A = rand(1.:9., 6, 4)
@show m,n = size(A);(m, n) = size(A) = (6, 4)
Here is a standard call:
Q,R = qr(A);
Q6×6 LinearAlgebra.QRCompactWYQ{Float64, Matrix{Float64}, Matrix{Float64}}R4×4 Matrix{Float64}:
-18.5203 -10.853 -7.93725 -15.9285
0.0 -4.38324 -0.19555 -2.53882
0.0 0.0 -5.09527 -5.70627
0.0 0.0 0.0 -5.85454If you look carefully, you see that we seemingly got a full but a thin . However, the above is not a standard matrix type. If you convert it to a true matrix, then it reverts to the thin form.
Tip
To enter the accented character Q̂, type Q\hat followed by Tab.
Q̂ = Matrix(Q)6×4 Matrix{Float64}:
-0.431959 0.841397 0.24808 -0.114664
-0.377964 -0.204862 -0.384659 0.125626
-0.10799 0.039243 -0.618325 -0.657812
-0.485954 -0.39376 0.379596 0.0980669
-0.431959 -0.299311 0.291859 -0.516703
-0.485954 0.0625228 -0.422958 0.51162We can test that is an orthogonal matrix:
opnorm(Q' * Q - I)4.617485372336832e-16The thin cannot be an orthogonal matrix, because it is not square, but it is still ONC:
Q̂' * Q̂ - I4×4 Matrix{Float64}:
-1.11022e-16 4.78643e-17 9.20989e-17 -3.82932e-18
4.78643e-17 -4.44089e-16 -1.9853e-16 -6.07525e-17
9.20989e-17 -1.9853e-16 -2.22045e-16 -3.01771e-17
-3.82932e-18 -6.07525e-17 -3.01771e-17 2.22045e-16Example 3.3.2
We’ll repeat the experiment of Example 3.2.1, which exposed instability in the normal equations.
t = range(0, 3, 400)
f = [ x -> sin(x)^2, x -> cos((1 + 1e-7) * x)^2, x -> 1. ]
A = [ f(t) for t in t, f in f ]
x = [1., 2, 1]
b = A * x;The error in the solution by Function 3.3.2 is similar to the bound predicted by the condition number.
observed_error = norm(FNC.lsqrfact(A, b) - x) / norm(x);
@show observed_error;
@show error_bound = cond(A) * eps();observed_error = 4.516380901323573e-10
error_bound = cond(A) * eps() = 4.053030228813602e-9
3.4 Computing QR factorizations¶
Example 3.4.1
We will use Householder reflections to produce a QR factorization of a random matrix.
Tip
The rand function can select randomly from within the interval , or from a vector or range that you specify.
A = rand(float(1:9), 6, 4)
m,n = size(A)(6, 4)Our first step is to introduce zeros below the diagonal in column 1 by using (3.4.4) and (3.4.1).
Tip
I can stand for an identity matrix of any size, inferred from the context when needed.
z = A[:, 1];
v = normalize(z - norm(z) * [1; zeros(m-1)])
P₁ = I - 2v * v' # reflector6×6 Matrix{Float64}:
0.306786 0.460179 0.460179 0.230089 0.613572 0.230089
0.460179 0.694518 -0.305482 -0.152741 -0.40731 -0.152741
0.460179 -0.305482 0.694518 -0.152741 -0.40731 -0.152741
0.230089 -0.152741 -0.152741 0.923629 -0.203655 -0.0763706
0.613572 -0.40731 -0.40731 -0.203655 0.45692 -0.203655
0.230089 -0.152741 -0.152741 -0.0763706 -0.203655 0.923629We check that this reflector introduces zeros as it should:
P₁ * z6-element Vector{Float64}:
13.038404810405298
1.1102230246251565e-15
1.1102230246251565e-15
3.3306690738754696e-16
-2.7755575615628914e-17
3.3306690738754696e-16Now we replace by .
A = P₁ * A6×4 Matrix{Float64}:
13.0384 10.5841 11.7346 4.14161
1.11022e-15 1.96541 5.85704 2.8975
1.11022e-15 -3.03459 2.85704 2.8975
3.33067e-16 5.4827 2.42852 1.94875
-2.77556e-17 -3.71279 -3.19062 3.52999
3.33067e-16 5.4827 3.42852 1.94875We are set to put zeros into column 2. We must not use row 1 in any way, lest it destroy the zeros we just introduced. So we leave it out of the next reflector.
z = A[2:m, 2]
v = normalize(z - norm(z) * [1; zeros(m-2)])
P₂ = I - 2v * v'5×5 Matrix{Float64}:
0.210742 -0.325387 0.587887 -0.398107 0.587887
-0.325387 0.865853 0.242368 -0.164127 0.242368
0.587887 0.242368 0.562106 0.296534 -0.437894
-0.398107 -0.164127 0.296534 0.799192 0.296534
0.587887 0.242368 -0.437894 0.296534 0.562106We now apply this reflector to rows 2 and below only.
A[2:m, :] = P₂ * A[2:m, :]
A6×4 Matrix{Float64}:
13.0384 10.5841 11.7346 4.14161
2.7538e-16 9.32612 5.01817 0.553789
7.66043e-16 -6.95149e-16 2.5112 1.93126
9.54908e-16 2.13289e-15 3.05336 3.69448
-4.48856e-16 -1.48753e-15 -3.61375 2.34781
9.54908e-16 2.13289e-15 4.05336 3.69448We need to iterate the process for the last two columns.
for j in 3:n
z = A[j:m, j]
v = normalize(z - norm(z) * [1; zeros(m-j)])
P = I - 2v * v'
A[j:m, :] = P * A[j:m, :]
endWe have now reduced the original to an upper triangular matrix using four orthogonal Householder reflections:
R = triu(A)6×4 Matrix{Float64}:
13.0384 10.5841 11.7346 4.14161
0.0 9.32612 5.01817 0.553789
0.0 0.0 6.71699 3.36773
0.0 0.0 0.0 5.01984
0.0 0.0 0.0 0.0
0.0 0.0 0.0 0.0