19 #ifndef EIGEN_LEVENBERGMARQUARDT_H
20 #define EIGEN_LEVENBERGMARQUARDT_H
24 namespace LevenbergMarquardtSpace {
28 ImproperInputParameters = 0,
29 RelativeReductionTooSmall = 1,
30 RelativeErrorTooSmall = 2,
31 RelativeErrorAndReductionTooSmall = 3,
33 TooManyFunctionEvaluation = 5,
41 template <
typename _Scalar,
int NX=Dynamic,
int NY=Dynamic>
44 typedef _Scalar Scalar;
46 InputsAtCompileTime = NX,
47 ValuesAtCompileTime = NY
49 typedef Matrix<Scalar,InputsAtCompileTime,1> InputType;
50 typedef Matrix<Scalar,ValuesAtCompileTime,1> ValueType;
51 typedef Matrix<Scalar,ValuesAtCompileTime,InputsAtCompileTime> JacobianType;
52 typedef ColPivHouseholderQR<JacobianType> QRSolver;
53 const int m_inputs, m_values;
55 DenseFunctor() : m_inputs(InputsAtCompileTime), m_values(ValuesAtCompileTime) {}
56 DenseFunctor(
int inputs,
int values) : m_inputs(inputs), m_values(values) {}
58 int inputs()
const {
return m_inputs; }
59 int values()
const {
return m_values; }
68 template <
typename _Scalar,
typename _Index>
71 typedef _Scalar Scalar;
73 typedef Matrix<Scalar,Dynamic,1> InputType;
74 typedef Matrix<Scalar,Dynamic,1> ValueType;
75 typedef SparseMatrix<Scalar, ColMajor, Index> JacobianType;
76 typedef SparseQR<JacobianType, COLAMDOrdering<int> > QRSolver;
78 InputsAtCompileTime = Dynamic,
79 ValuesAtCompileTime = Dynamic
82 SparseFunctor(
int inputs,
int values) : m_inputs(inputs), m_values(values) {}
84 int inputs()
const {
return m_inputs; }
85 int values()
const {
return m_values; }
87 const int m_inputs, m_values;
96 template <
typename QRSolver,
typename VectorType>
97 void lmpar2(
const QRSolver &qr,
const VectorType &diag,
const VectorType &qtb,
98 typename VectorType::Scalar m_delta,
typename VectorType::Scalar &par,
109 template<
typename _FunctorType>
113 typedef _FunctorType FunctorType;
114 typedef typename FunctorType::QRSolver QRSolver;
115 typedef typename FunctorType::JacobianType JacobianType;
116 typedef typename JacobianType::Scalar Scalar;
117 typedef typename JacobianType::RealScalar RealScalar;
118 typedef typename JacobianType::Index Index;
119 typedef typename QRSolver::Index PermIndex;
120 typedef Matrix<Scalar,Dynamic,1> FVectorType;
121 typedef PermutationMatrix<Dynamic,Dynamic> PermutationType;
124 : m_functor(functor),m_nfev(0),m_njev(0),m_fnorm(0.0),m_gnorm(0),
125 m_isInitialized(
false),m_info(InvalidInput)
128 m_useExternalScaling=
false;
131 LevenbergMarquardtSpace::Status minimize(FVectorType &x);
132 LevenbergMarquardtSpace::Status minimizeInit(FVectorType &x);
133 LevenbergMarquardtSpace::Status minimizeOneStep(FVectorType &x);
134 LevenbergMarquardtSpace::Status lmder1(
136 const Scalar tol = std::sqrt(NumTraits<Scalar>::epsilon())
138 static LevenbergMarquardtSpace::Status lmdif1(
139 FunctorType &functor,
142 const Scalar tol = std::sqrt(NumTraits<Scalar>::epsilon())
150 m_ftol = std::sqrt(NumTraits<RealScalar>::epsilon());
151 m_xtol = std::sqrt(NumTraits<RealScalar>::epsilon());
157 void setXtol(RealScalar xtol) { m_xtol = xtol; }
160 void setFtol(RealScalar ftol) { m_ftol = ftol; }
163 void setGtol(RealScalar gtol) { m_gtol = gtol; }
166 void setFactor(RealScalar factor) { m_factor = factor; }
178 FVectorType&
diag() {
return m_diag; }
184 Index
nfev() {
return m_nfev; }
187 Index
njev() {
return m_njev; }
190 RealScalar
fnorm() {
return m_fnorm; }
193 RealScalar
gnorm() {
return m_gnorm; }
200 FVectorType&
fvec() {
return m_fvec; }
231 JacobianType m_rfactor;
232 FunctorType &m_functor;
233 FVectorType m_fvec, m_qtf, m_diag;
248 bool m_useExternalScaling;
249 PermutationType m_permutation;
250 FVectorType m_wa1, m_wa2, m_wa3, m_wa4;
252 bool m_isInitialized;
253 ComputationInfo m_info;
256 template<
typename FunctorType>
257 LevenbergMarquardtSpace::Status
258 LevenbergMarquardt<FunctorType>::minimize(FVectorType &x)
260 LevenbergMarquardtSpace::Status status = minimizeInit(x);
261 if (status==LevenbergMarquardtSpace::ImproperInputParameters) {
262 m_isInitialized =
true;
267 status = minimizeOneStep(x);
268 }
while (status==LevenbergMarquardtSpace::Running);
269 m_isInitialized =
true;
273 template<
typename FunctorType>
274 LevenbergMarquardtSpace::Status
275 LevenbergMarquardt<FunctorType>::minimizeInit(FVectorType &x)
278 m = m_functor.values();
280 m_wa1.resize(n); m_wa2.resize(n); m_wa3.resize(n);
286 if (!m_useExternalScaling)
288 eigen_assert( (!m_useExternalScaling || m_diag.size()==n) ||
"When m_useExternalScaling is set, the caller must provide a valid 'm_diag'");
296 if (n <= 0 || m < n || m_ftol < 0. || m_xtol < 0. || m_gtol < 0. || m_maxfev <= 0 || m_factor <= 0.){
297 m_info = InvalidInput;
298 return LevenbergMarquardtSpace::ImproperInputParameters;
301 if (m_useExternalScaling)
302 for (Index j = 0; j < n; ++j)
305 m_info = InvalidInput;
306 return LevenbergMarquardtSpace::ImproperInputParameters;
312 if ( m_functor(x, m_fvec) < 0)
313 return LevenbergMarquardtSpace::UserAsked;
314 m_fnorm = m_fvec.stableNorm();
320 return LevenbergMarquardtSpace::NotStarted;
323 template<
typename FunctorType>
324 LevenbergMarquardtSpace::Status
325 LevenbergMarquardt<FunctorType>::lmder1(
331 m = m_functor.values();
334 if (n <= 0 || m < n || tol < 0.)
335 return LevenbergMarquardtSpace::ImproperInputParameters;
340 m_maxfev = 100*(n+1);
346 template<
typename FunctorType>
347 LevenbergMarquardtSpace::Status
348 LevenbergMarquardt<FunctorType>::lmdif1(
349 FunctorType &functor,
356 Index m = functor.values();
359 if (n <= 0 || m < n || tol < 0.)
360 return LevenbergMarquardtSpace::ImproperInputParameters;
362 NumericalDiff<FunctorType> numDiff(functor);
364 LevenbergMarquardt<NumericalDiff<FunctorType> > lm(numDiff);
367 lm.setMaxfev(200*(n+1));
369 LevenbergMarquardtSpace::Status info = LevenbergMarquardtSpace::Status(lm.minimize(x));
377 #endif // EIGEN_LEVENBERGMARQUARDT_H
FVectorType & fvec()
Definition: LevenbergMarquardt/LevenbergMarquardt.h:200
JacobianType & jacobian()
Definition: LevenbergMarquardt/LevenbergMarquardt.h:204
Performs non linear optimization over a non-linear function, using a variant of the Levenberg Marquar...
Definition: LevenbergMarquardt/LevenbergMarquardt.h:110
void setEpsilon(RealScalar epsfcn)
Definition: LevenbergMarquardt/LevenbergMarquardt.h:169
Index nfev()
Definition: LevenbergMarquardt/LevenbergMarquardt.h:184
FVectorType & diag()
Definition: LevenbergMarquardt/LevenbergMarquardt.h:178
ComputationInfo info() const
Reports whether the minimization was successful.
Definition: LevenbergMarquardt/LevenbergMarquardt.h:224
void setMaxfev(Index maxfev)
Definition: LevenbergMarquardt/LevenbergMarquardt.h:172
RealScalar lm_param(void)
Definition: LevenbergMarquardt/LevenbergMarquardt.h:196
void setFtol(RealScalar ftol)
Definition: LevenbergMarquardt/LevenbergMarquardt.h:160
void resetParameters()
Definition: LevenbergMarquardt/LevenbergMarquardt.h:146
RealScalar fnorm()
Definition: LevenbergMarquardt/LevenbergMarquardt.h:190
void setXtol(RealScalar xtol)
Definition: LevenbergMarquardt/LevenbergMarquardt.h:157
RealScalar gnorm()
Definition: LevenbergMarquardt/LevenbergMarquardt.h:193
Index njev()
Definition: LevenbergMarquardt/LevenbergMarquardt.h:187
JacobianType & matrixR()
Definition: LevenbergMarquardt/LevenbergMarquardt.h:209
PermutationType permutation()
Definition: LevenbergMarquardt/LevenbergMarquardt.h:213
void setFactor(RealScalar factor)
Definition: LevenbergMarquardt/LevenbergMarquardt.h:166
Index iterations()
Definition: LevenbergMarquardt/LevenbergMarquardt.h:181
void setGtol(RealScalar gtol)
Definition: LevenbergMarquardt/LevenbergMarquardt.h:163
void setExternalScaling(bool value)
Definition: LevenbergMarquardt/LevenbergMarquardt.h:175