19#ifndef EIGEN_LEVENBERGMARQUARDT_H
20#define EIGEN_LEVENBERGMARQUARDT_H
23#include "./InternalHeaderCheck.h"
26namespace LevenbergMarquardtSpace {
30 ImproperInputParameters = 0,
31 RelativeReductionTooSmall = 1,
32 RelativeErrorTooSmall = 2,
33 RelativeErrorAndReductionTooSmall = 3,
35 TooManyFunctionEvaluation = 5,
43template <
typename Scalar_,
int NX = Dynamic,
int NY = Dynamic>
45 typedef Scalar_ Scalar;
46 enum { InputsAtCompileTime = NX, ValuesAtCompileTime = NY };
47 typedef Matrix<Scalar, InputsAtCompileTime, 1> InputType;
48 typedef Matrix<Scalar, ValuesAtCompileTime, 1> ValueType;
49 typedef Matrix<Scalar, ValuesAtCompileTime, InputsAtCompileTime> JacobianType;
50 typedef ColPivHouseholderQR<JacobianType> QRSolver;
51 const int m_inputs, m_values;
53 DenseFunctor() : m_inputs(InputsAtCompileTime), m_values(ValuesAtCompileTime) {}
54 DenseFunctor(
int inputs,
int values) : m_inputs(inputs), m_values(values) {}
56 int inputs()
const {
return m_inputs; }
57 int values()
const {
return m_values; }
66template <
typename Scalar_,
typename Index_>
68 typedef Scalar_ Scalar;
70 typedef Matrix<Scalar, Dynamic, 1> InputType;
71 typedef Matrix<Scalar, Dynamic, 1> ValueType;
72 typedef SparseMatrix<Scalar, ColMajor, Index> JacobianType;
73 typedef SparseQR<JacobianType, COLAMDOrdering<int> > QRSolver;
74 enum { InputsAtCompileTime =
Dynamic, ValuesAtCompileTime =
Dynamic };
76 SparseFunctor(
int inputs,
int values) : m_inputs(inputs), m_values(values) {}
78 int inputs()
const {
return m_inputs; }
79 int values()
const {
return m_values; }
81 const int m_inputs, m_values;
89template <
typename QRSolver,
typename VectorType>
90void lmpar2(
const QRSolver &qr,
const VectorType &diag,
const VectorType &qtb,
typename VectorType::Scalar m_delta,
91 typename VectorType::Scalar &par, VectorType &x);
101template <
typename FunctorType_>
102class LevenbergMarquardt : internal::no_assignment_operator {
104 typedef FunctorType_ FunctorType;
105 typedef typename FunctorType::QRSolver QRSolver;
106 typedef typename FunctorType::JacobianType JacobianType;
107 typedef typename JacobianType::Scalar Scalar;
108 typedef typename JacobianType::RealScalar RealScalar;
109 typedef typename QRSolver::StorageIndex PermIndex;
114 LevenbergMarquardt(FunctorType &functor)
115 : m_functor(functor),
120 m_isInitialized(
false),
123 m_useExternalScaling =
false;
126 LevenbergMarquardtSpace::Status minimize(FVectorType &x);
127 LevenbergMarquardtSpace::Status minimizeInit(FVectorType &x);
128 LevenbergMarquardtSpace::Status minimizeOneStep(FVectorType &x);
130 static LevenbergMarquardtSpace::Status lmdif1(FunctorType &functor, FVectorType &x, Index *nfev,
167 RealScalar
xtol()
const {
return m_xtol; }
170 RealScalar
ftol()
const {
return m_ftol; }
173 RealScalar
gtol()
const {
return m_gtol; }
176 RealScalar
factor()
const {
return m_factor; }
179 RealScalar
epsilon()
const {
return m_epsfcn; }
182 Index
maxfev()
const {
return m_maxfev; }
185 FVectorType &
diag() {
return m_diag; }
191 Index
nfev() {
return m_nfev; }
194 Index
njev() {
return m_njev; }
197 RealScalar
fnorm() {
return m_fnorm; }
200 RealScalar
gnorm() {
return m_gnorm; }
207 FVectorType &
fvec() {
return m_fvec; }
235 JacobianType m_rfactor;
236 FunctorType &m_functor;
237 FVectorType m_fvec, m_qtf, m_diag;
252 bool m_useExternalScaling;
253 PermutationType m_permutation;
254 FVectorType m_wa1, m_wa2, m_wa3, m_wa4;
256 bool m_isInitialized;
260template <
typename FunctorType>
261LevenbergMarquardtSpace::Status LevenbergMarquardt<FunctorType>::minimize(FVectorType &x) {
262 LevenbergMarquardtSpace::Status status = minimizeInit(x);
263 if (status == LevenbergMarquardtSpace::ImproperInputParameters) {
264 m_isInitialized =
true;
269 status = minimizeOneStep(x);
270 }
while (status == LevenbergMarquardtSpace::Running);
271 m_isInitialized =
true;
275template <
typename FunctorType>
276LevenbergMarquardtSpace::Status LevenbergMarquardt<FunctorType>::minimizeInit(FVectorType &x) {
278 m = m_functor.values();
288 if (!m_useExternalScaling) m_diag.resize(n);
289 eigen_assert((!m_useExternalScaling || m_diag.size() == n) &&
290 "When m_useExternalScaling is set, the caller must provide a valid 'm_diag'");
298 if (n <= 0 || m < n || m_ftol < 0. || m_xtol < 0. || m_gtol < 0. || m_maxfev <= 0 || m_factor <= 0.) {
300 return LevenbergMarquardtSpace::ImproperInputParameters;
303 if (m_useExternalScaling)
304 for (
Index j = 0; j < n; ++j)
305 if (m_diag[j] <= 0.) {
307 return LevenbergMarquardtSpace::ImproperInputParameters;
313 if (m_functor(x, m_fvec) < 0)
return LevenbergMarquardtSpace::UserAsked;
314 m_fnorm = m_fvec.stableNorm();
320 return LevenbergMarquardtSpace::NotStarted;
323template <
typename FunctorType>
324LevenbergMarquardtSpace::Status LevenbergMarquardt<FunctorType>::lmder1(FVectorType &x,
const Scalar tol) {
326 m = m_functor.values();
329 if (n <= 0 || m < n || tol < 0.)
return LevenbergMarquardtSpace::ImproperInputParameters;
334 m_maxfev = 100 * (n + 1);
339template <
typename FunctorType>
340LevenbergMarquardtSpace::Status LevenbergMarquardt<FunctorType>::lmdif1(FunctorType &functor, FVectorType &x,
343 Index m = functor.values();
346 if (n <= 0 || m < n || tol < 0.)
return LevenbergMarquardtSpace::ImproperInputParameters;
353 lm.setMaxfev(200 * (n + 1));
355 LevenbergMarquardtSpace::Status info = LevenbergMarquardtSpace::Status(lm.minimize(x));
356 if (nfev) *nfev = lm.nfev();
Performs non linear optimization over a non-linear function, using a variant of the Levenberg Marquar...
Definition LevenbergMarquardt.h:102
RealScalar epsilon() const
Definition LevenbergMarquardt.h:179
ComputationInfo info() const
Reports whether the minimization was successful.
Definition LevenbergMarquardt.h:231
void setFtol(RealScalar ftol)
Definition LevenbergMarquardt.h:149
RealScalar gtol() const
Definition LevenbergMarquardt.h:173
void resetParameters()
Definition LevenbergMarquardt.h:134
RealScalar ftol() const
Definition LevenbergMarquardt.h:170
FVectorType & diag()
Definition LevenbergMarquardt.h:185
Index maxfev() const
Definition LevenbergMarquardt.h:182
void setExternalScaling(bool value)
Definition LevenbergMarquardt.h:164
RealScalar xtol() const
Definition LevenbergMarquardt.h:167
RealScalar gnorm()
Definition LevenbergMarquardt.h:200
Index iterations()
Definition LevenbergMarquardt.h:188
JacobianType & matrixR()
Definition LevenbergMarquardt.h:216
Index njev()
Definition LevenbergMarquardt.h:194
void setGtol(RealScalar gtol)
Definition LevenbergMarquardt.h:152
RealScalar lm_param(void)
Definition LevenbergMarquardt.h:203
RealScalar fnorm()
Definition LevenbergMarquardt.h:197
void setEpsilon(RealScalar epsfcn)
Definition LevenbergMarquardt.h:158
void setXtol(RealScalar xtol)
Definition LevenbergMarquardt.h:146
FVectorType & fvec()
Definition LevenbergMarquardt.h:207
void setFactor(RealScalar factor)
Definition LevenbergMarquardt.h:155
Index nfev()
Definition LevenbergMarquardt.h:191
PermutationType permutation()
Definition LevenbergMarquardt.h:220
RealScalar factor() const
Definition LevenbergMarquardt.h:176
void setMaxfev(Index maxfev)
Definition LevenbergMarquardt.h:161
JacobianType & jacobian()
Definition LevenbergMarquardt.h:211
Definition NumericalDiff.h:35
Namespace containing all symbols from the Eigen library.
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index