Eigen-unsupported  3.4.90 (git rev 67eeba6e720c5745abc77ae6c92ce0a44aa7b7ae)
Eigen::LevenbergMarquardt< FunctorType_ > Class Template Reference

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

Inherits internal::no_assignment_operator.

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 ()
 
JacobianType & jacobian ()
 
RealScalar lm_param (void)
 
JacobianType & matrixR ()
 
Index maxfev () const
 
Index nfev ()
 
Index njev ()
 
PermutationType permutation ()
 
void resetParameters ()
 
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
 

Member Function Documentation

◆ diag()

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

◆ epsilon()

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

◆ factor()

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

◆ fnorm()

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

◆ ftol()

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

◆ fvec()

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

◆ gnorm()

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

◆ gtol()

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

◆ 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

◆ iterations()

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

◆ jacobian()

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

◆ lm_param()

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

◆ 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()

◆ maxfev()

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

◆ nfev()

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

◆ njev()

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

◆ permutation()

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

the permutation used in the QR factorization

◆ resetParameters()

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

Sets the default parameters

◆ setEpsilon()

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

Sets the error precision

◆ 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()

◆ setFactor()

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

Sets the step bound for the diagonal shift

◆ setFtol()

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

Sets the tolerance for the norm of the vector function

◆ 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

◆ setMaxfev()

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

Sets the maximum number of function evaluation

◆ setXtol()

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

Sets the tolerance for the norm of the solution vector

◆ xtol()

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

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