compbio
Classes | Public Types | Public Member Functions | Static Public Member Functions | Public Attributes | List of all members
Eigen::LevenbergMarquardt< _FunctorType > Class Template Reference

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

#include <LevenbergMarquardt.h>

Inheritance diagram for Eigen::LevenbergMarquardt< _FunctorType >:
Eigen::internal::no_assignment_operator

Classes

struct  Parameters
 

Public Types

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

Public Member Functions

 LevenbergMarquardt (FunctorType &functor)
 
LevenbergMarquardtSpace::Status minimize (FVectorType &x)
 
LevenbergMarquardtSpace::Status minimizeInit (FVectorType &x)
 
LevenbergMarquardtSpace::Status minimizeOneStep (FVectorType &x)
 
LevenbergMarquardtSpace::Status lmder1 (FVectorType &x, const Scalar tol=std::sqrt(NumTraits< Scalar >::epsilon()))
 
void resetParameters ()
 Sets the default parameters.
 
void setXtol (RealScalar xtol)
 Sets the tolerance for the norm of the solution vector.
 
void setFtol (RealScalar ftol)
 Sets the tolerance for the norm of the vector function.
 
void setGtol (RealScalar gtol)
 Sets the tolerance for the norm of the gradient of the error vector.
 
void setFactor (RealScalar factor)
 Sets the step bound for the diagonal shift.
 
void setEpsilon (RealScalar epsfcn)
 Sets the error precision.
 
void setMaxfev (Index maxfev)
 Sets the maximum number of function evaluation.
 
void setExternalScaling (bool value)
 Use an external Scaling. More...
 
RealScalar xtol () const
 
RealScalar ftol () const
 
RealScalar gtol () const
 
RealScalar factor () const
 
RealScalar epsilon () const
 
Index maxfev () const
 
FVectorTypediag ()
 
Index iterations ()
 
Index nfev ()
 
Index njev ()
 
RealScalar fnorm ()
 
RealScalar gnorm ()
 
RealScalar lm_param (void)
 
FVectorTypefvec ()
 
JacobianType & jacobian ()
 
JacobianType & matrixR ()
 
PermutationType permutation ()
 the permutation used in the QR factorization
 
ComputationInfo info () const
 Reports whether the minimization was successful. More...
 
 LevenbergMarquardt (FunctorType &_functor)
 
LevenbergMarquardtSpace::Status lmder1 (FVectorType &x, const Scalar tol=sqrt_epsilon())
 
LevenbergMarquardtSpace::Status minimize (FVectorType &x)
 
LevenbergMarquardtSpace::Status minimizeInit (FVectorType &x)
 
LevenbergMarquardtSpace::Status minimizeOneStep (FVectorType &x)
 
LevenbergMarquardtSpace::Status lmstr1 (FVectorType &x, const Scalar tol=sqrt_epsilon())
 
LevenbergMarquardtSpace::Status minimizeOptimumStorage (FVectorType &x)
 
LevenbergMarquardtSpace::Status minimizeOptimumStorageInit (FVectorType &x)
 
LevenbergMarquardtSpace::Status minimizeOptimumStorageOneStep (FVectorType &x)
 
void resetParameters (void)
 
Scalar lm_param (void)
 

Static Public Member Functions

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

Public Attributes

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

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

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

§ 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

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

§ 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: