FastLevenbergMarquardt.lmsolve
— Functionlmsolve(
fun,
jac,
x0::StaticVector{N, <:AbstractFloat},
data = nothing,
lb::Union{Nothing, Real, AbstractVector{<:Real}} = nothing,
ub::Union{Nothing, Real, AbstractVector{<:Real}} = nothing;
kwargs...
) -> x, F, info, iter, nfev, njev
lmsolve!(
fun!,
jac!,
x0::AbstractVector{<:AbstractFloat},
m::Integer = length(x0),
data = nothing,
lb::Union{Nothing, Real, AbstractVector{<:Real}} = nothing,
ub::Union{Nothing, Real, AbstractVector{<:Real}} = nothing;
kwargs...,
) -> x, F, info, iter, nfev, njev, LM, solver
lmsolve!(
fun!,
jac!,
x0::AbstractVector{<:AbstractFloat},
f::AbstractVector{<:AbstractFloat},
J::AbstractMatrix{<:AbstractFloat},
data = nothing,
lb::Union{Nothing, Real, AbstractVector{<:Real}} = nothing,
ub::Union{Nothing, Real, AbstractVector{<:Real}} = nothing;
kwargs...,
) -> x, F, info, iter, nfev, njev, LM, solver
lmsolve!(
fun!,
jac!,
LM::LMWorkspace,
data = nothing,
lb::Union{Nothing, Real, AbstractVector{<:Real}} = nothing,
ub::Union{Nothing, Real, AbstractVector{<:Real}} = nothing;
kwargs...
) -> x, F, info, iter, nfev, njev, LM, solver
Minimize F(x) = ||f(x)||^2
using the Levenberg-Marquardt algorithm.
Arguments
fun/fun!
: function to be minimized||f||^2
,f = fun(x, data)
,f = fun!(f, x, data)
jac/jac!
: jacobian off
,J = jac(x, data)
,J = jac!(J, x, data)
x0::AbstractVector{<:AbstractFloat}
: initial guessm::Integer = length(x0)
: number of function valuesdata = nothing
: data passed tofun/fun!
andjac/jac!
f::AbstractVector{<:AbstractFloat}
: preallocated function vectorJ::AbstractMatrix{<:AbstractFloat}
: preallocated Jacobian matrixLM::LMWorkspace
: preallocated workspacelb::Union{Nothing, Real, AbstractVector{<:Real}} = nothing
: lower bounds forx
. Vectors must have same length asx
ub::Union{Nothing, Real, AbstractVector{<:Real}} = nothing
: upper bounds forx
. Vectors must have same length asx
Keywords
solver::Union{Nothing, Symbol} = nothing
: linear solver forlmsolve!
(lmsolve
always uses Cholesky factorization)nothing
: QR for dense Jacobian, Cholesky for sparse Jacobian:cholesky
: Cholesky factorization:qr
: QR factorization
ftol::Real = eps(eltype(x))
: relative tolerance for function: both actual and predicted reductions are less thanftol
xtol::Real = 1e-10
: relative tolerance for change inx
:all(abs(x - xk) < xtol * (xtol + abs(x)))
gtol::Real = eps(eltype(x))
: tolerance for gradient:norm(g, Inf) < gtol
maxit::Integer = 1000
: maximum number of iterationsfactor::Real = 1e-6
: initial factor for dampingfactoraccept::Real = 13
: factor for decreasing damping on good stepfactorreject::Real = 3
: factor for increasing damping on bad stepfactorupdate::Symbol = :marquardt
: factor update method∈ (:marquardt, :nielsen)
minscale::Real = 1e-12
: diagonal scaling lower boundmaxscale::Real = 1e16
: diagonal scaling upper boundminfactor::Real = 1e-28
: damping factor lower boundmaxfactor::Real = 1e32
: damping factor upper bound
Returns
x/LM.x
: solutionF
: final objectiveinfo::Int
: convergence status1
: both actual and predicted reductions are less thanftol
2
: relative difference between two consecutive iterates is less thanxtol
3
: inf norm of the gradient is less thangtol
-1
:maxit
reached
iter::Int
: number of iterationsnfev::Int
: number of function evaluationsnjev::Int
: number of Jacobian evaluationsLM::LMWorkspace
: workspacesolver::AbstractSolver
: solver
Notes
In the returned LMWorkspace
, only LM.x
and LM.f
are guaranteed to be updated. That is, LM.J
might not be the Jacobian at the returned x
.
FastLevenbergMarquardt.lmsolve!
— Functionlmsolve(
fun,
jac,
x0::StaticVector{N, <:AbstractFloat},
data = nothing,
lb::Union{Nothing, Real, AbstractVector{<:Real}} = nothing,
ub::Union{Nothing, Real, AbstractVector{<:Real}} = nothing;
kwargs...
) -> x, F, info, iter, nfev, njev
lmsolve!(
fun!,
jac!,
x0::AbstractVector{<:AbstractFloat},
m::Integer = length(x0),
data = nothing,
lb::Union{Nothing, Real, AbstractVector{<:Real}} = nothing,
ub::Union{Nothing, Real, AbstractVector{<:Real}} = nothing;
kwargs...,
) -> x, F, info, iter, nfev, njev, LM, solver
lmsolve!(
fun!,
jac!,
x0::AbstractVector{<:AbstractFloat},
f::AbstractVector{<:AbstractFloat},
J::AbstractMatrix{<:AbstractFloat},
data = nothing,
lb::Union{Nothing, Real, AbstractVector{<:Real}} = nothing,
ub::Union{Nothing, Real, AbstractVector{<:Real}} = nothing;
kwargs...,
) -> x, F, info, iter, nfev, njev, LM, solver
lmsolve!(
fun!,
jac!,
LM::LMWorkspace,
data = nothing,
lb::Union{Nothing, Real, AbstractVector{<:Real}} = nothing,
ub::Union{Nothing, Real, AbstractVector{<:Real}} = nothing;
kwargs...
) -> x, F, info, iter, nfev, njev, LM, solver
Minimize F(x) = ||f(x)||^2
using the Levenberg-Marquardt algorithm.
Arguments
fun/fun!
: function to be minimized||f||^2
,f = fun(x, data)
,f = fun!(f, x, data)
jac/jac!
: jacobian off
,J = jac(x, data)
,J = jac!(J, x, data)
x0::AbstractVector{<:AbstractFloat}
: initial guessm::Integer = length(x0)
: number of function valuesdata = nothing
: data passed tofun/fun!
andjac/jac!
f::AbstractVector{<:AbstractFloat}
: preallocated function vectorJ::AbstractMatrix{<:AbstractFloat}
: preallocated Jacobian matrixLM::LMWorkspace
: preallocated workspacelb::Union{Nothing, Real, AbstractVector{<:Real}} = nothing
: lower bounds forx
. Vectors must have same length asx
ub::Union{Nothing, Real, AbstractVector{<:Real}} = nothing
: upper bounds forx
. Vectors must have same length asx
Keywords
solver::Union{Nothing, Symbol} = nothing
: linear solver forlmsolve!
(lmsolve
always uses Cholesky factorization)nothing
: QR for dense Jacobian, Cholesky for sparse Jacobian:cholesky
: Cholesky factorization:qr
: QR factorization
ftol::Real = eps(eltype(x))
: relative tolerance for function: both actual and predicted reductions are less thanftol
xtol::Real = 1e-10
: relative tolerance for change inx
:all(abs(x - xk) < xtol * (xtol + abs(x)))
gtol::Real = eps(eltype(x))
: tolerance for gradient:norm(g, Inf) < gtol
maxit::Integer = 1000
: maximum number of iterationsfactor::Real = 1e-6
: initial factor for dampingfactoraccept::Real = 13
: factor for decreasing damping on good stepfactorreject::Real = 3
: factor for increasing damping on bad stepfactorupdate::Symbol = :marquardt
: factor update method∈ (:marquardt, :nielsen)
minscale::Real = 1e-12
: diagonal scaling lower boundmaxscale::Real = 1e16
: diagonal scaling upper boundminfactor::Real = 1e-28
: damping factor lower boundmaxfactor::Real = 1e32
: damping factor upper bound
Returns
x/LM.x
: solutionF
: final objectiveinfo::Int
: convergence status1
: both actual and predicted reductions are less thanftol
2
: relative difference between two consecutive iterates is less thanxtol
3
: inf norm of the gradient is less thangtol
-1
:maxit
reached
iter::Int
: number of iterationsnfev::Int
: number of function evaluationsnjev::Int
: number of Jacobian evaluationsLM::LMWorkspace
: workspacesolver::AbstractSolver
: solver
Notes
In the returned LMWorkspace
, only LM.x
and LM.f
are guaranteed to be updated. That is, LM.J
might not be the Jacobian at the returned x
.