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 no_assignment_operator.

Public Member Functions

FVectorType & diag ()
 
RealScalar fnorm ()
 
FVectorType & fvec ()
 
RealScalar gnorm ()
 
ComputationInfo info () const
 Reports whether the minimization was successful. More...
 
Index iterations ()
 
JacobianType & jacobian ()
 
RealScalar lm_param (void)
 
JacobianType & matrixR ()
 
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)
 

Member Function Documentation

§ diag()

FVectorType& diag ( )
inline
Returns
a reference to the diagonal of the jacobian

§ fnorm()

RealScalar fnorm ( )
inline
Returns
the norm of current vector function

§ fvec()

FVectorType& fvec ( )
inline
Returns
a reference to the current vector function

§ gnorm()

RealScalar gnorm ( )
inline
Returns
the norm of the gradient of the error

§ info()

ComputationInfo info ( ) const
inline

Reports whether the minimization was successful.

Returns
Success if the minimization was succesful, NumericalIssue if a numerical problem arises during the minimization process, for exemple 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

References LevenbergMarquardt< _FunctorType >::setFtol(), LevenbergMarquardt< _FunctorType >::setMaxfev(), and LevenbergMarquardt< _FunctorType >::setXtol().

§ iterations()

Index iterations ( )
inline
Returns
the number of iterations performed

§ jacobian()

JacobianType& jacobian ( )
inline
Returns
a reference to the matrix where the current Jacobian matrix is stored

§ lm_param()

RealScalar lm_param ( void  )
inline
Returns
the LevenbergMarquardt parameter

§ matrixR()

JacobianType& matrixR ( )
inline
Returns
a reference to the triangular matrix R from the QR of the jacobian matrix.
See also
jacobian()

§ nfev()

Index nfev ( )
inline
Returns
the number of functions evaluation

§ njev()

Index njev ( )
inline
Returns
the number of jacobian evaluation

§ permutation()

PermutationType permutation ( )
inline

the permutation used in the QR factorization

§ resetParameters()

void resetParameters ( )
inline

Sets the default parameters

§ setEpsilon()

void setEpsilon ( RealScalar  epsfcn)
inline

Sets the error precision

§ setExternalScaling()

void setExternalScaling ( bool  value)
inline

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

§ setFactor()

void setFactor ( RealScalar  factor)
inline

Sets the step bound for the diagonal shift

§ setFtol()

void setFtol ( RealScalar  ftol)
inline

Sets the tolerance for the norm of the vector function

Referenced by LevenbergMarquardt< _FunctorType >::info().

§ setGtol()

void setGtol ( RealScalar  gtol)
inline

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

§ setMaxfev()

void setMaxfev ( Index  maxfev)
inline

Sets the maximum number of function evaluation

Referenced by LevenbergMarquardt< _FunctorType >::info().

§ setXtol()

void setXtol ( RealScalar  xtol)
inline

Sets the tolerance for the norm of the solution vector

Referenced by LevenbergMarquardt< _FunctorType >::info().


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