19 #ifndef EIGEN_LEVENBERGMARQUARDT_H
20 #define EIGEN_LEVENBERGMARQUARDT_H
23 #include "./InternalHeaderCheck.h"
26 namespace LevenbergMarquardtSpace {
30 ImproperInputParameters = 0,
31 RelativeReductionTooSmall = 1,
32 RelativeErrorTooSmall = 2,
33 RelativeErrorAndReductionTooSmall = 3,
35 TooManyFunctionEvaluation = 5,
43 template <
typename Scalar_,
int NX=Dynamic,
int NY=Dynamic>
46 typedef Scalar_ Scalar;
48 InputsAtCompileTime = NX,
49 ValuesAtCompileTime = NY
51 typedef Matrix<Scalar,InputsAtCompileTime,1> InputType;
52 typedef Matrix<Scalar,ValuesAtCompileTime,1> ValueType;
53 typedef Matrix<Scalar,ValuesAtCompileTime,InputsAtCompileTime> JacobianType;
54 typedef ColPivHouseholderQR<JacobianType> QRSolver;
55 const int m_inputs, m_values;
57 DenseFunctor() : m_inputs(InputsAtCompileTime), m_values(ValuesAtCompileTime) {}
58 DenseFunctor(
int inputs,
int values) : m_inputs(inputs), m_values(values) {}
60 int inputs()
const {
return m_inputs; }
61 int values()
const {
return m_values; }
70 template <
typename Scalar_,
typename Index_>
73 typedef Scalar_ Scalar;
77 typedef SparseMatrix<Scalar, ColMajor, Index> JacobianType;
78 typedef SparseQR<JacobianType, COLAMDOrdering<int> > QRSolver;
84 SparseFunctor(
int inputs,
int values) : m_inputs(inputs), m_values(values) {}
86 int inputs()
const {
return m_inputs; }
87 int values()
const {
return m_values; }
89 const int m_inputs, m_values;
98 template <
typename QRSolver,
typename VectorType>
99 void lmpar2(
const QRSolver &qr,
const VectorType &diag,
const VectorType &qtb,
100 typename VectorType::Scalar m_delta,
typename VectorType::Scalar &par,
111 template<
typename FunctorType_>
115 typedef FunctorType_ FunctorType;
116 typedef typename FunctorType::QRSolver QRSolver;
117 typedef typename FunctorType::JacobianType JacobianType;
118 typedef typename JacobianType::Scalar Scalar;
119 typedef typename JacobianType::RealScalar RealScalar;
120 typedef typename QRSolver::StorageIndex PermIndex;
125 : m_functor(functor),m_nfev(0),m_njev(0),m_fnorm(0.0),m_gnorm(0),
129 m_useExternalScaling=
false;
132 LevenbergMarquardtSpace::Status minimize(
FVectorType &x);
133 LevenbergMarquardtSpace::Status minimizeInit(
FVectorType &x);
134 LevenbergMarquardtSpace::Status minimizeOneStep(
FVectorType &x);
135 LevenbergMarquardtSpace::Status lmder1(
139 static LevenbergMarquardtSpace::Status lmdif1(
140 FunctorType &functor,
181 RealScalar
xtol()
const {
return m_xtol; }
184 RealScalar
ftol()
const {
return m_ftol; }
187 RealScalar
gtol()
const {
return m_gtol; }
190 RealScalar
factor()
const {
return m_factor; }
193 RealScalar
epsilon()
const {
return m_epsfcn; }
205 Index
nfev() {
return m_nfev; }
208 Index
njev() {
return m_njev; }
211 RealScalar
fnorm() {
return m_fnorm; }
214 RealScalar
gnorm() {
return m_gnorm; }
252 JacobianType m_rfactor;
253 FunctorType &m_functor;
254 FVectorType m_fvec, m_qtf, m_diag;
269 bool m_useExternalScaling;
270 PermutationType m_permutation;
271 FVectorType m_wa1, m_wa2, m_wa3, m_wa4;
273 bool m_isInitialized;
277 template<
typename FunctorType>
278 LevenbergMarquardtSpace::Status
279 LevenbergMarquardt<FunctorType>::minimize(FVectorType &x)
281 LevenbergMarquardtSpace::Status status = minimizeInit(x);
282 if (status==LevenbergMarquardtSpace::ImproperInputParameters) {
283 m_isInitialized =
true;
288 status = minimizeOneStep(x);
289 }
while (status==LevenbergMarquardtSpace::Running);
290 m_isInitialized =
true;
294 template<
typename FunctorType>
295 LevenbergMarquardtSpace::Status
296 LevenbergMarquardt<FunctorType>::minimizeInit(FVectorType &x)
299 m = m_functor.values();
301 m_wa1.resize(n); m_wa2.resize(n); m_wa3.resize(n);
307 if (!m_useExternalScaling)
309 eigen_assert( (!m_useExternalScaling || m_diag.size()==n) &&
"When m_useExternalScaling is set, the caller must provide a valid 'm_diag'");
317 if (n <= 0 || m < n || m_ftol < 0. || m_xtol < 0. || m_gtol < 0. || m_maxfev <= 0 || m_factor <= 0.){
319 return LevenbergMarquardtSpace::ImproperInputParameters;
322 if (m_useExternalScaling)
323 for (
Index j = 0; j < n; ++j)
327 return LevenbergMarquardtSpace::ImproperInputParameters;
333 if ( m_functor(x, m_fvec) < 0)
334 return LevenbergMarquardtSpace::UserAsked;
335 m_fnorm = m_fvec.stableNorm();
341 return LevenbergMarquardtSpace::NotStarted;
344 template<
typename FunctorType>
345 LevenbergMarquardtSpace::Status
346 LevenbergMarquardt<FunctorType>::lmder1(
352 m = m_functor.values();
355 if (n <= 0 || m < n || tol < 0.)
356 return LevenbergMarquardtSpace::ImproperInputParameters;
361 m_maxfev = 100*(n+1);
367 template<
typename FunctorType>
368 LevenbergMarquardtSpace::Status
369 LevenbergMarquardt<FunctorType>::lmdif1(
370 FunctorType &functor,
377 Index m = functor.values();
380 if (n <= 0 || m < n || tol < 0.)
381 return LevenbergMarquardtSpace::ImproperInputParameters;
383 NumericalDiff<FunctorType> numDiff(functor);
385 LevenbergMarquardt<NumericalDiff<FunctorType> > lm(numDiff);
388 lm.setMaxfev(200*(n+1));
390 LevenbergMarquardtSpace::Status info = LevenbergMarquardtSpace::Status(lm.minimize(x));
Performs non linear optimization over a non-linear function, using a variant of the Levenberg Marquar...
Definition: LevenbergMarquardt/LevenbergMarquardt.h:113
RealScalar epsilon() const
Definition: LevenbergMarquardt/LevenbergMarquardt.h:193
ComputationInfo info() const
Reports whether the minimization was successful.
Definition: LevenbergMarquardt/LevenbergMarquardt.h:245
void setFtol(RealScalar ftol)
Definition: LevenbergMarquardt/LevenbergMarquardt.h:163
RealScalar gtol() const
Definition: LevenbergMarquardt/LevenbergMarquardt.h:187
void resetParameters()
Definition: LevenbergMarquardt/LevenbergMarquardt.h:147
RealScalar ftol() const
Definition: LevenbergMarquardt/LevenbergMarquardt.h:184
Index maxfev() const
Definition: LevenbergMarquardt/LevenbergMarquardt.h:196
void setExternalScaling(bool value)
Definition: LevenbergMarquardt/LevenbergMarquardt.h:178
RealScalar xtol() const
Definition: LevenbergMarquardt/LevenbergMarquardt.h:181
RealScalar gnorm()
Definition: LevenbergMarquardt/LevenbergMarquardt.h:214
FVectorType & fvec()
Definition: LevenbergMarquardt/LevenbergMarquardt.h:221
Index iterations()
Definition: LevenbergMarquardt/LevenbergMarquardt.h:202
JacobianType & jacobian()
Definition: LevenbergMarquardt/LevenbergMarquardt.h:225
FVectorType & diag()
Definition: LevenbergMarquardt/LevenbergMarquardt.h:199
Index njev()
Definition: LevenbergMarquardt/LevenbergMarquardt.h:208
void setGtol(RealScalar gtol)
Definition: LevenbergMarquardt/LevenbergMarquardt.h:166
RealScalar lm_param(void)
Definition: LevenbergMarquardt/LevenbergMarquardt.h:217
RealScalar fnorm()
Definition: LevenbergMarquardt/LevenbergMarquardt.h:211
void setEpsilon(RealScalar epsfcn)
Definition: LevenbergMarquardt/LevenbergMarquardt.h:172
void setXtol(RealScalar xtol)
Definition: LevenbergMarquardt/LevenbergMarquardt.h:160
void setFactor(RealScalar factor)
Definition: LevenbergMarquardt/LevenbergMarquardt.h:169
Index nfev()
Definition: LevenbergMarquardt/LevenbergMarquardt.h:205
PermutationType permutation()
Definition: LevenbergMarquardt/LevenbergMarquardt.h:234
RealScalar factor() const
Definition: LevenbergMarquardt/LevenbergMarquardt.h:190
JacobianType & matrixR()
Definition: LevenbergMarquardt/LevenbergMarquardt.h:230
void setMaxfev(Index maxfev)
Definition: LevenbergMarquardt/LevenbergMarquardt.h:175
Namespace containing all symbols from the Eigen library.
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)