ITK  4.13.0
Insight Segmentation and Registration Toolkit
Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | Static Protected Member Functions | Private Attributes | List of all members
itk::LBFGS2Optimizerv4 Class Reference

#include <itkLBFGS2Optimizerv4.h>

+ Inheritance diagram for itk::LBFGS2Optimizerv4:
+ Collaboration diagram for itk::LBFGS2Optimizerv4:

Detailed Description

Wrap of the libLBFGS[1] algorithm for use in ITKv4 registration framework. LibLBFGS is a translation of LBFGS code by Nocedal [2] and adds the orthantwise limited-memmory Quais-Newton method [3] for optimization with L1-norm on the parameters.

LBFGS is a quasi-Newton method uses an approximate estimate of the inverse Hessian $ (\nabla^2 f(x) )^-1 $ to scale the gradient step:

\[ x_{n+1} = x_n - s (\nabla^2 f(x_n) )^-1 \nabla f(x) \]

with $ s $ the step size.

The inverse Hessian is approximated from the gradients of previous iteration and thus only the gradient of the objective function is required.

The step size $ s $ is determined through line search which defaults to the approach by More and Thuente [4]. This line search approach finds a step size such that

\[ \lVert \nabla f(x + s (\nabla^2 f(x_n) )^{-1} \nabla f(x) ) \rVert \le \nu \lVert \nabla f(x) \rVert \]

The parameter $\nu$ is set through SetLineSearchAccuracy() (default 0.9) and SetGradientLineSearchAccuracy()

Instead of the More-Tunete method, backtracking with three different conditions [7] are availabe and can be set through SetLineSearch():

The optimization stops when either the gradient satisfies the condition

\[ \lVert \nabla f(x) \rVert \le \epsilon \max(1, \lVert X \rVert) \]

or a maximum number of function evaluations has been reached. The tolerance $\epsilon$ is set through SetSolutionAccuracy() (default 1e-5) and the maximum number of function evaluations is set through SetMaximumIterations() (default 0 = no maximum).

References:

[1] libLBFGS

[2] NETLIB lbfgs

[3] Galen Andrew and Jianfeng Gao. Scalable training of L1-regularized log-linear models. 24th International Conference on Machine Learning, pp. 33-40, 2007.

[4] Jorge Nocedal. Updating Quasi-Newton Matrices with Limited Storage. Mathematics of Computation, Vol. 35, No. 151, pp. 773–782, 1980.

[5] Dong C. Liu and Jorge Nocedal. On the limited memory BFGS method for large scale optimization. Mathematical Programming B, Vol. 45, No. 3, pp. 503-528, 1989.

[6] More, J. J. and D. J. Thuente. Line Search Algorithms with Guaranteed Sufficient Decrease. ACM Transactions on Mathematical Software 20, no. 3 (1994): 286–307.

[7] John E. Dennis and Robert B. Schnabel. Numerical Methods for Unconstrained Optimization and Nonlinear Equations, Englewood Cliffs, 1983.

Definition at line 102 of file itkLBFGS2Optimizerv4.h.

Public Types

typedef SmartPointer< const SelfConstPointer
 
enum  LineSearchMethod {
  LINESEARCH_DEFAULT = 0,
  LINESEARCH_MORETHUENTE = 0,
  LINESEARCH_BACKTRACKING_ARMIJO = 1,
  LINESEARCH_BACKTRACKING = 2,
  LINESEARCH_BACKTRACKING_WOLFE = 2,
  LINESEARCH_BACKTRACKING_STRONG_WOLFE = 3
}
 
typedef Superclass::MetricType MetricType
 
typedef Superclass::ParametersType ParametersType
 
typedef SmartPointer< SelfPointer
 
typedef double PrecisionType
 
typedef Superclass::ScalesType ScalesType
 
typedef LBFGS2Optimizerv4 Self
 
typedef
ObjectToObjectOptimizerBaseTemplate
< double > 
Superclass
 
- Public Types inherited from itk::ObjectToObjectOptimizerBaseTemplate< double >
typedef SmartPointer< const SelfConstPointer
 
typedef MetricType::DerivativeType DerivativeType
 
typedef MetricType::MeasureType MeasureType
 
typedef
ObjectToObjectMetricBaseTemplate
< double > 
MetricType
 
typedef MetricType::Pointer MetricTypePointer
 
typedef
MetricType::NumberOfParametersType 
NumberOfParametersType
 
typedef OptimizerParameters
< double > 
ParametersType
 
typedef SmartPointer< SelfPointer
 
typedef
OptimizerParameterScalesEstimatorTemplate
< double > 
ScalesEstimatorType
 
typedef OptimizerParameters
< double > 
ScalesType
 
typedef
ObjectToObjectOptimizerBaseTemplate 
Self
 
typedef std::ostringstream StopConditionDescriptionType
 
typedef std::string StopConditionReturnStringType
 
typedef Object Superclass
 
- Public Types inherited from itk::Object
typedef SmartPointer< const SelfConstPointer
 
typedef SmartPointer< SelfPointer
 
typedef Object Self
 
typedef LightObject Superclass
 
- Public Types inherited from itk::LightObject
typedef SmartPointer< const SelfConstPointer
 
typedef SmartPointer< SelfPointer
 
typedef LightObject Self
 

Public Member Functions

virtual ::itk::LightObject::Pointer CreateAnother () const
 
virtual PrecisionType GetCurrentGradientNorm () const
 
virtual PrecisionType GetCurrentNumberOfEvaluations () const
 
virtual PrecisionType GetCurrentParameterNorm () const
 
virtual PrecisionType GetCurrentStepSize () const
 
virtual const char * GetNameOfClass () const
 
virtual const
StopConditionReturnStringType 
GetStopConditionDescription () const override
 
virtual void SetScales (const ScalesType &) override
 
virtual void SetWeights (const ScalesType) override
 
virtual void StartOptimization (bool doOnlyInitialization=false) override
 
void SetHessianApproximationAccuracy (int m)
 
int GetHessianApproximationAccuracy () const
 
void SetSolutionAccuracy (PrecisionType epsilon)
 
PrecisionType GetSolutionAccuracy () const
 
void SetDeltaConvergenceDistance (int nPast)
 
int GetDeltaConvergenceDistance () const
 
void SetDeltaConvergenceTolerance (PrecisionType tol)
 
PrecisionType GetDeltaConvergenceTolerance () const
 
void SetMaximumIterations (int maxIterations)
 
int GetMaximumIterations () const
 
virtual SizeValueType GetNumberOfIterations () const override
 
virtual void SetNumberOfIterations (const SizeValueType _arg) override
 
void SetLineSearch (const LineSearchMethod &linesearch)
 
LineSearchMethod GetLineSearch () const
 
void SetMaximumLineSearchEvaluations (int n)
 
int GetMaximumLineSearchEvaluations () const
 
void SetMinimumLineSearchStep (PrecisionType step)
 
PrecisionType GetMinimumLineSearchStep () const
 
void SetMaximumLineSearchStep (PrecisionType step)
 
PrecisionType GetMaximumLineSearchStep () const
 
void SetLineSearchAccuracy (PrecisionType ftol)
 
PrecisionType GetLineSearchAccuracy () const
 
void SetWolfeCoefficient (PrecisionType wc)
 
PrecisionType GetWolfeCoefficient () const
 
void SetLineSearchGradientAccuracy (PrecisionType gtol)
 
PrecisionType GetLineSearchGradientAccuracy () const
 
void SetMachinePrecisionTolerance (PrecisionType xtol)
 
PrecisionType GetMachinePrecisionTolerance () const
 
void SetOrthantwiseCoefficient (PrecisionType orthant_c)
 
PrecisionType GetOrthantwiseCoefficient () const
 
void SetOrthantwiseStart (int start)
 
int GetOrthantwiseStart () const
 
void SetOrthantwiseEnd (int end)
 
int GetOrthantwiseEnd () const
 
- Public Member Functions inherited from itk::ObjectToObjectOptimizerBaseTemplate< double >
virtual SizeValueType GetCurrentIteration () const
 
virtual const MeasureTypeGetCurrentMetricValue () const
 
virtual const ParametersTypeGetCurrentPosition () const
 
virtual const ThreadIdTypeGetNumberOfThreads () const
 
virtual const ScalesTypeGetScales () const
 
virtual const bool & GetScalesAreIdentity () const
 
bool GetScalesInitialized () const
 
virtual const MeasureTypeGetValue () const
 
virtual const ScalesTypeGetWeights () const
 
virtual const bool & GetWeightsAreIdentity () const
 
virtual void SetNumberOfThreads (ThreadIdType number)
 
virtual void SetScalesEstimator (ScalesEstimatorType *_arg)
 
virtual void SetWeights (ScalesType _arg)
 
virtual void SetMetric (MetricType *_arg)
 
virtual MetricTypeGetModifiableMetric ()
 
virtual const MetricTypeGetMetric () const
 
virtual void SetScales (const ScalesType &scales)
 
virtual void SetDoEstimateScales (bool _arg)
 
virtual const bool & GetDoEstimateScales () const
 
virtual void DoEstimateScalesOn ()
 
virtual void DoEstimateScalesOff ()
 
- Public Member Functions inherited from itk::Object
unsigned long AddObserver (const EventObject &event, Command *)
 
unsigned long AddObserver (const EventObject &event, Command *) const
 
virtual void DebugOff () const
 
virtual void DebugOn () const
 
CommandGetCommand (unsigned long tag)
 
bool GetDebug () const
 
MetaDataDictionaryGetMetaDataDictionary ()
 
const MetaDataDictionaryGetMetaDataDictionary () const
 
virtual ModifiedTimeType GetMTime () const
 
virtual const TimeStampGetTimeStamp () const
 
bool HasObserver (const EventObject &event) const
 
void InvokeEvent (const EventObject &)
 
void InvokeEvent (const EventObject &) const
 
virtual void Modified () const
 
virtual void Register () const override
 
void RemoveAllObservers ()
 
void RemoveObserver (unsigned long tag)
 
void SetDebug (bool debugFlag) const
 
void SetMetaDataDictionary (const MetaDataDictionary &rhs)
 
virtual void SetReferenceCount (int) override
 
virtual void UnRegister () const noexceptoverride
 
virtual void SetObjectName (std::string _arg)
 
virtual const std::string & GetObjectName () const
 
- Public Member Functions inherited from itk::LightObject
virtual void Delete ()
 
virtual int GetReferenceCount () const
 
 itkCloneMacro (Self)
 
void Print (std::ostream &os, Indent indent=0) const
 

Static Public Member Functions

static Pointer New ()
 
- Static Public Member Functions inherited from itk::Object
static bool GetGlobalWarningDisplay ()
 
static void GlobalWarningDisplayOff ()
 
static void GlobalWarningDisplayOn ()
 
static Pointer New ()
 
static void SetGlobalWarningDisplay (bool flag)
 
- Static Public Member Functions inherited from itk::LightObject
static void BreakOnError ()
 
static Pointer New ()
 

Protected Member Functions

PrecisionType EvaluateCost (const PrecisionType *x, PrecisionType *g, const int n, const PrecisionType step)
 
 LBFGS2Optimizerv4 ()
 
virtual void PrintSelf (std::ostream &os, Indent indent) const override
 
int UpdateProgress (const PrecisionType *x, const PrecisionType *g, const PrecisionType fx, const PrecisionType xnorm, const PrecisionType gnorm, const PrecisionType step, int n, int k, int ls)
 
virtual ~LBFGS2Optimizerv4 () override
 
- Protected Member Functions inherited from itk::ObjectToObjectOptimizerBaseTemplate< double >
 ObjectToObjectOptimizerBaseTemplate ()
 
virtual ~ObjectToObjectOptimizerBaseTemplate () override
 
- Protected Member Functions inherited from itk::Object
 Object ()
 
bool PrintObservers (std::ostream &os, Indent indent) const
 
virtual void SetTimeStamp (const TimeStamp &time)
 
virtual ~Object () override
 
- Protected Member Functions inherited from itk::LightObject
virtual LightObject::Pointer InternalClone () const
 
 LightObject ()
 
virtual void PrintHeader (std::ostream &os, Indent indent) const
 
virtual void PrintTrailer (std::ostream &os, Indent indent) const
 
virtual ~LightObject ()
 

Static Protected Member Functions

static PrecisionType EvaluateCostCallback (void *instance, const PrecisionType *x, PrecisionType *g, const int n, const PrecisionType step)
 
static int UpdateProgressCallback (void *Instance, const PrecisionType *x, const PrecisionType *g, const PrecisionType fx, const PrecisionType xnorm, const PrecisionType gnorm, const PrecisionType step, int n, int k, int ls)
 

Private Attributes

const double * m_CurrentGradient
 
double m_CurrentGradientNorm
 
int m_CurrentNumberOfEvaluations
 
const double * m_CurrentParameter
 
double m_CurrentParameterNorm
 
double m_CurrentStepSize
 
AutoPointer
< PrivateImplementationHolder > 
m_Pimpl
 
int m_StatusCode
 

Additional Inherited Members

- Protected Attributes inherited from itk::ObjectToObjectOptimizerBaseTemplate< double >
SizeValueType m_CurrentIteration
 
MeasureType m_CurrentMetricValue
 
bool m_DoEstimateScales
 
MetricTypePointer m_Metric
 
SizeValueType m_NumberOfIterations
 
ThreadIdType m_NumberOfThreads
 
ScalesType m_Scales
 
bool m_ScalesAreIdentity
 
ScalesEstimatorType::Pointer m_ScalesEstimator
 
ScalesType m_Weights
 
bool m_WeightsAreIdentity
 
- Protected Attributes inherited from itk::LightObject
AtomicInt< int > m_ReferenceCount
 

Member Typedef Documentation

Definition at line 166 of file itkLBFGS2Optimizerv4.h.

typedef Superclass::MetricType itk::LBFGS2Optimizerv4::MetricType

Definition at line 168 of file itkLBFGS2Optimizerv4.h.

typedef Superclass::ParametersType itk::LBFGS2Optimizerv4::ParametersType

Definition at line 169 of file itkLBFGS2Optimizerv4.h.

Definition at line 165 of file itkLBFGS2Optimizerv4.h.

currently only double is used in lbfgs need to figure out how to make it a template parameter and set the required define so lbfgs.h uses the corrcet version

Definition at line 159 of file itkLBFGS2Optimizerv4.h.

typedef Superclass::ScalesType itk::LBFGS2Optimizerv4::ScalesType

Definition at line 170 of file itkLBFGS2Optimizerv4.h.

Standard "Self" typedef.

Definition at line 163 of file itkLBFGS2Optimizerv4.h.

Definition at line 164 of file itkLBFGS2Optimizerv4.h.

Member Enumeration Documentation

Enumerator
LINESEARCH_DEFAULT 

The default algorithm (MoreThuente method).

LINESEARCH_MORETHUENTE 

MoreThuente method proposd by More and Thuente.

LINESEARCH_BACKTRACKING_ARMIJO 

Backtracking method with the Armijo condition. The backtracking method finds the step length such that it satisfies the sufficient decrease (Armijo) condition,

  • f(x + a * d) <= f(x) + lbfgs_parameter_t::ftol * a * g(x)^T d,

where x is the current point, d is the current search direction, and a is the step length.

LINESEARCH_BACKTRACKING 

The backtracking method with the defualt (regular Wolfe) condition.

LINESEARCH_BACKTRACKING_WOLFE 

Backtracking method with regular Wolfe condition. The backtracking method finds the step length such that it satisfies both the Armijo condition (LINESEARCH_BACKTRACKING_ARMIJO) and the curvature condition,

  • g(x + a * d)^T d >= lbfgs_parameter_t::wolfe * g(x)^T d,

where x is the current point, d is the current search direction, and a is the step length.

LINESEARCH_BACKTRACKING_STRONG_WOLFE 

Backtracking method with strong Wolfe condition. The backtracking method finds the step length such that it satisfies both the Armijo condition (LINESEARCH_BACKTRACKING_ARMIJO) and the following condition,

  • |g(x + a * d)^T d| <= lbfgs_parameter_t::wolfe * |g(x)^T d|,

where x is the current point, d is the current search direction, and a is the step length.

Definition at line 107 of file itkLBFGS2Optimizerv4.h.

Constructor & Destructor Documentation

itk::LBFGS2Optimizerv4::LBFGS2Optimizerv4 ( )
protected
virtual itk::LBFGS2Optimizerv4::~LBFGS2Optimizerv4 ( )
overrideprotectedvirtual

Member Function Documentation

virtual::itk::LightObject::Pointer itk::LBFGS2Optimizerv4::CreateAnother ( ) const
virtual

Create an object from an instance, potentially deferring to a factory. This method allows you to create an instance of an object that is exactly the same type as the referring object. This is useful in cases where an object has been cast back to a base class.

Reimplemented from itk::Object.

PrecisionType itk::LBFGS2Optimizerv4::EvaluateCost ( const PrecisionType x,
PrecisionType g,
const int  n,
const PrecisionType  step 
)
protected
static PrecisionType itk::LBFGS2Optimizerv4::EvaluateCostCallback ( void *  instance,
const PrecisionType x,
PrecisionType g,
const int  n,
const PrecisionType  step 
)
staticprotected

Function evluation callback from libLBFGS frowrad to instance

virtual PrecisionType itk::LBFGS2Optimizerv4::GetCurrentGradientNorm ( ) const
virtual

Get gradient norm of current iteration

virtual PrecisionType itk::LBFGS2Optimizerv4::GetCurrentNumberOfEvaluations ( ) const
virtual

Get number of evaluations for current iteration

virtual PrecisionType itk::LBFGS2Optimizerv4::GetCurrentParameterNorm ( ) const
virtual

Get paramater norm of current iteration

virtual PrecisionType itk::LBFGS2Optimizerv4::GetCurrentStepSize ( ) const
virtual

Get stepsize of current iteration

int itk::LBFGS2Optimizerv4::GetDeltaConvergenceDistance ( ) const

Set/Ger distance for delta-based convergence test. This parameter determines the distance, in iterations, to compute the rate of decrease of the objective function. If the value of this parameter is zero, the library does not perform the delta-based convergence test. The default value is 0.

PrecisionType itk::LBFGS2Optimizerv4::GetDeltaConvergenceTolerance ( ) const

Delta for convergence test. This parameter determines the minimum rate of decrease of the objective function. The library stops iterations when the following condition is met: $(f' - f) / f < \delta$, where f' is the objective value of past iterations ago, and f is the objective value of the current iteration. The default value is 0.

int itk::LBFGS2Optimizerv4::GetHessianApproximationAccuracy ( ) const

Set/Get the number of corrections to approximate the inverse hessian matrix. The L-BFGS routine stores the computation results of previous m iterations to approximate the inverse hessian matrix of the current iteration. This parameter controls the size of the limited memories (corrections). The default value is 6. Values less than 3 are not recommended. Large values will result in excessive computing time.

LineSearchMethod itk::LBFGS2Optimizerv4::GetLineSearch ( ) const

The line search algorithm. This parameter specifies a line search algorithm to be used by the L-BFGS routine. See lbfgs.h for enumeration of line search type. Defaults to More-Thuente's method.

PrecisionType itk::LBFGS2Optimizerv4::GetLineSearchAccuracy ( ) const

A parameter to control the accuracy of the line search routine. The default value is 1e-4. This parameter should be greater than zero and smaller than 0.5.

PrecisionType itk::LBFGS2Optimizerv4::GetLineSearchGradientAccuracy ( ) const

A parameter to control the gradient accuracy of the More-Thuente line search routine. The default value is 0.9. If the function and gradient evaluations are inexpensive with respect to the cost of the iteration (which is sometimes the case when solving very large problems) it may be advantageous to set this parameter to a small value. A typical small value is 0.1. This parameter shuold be greater than the ftol parameter (1e-4) and smaller than 1.0.

PrecisionType itk::LBFGS2Optimizerv4::GetMachinePrecisionTolerance ( ) const

The machine precision for floating-point values. This parameter must be a positive value set by a client program to estimate the machine precision. The line search routine will terminate with the status code (LBFGSERR_ROUNDING_ERROR) if the relative width of the interval of uncertainty is less than this parameter.

int itk::LBFGS2Optimizerv4::GetMaximumIterations ( ) const

The maximum number of iterations. The lbfgs() function terminates an optimization process with LBFGSERR_MAXIMUMITERATION status code when the iteration count exceedes this parameter. Setting this parameter to zero continues an optimization process until a convergence or error. The default value is 0.

int itk::LBFGS2Optimizerv4::GetMaximumLineSearchEvaluations ( ) const

The maximum number of trials for the line search. This parameter controls the number of function and gradients evaluations per iteration for the line search routine. The default value is 20.

PrecisionType itk::LBFGS2Optimizerv4::GetMaximumLineSearchStep ( ) const

The maximum step of the line search. The default value is 1e+20. This value need not be modified unless the exponents are too large for the machine being used, or unless the problem is extremely badly scaled (in which case the exponents should be increased).

PrecisionType itk::LBFGS2Optimizerv4::GetMinimumLineSearchStep ( ) const

The minimum step of the line search routine. The default value is 1e-20. This value need not be modified unless the exponents are too large for the machine being used, or unless the problem is extremely badly scaled (in which case the exponents should be increased).

virtual const char* itk::LBFGS2Optimizerv4::GetNameOfClass ( ) const
virtual

Run-time type information (and related methods).

Reimplemented from itk::ObjectToObjectOptimizerBaseTemplate< double >.

virtual SizeValueType itk::LBFGS2Optimizerv4::GetNumberOfIterations ( ) const
inlineoverridevirtual

Aliased to Set/Get MaximumIterations to match base class interface.

Reimplemented from itk::ObjectToObjectOptimizerBaseTemplate< double >.

Definition at line 253 of file itkLBFGS2Optimizerv4.h.

PrecisionType itk::LBFGS2Optimizerv4::GetOrthantwiseCoefficient ( ) const

Coeefficient for the L1 norm of variables. This parameter should be set to zero for standard minimization problems. Setting this parameter to a positive value activates Orthant-Wise Limited-memory Quasi-Newton (OWL-QN) method, which minimizes the objective function F(x) combined with the L1 norm |x| of the variables, $F(x) + C |x|}. $. This parameter is the coefficient for the |x|, i.e., C. As the L1 norm |x| is not differentiable at zero, the library modifies function and gradient evaluations from a client program suitably; a client program thus have only to return the function value F(x) and gradients G(x) as usual. The default value is zero.

int itk::LBFGS2Optimizerv4::GetOrthantwiseEnd ( ) const

End index for computing L1 norm of the variables. This parameter is valid only for OWL-QN method (i.e., $ orthantwise_c != 0 $). This parameter e, (0 < e <= N) specifies the index number at which the library stops computing the L1 norm of the variables x,

int itk::LBFGS2Optimizerv4::GetOrthantwiseStart ( ) const

Start index for computing L1 norm of the variables. This parameter is valid only for OWL-QN method (i.e., $ orthantwise_c != 0 $). This parameter b (0 <= b < N) specifies the index number from which the library computes the L1 norm of the variables x,

\[ |x| := |x_{b}| + |x_{b+1}| + ... + |x_{N}| . \]

In other words, variables $x_1, ..., x_{b-1}$ are not used for computing the L1 norm. Setting b, (0 < b < N), one can protect variables, $x_1, ..., x_{b-1}$ (e.g., a bias term of logistic regression) from being regularized. The default value is zero.

PrecisionType itk::LBFGS2Optimizerv4::GetSolutionAccuracy ( ) const

Set/Get epsilon for convergence test. This parameter determines the accuracy with which the solution is to be found. A minimization terminates when $||g|| < \epsilon * max(1, ||x||)$, where ||.|| denotes the Euclidean (L2) norm. The default value is 1e-5.

virtual const StopConditionReturnStringType itk::LBFGS2Optimizerv4::GetStopConditionDescription ( ) const
overridevirtual

Stop condition return string type

Implements itk::ObjectToObjectOptimizerBaseTemplate< double >.

PrecisionType itk::LBFGS2Optimizerv4::GetWolfeCoefficient ( ) const

A coefficient for the Wolfe condition. This parameter is valid only when the backtracking line-search algorithm is used with the Wolfe condition, LINESEARCH_BACKTRACKING_STRONG_WOLFE or LINESEARCH_BACKTRACKING_WOLFE . The default value is 0.9. This parameter should be greater than the ftol parameter and smaller than 1.0.

static Pointer itk::LBFGS2Optimizerv4::New ( )
static

Method for creation through the object factory.

virtual void itk::LBFGS2Optimizerv4::PrintSelf ( std::ostream &  os,
Indent  indent 
) const
overrideprotectedvirtual

Methods invoked by Print() to print information about the object including superclasses. Typically not called by the user (use Print() instead) but used in the hierarchical print process to combine the output of several classes.

Reimplemented from itk::ObjectToObjectOptimizerBaseTemplate< double >.

void itk::LBFGS2Optimizerv4::SetDeltaConvergenceDistance ( int  nPast)

Set/Ger distance for delta-based convergence test. This parameter determines the distance, in iterations, to compute the rate of decrease of the objective function. If the value of this parameter is zero, the library does not perform the delta-based convergence test. The default value is 0.

void itk::LBFGS2Optimizerv4::SetDeltaConvergenceTolerance ( PrecisionType  tol)

Delta for convergence test. This parameter determines the minimum rate of decrease of the objective function. The library stops iterations when the following condition is met: $(f' - f) / f < \delta$, where f' is the objective value of past iterations ago, and f is the objective value of the current iteration. The default value is 0.

void itk::LBFGS2Optimizerv4::SetHessianApproximationAccuracy ( int  m)

Set/Get the number of corrections to approximate the inverse hessian matrix. The L-BFGS routine stores the computation results of previous m iterations to approximate the inverse hessian matrix of the current iteration. This parameter controls the size of the limited memories (corrections). The default value is 6. Values less than 3 are not recommended. Large values will result in excessive computing time.

void itk::LBFGS2Optimizerv4::SetLineSearch ( const LineSearchMethod linesearch)

The line search algorithm. This parameter specifies a line search algorithm to be used by the L-BFGS routine. See lbfgs.h for enumeration of line search type. Defaults to More-Thuente's method.

void itk::LBFGS2Optimizerv4::SetLineSearchAccuracy ( PrecisionType  ftol)

A parameter to control the accuracy of the line search routine. The default value is 1e-4. This parameter should be greater than zero and smaller than 0.5.

void itk::LBFGS2Optimizerv4::SetLineSearchGradientAccuracy ( PrecisionType  gtol)

A parameter to control the gradient accuracy of the More-Thuente line search routine. The default value is 0.9. If the function and gradient evaluations are inexpensive with respect to the cost of the iteration (which is sometimes the case when solving very large problems) it may be advantageous to set this parameter to a small value. A typical small value is 0.1. This parameter shuold be greater than the ftol parameter (1e-4) and smaller than 1.0.

void itk::LBFGS2Optimizerv4::SetMachinePrecisionTolerance ( PrecisionType  xtol)

The machine precision for floating-point values. This parameter must be a positive value set by a client program to estimate the machine precision. The line search routine will terminate with the status code (LBFGSERR_ROUNDING_ERROR) if the relative width of the interval of uncertainty is less than this parameter.

void itk::LBFGS2Optimizerv4::SetMaximumIterations ( int  maxIterations)

The maximum number of iterations. The lbfgs() function terminates an optimization process with LBFGSERR_MAXIMUMITERATION status code when the iteration count exceedes this parameter. Setting this parameter to zero continues an optimization process until a convergence or error. The default value is 0.

void itk::LBFGS2Optimizerv4::SetMaximumLineSearchEvaluations ( int  n)

The maximum number of trials for the line search. This parameter controls the number of function and gradients evaluations per iteration for the line search routine. The default value is 20.

void itk::LBFGS2Optimizerv4::SetMaximumLineSearchStep ( PrecisionType  step)

The maximum step of the line search. The default value is 1e+20. This value need not be modified unless the exponents are too large for the machine being used, or unless the problem is extremely badly scaled (in which case the exponents should be increased).

void itk::LBFGS2Optimizerv4::SetMinimumLineSearchStep ( PrecisionType  step)

The minimum step of the line search routine. The default value is 1e-20. This value need not be modified unless the exponents are too large for the machine being used, or unless the problem is extremely badly scaled (in which case the exponents should be increased).

virtual void itk::LBFGS2Optimizerv4::SetNumberOfIterations ( const SizeValueType  _arg)
inlineoverridevirtual

Aliased to Set/Get MaximumIterations to match base class interface.

Reimplemented from itk::ObjectToObjectOptimizerBaseTemplate< double >.

Definition at line 254 of file itkLBFGS2Optimizerv4.h.

void itk::LBFGS2Optimizerv4::SetOrthantwiseCoefficient ( PrecisionType  orthant_c)

Coeefficient for the L1 norm of variables. This parameter should be set to zero for standard minimization problems. Setting this parameter to a positive value activates Orthant-Wise Limited-memory Quasi-Newton (OWL-QN) method, which minimizes the objective function F(x) combined with the L1 norm |x| of the variables, $F(x) + C |x|}. $. This parameter is the coefficient for the |x|, i.e., C. As the L1 norm |x| is not differentiable at zero, the library modifies function and gradient evaluations from a client program suitably; a client program thus have only to return the function value F(x) and gradients G(x) as usual. The default value is zero.

void itk::LBFGS2Optimizerv4::SetOrthantwiseEnd ( int  end)

End index for computing L1 norm of the variables. This parameter is valid only for OWL-QN method (i.e., $ orthantwise_c != 0 $). This parameter e, (0 < e <= N) specifies the index number at which the library stops computing the L1 norm of the variables x,

void itk::LBFGS2Optimizerv4::SetOrthantwiseStart ( int  start)

Start index for computing L1 norm of the variables. This parameter is valid only for OWL-QN method (i.e., $ orthantwise_c != 0 $). This parameter b (0 <= b < N) specifies the index number from which the library computes the L1 norm of the variables x,

\[ |x| := |x_{b}| + |x_{b+1}| + ... + |x_{N}| . \]

In other words, variables $x_1, ..., x_{b-1}$ are not used for computing the L1 norm. Setting b, (0 < b < N), one can protect variables, $x_1, ..., x_{b-1}$ (e.g., a bias term of logistic regression) from being regularized. The default value is zero.

virtual void itk::LBFGS2Optimizerv4::SetScales ( const ScalesType )
overridevirtual

This optimizer does not support scaling of the derivatives.

void itk::LBFGS2Optimizerv4::SetSolutionAccuracy ( PrecisionType  epsilon)

Set/Get epsilon for convergence test. This parameter determines the accuracy with which the solution is to be found. A minimization terminates when $||g|| < \epsilon * max(1, ||x||)$, where ||.|| denotes the Euclidean (L2) norm. The default value is 1e-5.

virtual void itk::LBFGS2Optimizerv4::SetWeights ( const ScalesType  )
overridevirtual

This optimizer does not support weighting of the derivatives.

void itk::LBFGS2Optimizerv4::SetWolfeCoefficient ( PrecisionType  wc)

A coefficient for the Wolfe condition. This parameter is valid only when the backtracking line-search algorithm is used with the Wolfe condition, LINESEARCH_BACKTRACKING_STRONG_WOLFE or LINESEARCH_BACKTRACKING_WOLFE . The default value is 0.9. This parameter should be greater than the ftol parameter and smaller than 1.0.

virtual void itk::LBFGS2Optimizerv4::StartOptimization ( bool  doOnlyInitialization = false)
overridevirtual

Start optimization with an initial value.

Reimplemented from itk::ObjectToObjectOptimizerBaseTemplate< double >.

int itk::LBFGS2Optimizerv4::UpdateProgress ( const PrecisionType x,
const PrecisionType g,
const PrecisionType  fx,
const PrecisionType  xnorm,
const PrecisionType  gnorm,
const PrecisionType  step,
int  n,
int  k,
int  ls 
)
protected

Update the progress as reported from libLBFSG and notify itkObject

static int itk::LBFGS2Optimizerv4::UpdateProgressCallback ( void *  Instance,
const PrecisionType x,
const PrecisionType g,
const PrecisionType  fx,
const PrecisionType  xnorm,
const PrecisionType  gnorm,
const PrecisionType  step,
int  n,
int  k,
int  ls 
)
staticprotected

Progress callback from libLBFGS forwards it to the specific instance

Member Data Documentation

const double* itk::LBFGS2Optimizerv4::m_CurrentGradient
private

Progress update variables

Definition at line 462 of file itkLBFGS2Optimizerv4.h.

double itk::LBFGS2Optimizerv4::m_CurrentGradientNorm
private

Definition at line 467 of file itkLBFGS2Optimizerv4.h.

int itk::LBFGS2Optimizerv4::m_CurrentNumberOfEvaluations
private

Definition at line 468 of file itkLBFGS2Optimizerv4.h.

const double* itk::LBFGS2Optimizerv4::m_CurrentParameter
private

Definition at line 463 of file itkLBFGS2Optimizerv4.h.

double itk::LBFGS2Optimizerv4::m_CurrentParameterNorm
private

Definition at line 466 of file itkLBFGS2Optimizerv4.h.

double itk::LBFGS2Optimizerv4::m_CurrentStepSize
private

Definition at line 465 of file itkLBFGS2Optimizerv4.h.

AutoPointer<PrivateImplementationHolder> itk::LBFGS2Optimizerv4::m_Pimpl
private

Definition at line 457 of file itkLBFGS2Optimizerv4.h.

int itk::LBFGS2Optimizerv4::m_StatusCode
private

Definition at line 470 of file itkLBFGS2Optimizerv4.h.


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