CGLS
Krylov.cgls — Function(x, stats) = cgls(A, b::AbstractVector{FC};
M=I, λ::T=zero(T), atol::T=√eps(T), rtol::T=√eps(T),
radius::T=zero(T), itmax::Int=0, verbose::Int=0, history::Bool=false,
ldiv::Bool=false, callback=solver->false)T is an AbstractFloat such as Float32, Float64 or BigFloat. FC is T or Complex{T}.
Solve the regularized linear least-squares problem
minimize ‖b - Ax‖₂² + λ‖x‖₂²using the Conjugate Gradient (CG) method, where λ ≥ 0 is a regularization parameter. This method is equivalent to applying CG to the normal equations
(AᵀA + λI) x = Aᵀbbut is more stable.
CGLS produces monotonic residuals ‖r‖₂ but not optimality residuals ‖Aᵀr‖₂. It is formally equivalent to LSQR, though can be slightly less accurate, but simpler to implement.
The callback is called as callback(solver) and should return true if the main loop should terminate, and false otherwise.
References
- M. R. Hestenes and E. Stiefel. Methods of conjugate gradients for solving linear systems, Journal of Research of the National Bureau of Standards, 49(6), pp. 409–436, 1952.
- A. Björck, T. Elfving and Z. Strakos, Stability of Conjugate Gradient and Lanczos Methods for Linear Least Squares Problems, SIAM Journal on Matrix Analysis and Applications, 19(3), pp. 720–736, 1998.
Krylov.cgls! — Functionsolver = cgls!(solver::CglsSolver, A, b; kwargs...)where kwargs are keyword arguments of cgls.
See CglsSolver for more details about the solver.
CRLS
Krylov.crls — Function(x, stats) = crls(A, b::AbstractVector{FC};
M=I, λ::T=zero(T), atol::T=√eps(T), rtol::T=√eps(T),
radius::T=zero(T), itmax::Int=0, verbose::Int=0, history::Bool=false,
ldiv::Bool=false, callback=solver->false)T is an AbstractFloat such as Float32, Float64 or BigFloat. FC is T or Complex{T}.
Solve the linear least-squares problem
minimize ‖b - Ax‖₂² + λ‖x‖₂²using the Conjugate Residuals (CR) method. This method is equivalent to applying MINRES to the normal equations
(AᵀA + λI) x = Aᵀb.This implementation recurs the residual r := b - Ax.
CRLS produces monotonic residuals ‖r‖₂ and optimality residuals ‖Aᵀr‖₂. It is formally equivalent to LSMR, though can be substantially less accurate, but simpler to implement.
The callback is called as callback(solver) and should return true if the main loop should terminate, and false otherwise.
Reference
- D. C.-L. Fong, Minimum-Residual Methods for Sparse, Least-Squares using Golubg-Kahan Bidiagonalization, Ph.D. Thesis, Stanford University, 2011.
Krylov.crls! — Functionsolver = crls!(solver::CrlsSolver, A, b; kwargs...)where kwargs are keyword arguments of crls.
See CrlsSolver for more details about the solver.
LSLQ
Krylov.lslq — Function(x, stats) = lslq(A, b::AbstractVector{FC};
M=I, N=I, sqd::Bool=false, λ::T=zero(T),
atol::T=√eps(T), btol::T=√eps(T), etol::T=√eps(T),
window::Int=5, utol::T=√eps(T), itmax::Int=0,
σ::T=zero(T), transfer_to_lsqr::Bool=false,
conlim::T=1/√eps(T), verbose::Int=0, history::Bool=false,
ldiv::Bool=false, callback=solver->false)T is an AbstractFloat such as Float32, Float64 or BigFloat. FC is T or Complex{T}.
Solve the regularized linear least-squares problem
minimize ‖b - Ax‖₂² + λ²‖x‖₂²using the LSLQ method, where λ ≥ 0 is a regularization parameter. LSLQ is formally equivalent to applying SYMMLQ to the normal equations
(AᵀA + λ²I) x = Aᵀbbut is more stable.
Main features
- the solution estimate is updated along orthogonal directions
- the norm of the solution estimate ‖xᴸₖ‖₂ is increasing
- the error ‖eₖ‖₂ := ‖xᴸₖ - x*‖₂ is decreasing
- it is possible to transition cheaply from the LSLQ iterate to the LSQR iterate if there is an advantage (there always is in terms of error)
- if
Ais rank deficient, identify the minimum least-squares solution
Optional arguments
M: a symmetric and positive definite dual preconditionerN: a symmetric and positive definite primal preconditionersqdindicates that we are solving a symmetric and quasi-definite system withλ=1
If λ > 0, we solve the symmetric and quasi-definite system
[ E A ] [ r ] [ b ]
[ Aᵀ -λ²F ] [ x ] = [ 0 ],where E and F are symmetric and positive definite. Preconditioners M = E⁻¹ ≻ 0 and N = F⁻¹ ≻ 0 may be provided in the form of linear operators. If sqd=true, λ is set to the common value 1.
The system above represents the optimality conditions of
minimize ‖b - Ax‖²_E⁻¹ + λ²‖x‖²_F.For a symmetric and positive definite matrix K, the K-norm of a vector x is ‖x‖²_K = xᵀKx. LSLQ is then equivalent to applying SYMMLQ to (AᵀE⁻¹A + λ²F)x = AᵀE⁻¹b with r = E⁻¹(b - Ax).
If λ = 0, we solve the symmetric and indefinite system
[ E A ] [ r ] [ b ]
[ Aᵀ 0 ] [ x ] = [ 0 ].The system above represents the optimality conditions of
minimize ‖b - Ax‖²_E⁻¹.In this case, N can still be specified and indicates the weighted norm in which x and Aᵀr should be measured. r can be recovered by computing E⁻¹(b - Ax).
λis a regularization parameter (see the problem statement above)σis an underestimate of the smallest nonzero singular value ofA–-settingσtoo large will result in an error in the course of the iterationsatolis a stopping tolerance based on the residualbtolis a stopping tolerance used to detect zero-residual problemsetolis a stopping tolerance based on the lower bound on the errorwindowis the number of iterations used to accumulate a lower bound on the errorutolis a stopping tolerance based on the upper bound on the errortransfer_to_lsqrreturn the CG solution estimate (i.e., the LSQR point) instead of the LQ estimateitmaxis the maximum number of iterations (0 means no imposed limit)conlimis the limit on the estimated condition number ofAbeyond which the solution will be abandonedverbosedetermines verbosity.
Return values
lslq returns the tuple (x, stats) where
xis the LQ solution estimatestatscollects other statistics on the run in a LSLQStatsstats.err_lbndsis a vector of lower bounds on the LQ error–-the vector is empty ifwindowis set to zerostats.err_ubnds_lqis a vector of upper bounds on the LQ error–-the vector is empty ifσ == 0is left at zerostats.err_ubnds_cgis a vector of upper bounds on the CG error–-the vector is empty ifσ == 0is left at zerostats.error_with_bndis a boolean indicating whether there was an error in the upper bounds computation (cancellation errors, too large σ ...)
Stopping conditions
The iterations stop as soon as one of the following conditions holds true:
- the optimality residual is sufficiently small (
stats.status = "found approximate minimum least-squares solution") in the sense that either- ‖Aᵀr‖ / (‖A‖ ‖r‖) ≤ atol, or
- 1 + ‖Aᵀr‖ / (‖A‖ ‖r‖) ≤ 1
- an approximate zero-residual solution has been found (
stats.status = "found approximate zero-residual solution") in the sense that either- ‖r‖ / ‖b‖ ≤ btol + atol ‖A‖ * ‖xᴸ‖ / ‖b‖, or
- 1 + ‖r‖ / ‖b‖ ≤ 1
- the estimated condition number of
Ais too large in the sense that either- 1/cond(A) ≤ 1/conlim (
stats.status = "condition number exceeds tolerance"), or - 1 + 1/cond(A) ≤ 1 (
stats.status = "condition number seems too large for this machine")
- 1/cond(A) ≤ 1/conlim (
- the lower bound on the LQ forward error is less than etol * ‖xᴸ‖
- the upper bound on the CG forward error is less than utol * ‖xᶜ‖
The callback is called as callback(solver) and should return true if the main loop should terminate, and false otherwise.
References
- R. Estrin, D. Orban and M. A. Saunders, Euclidean-norm error bounds for SYMMLQ and CG, SIAM Journal on Matrix Analysis and Applications, 40(1), pp. 235–253, 2019.
- R. Estrin, D. Orban and M. A. Saunders, LSLQ: An Iterative Method for Linear Least-Squares with an Error Minimization Property, SIAM Journal on Matrix Analysis and Applications, 40(1), pp. 254–275, 2019.
Krylov.lslq! — Functionsolver = lslq!(solver::LslqSolver, A, b; kwargs...)where kwargs are keyword arguments of lslq.
See LslqSolver for more details about the solver.
LSQR
Krylov.lsqr — Function(x, stats) = lsqr(A, b::AbstractVector{FC};
M=I, N=I, sqd::Bool=false, λ::T=zero(T),
axtol::T=√eps(T), btol::T=√eps(T),
atol::T=zero(T), rtol::T=zero(T),
etol::T=√eps(T), window::Int=5,
itmax::Int=0, conlim::T=1/√eps(T),
radius::T=zero(T), verbose::Int=0, history::Bool=false,
ldiv::Bool=false, callback=solver->false)T is an AbstractFloat such as Float32, Float64 or BigFloat. FC is T or Complex{T}.
Solve the regularized linear least-squares problem
minimize ‖b - Ax‖₂² + λ²‖x‖₂²using the LSQR method, where λ ≥ 0 is a regularization parameter. LSQR is formally equivalent to applying CG to the normal equations
(AᵀA + λ²I) x = Aᵀb(and therefore to CGLS) but is more stable.
LSQR produces monotonic residuals ‖r‖₂ but not optimality residuals ‖Aᵀr‖₂. It is formally equivalent to CGLS, though can be slightly more accurate.
If λ > 0, LSQR solves the symmetric and quasi-definite system
[ E A ] [ r ] [ b ]
[ Aᵀ -λ²F ] [ x ] = [ 0 ],where E and F are symmetric and positive definite. Preconditioners M = E⁻¹ ≻ 0 and N = F⁻¹ ≻ 0 may be provided in the form of linear operators. If sqd=true, λ is set to the common value 1.
The system above represents the optimality conditions of
minimize ‖b - Ax‖²_E⁻¹ + λ²‖x‖²_F.For a symmetric and positive definite matrix K, the K-norm of a vector x is ‖x‖²_K = xᵀKx. LSQR is then equivalent to applying CG to (AᵀE⁻¹A + λ²F)x = AᵀE⁻¹b with r = E⁻¹(b - Ax).
If λ = 0, we solve the symmetric and indefinite system
[ E A ] [ r ] [ b ]
[ Aᵀ 0 ] [ x ] = [ 0 ].The system above represents the optimality conditions of
minimize ‖b - Ax‖²_E⁻¹.In this case, N can still be specified and indicates the weighted norm in which x and Aᵀr should be measured. r can be recovered by computing E⁻¹(b - Ax).
The callback is called as callback(solver) and should return true if the main loop should terminate, and false otherwise.
Reference
- C. C. Paige and M. A. Saunders, LSQR: An Algorithm for Sparse Linear Equations and Sparse Least Squares, ACM Transactions on Mathematical Software, 8(1), pp. 43–71, 1982.
Krylov.lsqr! — Functionsolver = lsqr!(solver::LsqrSolver, A, b; kwargs...)where kwargs are keyword arguments of lsqr.
See LsqrSolver for more details about the solver.
LSMR
Krylov.lsmr — Function(x, stats) = lsmr(A, b::AbstractVector{FC};
M=I, N=I, sqd::Bool=false, λ::T=zero(T),
axtol::T=√eps(T), btol::T=√eps(T),
atol::T=zero(T), rtol::T=zero(T),
etol::T=√eps(T), window::Int=5,
itmax::Int=0, conlim::T=1/√eps(T),
radius::T=zero(T), verbose::Int=0,
history::Bool=false, ldiv::Bool=false,
callback=solver->false)T is an AbstractFloat such as Float32, Float64 or BigFloat. FC is T or Complex{T}.
Solve the regularized linear least-squares problem
minimize ‖b - Ax‖₂² + λ²‖x‖₂²using the LSMR method, where λ ≥ 0 is a regularization parameter. LSMR is formally equivalent to applying MINRES to the normal equations
(AᵀA + λ²I) x = Aᵀb(and therefore to CRLS) but is more stable.
LSMR produces monotonic residuals ‖r‖₂ and optimality residuals ‖Aᵀr‖₂. It is formally equivalent to CRLS, though can be substantially more accurate.
LSMR can be also used to find a null vector of a singular matrix A by solving the problem min ‖Aᵀx - b‖ with any nonzero vector b. At a minimizer, the residual vector r = b - Aᵀx will satisfy Ar = 0.
If λ > 0, we solve the symmetric and quasi-definite system
[ E A ] [ r ] [ b ]
[ Aᵀ -λ²F ] [ x ] = [ 0 ],where E and F are symmetric and positive definite. Preconditioners M = E⁻¹ ≻ 0 and N = F⁻¹ ≻ 0 may be provided in the form of linear operators. If sqd=true, λ is set to the common value 1.
The system above represents the optimality conditions of
minimize ‖b - Ax‖²_E⁻¹ + λ²‖x‖²_F.For a symmetric and positive definite matrix K, the K-norm of a vector x is ‖x‖²_K = xᵀKx. LSMR is then equivalent to applying MINRES to (AᵀE⁻¹A + λ²F)x = AᵀE⁻¹b with r = E⁻¹(b - Ax).
If λ = 0, we solve the symmetric and indefinite system
[ E A ] [ r ] [ b ]
[ Aᵀ 0 ] [ x ] = [ 0 ].The system above represents the optimality conditions of
minimize ‖b - Ax‖²_E⁻¹.In this case, N can still be specified and indicates the weighted norm in which x and Aᵀr should be measured. r can be recovered by computing E⁻¹(b - Ax).
The callback is called as callback(solver) and should return true if the main loop should terminate, and false otherwise.
Reference
- D. C.-L. Fong and M. A. Saunders, LSMR: An Iterative Algorithm for Sparse Least Squares Problems, SIAM Journal on Scientific Computing, 33(5), pp. 2950–2971, 2011.
Krylov.lsmr! — Functionsolver = lsmr!(solver::LsmrSolver, A, b; kwargs...)where kwargs are keyword arguments of lsmr.
See LsmrSolver for more details about the solver.