Note the packages needed in block [1] if not already installed. Most works except for flop counting in Ex. 2.5.5
#import Pkg; Pkg.add("Plots") # Add on first run
#import Pkg; Pkg.add("DataFrames") # Add on first run
using Plots
using Polynomials
using LinearAlgebra,SparseArrays
using DataFrames
include("FNC.jl")
Main.FNC
We create two vectors for data about the population of China. The first has the years of census data, the other has the numbers of millions of people.
year = 1980:10:2010
pop = [984.736, 1148.364, 1263.638, 1330.141];
It's convenient to measure time in years since 1980. We use .-
to subtract a scalar from a vector elementwise.
t = year .- 1980
y = pop;
Now we have four data points $(t_1,y_1),\dots,(t_4,y_4)$, so $n=4$ and we seek an interpolating cubic polynomial. We construct the associated Vandermonde matrix:
V = [ t[i]^j for i=1:4, j=0:3 ]
4×4 Matrix{Int64}: 1 0 0 0 1 10 100 1000 1 20 400 8000 1 30 900 27000
To solve for the vector of polynomial coefficients, we use a backslash:
c = V \ y
4-element Vector{Float64}: 984.736 18.76660000000003 -0.2396850000000029 -6.949999999993039e-5
The algorithms used by the backslash operator are the main topic of this chapter. For now, observe that the coefficients of the cubic polynomial vary over several orders of magnitude, which is typical in this context. By our definitions, these coefficients are given in ascending order of power in $t$.
We can use the resulting polynomial to estimate the population of China in 2005:
p = Polynomial(c) # construct a polynomial from coeffs - changed from poly()
p(2005-1980) # apply the 1980 time shift
1303.0119375
The official figure is 1297.8, so our result is not bad.
We can visualize the interpolation process. First, we plot the data as points. We'll shift the $t$ variable back to actual years.
scatter(1980 .+ t,y,label="actual",
legend=:topleft,xlabel="year",ylabel="population (millions)",title="Population of China")
We want to superimpose a plot of the polynomial. We do that by evaluating it at a lot of points in the interval. Note the use of plot!
to add to the current plot, rather than replacing it.
tt = LinRange(0,30,500) # 500 times from 0 to 30 years
yy = p.(tt) # use dot to apply to all elements of the vector
plot!(1980 .+ tt,yy, label="interpolant")
Let's redo it, this time continuing the curve outside of the original date range.
scatter(1980 .+ t,y,label="actual",
legend=:topleft,xlabel="year",ylabel="population (millions)",title="Population of China")
tt = LinRange(-10,50,300)
plot!(1980 .+ tt,p.(tt),label="interpolant")
While the interpolation is plausible, the extrapolation to the future is highly questionable! As a rule, extrapolation more than a short distance beyond the original interval is not reliable.
Square brackets are used to enclose elements of a matrix or vector. Use spaces for horizontal concatenation, and semicolons or new lines to indicate vertical concatenation.
A = [ 1 2 3 4 5; 50 40 30 20 10
pi sqrt(2) exp(1) (1+sqrt(5))/2 log(3) ]
3×5 Matrix{Float64}: 1.0 2.0 3.0 4.0 5.0 50.0 40.0 30.0 20.0 10.0 3.14159 1.41421 2.71828 1.61803 1.09861
m,n = size(A)
(3, 5)
A vector is not quite the same thing as a matrix. It has only one dimension, not two. Separate its elements by commas.
x = [ 3, 3, 0, 1, 0 ]
size(x)
(5,)
For many purposes, though, an $n$-vector in Julia is a lot like an $n\times 1$ column vector.
size( [3;3;0;1;0] )
(5,)
Concatenated elements within brackets may be matrices for a block representation, as long as all the block sizes are compatible.
AA = [ A; A ]
6×5 Matrix{Float64}: 1.0 2.0 3.0 4.0 5.0 50.0 40.0 30.0 20.0 10.0 3.14159 1.41421 2.71828 1.61803 1.09861 1.0 2.0 3.0 4.0 5.0 50.0 40.0 30.0 20.0 10.0 3.14159 1.41421 2.71828 1.61803 1.09861
B = [ zeros(3,2) ones(3,1) ]
3×3 Matrix{Float64}: 0.0 0.0 1.0 0.0 0.0 1.0 0.0 0.0 1.0
The dot-quote .'
transposes a matrix. A single quote '
on its own performs the hermitian (transpose and complex conjugation). For a real matrix, the two operations are the same.
A'
5×3 adjoint(::Matrix{Float64}) with eltype Float64: 1.0 50.0 3.14159 2.0 40.0 1.41421 3.0 30.0 2.71828 4.0 20.0 1.61803 5.0 10.0 1.09861
If x
is simply a vector, then its transpose has a row shape.
x'
1×5 adjoint(::Vector{Int64}) with eltype Int64: 3 3 0 1 0
There are many convenient shorthand ways of building vectors and matrices other than entering all of their entries directly or in a loop. To get a vector with evenly spaced entries between two endpoints, you have two options.
y = 1:4 # start:stop
1:4
z = ( 0:3:12 )' # start:step:stop
1×5 adjoint(::StepRange{Int64, Int64}) with eltype Int64: 0 3 6 9 12
(Technically, y
above is not a vector but a range. It behaves identically in most circumstances.)
s = LinRange(-1,1,5)
5-element LinRange{Float64, Int64}: -1.0,-0.5,0.0,0.5,1.0
Accessing an element is done by giving one (for a vector) or two index values in square brackets. The keyword end
as an index refers to the last position in the corresponding dimension.
a = A[2,end-1]
20.0
x[2]
3
The indices can be vectors or ranges, in which case a block of the matrix is accessed.
A[1:2,end-2:end] # first two rows, last three columns
2×3 Matrix{Float64}: 3.0 4.0 5.0 30.0 20.0 10.0
If a dimension has only the index :
(a colon), then it refers to all the entries in that dimension of the matrix.
A[:,1:2:end] # all of the odd columns
3×3 Matrix{Float64}: 1.0 3.0 5.0 50.0 30.0 10.0 3.14159 2.71828 1.09861
The matrix and vector senses of addition, subtraction, scalar multiplication, multiplication, and power are all handled by the usual symbols. If matrix sizes are such that the operation is not defined, an error message will result.
B = diagm( 0=>[-1,0,-5] ) # create a diagonal matrix
3×3 Matrix{Int64}: -1 0 0 0 0 0 0 0 -5
BA = B*A # matrix product
3×5 Matrix{Float64}: -1.0 -2.0 -3.0 -4.0 -5.0 0.0 0.0 0.0 0.0 0.0 -15.708 -7.07107 -13.5914 -8.09017 -5.49306
A*B
causes an error.
A*B
DimensionMismatch("matrix A has dimensions (3,5), matrix B has dimensions (3,3)") Stacktrace: [1] _generic_matmatmul!(C::Matrix{Float64}, tA::Char, tB::Char, A::Matrix{Float64}, B::Matrix{Int64}, _add::LinearAlgebra.MulAddMul{true, true, Bool, Bool}) @ LinearAlgebra /opt/julia-1.7.0/share/julia/stdlib/v1.7/LinearAlgebra/src/matmul.jl:810 [2] generic_matmatmul!(C::Matrix{Float64}, tA::Char, tB::Char, A::Matrix{Float64}, B::Matrix{Int64}, _add::LinearAlgebra.MulAddMul{true, true, Bool, Bool}) @ LinearAlgebra /opt/julia-1.7.0/share/julia/stdlib/v1.7/LinearAlgebra/src/matmul.jl:798 [3] mul! @ /opt/julia-1.7.0/share/julia/stdlib/v1.7/LinearAlgebra/src/matmul.jl:302 [inlined] [4] mul! @ /opt/julia-1.7.0/share/julia/stdlib/v1.7/LinearAlgebra/src/matmul.jl:275 [inlined] [5] *(A::Matrix{Float64}, B::Matrix{Int64}) @ LinearAlgebra /opt/julia-1.7.0/share/julia/stdlib/v1.7/LinearAlgebra/src/matmul.jl:153 [6] top-level scope @ In[27]:1 [7] eval @ ./boot.jl:373 [inlined] [8] include_string(mapexpr::typeof(REPL.softscope), mod::Module, code::String, filename::String) @ Base ./loading.jl:1196
A square matrix raised to an integer power is the same as repeated matrix multiplication.
B^3 # same as B*B*B
3×3 Matrix{Int64}: -1 0 0 0 0 0 0 0 -125
In many cases, one instead wants to treat a matrix or vector as a mere array and simply apply a single operation to each element of it. For multiplication, division, and power, the corresponding operators start with a dot.
C = -A;
A*C
would be an error.
elementwise = A.*C
3×5 Matrix{Float64}: -1.0 -4.0 -9.0 -16.0 -25.0 -2500.0 -1600.0 -900.0 -400.0 -100.0 -9.8696 -2.0 -7.38906 -2.61803 -1.20695
The two operands of a dot operator have to have the same size—unless one is a scalar, in which case it is expanded or "broadcast" to be the same size as the other operand.
xtotwo = x.^2
5-element Vector{Int64}: 9 9 0 1 0
twotox = 2.0.^x
5-element Vector{Float64}: 8.0 8.0 1.0 2.0 1.0
Most of the mathematical functions, such as cos, sin, log, exp and sqrt, expect scalars as operands. However, you can broadcast any function across a vector or array by using a special dot syntax.
println(cos.(pi*x)); # vectorize a single function
println(@. cos(pi*x)); # vectorize an entire expression
[-1.0, -1.0, 1.0, -1.0, 1.0] [-1.0, -1.0, 1.0, -1.0, 1.0]
For a square matrix $A$, the command A\b
is mathematically equivalent to $A^{-1}b$. This command is not part of the core Julia, though, so it has to be explicitly loaded before the first use in a session.
A = [1 0 -1; 2 2 1; -1 -3 0]
3×3 Matrix{Int64}: 1 0 -1 2 2 1 -1 -3 0
b = [1,2,3]
3-element Vector{Int64}: 1 2 3
x = A\b
3-element Vector{Float64}: 2.1428571428571432 -1.7142857142857144 1.1428571428571428
One way to check the answer is to compute a quantity known as the residual. It is (hopefully) close to machine precision, scaled by the size of the entries of the data.
residual = b - A*x
3-element Vector{Float64}: -4.440892098500626e-16 -4.440892098500626e-16 0.0
If the matrix $A$ is singular, you may get an error ("exception" in Julia-speak).
A = [0 1; 0 0]
b = [1,-1]
x = A\b
LAPACKException(1) Stacktrace: [1] chklapackerror(ret::Int64) @ LinearAlgebra.LAPACK /opt/julia-1.7.0/share/julia/stdlib/v1.7/LinearAlgebra/src/lapack.jl:43 [2] trtrs!(uplo::Char, trans::Char, diag::Char, A::Matrix{Float64}, B::Vector{Float64}) @ LinearAlgebra.LAPACK /opt/julia-1.7.0/share/julia/stdlib/v1.7/LinearAlgebra/src/lapack.jl:3431 [3] ldiv! @ /opt/julia-1.7.0/share/julia/stdlib/v1.7/LinearAlgebra/src/triangular.jl:727 [inlined] [4] \(A::UpperTriangular{Int64, Matrix{Int64}}, B::Vector{Int64}) @ LinearAlgebra /opt/julia-1.7.0/share/julia/stdlib/v1.7/LinearAlgebra/src/triangular.jl:1661 [5] \(A::Matrix{Int64}, B::Vector{Int64}) @ LinearAlgebra /opt/julia-1.7.0/share/julia/stdlib/v1.7/LinearAlgebra/src/generic.jl:1140 [6] top-level scope @ In[38]:3 [7] eval @ ./boot.jl:373 [inlined] [8] include_string(mapexpr::typeof(REPL.softscope), mod::Module, code::String, filename::String) @ Base ./loading.jl:1196
It's not exactly user-friendly here. Moreover, detecting singularity is a lot like checking whether two floating point numbers are exactly equal: because of roundoff, it could be missed. We're headed toward a more robust way to fully describe the situation.
It's easy to get just the lower triangular part of any matrix using the tril
command.
A = rand(1.:9.,5,5)
L = tril(A)
5×5 Matrix{Float64}: 8.0 0.0 0.0 0.0 0.0 9.0 7.0 0.0 0.0 0.0 4.0 5.0 2.0 0.0 0.0 2.0 3.0 7.0 3.0 0.0 7.0 9.0 2.0 8.0 6.0
We'll set up and solve a linear system with this matrix.
b = ones(5)
x = FNC.forwardsub(L,b)
5-element Vector{Float64}: 0.125 -0.017857142857142856 0.29464285714285715 -0.4196428571428572 0.5089285714285715
It's not clear what the error in this answer is. However, the residual, while not zero, is comparable to $\varepsilon_\text{mach}$ in size.
b - L*x
5-element Vector{Float64}: 0.0 0.0 -2.220446049250313e-16 2.220446049250313e-16 -2.220446049250313e-16
Next we'll engineer a problem to which we know the exact answer.
alpha = 0.3;
beta = 2.2;
U = diagm(0=>ones(5),1=>[-1,-1,-1,-1])
U[1,[4,5]] = [ alpha-beta, beta ]
U
5×5 Matrix{Float64}: 1.0 -1.0 0.0 -1.9 2.2 0.0 1.0 -1.0 0.0 0.0 0.0 0.0 1.0 -1.0 0.0 0.0 0.0 0.0 1.0 -1.0 0.0 0.0 0.0 0.0 1.0
x_exact = ones(5)
b = [alpha,0,0,0,1]
x = FNC.backsub(U,b)
err = x - x_exact
5-element Vector{Float64}: 2.220446049250313e-16 0.0 0.0 0.0 0.0
Everything seems OK here. But another example, with a different value for $\beta$, is more troubling.
alpha = 0.3;
beta = 1e12;
U = diagm(0=>ones(5),1=>[-1,-1,-1,-1])
U[1,[4,5]] = [ alpha-beta, beta ]
b = [alpha,0,0,0,1]
x = FNC.backsub(U,b)
err = x - x_exact
5-element Vector{Float64}: -4.882812499995559e-5 0.0 0.0 0.0 0.0
It's not so good to get four digits of accuracy after starting with sixteen! But the source of the error is not hard to track down. Solving for $x_1$ performs $(\alpha-\beta)+\beta$ in the first row. Since $|\alpha|$ is so much smaller than $|\beta|$, this a recipe for losing digits to subtractive cancellation.
We create a 4-by-4 linear system with the matrix
A = [
2 0 4 3
-4 5 -7 -10
1 15 2 -4.5
-2 0 2 -13
];
and with the right-hand side
b = [ 4, 9, 29, 40 ];
We define an augmented matrix by tacking $b$ on the end as a new column.
S = [A b]
4×5 Matrix{Float64}: 2.0 0.0 4.0 3.0 4.0 -4.0 5.0 -7.0 -10.0 9.0 1.0 15.0 2.0 -4.5 29.0 -2.0 0.0 2.0 -13.0 40.0
The goal is to introduce zeros into the lower triangle of this matrix. By using only elementary row operations, we ensure that the matrix $S$ always represents a linear system that is equivalent to the original. We proceed from left to right and top to bottom. The first step is to put a zero in the (2,1) location using a multiple of row 1:
@show mult21 = S[2,1]/S[1,1];
S[2,:] -= mult21*S[1,:]; # -= means "subtract the following from"
S
mult21 = S[2, 1] / S[1, 1] = -2.0
4×5 Matrix{Float64}: 2.0 0.0 4.0 3.0 4.0 0.0 5.0 1.0 -4.0 17.0 1.0 15.0 2.0 -4.5 29.0 -2.0 0.0 2.0 -13.0 40.0
We repeat the process for the (3,1) and (4,1) entries.
@show mult31 = S[3,1]/S[1,1];
S[3,:] -= mult31*S[1,:];
@show mult41 = S[4,1]/S[1,1];
S[4,:] -= mult41*S[1,:];
S
mult31 = S[3, 1] / S[1, 1] = 0.5 mult41 = S[4, 1] / S[1, 1] = -1.0
4×5 Matrix{Float64}: 2.0 0.0 4.0 3.0 4.0 0.0 5.0 1.0 -4.0 17.0 0.0 15.0 0.0 -6.0 27.0 0.0 0.0 6.0 -10.0 44.0
The first column has the zero structure we want. To avoid interfering with that, we no longer add multiples of row 1 to anything. Instead, to handle column 2, we use multiples of row 2. We'll also exploit the highly repetitive nature of the operations to write them as a loop.
for i = 3:4
mult = S[i,2]/S[2,2]
S[i,:] -= mult*S[2,:]
end
S
4×5 Matrix{Float64}: 2.0 0.0 4.0 3.0 4.0 0.0 5.0 1.0 -4.0 17.0 0.0 0.0 -3.0 6.0 -24.0 0.0 0.0 6.0 -10.0 44.0
We finish out the triangularization with a zero in the (4,3) place. It's a little silly to use a loop for just one iteration, but the point is to establish a pattern.
for i = 4
mult = S[i,3]/S[3,3]
S[i,:] -= mult*S[3,:]
end
S
4×5 Matrix{Float64}: 2.0 0.0 4.0 3.0 4.0 0.0 5.0 1.0 -4.0 17.0 0.0 0.0 -3.0 6.0 -24.0 0.0 0.0 0.0 2.0 -4.0
Recall that $S$ is an augmented matrix: it represents the system $Ux=z$, where
U = S[:,1:4]
4×4 Matrix{Float64}: 2.0 0.0 4.0 3.0 0.0 5.0 1.0 -4.0 0.0 0.0 -3.0 6.0 0.0 0.0 0.0 2.0
z = S[:,5]
4-element Vector{Float64}: 4.0 17.0 -24.0 -4.0
The solutions to this system are identical to those of the original system, but this one can be solved by backward substitution.
x = FNC.backsub(U,z)
4-element Vector{Float64}: -3.0 1.0 4.0 -2.0
b - A*x
4-element Vector{Float64}: 0.0 0.0 0.0 0.0
We revisit the previous example using algebra to express the row operations on $A$.
A = [2 0 4 3 ; -4 5 -7 -10 ; 1 15 2 -4.5 ; -2 0 2 -13];
We use the identity and its columns heavily.
I = diagm(0=>ones(4))
4×4 Matrix{Float64}: 1.0 0.0 0.0 0.0 0.0 1.0 0.0 0.0 0.0 0.0 1.0 0.0 0.0 0.0 0.0 1.0
The first step is to put a zero in the (2,1) location using a multiple of row 1:
mult21 = A[2,1]/A[1,1]
L21 = I - mult21*I[:,2]*I[:,1]'
A = L21*A
4×4 Matrix{Float64}: 2.0 0.0 4.0 3.0 0.0 5.0 1.0 -4.0 1.0 15.0 2.0 -4.5 -2.0 0.0 2.0 -13.0
We repeat the process for the (3,1) and (4,1) entries.
mult31 = A[3,1]/A[1,1];
L31 = I - mult31*I[:,3]*I[:,1]';
A = L31*A;
mult41 = A[4,1]/A[1,1];
L41 = I - mult41*I[:,4]*I[:,1]';
A = L41*A
4×4 Matrix{Float64}: 2.0 0.0 4.0 3.0 0.0 5.0 1.0 -4.0 0.0 15.0 0.0 -6.0 0.0 0.0 6.0 -10.0
And so on, following the pattern as before.
A = [2 0 4 3; -4 5 -7 -10; 1 15 2 -4.5; -2 0 2 -13];
L,U = FNC.lufact(A)
([1.0 0.0 0.0 0.0; -2.0 1.0 0.0 0.0; 0.5 3.0 1.0 0.0; -1.0 0.0 -2.0 1.0], [2.0 0.0 4.0 3.0; 0.0 5.0 1.0 -4.0; 0.0 0.0 -3.0 6.0; 0.0 0.0 0.0 2.0])
LtimesU = L*U
4×4 Matrix{Float64}: 2.0 0.0 4.0 3.0 -4.0 5.0 -7.0 -10.0 1.0 15.0 2.0 -4.5 -2.0 0.0 2.0 -13.0
It's best to compare two floating-point quantities by taking their difference.
A - LtimesU
4×4 Matrix{Float64}: 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
(Usually we can expect "zero" only up to machine precision. However, all the exact numbers in this example are also floating-point numbers.)
To solve a linear system, we no longer need the matrix $A$.
b = [4,9,29,40]
z = FNC.forwardsub(L,b)
x = FNC.backsub(U,z)
4-element Vector{Float64}: -3.0 1.0 4.0 -2.0
b - A*x
4-element Vector{Float64}: 0.0 0.0 0.0 0.0
n = 6
A = randn(n,n)
x = ones(n)
y = zeros(n)
for i = 1:n
for j = 1:n
y[i] = y[i] + A[i,j]*x[j] # 2 flops
end
end
Each of the loops implies a summation of flops. The total flop count for this algorithm is [ \sum{i=1}^n \sum{j=1}^n 2 = \sum_{i=1}^n 2n = 2n^2. ] Since the matrix $A$ has $n^2$ elements, all of which have to be involved in the product, it seems unlikely that we could get a flop count that is smaller than $O(n^2)$.
Let's run an experiment with the built-in matrix-vector multiplication. We assume that flops dominate the computation time and thus measure elapsed time.
n = 400:400:4000
t = zeros(size(n))
for (i,n) in enumerate(n)
A = randn(n,n)
x = randn(n)
t[i] = @elapsed for j = 1:10; A*x; end
end
The reason for doing multiple repetitions at each value of $n$ is to avoid having times so short that the resolution of the timer is a factor.
DataFrame(size=n,time=t)
10 rows × 2 columns
size | time | |
---|---|---|
Int64 | Float64 | |
1 | 400 | 0.000509145 |
2 | 800 | 0.00102902 |
3 | 1200 | 0.00260645 |
4 | 1600 | 0.00545276 |
5 | 2000 | 0.0276806 |
6 | 2400 | 0.0519254 |
7 | 2800 | 0.0617959 |
8 | 3200 | 0.0664739 |
9 | 3600 | 0.0871947 |
10 | 4000 | 0.104178 |
Let's repeat the experiment of the previous example for more, and larger, values of $n$.
n = 400:200:6000
t = zeros(size(n))
for (i,n) in enumerate(n)
A = randn(n,n)
x = randn(n)
t[i] = @elapsed for j = 1:10; A*x; end
end
Plotting the time as a function of $n$ on log-log scales is equivalent to plotting the logs of the variables, but is formatted more neatly.
plot(n,t,m=:o,
xaxis=(:log10,"\$n\$"),yaxis=(:log10,"elapsed time (sec)"),
title="Timing of matrix-vector multiplications",label="data",leg=false)
You can see that while the full story is complicated, the graph is trending to a straight line of positive slope. For comparison, we can plot a line that represents $O(n^2)$ growth exactly. (All such lines have slope equal to 2.)
plot!(n,(n/n[end]).^2*t[end],l=:dash,
label="\$O(n^3)\$",legend=:topleft)
We'll test the conclusion of $O(n^3)$ flops experimentally, using the built-in lu
function instead of the purely instructive lufact
.
# Code not working See https://fncbook.github.io/fnc/linsys/demos/flops-lufact.html, which doesn't work
n = 200:100:2400
t = zeros(size(n))
for (i,n) in enumerate(n)
A = randn(n,n)
t[i] = @elapsed for j = 1:6; lu(A); end
end
ReadOnlyMemoryError() Stacktrace: [1] getrf!(A::Matrix{Float64}) @ LinearAlgebra.LAPACK /opt/julia-1.7.0/share/julia/stdlib/v1.7/LinearAlgebra/src/lapack.jl:575 [2] #lu!#146 @ /opt/julia-1.7.0/share/julia/stdlib/v1.7/LinearAlgebra/src/lu.jl:81 [inlined] [3] #lu#153 @ /opt/julia-1.7.0/share/julia/stdlib/v1.7/LinearAlgebra/src/lu.jl:279 [inlined] [4] lu (repeats 2 times) @ /opt/julia-1.7.0/share/julia/stdlib/v1.7/LinearAlgebra/src/lu.jl:278 [inlined] [5] macro expansion @ ./In[134]:5 [inlined] [6] macro expansion @ ./timing.jl:299 [inlined] [7] top-level scope @ ./In[134]:5 [8] eval @ ./boot.jl:373 [inlined] [9] include_string(mapexpr::typeof(REPL.softscope), mod::Module, code::String, filename::String) @ Base ./loading.jl:1196
We plot the timings on a log-log graph and compare it to $O(n^3)$. The result could vary significantly from machine to machine.
plot(n,t,m=:o,
xaxis=(:log10,"\$n\$"),yaxis=(:log10,"elapsed time"),label="data")
plot!(n,(n/n[end]).^3*t[end],l=:dash,
label="\$O(n^3)\$")
┌ Warning: No strict ticks found └ @ PlotUtils /home/hermanr/.julia/packages/PlotUtils/xekml/src/ticks.jl:333 ┌ Warning: No strict ticks found └ @ PlotUtils /home/hermanr/.julia/packages/PlotUtils/xekml/src/ticks.jl:333 ┌ Warning: Invalid negative or zero value 0.0 found at series index 1 for log10 based yscale └ @ Plots /home/hermanr/.julia/packages/Plots/FI0vT/src/utils.jl:95
Here is the previously solved system.
A = [2 0 4 3; -4 5 -7 -10; 1 15 2 -4.5; -2 0 2 -13];
b = [4,9,29,40];
It has a perfectly good solution, obtainable through LU factorization.
L,U = FNC.lufact(A)
x = FNC.backsub( U, FNC.forwardsub(L,b) )
4-element Vector{Float64}: -3.0 1.0 4.0 -2.0
If we swap the second and fourth equations, nothing essential is changed, and Julia still finds the solution.
A[[2,4],:] = A[[4,2],:]
b[[2,4]] = b[[4,2]]
x = A\b
4-element Vector{Float64}: -3.0000000000000235 1.0000000000000018 4.000000000000008 -1.9999999999999953
However, LU factorization fails.
L,U = FNC.lufact(A)
L
4×4 Matrix{Float64}: 1.0 0.0 0.0 0.0 -1.0 1.0 0.0 0.0 0.5 Inf 1.0 0.0 -2.0 Inf NaN 1.0
Here is the system that "broke" LU factorization for us.
A = [2 0 4 3; -4 5 -7 -10; 1 15 2 -4.5; -2 0 2 -13];
b = [4,9,29,40];
When we use the lu
function (from LinearAlgebra
) with three outputs, we get the elements of the PLU factorization.
L,U,p = lu(A);
L
4×4 Matrix{Float64}: 1.0 0.0 0.0 0.0 -0.25 1.0 0.0 0.0 0.5 -0.153846 1.0 0.0 -0.5 0.153846 0.0833333 1.0
U
4×4 Matrix{Float64}: -4.0 5.0 -7.0 -10.0 0.0 16.25 0.25 -7.0 0.0 0.0 5.53846 -9.07692 0.0 0.0 0.0 -0.166667
p
4-element Vector{Int64}: 2 3 4 1
As you see above, the p
return is a vector permutation of 1:n
, rather than the permutation matrix P
. We can recover the latter as follows:
I = diagm(0=>ones(4))
P = I[p,:]
4×4 Matrix{Float64}: 0.0 1.0 0.0 0.0 0.0 0.0 1.0 0.0 0.0 0.0 0.0 1.0 1.0 0.0 0.0 0.0
However, this is rarely necessary in practice (and the vector requires a lot less storage). We can the linear system, for example, using only p
.
x = FNC.backsub( U, FNC.forwardsub(L,b[p,:]) )
4-element Vector{Float64}: -3.0000000000000364 1.0000000000000033 4.000000000000012 -1.999999999999992
If you call lu
with just one output, it is a "factorization object". You can access the individual parts of it using a dot syntax.
fact = lu(A)
fact.L
4×4 Matrix{Float64}: 1.0 0.0 0.0 0.0 -0.25 1.0 0.0 0.0 0.5 -0.153846 1.0 0.0 -0.5 0.153846 0.0833333 1.0
The factorization object can be used efficiently to solve linear systems by the backslash.
x = fact\b
4-element Vector{Float64}: -3.0000000000000235 1.0000000000000018 4.000000000000008 -1.9999999999999953
The idea here is that if you have to solve many different linear systems for the same matrix, you can perform the computationally expensive factorization just once, and repeat only the much faster triangular solves for the different right-hand sides.
In Julia the standard LinearAlgebra
package has a norm
command for vector norms.
x = [2,-3,1,-1]
twonorm = norm(x) # or norm(x,2)
3.872983346207417
infnorm = norm(x,Inf)
3.0
onenorm = norm(x,1)
7.0
A = [ 2 0; 1 -1 ]
2×2 Matrix{Int64}: 2 0 1 -1
In Julia one uses norm
for vector norms and opnorm
for induced matrix norms. The default matrix norm is the 2-norm.
twonorm = opnorm(A)
2.2882456112707374
(A potential snag is that norm
does work on a matrix but treats it like a vector of stacked columns, giving a different result.)
You can get the 1-norm as well.
onenorm = opnorm(A,1)
3.0
The 1-norm is equivalent to
maximum( sum(abs.(A),dims=1) ) # sum down the rows (1st matrix dimension)
3
Similarly, we can get the $\infty$-norm and check our formula for it.
infnorm = opnorm(A,Inf)
2.0
maximum( sum(abs.(A),dims=2) ) # sum across columns (2nd matrix dimension)
2
Here we illustrate the geometric interpretation of the 2-norm. First, we will sample a lot of vectors on the unit circle in $\mathbb{R}^2$.
theta = 2*pi*(0:1/600:1)
x = @. [ cos(theta) sin(theta) ]' # 601 unit columns
plot(x[1,:],x[2,:],aspect_ratio=1,layout=(1,2),subplot=1,
title="Unit circle",leg=:none,xlabel="\$x_1\$",ylabel="\$x_2\$")
We can apply A
to every column of x
simply by using a matrix multiplication.
Ax = A*x;
We superimpose the image of the unit circle with the circle whose radius is $\|A\|_2$, and display the plots side by side.
plot!(Ax[1,:],Ax[2,:],subplot=2,aspect_ratio=1,
title="Image under map",leg=:none,xlabel="\$x_1\$",ylabel="\$x_2\$")
plot!(twonorm*x[1,:],twonorm*x[2,:],subplot=2,l=:dash)
Julia has a function cond
to compute matrix condition numbers. By default, the 2-norm is used. As an example, the family of Hilbert matrices is famously badly conditioned. Here is the $7\times 7$ case.
A = [ 1/(i+j) for i=1:7, j=1:7 ]
7×7 Matrix{Float64}: 0.5 0.333333 0.25 0.2 0.166667 0.142857 0.125 0.333333 0.25 0.2 0.166667 0.142857 0.125 0.111111 0.25 0.2 0.166667 0.142857 0.125 0.111111 0.1 0.2 0.166667 0.142857 0.125 0.111111 0.1 0.0909091 0.166667 0.142857 0.125 0.111111 0.1 0.0909091 0.0833333 0.142857 0.125 0.111111 0.1 0.0909091 0.0833333 0.0769231 0.125 0.111111 0.1 0.0909091 0.0833333 0.0769231 0.0714286
kappa = cond(A)
1.697836309221093e9
Next we engineer a linear system problem to which we know the exact answer.
x_exact = 1.:7.
b = A*x_exact
7-element Vector{Float64}: 5.282142857142857 4.342063492063492 3.7130952380952382 3.2538239538239537 2.900613275613275 2.619197469197469 2.389063714063714
Now we perturb the data randomly with a vector of norm $10^{-12}$.
dA = randn(size(A)); dA = 1e-12*(dA/norm(dA));
db = randn(size(b)); db = 1e-12*(db/norm(db));
We solve the perturbed problem using built-in pivoted LU and see how the solution was changed.
x = (A+dA) \ (b+db);
dx = x - x_exact;
Here is the relative error in the solution.
rel_error = norm(dx) / norm(x_exact)
7.801177505809046e-5
And here are upper bounds predicted using the condition number of the original matrix.
@show b_bound = kappa * 1e-12/norm(b);
@show A_bound = kappa * 1e-12/norm(A);
b_bound = (kappa * 1.0e-12) / norm(b) = 0.00017690558572904666 A_bound = (kappa * 1.0e-12) / norm(A) = 0.001442507718610632
Even if we don't make any manual perturbations to the data, machine epsilon does when we solve the linear system numerically.
x = A\b;
@show rel_error = norm(x - x_exact) / norm(x_exact);
@show rounding_bound = kappa*eps();
rel_error = norm(x - x_exact) / norm(x_exact) = 3.619160875001057e-8 rounding_bound = kappa * eps() = 3.7699539250837087e-7
Because $\kappa\approx 10^8$, it's possible to lose 8 digits of accuracy in the process of passing from $A$ and $b$ to $x$. That's independent of the algorithm; it's inevitable once the data are expressed in double precision.
Larger Hilbert matrices are even more poorly conditioned.
A = [ 1/(i+j) for i=1:14, j=1:14 ];
kappa = cond(A)
4.756720805226143e17
Before we compute the solution, note that $\kappa$ exceeds 1/eps()
. In principle we therefore might end up with an answer that is completely wrong (i.e., a relative error greater than 100%).
rounding_bound = kappa*eps()
105.62041919351157
x_exact = 1.:14.;
b = A*x_exact;
x = A\b;
We got an answer. But in fact the error does exceed 100%.
relative_error = norm(x_exact - x) / norm(x_exact)
5.839548471598157
Here is a matrix with both lower and upper bandwidth equal to one. Such a matrix is called tridiagonal.
A = [ 2 -1 0 0 0 0
4 2 -1 0 0 0
0 3 0 -1 0 0
0 0 2 2 -1 0
0 0 0 1 1 -1
0 0 0 0 0 2 ]
6×6 Matrix{Int64}: 2 -1 0 0 0 0 4 2 -1 0 0 0 0 3 0 -1 0 0 0 0 2 2 -1 0 0 0 0 1 1 -1 0 0 0 0 0 2
We can extract the elements on any diagonal using the diag
command. The "main" or central diagonal is numbered zero, above and to the right of that is positive, and below and to the left is negative.
diag_main = diag(A,0)
6-element Vector{Int64}: 2 2 0 2 1 2
diag_plusone = diag(A,1)
5-element Vector{Int64}: -1 -1 -1 -1 -1
diag_minusone = diag(A,-1)
5-element Vector{Int64}: 4 3 2 1 0
We can also put whatever numbers we like onto any diagonal with the diagm
command.
A = A + diagm(2=>[pi,8,6,7])
6×6 Matrix{Float64}: 2.0 -1.0 3.14159 0.0 0.0 0.0 4.0 2.0 -1.0 8.0 0.0 0.0 0.0 3.0 0.0 -1.0 6.0 0.0 0.0 0.0 2.0 2.0 -1.0 7.0 0.0 0.0 0.0 1.0 1.0 -1.0 0.0 0.0 0.0 0.0 0.0 2.0
L,U = FNC.lufact(A)
L
6×6 Matrix{Float64}: 1.0 0.0 0.0 0.0 0.0 0.0 2.0 1.0 0.0 0.0 0.0 0.0 0.0 0.75 1.0 0.0 0.0 0.0 0.0 0.0 0.36614 1.0 0.0 0.0 0.0 0.0 0.0 0.219155 1.0 0.0 0.0 0.0 0.0 0.0 0.0 1.0
U
6×6 Matrix{Float64}: 2.0 -1.0 3.14159 0.0 0.0 0.0 0.0 4.0 -7.28319 8.0 0.0 0.0 0.0 0.0 5.46239 -7.0 6.0 0.0 0.0 0.0 0.0 4.56298 -3.19684 7.0 0.0 0.0 0.0 0.0 1.7006 -2.53408 0.0 0.0 0.0 0.0 0.0 2.0
Observe above that the lower and upper bandwidths of $A$ are preserved in the $L$ and $U$ results.
We'll use a large banded matrix to observe the speedup possible in LU factorization. We'll also need to load in a (standard) package for sparse matrices.
If we use an ordinary "dense" matrix, then there's no way to exploit a banded structure such as tridiagonality.
n = 10000
A = diagm(0=>1:n,1=>n-1:-1:1,-1=>ones(n-1))
10000×10000 Matrix{Float64}: 1.0 9999.0 0.0 0.0 0.0 … 0.0 0.0 0.0 0.0 1.0 2.0 9998.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 1.0 3.0 9997.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 1.0 4.0 9996.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 1.0 5.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 1.0 … 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 … 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 ⋮ ⋱ 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 … 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 … 4.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 9997.0 3.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 1.0 9998.0 2.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 1.0 9999.0 1.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 1.0 10000.0
@time lu(A);
ReadOnlyMemoryError() Stacktrace: [1] getrf!(A::Matrix{Float64}) @ LinearAlgebra.LAPACK /opt/julia-1.7.0/share/julia/stdlib/v1.7/LinearAlgebra/src/lapack.jl:575 [2] #lu!#146 @ /opt/julia-1.7.0/share/julia/stdlib/v1.7/LinearAlgebra/src/lu.jl:81 [inlined] [3] #lu#153 @ /opt/julia-1.7.0/share/julia/stdlib/v1.7/LinearAlgebra/src/lu.jl:279 [inlined] [4] lu (repeats 2 times) @ /opt/julia-1.7.0/share/julia/stdlib/v1.7/LinearAlgebra/src/lu.jl:278 [inlined] [5] top-level scope @ ./timing.jl:220 [inlined] [6] top-level scope @ ./In[118]:0 [7] eval @ ./boot.jl:373 [inlined] [8] include_string(mapexpr::typeof(REPL.softscope), mod::Module, code::String, filename::String) @ Base ./loading.jl:1196
If instead we construct a proper "sparse" matrix, though, the speedup can be dramatic.
A = spdiagm(0=>1:n,1=>n-1:-1:1,-1=>ones(n-1))
10000×10000 SparseMatrixCSC{Float64, Int64} with 29998 stored entries: ⠻⣦⡀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀ ⠀⠈⠻⣦⡀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀ ⠀⠀⠀⠈⠻⣦⡀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀ ⠀⠀⠀⠀⠀⠈⠻⣦⡀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀ ⠀⠀⠀⠀⠀⠀⠀⠈⠻⣦⡀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀ ⠀⠀⠀⠀⠀⠀⠀⠀⠀⠈⠻⣦⡀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀ ⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠈⠻⣦⡀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀ ⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠈⠻⣦⡀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀ ⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠈⠻⣦⡀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀ ⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠈⠻⣦⡀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀ ⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠈⠻⣦⡀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀ ⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠈⠻⣦⡀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀ ⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠈⠻⣦⡀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀ ⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠈⠻⣦⡀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀ ⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠈⠻⣦⡀⠀⠀⠀⠀⠀⠀⠀⠀⠀ ⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠈⠻⣦⡀⠀⠀⠀⠀⠀⠀⠀ ⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠈⠻⣦⡀⠀⠀⠀⠀⠀ ⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠈⠻⣦⡀⠀⠀⠀ ⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠈⠻⣦⡀⠀ ⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠈⠻⣦
@time lu(A);
0.163473 seconds (420.51 k allocations: 33.359 MiB, 93.40% compilation time)
We begin with a symmetric $A$.
A = [ 2 4 4 2
4 5 8 -5
4 8 6 2
2 -5 2 -26 ];
Carrying out our usual elimination in the first column leads us to
L1 = diagm(0=>ones(4))
L1[2:4,1] = [-2,-2,-1]
A1 = L1*A
4×4 Matrix{Float64}: 2.0 4.0 4.0 2.0 0.0 -3.0 0.0 -9.0 0.0 0.0 -2.0 -2.0 0.0 -9.0 -2.0 -28.0
But now let's note that if we transpose this result, we have the same first column as before! So we could apply again and then transpose back.
A2 = (L1*A1')'
4×4 adjoint(::Matrix{Float64}) with eltype Float64: 2.0 0.0 0.0 0.0 0.0 -3.0 0.0 -9.0 0.0 0.0 -2.0 -2.0 0.0 -9.0 -2.0 -28.0
Using transpose identities, this is just
A2 = A1*L1'
4×4 Matrix{Float64}: 2.0 0.0 0.0 0.0 0.0 -3.0 0.0 -9.0 0.0 0.0 -2.0 -2.0 0.0 -9.0 -2.0 -28.0
Now you can see how we proceed down and to the right, eliminating in a column and then symmetrically in the corresponding row.
L2 = diagm(0=>ones(4))
L2[3:4,2] = [0,-3]
A3 = L2*A2*L2'
4×4 Matrix{Float64}: 2.0 0.0 0.0 0.0 0.0 -3.0 0.0 0.0 0.0 0.0 -2.0 -2.0 0.0 0.0 -2.0 -1.0
Finally, we arrive at a diagonal matrix.
L3 = diagm(0=>ones(4))
L3[4,3] = -1
D = L3*A3*L3'
4×4 Matrix{Float64}: 2.0 0.0 0.0 0.0 0.0 -3.0 0.0 0.0 0.0 0.0 -2.0 0.0 0.0 0.0 0.0 1.0
A randomly chosen matrix is extremely unlikely to be symmetric. However, there is a simple way to symmetrize one.
A = rand(1.:9.,4,4)
B = A + A'
4×4 Matrix{Float64}: 16.0 6.0 11.0 16.0 6.0 14.0 13.0 12.0 11.0 13.0 16.0 7.0 16.0 12.0 7.0 2.0
Similarly, a random symmetric matrix is unlikely to be positive definite. The Cholesky algorithm always detects a non-PD matrix by quitting with an error.
cholesky(B)
PosDefException: matrix is not positive definite; Cholesky factorization failed. Stacktrace: [1] checkpositivedefinite @ /opt/julia-1.7.0/share/julia/stdlib/v1.7/LinearAlgebra/src/factorization.jl:18 [inlined] [2] cholesky!(A::Hermitian{Float64, Matrix{Float64}}, ::Val{false}; check::Bool) @ LinearAlgebra /opt/julia-1.7.0/share/julia/stdlib/v1.7/LinearAlgebra/src/cholesky.jl:266 [3] cholesky!(A::Matrix{Float64}, ::Val{false}; check::Bool) @ LinearAlgebra /opt/julia-1.7.0/share/julia/stdlib/v1.7/LinearAlgebra/src/cholesky.jl:298 [4] #cholesky#143 @ /opt/julia-1.7.0/share/julia/stdlib/v1.7/LinearAlgebra/src/cholesky.jl:394 [inlined] [5] cholesky (repeats 2 times) @ /opt/julia-1.7.0/share/julia/stdlib/v1.7/LinearAlgebra/src/cholesky.jl:394 [inlined] [6] top-level scope @ In[128]:1 [7] eval @ ./boot.jl:373 [inlined] [8] include_string(mapexpr::typeof(REPL.softscope), mod::Module, code::String, filename::String) @ Base ./loading.jl:1196
It's not hard to manufacture an SPD matrix to try out the Cholesky factorization.
B = A'*A
cf = cholesky(B)
Cholesky{Float64, Matrix{Float64}} U factor: 4×4 UpperTriangular{Float64, Matrix{Float64}}: 13.0 10.7692 9.53846 9.0 ⋅ 7.07274 3.14985 0.152264 ⋅ ⋅ 5.30058 5.40964 ⋅ ⋅ ⋅ 2.17086
What's returned is a "factorization object." (This allows it to be used efficiently in various contexts.) Another step is required to extract the factor as a matrix.
R = Matrix(cf.U)
4×4 Matrix{Float64}: 13.0 10.7692 9.53846 9.0 0.0 7.07274 3.14985 0.152264 0.0 0.0 5.30058 5.40964 0.0 0.0 0.0 2.17086
norm(R'*R - B)
0.0