Eigen::LevenbergMarquardt< FunctorType_ > Class Template Reference

Performs non linear optimization over a non-linear function, using a variant of the Levenberg Marquardt algorithm. More...

Inherits internal::no_assignment_operator.

Classes

struct  Parameters
 

Public Types

typedef FunctorType_ FunctorType
 
typedef Matrix< Scalar, Dynamic, 1 > FVectorType
 
typedef Matrix< Scalar, Dynamic, 1 > FVectorType
 
typedef DenseIndex Index
 
typedef FunctorType::JacobianType JacobianType
 
typedef Matrix< Scalar, Dynamic, DynamicJacobianType
 
typedef QRSolver::StorageIndex PermIndex
 
typedef PermutationMatrix< Dynamic, Dynamic, int > PermutationType
 
typedef FunctorType::QRSolver QRSolver
 
typedef JacobianType::RealScalar RealScalar
 
typedef JacobianType::Scalar Scalar
 

Public Member Functions

FVectorTypediag ()
 
RealScalar epsilon () const
 
RealScalar factor () const
 
RealScalar fnorm ()
 
RealScalar ftol () const
 
FVectorTypefvec ()
 
RealScalar gnorm ()
 
RealScalar gtol () const
 
ComputationInfo info () const
 Reports whether the minimization was successful. More...
 
Index iterations ()
 
JacobianTypejacobian ()
 
 LevenbergMarquardt (FunctorType &_functor)
 
 LevenbergMarquardt (FunctorType &functor)
 
RealScalar lm_param (void)
 
Scalar lm_param (void)
 
LevenbergMarquardtSpace::Status lmder1 (FVectorType &x, const Scalar tol=sqrt_epsilon())
 
LevenbergMarquardtSpace::Status lmder1 (FVectorType &x, const Scalar tol=std::sqrt(NumTraits< Scalar >::epsilon()))
 
LevenbergMarquardtSpace::Status lmstr1 (FVectorType &x, const Scalar tol=sqrt_epsilon())
 
JacobianTypematrixR ()
 
Index maxfev () const
 
LevenbergMarquardtSpace::Status minimize (FVectorType &x)
 
LevenbergMarquardtSpace::Status minimize (FVectorType &x)
 
LevenbergMarquardtSpace::Status minimizeInit (FVectorType &x)
 
LevenbergMarquardtSpace::Status minimizeInit (FVectorType &x)
 
LevenbergMarquardtSpace::Status minimizeOneStep (FVectorType &x)
 
LevenbergMarquardtSpace::Status minimizeOneStep (FVectorType &x)
 
LevenbergMarquardtSpace::Status minimizeOptimumStorage (FVectorType &x)
 
LevenbergMarquardtSpace::Status minimizeOptimumStorageInit (FVectorType &x)
 
LevenbergMarquardtSpace::Status minimizeOptimumStorageOneStep (FVectorType &x)
 
Index nfev ()
 
Index njev ()
 
PermutationType permutation ()
 
void resetParameters ()
 
void resetParameters (void)
 
void setEpsilon (RealScalar epsfcn)
 
void setExternalScaling (bool value)
 
void setFactor (RealScalar factor)
 
void setFtol (RealScalar ftol)
 
void setGtol (RealScalar gtol)
 
void setMaxfev (Index maxfev)
 
void setXtol (RealScalar xtol)
 
RealScalar xtol () const
 

Static Public Member Functions

static LevenbergMarquardtSpace::Status lmdif1 (FunctorType &functor, FVectorType &x, Index *nfev, const Scalar tol=sqrt_epsilon())
 
static LevenbergMarquardtSpace::Status lmdif1 (FunctorType &functor, FVectorType &x, Index *nfev, const Scalar tol=std::sqrt(NumTraits< Scalar >::epsilon()))
 

Public Attributes

FVectorType diag
 
JacobianType fjac
 
Scalar fnorm
 
FVectorType fvec
 
Scalar gnorm
 
Index iter
 
Index nfev
 
Index njev
 
Parameters parameters
 
PermutationMatrix< Dynamic, Dynamicpermutation
 
FVectorType qtf
 
bool useExternalScaling
 

Private Member Functions

LevenbergMarquardtoperator= (const LevenbergMarquardt &)
 

Static Private Member Functions

static Scalar sqrt_epsilon ()
 

Private Attributes

Scalar actred
 
Scalar delta
 
Scalar dirder
 
Scalar fnorm1
 
FunctorTypefunctor
 
Index m
 
RealScalar m_delta
 
FVectorType m_diag
 
RealScalar m_epsfcn
 
RealScalar m_factor
 
JacobianType m_fjac
 
RealScalar m_fnorm
 
RealScalar m_ftol
 
FunctorTypem_functor
 
FVectorType m_fvec
 
RealScalar m_gnorm
 
RealScalar m_gtol
 
ComputationInfo m_info
 
bool m_isInitialized
 
Index m_iter
 
Index m_maxfev
 
Index m_nfev
 
Index m_njev
 
RealScalar m_par
 
PermutationType m_permutation
 
FVectorType m_qtf
 
JacobianType m_rfactor
 
bool m_useExternalScaling
 
FVectorType m_wa1
 
FVectorType m_wa2
 
FVectorType m_wa3
 
FVectorType m_wa4
 
RealScalar m_xtol
 
Index n
 
Scalar par
 
Scalar pnorm
 
Scalar prered
 
Scalar ratio
 
Scalar sum
 
Scalar temp
 
Scalar temp1
 
Scalar temp2
 
FVectorType wa1
 
FVectorType wa2
 
FVectorType wa3
 
FVectorType wa4
 
Scalar xnorm
 

Detailed Description

template<typename FunctorType_>
class Eigen::LevenbergMarquardt< FunctorType_ >

Performs non linear optimization over a non-linear function, using a variant of the Levenberg Marquardt algorithm.

Check wikipedia for more information. http://en.wikipedia.org/wiki/Levenberg%E2%80%93Marquardt_algorithm

Definition at line 112 of file LevenbergMarquardt/LevenbergMarquardt.h.

Member Typedef Documentation

◆ FunctorType

template<typename FunctorType_ >
typedef FunctorType_ Eigen::LevenbergMarquardt< FunctorType_ >::FunctorType

Definition at line 115 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ FVectorType [1/2]

template<typename FunctorType_ >
typedef Matrix<Scalar,Dynamic,1> Eigen::LevenbergMarquardt< FunctorType_ >::FVectorType

Definition at line 121 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ FVectorType [2/2]

template<typename FunctorType_ >
typedef Matrix< Scalar, Dynamic, 1 > Eigen::LevenbergMarquardt< FunctorType_ >::FVectorType

Definition at line 78 of file NonLinearOptimization/LevenbergMarquardt.h.

◆ Index

template<typename FunctorType_ >
typedef DenseIndex Eigen::LevenbergMarquardt< FunctorType_ >::Index

Definition at line 60 of file NonLinearOptimization/LevenbergMarquardt.h.

◆ JacobianType [1/2]

template<typename FunctorType_ >
typedef FunctorType::JacobianType Eigen::LevenbergMarquardt< FunctorType_ >::JacobianType

Definition at line 117 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ JacobianType [2/2]

template<typename FunctorType_ >
typedef Matrix< Scalar, Dynamic, Dynamic > Eigen::LevenbergMarquardt< FunctorType_ >::JacobianType

Definition at line 79 of file NonLinearOptimization/LevenbergMarquardt.h.

◆ PermIndex

template<typename FunctorType_ >
typedef QRSolver::StorageIndex Eigen::LevenbergMarquardt< FunctorType_ >::PermIndex

Definition at line 120 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ PermutationType

template<typename FunctorType_ >
typedef PermutationMatrix<Dynamic,Dynamic,int> Eigen::LevenbergMarquardt< FunctorType_ >::PermutationType

Definition at line 122 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ QRSolver

template<typename FunctorType_ >
typedef FunctorType::QRSolver Eigen::LevenbergMarquardt< FunctorType_ >::QRSolver

Definition at line 116 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ RealScalar

template<typename FunctorType_ >
typedef JacobianType::RealScalar Eigen::LevenbergMarquardt< FunctorType_ >::RealScalar

Definition at line 119 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ Scalar

template<typename FunctorType_ >
typedef JacobianType::Scalar Eigen::LevenbergMarquardt< FunctorType_ >::Scalar

Definition at line 118 of file LevenbergMarquardt/LevenbergMarquardt.h.

Constructor & Destructor Documentation

◆ LevenbergMarquardt() [1/2]

template<typename FunctorType_ >
Eigen::LevenbergMarquardt< FunctorType_ >::LevenbergMarquardt ( FunctorType functor)
inline

◆ LevenbergMarquardt() [2/2]

Member Function Documentation

◆ diag()

template<typename FunctorType_ >
FVectorType& Eigen::LevenbergMarquardt< FunctorType_ >::diag ( )
inline
Returns
a reference to the diagonal of the jacobian

Definition at line 199 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ epsilon()

template<typename FunctorType_ >
RealScalar Eigen::LevenbergMarquardt< FunctorType_ >::epsilon ( ) const
inline
Returns
the error precision

Definition at line 193 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ factor()

template<typename FunctorType_ >
RealScalar Eigen::LevenbergMarquardt< FunctorType_ >::factor ( ) const
inline
Returns
the step bound for the diagonal shift

Definition at line 190 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ fnorm()

template<typename FunctorType_ >
RealScalar Eigen::LevenbergMarquardt< FunctorType_ >::fnorm ( )
inline
Returns
the norm of current vector function

Definition at line 211 of file LevenbergMarquardt/LevenbergMarquardt.h.

211 {return m_fnorm; }

◆ ftol()

template<typename FunctorType_ >
RealScalar Eigen::LevenbergMarquardt< FunctorType_ >::ftol ( ) const
inline
Returns
the tolerance for the norm of the vector function

Definition at line 184 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ fvec()

template<typename FunctorType_ >
FVectorType& Eigen::LevenbergMarquardt< FunctorType_ >::fvec ( )
inline
Returns
a reference to the current vector function

Definition at line 221 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ gnorm()

template<typename FunctorType_ >
RealScalar Eigen::LevenbergMarquardt< FunctorType_ >::gnorm ( )
inline
Returns
the norm of the gradient of the error

Definition at line 214 of file LevenbergMarquardt/LevenbergMarquardt.h.

214 {return m_gnorm; }

◆ gtol()

template<typename FunctorType_ >
RealScalar Eigen::LevenbergMarquardt< FunctorType_ >::gtol ( ) const
inline
Returns
the tolerance for the norm of the gradient of the error vector

Definition at line 187 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ info()

template<typename FunctorType_ >
ComputationInfo Eigen::LevenbergMarquardt< FunctorType_ >::info ( ) const
inline

Reports whether the minimization was successful.

Returns
Success if the minimization was successful, NumericalIssue if a numerical problem arises during the minimization process, for example during the QR factorization NoConvergence if the minimization did not converge after the maximum number of function evaluation allowed InvalidInput if the input matrix is invalid

Definition at line 245 of file LevenbergMarquardt/LevenbergMarquardt.h.

246  {
247 
248  return m_info;
249  }

◆ iterations()

template<typename FunctorType_ >
Index Eigen::LevenbergMarquardt< FunctorType_ >::iterations ( )
inline
Returns
the number of iterations performed

Definition at line 202 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ jacobian()

template<typename FunctorType_ >
JacobianType& Eigen::LevenbergMarquardt< FunctorType_ >::jacobian ( )
inline
Returns
a reference to the matrix where the current Jacobian matrix is stored

Definition at line 225 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ lm_param() [1/2]

template<typename FunctorType_ >
RealScalar Eigen::LevenbergMarquardt< FunctorType_ >::lm_param ( void  )
inline
Returns
the LevenbergMarquardt parameter

Definition at line 217 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ lm_param() [2/2]

template<typename FunctorType_ >
Scalar Eigen::LevenbergMarquardt< FunctorType_ >::lm_param ( void  )
inline

◆ lmder1() [1/2]

template<typename FunctorType_ >
LevenbergMarquardtSpace::Status Eigen::LevenbergMarquardt< FunctorType_ >::lmder1 ( FVectorType x,
const Scalar  tol = sqrt_epsilon() 
)

◆ lmder1() [2/2]

template<typename FunctorType , typename Scalar >
LevenbergMarquardtSpace::Status Eigen::LevenbergMarquardt< FunctorType, Scalar >::lmder1 ( FVectorType x,
const Scalar  tol = std::sqrt(NumTraits<Scalar>::epsilon()) 
)

Definition at line 346 of file LevenbergMarquardt/LevenbergMarquardt.h.

350 {
351  n = x.size();
352  m = m_functor.values();
353 
354  /* check the input parameters for errors. */
355  if (n <= 0 || m < n || tol < 0.)
357 
358  resetParameters();
359  m_ftol = tol;
360  m_xtol = tol;
361  m_maxfev = 100*(n+1);
362 
363  return minimize(x);
364 }
LevenbergMarquardtSpace::Status minimize(FVectorType &x)

◆ lmdif1() [1/2]

template<typename FunctorType_ >
static LevenbergMarquardtSpace::Status Eigen::LevenbergMarquardt< FunctorType_ >::lmdif1 ( FunctorType functor,
FVectorType x,
Index nfev,
const Scalar  tol = sqrt_epsilon() 
)
static

◆ lmdif1() [2/2]

template<typename FunctorType , typename Scalar >
LevenbergMarquardtSpace::Status Eigen::LevenbergMarquardt< FunctorType, Scalar >::lmdif1 ( FunctorType functor,
FVectorType x,
Index nfev,
const Scalar  tol = std::sqrt(NumTraits<Scalar>::epsilon()) 
)
static

Definition at line 369 of file LevenbergMarquardt/LevenbergMarquardt.h.

375 {
376  Index n = x.size();
377  Index m = functor.values();
378 
379  /* check the input parameters for errors. */
380  if (n <= 0 || m < n || tol < 0.)
382 
383  NumericalDiff<FunctorType> numDiff(functor);
384  // embedded LevenbergMarquardt
385  LevenbergMarquardt<NumericalDiff<FunctorType> > lm(numDiff);
386  lm.setFtol(tol);
387  lm.setXtol(tol);
388  lm.setMaxfev(200*(n+1));
389 
391  if (nfev)
392  * nfev = lm.nfev();
393  return info;
394 }
ComputationInfo info() const
Reports whether the minimization was successful.
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index

◆ lmstr1()

template<typename FunctorType , typename Scalar >
LevenbergMarquardtSpace::Status Eigen::LevenbergMarquardt< FunctorType, Scalar >::lmstr1 ( FVectorType x,
const Scalar  tol = sqrt_epsilon() 
)

Definition at line 365 of file NonLinearOptimization/LevenbergMarquardt.h.

369 {
370  n = x.size();
371  m = functor.values();
372 
373  /* check the input parameters for errors. */
374  if (n <= 0 || m < n || tol < 0.)
376 
377  resetParameters();
378  parameters.ftol = tol;
379  parameters.xtol = tol;
380  parameters.maxfev = 100*(n+1);
381 
382  return minimizeOptimumStorage(x);
383 }
LevenbergMarquardtSpace::Status minimizeOptimumStorage(FVectorType &x)

◆ matrixR()

template<typename FunctorType_ >
JacobianType& Eigen::LevenbergMarquardt< FunctorType_ >::matrixR ( )
inline
Returns
a reference to the triangular matrix R from the QR of the jacobian matrix.
See also
jacobian()

Definition at line 230 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ maxfev()

template<typename FunctorType_ >
Index Eigen::LevenbergMarquardt< FunctorType_ >::maxfev ( ) const
inline
Returns
the maximum number of function evaluation

Definition at line 196 of file LevenbergMarquardt/LevenbergMarquardt.h.

196 {return m_maxfev; }

◆ minimize() [1/2]

template<typename FunctorType , typename Scalar >
LevenbergMarquardtSpace::Status Eigen::LevenbergMarquardt< FunctorType, Scalar >::minimize ( FVectorType x)

Definition at line 279 of file LevenbergMarquardt/LevenbergMarquardt.h.

280 {
283  m_isInitialized = true;
284  return status;
285  }
286  do {
287 // std::cout << " uv " << x.transpose() << "\n";
288  status = minimizeOneStep(x);
289  } while (status==LevenbergMarquardtSpace::Running);
290  m_isInitialized = true;
291  return status;
292 }
LevenbergMarquardtSpace::Status minimizeOneStep(FVectorType &x)
Definition: LMonestep.h:23
LevenbergMarquardtSpace::Status minimizeInit(FVectorType &x)

◆ minimize() [2/2]

template<typename FunctorType_ >
LevenbergMarquardtSpace::Status Eigen::LevenbergMarquardt< FunctorType_ >::minimize ( FVectorType x)

◆ minimizeInit() [1/2]

template<typename FunctorType , typename Scalar >
LevenbergMarquardtSpace::Status Eigen::LevenbergMarquardt< FunctorType, Scalar >::minimizeInit ( FVectorType x)

Definition at line 296 of file LevenbergMarquardt/LevenbergMarquardt.h.

297 {
298  n = x.size();
299  m = m_functor.values();
300 
302  m_wa4.resize(m);
303  m_fvec.resize(m);
304  //FIXME Sparse Case : Allocate space for the jacobian
305  m_fjac.resize(m, n);
306 // m_fjac.reserve(VectorXi::Constant(n,5)); // FIXME Find a better alternative
308  m_diag.resize(n);
309  eigen_assert( (!m_useExternalScaling || m_diag.size()==n) && "When m_useExternalScaling is set, the caller must provide a valid 'm_diag'");
310  m_qtf.resize(n);
311 
312  /* Function Body */
313  m_nfev = 0;
314  m_njev = 0;
315 
316  /* check the input parameters for errors. */
317  if (n <= 0 || m < n || m_ftol < 0. || m_xtol < 0. || m_gtol < 0. || m_maxfev <= 0 || m_factor <= 0.){
320  }
321 
323  for (Index j = 0; j < n; ++j)
324  if (m_diag[j] <= 0.)
325  {
328  }
329 
330  /* evaluate the function at the starting point */
331  /* and calculate its norm. */
332  m_nfev = 1;
333  if ( m_functor(x, m_fvec) < 0)
335  m_fnorm = m_fvec.stableNorm();
336 
337  /* initialize levenberg-marquardt parameter and iteration counter. */
338  m_par = 0.;
339  m_iter = 1;
340 
342 }
#define eigen_assert(x)
constexpr void resize(Index rows, Index cols)
std::ptrdiff_t j

◆ minimizeInit() [2/2]

template<typename FunctorType_ >
LevenbergMarquardtSpace::Status Eigen::LevenbergMarquardt< FunctorType_ >::minimizeInit ( FVectorType x)

◆ minimizeOneStep() [1/2]

template<typename FunctorType , typename Scalar >
LevenbergMarquardtSpace::Status Eigen::LevenbergMarquardt< FunctorType, Scalar >::minimizeOneStep ( FVectorType x)

Definition at line 23 of file LMonestep.h.

24 {
25  using std::abs;
26  using std::sqrt;
30  eigen_assert(x.size()==n); // check the caller is not cheating us
31 
32  temp = 0.0; xnorm = 0.0;
33  /* calculate the jacobian matrix. */
34  Index df_ret = m_functor.df(x, m_fjac);
35  if (df_ret<0)
37  if (df_ret>0)
38  // numerical diff, we evaluated the function df_ret times
39  m_nfev += df_ret;
40  else m_njev++;
41 
42  /* compute the qr factorization of the jacobian. */
43  for (int j = 0; j < x.size(); ++j)
44  m_wa2(j) = m_fjac.col(j).blueNorm();
45  QRSolver qrfac(m_fjac);
46  if(qrfac.info() != Success) {
49  }
50  // Make a copy of the first factor with the associated permutation
51  m_rfactor = qrfac.matrixR();
52  m_permutation = (qrfac.colsPermutation());
53 
54  /* on the first iteration and if external scaling is not used, scale according */
55  /* to the norms of the columns of the initial jacobian. */
56  if (m_iter == 1) {
58  for (Index j = 0; j < n; ++j)
59  m_diag[j] = (m_wa2[j]==0.)? 1. : m_wa2[j];
60 
61  /* on the first iteration, calculate the norm of the scaled x */
62  /* and initialize the step bound m_delta. */
63  xnorm = m_diag.cwiseProduct(x).stableNorm();
65  if (m_delta == 0.)
66  m_delta = m_factor;
67  }
68 
69  /* form (q transpose)*m_fvec and store the first n components in */
70  /* m_qtf. */
71  m_wa4 = m_fvec;
72  m_wa4 = qrfac.matrixQ().adjoint() * m_fvec;
73  m_qtf = m_wa4.head(n);
74 
75  /* compute the norm of the scaled gradient. */
76  m_gnorm = 0.;
77  if (m_fnorm != 0.)
78  for (Index j = 0; j < n; ++j)
79  if (m_wa2[m_permutation.indices()[j]] != 0.)
80  m_gnorm = (std::max)(m_gnorm, abs( m_rfactor.col(j).head(j+1).dot(m_qtf.head(j+1)/m_fnorm) / m_wa2[m_permutation.indices()[j]]));
81 
82  /* test for convergence of the gradient norm. */
83  if (m_gnorm <= m_gtol) {
84  m_info = Success;
86  }
87 
88  /* rescale if necessary. */
90  m_diag = m_diag.cwiseMax(m_wa2);
91 
92  do {
93  /* determine the levenberg-marquardt parameter. */
95 
96  /* store the direction p and x + p. calculate the norm of p. */
97  m_wa1 = -m_wa1;
98  m_wa2 = x + m_wa1;
99  pnorm = m_diag.cwiseProduct(m_wa1).stableNorm();
100 
101  /* on the first iteration, adjust the initial step bound. */
102  if (m_iter == 1)
104 
105  /* evaluate the function at x + p and calculate its norm. */
106  if ( m_functor(m_wa2, m_wa4) < 0)
108  ++m_nfev;
109  fnorm1 = m_wa4.stableNorm();
110 
111  /* compute the scaled actual reduction. */
112  actred = -1.;
113  if (Scalar(.1) * fnorm1 < m_fnorm)
114  actred = 1. - numext::abs2(fnorm1 / m_fnorm);
115 
116  /* compute the scaled predicted reduction and */
117  /* the scaled directional derivative. */
118  m_wa3 = m_rfactor.template triangularView<Upper>() * (m_permutation.inverse() *m_wa1);
119  temp1 = numext::abs2(m_wa3.stableNorm() / m_fnorm);
121  prered = temp1 + temp2 / Scalar(.5);
122  dirder = -(temp1 + temp2);
123 
124  /* compute the ratio of the actual to the predicted */
125  /* reduction. */
126  ratio = 0.;
127  if (prered != 0.)
128  ratio = actred / prered;
129 
130  /* update the step bound. */
131  if (ratio <= Scalar(.25)) {
132  if (actred >= 0.)
133  temp = RealScalar(.5);
134  if (actred < 0.)
135  temp = RealScalar(.5) * dirder / (dirder + RealScalar(.5) * actred);
136  if (RealScalar(.1) * fnorm1 >= m_fnorm || temp < RealScalar(.1))
137  temp = Scalar(.1);
138  /* Computing MIN */
140  m_par /= temp;
141  } else if (!(m_par != 0. && ratio < RealScalar(.75))) {
142  m_delta = pnorm / RealScalar(.5);
143  m_par = RealScalar(.5) * m_par;
144  }
145 
146  /* test for successful iteration. */
147  if (ratio >= RealScalar(1e-4)) {
148  /* successful iteration. update x, m_fvec, and their norms. */
149  x = m_wa2;
150  m_wa2 = m_diag.cwiseProduct(x);
151  m_fvec = m_wa4;
152  xnorm = m_wa2.stableNorm();
153  m_fnorm = fnorm1;
154  ++m_iter;
155  }
156 
157  /* tests for convergence. */
158  if (abs(actred) <= m_ftol && prered <= m_ftol && Scalar(.5) * ratio <= 1. && m_delta <= m_xtol * xnorm)
159  {
160  m_info = Success;
162  }
163  if (abs(actred) <= m_ftol && prered <= m_ftol && Scalar(.5) * ratio <= 1.)
164  {
165  m_info = Success;
167  }
168  if (m_delta <= m_xtol * xnorm)
169  {
170  m_info = Success;
172  }
173 
174  /* tests for termination and stringent tolerances. */
175  if (m_nfev >= m_maxfev)
176  {
179  }
180  if (abs(actred) <= NumTraits<Scalar>::epsilon() && prered <= NumTraits<Scalar>::epsilon() && Scalar(.5) * ratio <= 1.)
181  {
182  m_info = Success;
184  }
185  if (m_delta <= NumTraits<Scalar>::epsilon() * xnorm)
186  {
187  m_info = Success;
189  }
190  if (m_gnorm <= NumTraits<Scalar>::epsilon())
191  {
192  m_info = Success;
194  }
195 
196  } while (ratio < Scalar(1e-4));
197 
199 }
InverseReturnType inverse() const
IndicesType & indices()
NumericalIssue
void lmpar2(const QRSolver &qr, const VectorType &diag, const VectorType &qtb, typename VectorType::Scalar m_delta, typename VectorType::Scalar &par, VectorType &x)
Definition: LMpar.h:22
bool abs2(bool x)
Eigen::AutoDiffScalar< EIGEN_EXPR_BINARYOP_SCALAR_RETURN_TYPE(Eigen::internal::remove_all_t< DerType >, typename Eigen::internal::traits< Eigen::internal::remove_all_t< DerType >>::Scalar, product) > sqrt(const Eigen::AutoDiffScalar< DerType > &x)
CleanedUpDerType< DerType >::type() min(const AutoDiffScalar< DerType > &x, const T &y)
CleanedUpDerType< DerType >::type() max(const AutoDiffScalar< DerType > &x, const T &y)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
adouble abs(const adouble &x)
Definition: AdolcForward:74

◆ minimizeOneStep() [2/2]

template<typename FunctorType_ >
LevenbergMarquardtSpace::Status Eigen::LevenbergMarquardt< FunctorType_ >::minimizeOneStep ( FVectorType x)

◆ minimizeOptimumStorage()

template<typename FunctorType , typename Scalar >
LevenbergMarquardtSpace::Status Eigen::LevenbergMarquardt< FunctorType, Scalar >::minimizeOptimumStorage ( FVectorType x)

Definition at line 615 of file NonLinearOptimization/LevenbergMarquardt.h.

616 {
619  return status;
620  do {
622  } while (status==LevenbergMarquardtSpace::Running);
623  return status;
624 }
LevenbergMarquardtSpace::Status minimizeOptimumStorageOneStep(FVectorType &x)
LevenbergMarquardtSpace::Status minimizeOptimumStorageInit(FVectorType &x)

◆ minimizeOptimumStorageInit()

template<typename FunctorType , typename Scalar >
LevenbergMarquardtSpace::Status Eigen::LevenbergMarquardt< FunctorType, Scalar >::minimizeOptimumStorageInit ( FVectorType x)

Definition at line 387 of file NonLinearOptimization/LevenbergMarquardt.h.

388 {
389  n = x.size();
390  m = functor.values();
391 
392  wa1.resize(n); wa2.resize(n); wa3.resize(n);
393  wa4.resize(m);
394  fvec.resize(m);
395  // Only R is stored in fjac. Q is only used to compute 'qtf', which is
396  // Q.transpose()*rhs. qtf will be updated using givens rotation,
397  // instead of storing them in Q.
398  // The purpose it to only use a nxn matrix, instead of mxn here, so
399  // that we can handle cases where m>>n :
400  fjac.resize(n, n);
401  if (!useExternalScaling)
402  diag.resize(n);
403  eigen_assert( (!useExternalScaling || diag.size()==n) && "When useExternalScaling is set, the caller must provide a valid 'diag'");
404  qtf.resize(n);
405 
406  /* Function Body */
407  nfev = 0;
408  njev = 0;
409 
410  /* check the input parameters for errors. */
411  if (n <= 0 || m < n || parameters.ftol < 0. || parameters.xtol < 0. || parameters.gtol < 0. || parameters.maxfev <= 0 || parameters.factor <= 0.)
413 
414  if (useExternalScaling)
415  for (Index j = 0; j < n; ++j)
416  if (diag[j] <= 0.)
418 
419  /* evaluate the function at the starting point */
420  /* and calculate its norm. */
421  nfev = 1;
422  if ( functor(x, fvec) < 0)
424  fnorm = fvec.stableNorm();
425 
426  /* initialize levenberg-marquardt parameter and iteration counter. */
427  par = 0.;
428  iter = 1;
429 
431 }

◆ minimizeOptimumStorageOneStep()

template<typename FunctorType , typename Scalar >
LevenbergMarquardtSpace::Status Eigen::LevenbergMarquardt< FunctorType, Scalar >::minimizeOptimumStorageOneStep ( FVectorType x)

Definition at line 436 of file NonLinearOptimization/LevenbergMarquardt.h.

437 {
438  using std::abs;
439  using std::sqrt;
440 
441  eigen_assert(x.size()==n); // check the caller is not cheating us
442 
443  Index i, j;
444  bool sing;
445 
446  /* compute the qr factorization of the jacobian matrix */
447  /* calculated one row at a time, while simultaneously */
448  /* forming (q transpose)*fvec and storing the first */
449  /* n components in qtf. */
450  qtf.fill(0.);
451  fjac.fill(0.);
452  Index rownb = 2;
453  for (i = 0; i < m; ++i) {
454  if (functor.df(x, wa3, rownb) < 0) return LevenbergMarquardtSpace::UserAsked;
455  internal::rwupdt<Scalar>(fjac, wa3, qtf, fvec[i]);
456  ++rownb;
457  }
458  ++njev;
459 
460  /* if the jacobian is rank deficient, call qrfac to */
461  /* reorder its columns and update the components of qtf. */
462  sing = false;
463  for (j = 0; j < n; ++j) {
464  if (fjac(j,j) == 0.)
465  sing = true;
466  wa2[j] = fjac.col(j).head(j).stableNorm();
467  }
469  if (sing) {
470  wa2 = fjac.colwise().blueNorm();
471  // TODO We have no unit test covering this code path, do not modify
472  // until it is carefully tested
473  ColPivHouseholderQR<JacobianType> qrfac(fjac);
474  fjac = qrfac.matrixQR();
475  wa1 = fjac.diagonal();
476  fjac.diagonal() = qrfac.hCoeffs();
477  permutation = qrfac.colsPermutation();
478  // TODO : avoid this:
479  for(Index ii=0; ii< fjac.cols(); ii++) fjac.col(ii).segment(ii+1, fjac.rows()-ii-1) *= fjac(ii,ii); // rescale vectors
480 
481  for (j = 0; j < n; ++j) {
482  if (fjac(j,j) != 0.) {
483  sum = 0.;
484  for (i = j; i < n; ++i)
485  sum += fjac(i,j) * qtf[i];
486  temp = -sum / fjac(j,j);
487  for (i = j; i < n; ++i)
488  qtf[i] += fjac(i,j) * temp;
489  }
490  fjac(j,j) = wa1[j];
491  }
492  }
493 
494  /* on the first iteration and if external scaling is not used, scale according */
495  /* to the norms of the columns of the initial jacobian. */
496  if (iter == 1) {
497  if (!useExternalScaling)
498  for (j = 0; j < n; ++j)
499  diag[j] = (wa2[j]==0.)? 1. : wa2[j];
500 
501  /* on the first iteration, calculate the norm of the scaled x */
502  /* and initialize the step bound delta. */
503  xnorm = diag.cwiseProduct(x).stableNorm();
505  if (delta == 0.)
507  }
508 
509  /* compute the norm of the scaled gradient. */
510  gnorm = 0.;
511  if (fnorm != 0.)
512  for (j = 0; j < n; ++j)
513  if (wa2[permutation.indices()[j]] != 0.)
514  gnorm = (std::max)(gnorm, abs( fjac.col(j).head(j+1).dot(qtf.head(j+1)/fnorm) / wa2[permutation.indices()[j]]));
515 
516  /* test for convergence of the gradient norm. */
517  if (gnorm <= parameters.gtol)
519 
520  /* rescale if necessary. */
521  if (!useExternalScaling)
522  diag = diag.cwiseMax(wa2);
523 
524  do {
525 
526  /* determine the levenberg-marquardt parameter. */
527  internal::lmpar<Scalar>(fjac, permutation.indices(), diag, qtf, delta, par, wa1);
528 
529  /* store the direction p and x + p. calculate the norm of p. */
530  wa1 = -wa1;
531  wa2 = x + wa1;
532  pnorm = diag.cwiseProduct(wa1).stableNorm();
533 
534  /* on the first iteration, adjust the initial step bound. */
535  if (iter == 1)
536  delta = (std::min)(delta,pnorm);
537 
538  /* evaluate the function at x + p and calculate its norm. */
539  if ( functor(wa2, wa4) < 0)
541  ++nfev;
542  fnorm1 = wa4.stableNorm();
543 
544  /* compute the scaled actual reduction. */
545  actred = -1.;
546  if (Scalar(.1) * fnorm1 < fnorm)
547  actred = 1. - numext::abs2(fnorm1 / fnorm);
548 
549  /* compute the scaled predicted reduction and */
550  /* the scaled directional derivative. */
551  wa3 = fjac.topLeftCorner(n,n).template triangularView<Upper>() * (permutation.inverse() * wa1);
552  temp1 = numext::abs2(wa3.stableNorm() / fnorm);
554  prered = temp1 + temp2 / Scalar(.5);
555  dirder = -(temp1 + temp2);
556 
557  /* compute the ratio of the actual to the predicted */
558  /* reduction. */
559  ratio = 0.;
560  if (prered != 0.)
561  ratio = actred / prered;
562 
563  /* update the step bound. */
564  if (ratio <= Scalar(.25)) {
565  if (actred >= 0.)
566  temp = Scalar(.5);
567  if (actred < 0.)
568  temp = Scalar(.5) * dirder / (dirder + Scalar(.5) * actred);
569  if (Scalar(.1) * fnorm1 >= fnorm || temp < Scalar(.1))
570  temp = Scalar(.1);
571  /* Computing MIN */
572  delta = temp * (std::min)(delta, pnorm / Scalar(.1));
573  par /= temp;
574  } else if (!(par != 0. && ratio < Scalar(.75))) {
575  delta = pnorm / Scalar(.5);
576  par = Scalar(.5) * par;
577  }
578 
579  /* test for successful iteration. */
580  if (ratio >= Scalar(1e-4)) {
581  /* successful iteration. update x, fvec, and their norms. */
582  x = wa2;
583  wa2 = diag.cwiseProduct(x);
584  fvec = wa4;
585  xnorm = wa2.stableNorm();
586  fnorm = fnorm1;
587  ++iter;
588  }
589 
590  /* tests for convergence. */
591  if (abs(actred) <= parameters.ftol && prered <= parameters.ftol && Scalar(.5) * ratio <= 1. && delta <= parameters.xtol * xnorm)
593  if (abs(actred) <= parameters.ftol && prered <= parameters.ftol && Scalar(.5) * ratio <= 1.)
595  if (delta <= parameters.xtol * xnorm)
597 
598  /* tests for termination and stringent tolerances. */
599  if (nfev >= parameters.maxfev)
601  if (abs(actred) <= NumTraits<Scalar>::epsilon() && prered <= NumTraits<Scalar>::epsilon() && Scalar(.5) * ratio <= 1.)
603  if (delta <= NumTraits<Scalar>::epsilon() * xnorm)
605  if (gnorm <= NumTraits<Scalar>::epsilon())
607 
608  } while (ratio < Scalar(1e-4));
609 
611 }
int i
PermutationMatrix< Dynamic, Dynamic > permutation

◆ nfev()

template<typename FunctorType_ >
Index Eigen::LevenbergMarquardt< FunctorType_ >::nfev ( )
inline
Returns
the number of functions evaluation

Definition at line 205 of file LevenbergMarquardt/LevenbergMarquardt.h.

205 { return m_nfev; }

◆ njev()

template<typename FunctorType_ >
Index Eigen::LevenbergMarquardt< FunctorType_ >::njev ( )
inline
Returns
the number of jacobian evaluation

Definition at line 208 of file LevenbergMarquardt/LevenbergMarquardt.h.

208 { return m_njev; }

◆ operator=()

template<typename FunctorType_ >
LevenbergMarquardt& Eigen::LevenbergMarquardt< FunctorType_ >::operator= ( const LevenbergMarquardt< FunctorType_ > &  )
private

◆ permutation()

template<typename FunctorType_ >
PermutationType Eigen::LevenbergMarquardt< FunctorType_ >::permutation ( )
inline

the permutation used in the QR factorization

Definition at line 234 of file LevenbergMarquardt/LevenbergMarquardt.h.

234 {return m_permutation; }

◆ resetParameters() [1/2]

template<typename FunctorType_ >
void Eigen::LevenbergMarquardt< FunctorType_ >::resetParameters ( )
inline

Sets the default parameters

Definition at line 147 of file LevenbergMarquardt/LevenbergMarquardt.h.

148  {
149  using std::sqrt;
150 
151  m_factor = 100.;
152  m_maxfev = 400;
153  m_ftol = sqrt(NumTraits<RealScalar>::epsilon());
154  m_xtol = sqrt(NumTraits<RealScalar>::epsilon());
155  m_gtol = 0. ;
156  m_epsfcn = 0. ;
157  }

◆ resetParameters() [2/2]

template<typename FunctorType_ >
void Eigen::LevenbergMarquardt< FunctorType_ >::resetParameters ( void  )
inline

Definition at line 106 of file NonLinearOptimization/LevenbergMarquardt.h.

106 { parameters = Parameters(); }

◆ setEpsilon()

template<typename FunctorType_ >
void Eigen::LevenbergMarquardt< FunctorType_ >::setEpsilon ( RealScalar  epsfcn)
inline

Sets the error precision

Definition at line 172 of file LevenbergMarquardt/LevenbergMarquardt.h.

172 { m_epsfcn = epsfcn; }

◆ setExternalScaling()

template<typename FunctorType_ >
void Eigen::LevenbergMarquardt< FunctorType_ >::setExternalScaling ( bool  value)
inline

Use an external Scaling. If set to true, pass a nonzero diagonal to diag()

Definition at line 178 of file LevenbergMarquardt/LevenbergMarquardt.h.

178 {m_useExternalScaling = value; }

◆ setFactor()

template<typename FunctorType_ >
void Eigen::LevenbergMarquardt< FunctorType_ >::setFactor ( RealScalar  factor)
inline

Sets the step bound for the diagonal shift

Definition at line 169 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ setFtol()

template<typename FunctorType_ >
void Eigen::LevenbergMarquardt< FunctorType_ >::setFtol ( RealScalar  ftol)
inline

Sets the tolerance for the norm of the vector function

Definition at line 163 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ setGtol()

template<typename FunctorType_ >
void Eigen::LevenbergMarquardt< FunctorType_ >::setGtol ( RealScalar  gtol)
inline

Sets the tolerance for the norm of the gradient of the error vector

Definition at line 166 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ setMaxfev()

template<typename FunctorType_ >
void Eigen::LevenbergMarquardt< FunctorType_ >::setMaxfev ( Index  maxfev)
inline

Sets the maximum number of function evaluation

Definition at line 175 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ setXtol()

template<typename FunctorType_ >
void Eigen::LevenbergMarquardt< FunctorType_ >::setXtol ( RealScalar  xtol)
inline

Sets the tolerance for the norm of the solution vector

Definition at line 160 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ sqrt_epsilon()

template<typename FunctorType_ >
static Scalar Eigen::LevenbergMarquardt< FunctorType_ >::sqrt_epsilon ( )
inlinestaticprivate

Definition at line 50 of file NonLinearOptimization/LevenbergMarquardt.h.

51  {
52  using std::sqrt;
53  return sqrt(NumTraits<Scalar>::epsilon());
54  }

◆ xtol()

template<typename FunctorType_ >
RealScalar Eigen::LevenbergMarquardt< FunctorType_ >::xtol ( ) const
inline
Returns
the tolerance for the norm of the solution vector

Definition at line 181 of file LevenbergMarquardt/LevenbergMarquardt.h.

181 {return m_xtol; }

Member Data Documentation

◆ actred

template<typename FunctorType_ >
Scalar Eigen::LevenbergMarquardt< FunctorType_ >::actred
private

Definition at line 130 of file NonLinearOptimization/LevenbergMarquardt.h.

◆ delta

template<typename FunctorType_ >
Scalar Eigen::LevenbergMarquardt< FunctorType_ >::delta
private

Definition at line 128 of file NonLinearOptimization/LevenbergMarquardt.h.

◆ diag

template<typename FunctorType_ >
FVectorType Eigen::LevenbergMarquardt< FunctorType_ >::diag

Definition at line 109 of file NonLinearOptimization/LevenbergMarquardt.h.

◆ dirder

template<typename FunctorType_ >
Scalar Eigen::LevenbergMarquardt< FunctorType_ >::dirder
private

Definition at line 130 of file NonLinearOptimization/LevenbergMarquardt.h.

◆ fjac

template<typename FunctorType_ >
JacobianType Eigen::LevenbergMarquardt< FunctorType_ >::fjac

Definition at line 110 of file NonLinearOptimization/LevenbergMarquardt.h.

◆ fnorm

template<typename FunctorType_ >
Scalar Eigen::LevenbergMarquardt< FunctorType_ >::fnorm

Definition at line 115 of file NonLinearOptimization/LevenbergMarquardt.h.

◆ fnorm1

template<typename FunctorType_ >
Scalar Eigen::LevenbergMarquardt< FunctorType_ >::fnorm1
private

Definition at line 130 of file NonLinearOptimization/LevenbergMarquardt.h.

◆ functor

template<typename FunctorType_ >
FunctorType& Eigen::LevenbergMarquardt< FunctorType_ >::functor
private

Definition at line 121 of file NonLinearOptimization/LevenbergMarquardt.h.

◆ fvec

template<typename FunctorType_ >
FVectorType Eigen::LevenbergMarquardt< FunctorType_ >::fvec

Definition at line 109 of file NonLinearOptimization/LevenbergMarquardt.h.

◆ gnorm

template<typename FunctorType_ >
Scalar Eigen::LevenbergMarquardt< FunctorType_ >::gnorm

Definition at line 115 of file NonLinearOptimization/LevenbergMarquardt.h.

◆ iter

template<typename FunctorType_ >
Index Eigen::LevenbergMarquardt< FunctorType_ >::iter

Definition at line 114 of file NonLinearOptimization/LevenbergMarquardt.h.

◆ m

template<typename FunctorType_ >
Index Eigen::LevenbergMarquardt< FunctorType_ >::m
private

Definition at line 256 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ m_delta

template<typename FunctorType_ >
RealScalar Eigen::LevenbergMarquardt< FunctorType_ >::m_delta
private

Definition at line 268 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ m_diag

template<typename FunctorType_ >
FVectorType Eigen::LevenbergMarquardt< FunctorType_ >::m_diag
private

Definition at line 254 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ m_epsfcn

template<typename FunctorType_ >
RealScalar Eigen::LevenbergMarquardt< FunctorType_ >::m_epsfcn
private

Definition at line 266 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ m_factor

template<typename FunctorType_ >
RealScalar Eigen::LevenbergMarquardt< FunctorType_ >::m_factor
private

Definition at line 261 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ m_fjac

template<typename FunctorType_ >
JacobianType Eigen::LevenbergMarquardt< FunctorType_ >::m_fjac
private

Definition at line 251 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ m_fnorm

template<typename FunctorType_ >
RealScalar Eigen::LevenbergMarquardt< FunctorType_ >::m_fnorm
private

Definition at line 259 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ m_ftol

template<typename FunctorType_ >
RealScalar Eigen::LevenbergMarquardt< FunctorType_ >::m_ftol
private

Definition at line 263 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ m_functor

template<typename FunctorType_ >
FunctorType& Eigen::LevenbergMarquardt< FunctorType_ >::m_functor
private

Definition at line 253 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ m_fvec

template<typename FunctorType_ >
FVectorType Eigen::LevenbergMarquardt< FunctorType_ >::m_fvec
private

Definition at line 254 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ m_gnorm

template<typename FunctorType_ >
RealScalar Eigen::LevenbergMarquardt< FunctorType_ >::m_gnorm
private

Definition at line 260 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ m_gtol

template<typename FunctorType_ >
RealScalar Eigen::LevenbergMarquardt< FunctorType_ >::m_gtol
private

Definition at line 265 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ m_info

template<typename FunctorType_ >
ComputationInfo Eigen::LevenbergMarquardt< FunctorType_ >::m_info
private

Definition at line 274 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ m_isInitialized

template<typename FunctorType_ >
bool Eigen::LevenbergMarquardt< FunctorType_ >::m_isInitialized
private

Definition at line 273 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ m_iter

template<typename FunctorType_ >
Index Eigen::LevenbergMarquardt< FunctorType_ >::m_iter
private

Definition at line 267 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ m_maxfev

template<typename FunctorType_ >
Index Eigen::LevenbergMarquardt< FunctorType_ >::m_maxfev
private

Definition at line 262 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ m_nfev

template<typename FunctorType_ >
Index Eigen::LevenbergMarquardt< FunctorType_ >::m_nfev
private

Definition at line 257 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ m_njev

template<typename FunctorType_ >
Index Eigen::LevenbergMarquardt< FunctorType_ >::m_njev
private

Definition at line 258 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ m_par

template<typename FunctorType_ >
RealScalar Eigen::LevenbergMarquardt< FunctorType_ >::m_par
private

Definition at line 272 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ m_permutation

template<typename FunctorType_ >
PermutationType Eigen::LevenbergMarquardt< FunctorType_ >::m_permutation
private

Definition at line 270 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ m_qtf

template<typename FunctorType_ >
FVectorType Eigen::LevenbergMarquardt< FunctorType_ >::m_qtf
private

Definition at line 254 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ m_rfactor

template<typename FunctorType_ >
JacobianType Eigen::LevenbergMarquardt< FunctorType_ >::m_rfactor
private

Definition at line 252 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ m_useExternalScaling

template<typename FunctorType_ >
bool Eigen::LevenbergMarquardt< FunctorType_ >::m_useExternalScaling
private

Definition at line 269 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ m_wa1

template<typename FunctorType_ >
FVectorType Eigen::LevenbergMarquardt< FunctorType_ >::m_wa1
private

Definition at line 271 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ m_wa2

template<typename FunctorType_ >
FVectorType Eigen::LevenbergMarquardt< FunctorType_ >::m_wa2
private

Definition at line 271 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ m_wa3

template<typename FunctorType_ >
FVectorType Eigen::LevenbergMarquardt< FunctorType_ >::m_wa3
private

Definition at line 271 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ m_wa4

template<typename FunctorType_ >
FVectorType Eigen::LevenbergMarquardt< FunctorType_ >::m_wa4
private

Definition at line 271 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ m_xtol

template<typename FunctorType_ >
RealScalar Eigen::LevenbergMarquardt< FunctorType_ >::m_xtol
private

Definition at line 264 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ n

template<typename FunctorType_ >
Index Eigen::LevenbergMarquardt< FunctorType_ >::n
private

Definition at line 255 of file LevenbergMarquardt/LevenbergMarquardt.h.

◆ nfev

template<typename FunctorType_ >
Index Eigen::LevenbergMarquardt< FunctorType_ >::nfev

Definition at line 112 of file NonLinearOptimization/LevenbergMarquardt.h.

◆ njev

template<typename FunctorType_ >
Index Eigen::LevenbergMarquardt< FunctorType_ >::njev

Definition at line 113 of file NonLinearOptimization/LevenbergMarquardt.h.

◆ par

template<typename FunctorType_ >
Scalar Eigen::LevenbergMarquardt< FunctorType_ >::par
private

Definition at line 126 of file NonLinearOptimization/LevenbergMarquardt.h.

◆ parameters

template<typename FunctorType_ >
Parameters Eigen::LevenbergMarquardt< FunctorType_ >::parameters

Definition at line 108 of file NonLinearOptimization/LevenbergMarquardt.h.

◆ permutation

template<typename FunctorType_ >
PermutationMatrix<Dynamic,Dynamic> Eigen::LevenbergMarquardt< FunctorType_ >::permutation

Definition at line 111 of file NonLinearOptimization/LevenbergMarquardt.h.

◆ pnorm

template<typename FunctorType_ >
Scalar Eigen::LevenbergMarquardt< FunctorType_ >::pnorm
private

Definition at line 130 of file NonLinearOptimization/LevenbergMarquardt.h.

◆ prered

template<typename FunctorType_ >
Scalar Eigen::LevenbergMarquardt< FunctorType_ >::prered
private

Definition at line 130 of file NonLinearOptimization/LevenbergMarquardt.h.

◆ qtf

template<typename FunctorType_ >
FVectorType Eigen::LevenbergMarquardt< FunctorType_ >::qtf

Definition at line 109 of file NonLinearOptimization/LevenbergMarquardt.h.

◆ ratio

template<typename FunctorType_ >
Scalar Eigen::LevenbergMarquardt< FunctorType_ >::ratio
private

Definition at line 129 of file NonLinearOptimization/LevenbergMarquardt.h.

◆ sum

template<typename FunctorType_ >
Scalar Eigen::LevenbergMarquardt< FunctorType_ >::sum
private

Definition at line 126 of file NonLinearOptimization/LevenbergMarquardt.h.

◆ temp

template<typename FunctorType_ >
Scalar Eigen::LevenbergMarquardt< FunctorType_ >::temp
private

Definition at line 127 of file NonLinearOptimization/LevenbergMarquardt.h.

◆ temp1

template<typename FunctorType_ >
Scalar Eigen::LevenbergMarquardt< FunctorType_ >::temp1
private

Definition at line 127 of file NonLinearOptimization/LevenbergMarquardt.h.

◆ temp2

template<typename FunctorType_ >
Scalar Eigen::LevenbergMarquardt< FunctorType_ >::temp2
private

Definition at line 127 of file NonLinearOptimization/LevenbergMarquardt.h.

◆ useExternalScaling

template<typename FunctorType_ >
bool Eigen::LevenbergMarquardt< FunctorType_ >::useExternalScaling

Definition at line 116 of file NonLinearOptimization/LevenbergMarquardt.h.

◆ wa1

template<typename FunctorType_ >
FVectorType Eigen::LevenbergMarquardt< FunctorType_ >::wa1
private

Definition at line 124 of file NonLinearOptimization/LevenbergMarquardt.h.

◆ wa2

template<typename FunctorType_ >
FVectorType Eigen::LevenbergMarquardt< FunctorType_ >::wa2
private

Definition at line 124 of file NonLinearOptimization/LevenbergMarquardt.h.

◆ wa3

template<typename FunctorType_ >
FVectorType Eigen::LevenbergMarquardt< FunctorType_ >::wa3
private

Definition at line 124 of file NonLinearOptimization/LevenbergMarquardt.h.

◆ wa4

template<typename FunctorType_ >
FVectorType Eigen::LevenbergMarquardt< FunctorType_ >::wa4
private

Definition at line 124 of file NonLinearOptimization/LevenbergMarquardt.h.

◆ xnorm

template<typename FunctorType_ >
Scalar Eigen::LevenbergMarquardt< FunctorType_ >::xnorm
private

Definition at line 130 of file NonLinearOptimization/LevenbergMarquardt.h.


The documentation for this class was generated from the following files: