ITK
6.0.0
Insight Toolkit
|
#include <itkLBFGS2Optimizerv4.h>
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-memory Quasi-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 available 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 165 of file itkLBFGS2Optimizerv4.h.
Static Public Member Functions | |
static Pointer | New () |
Static Public Member Functions inherited from itk::GradientDescentOptimizerv4Template< TInternalComputationValueType > | |
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 val) |
Static Public Member Functions inherited from itk::LightObject | |
static void | BreakOnError () |
static Pointer | New () |
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 Member Functions | |
void | AdvanceOneStep () override |
void | SetMinimumConvergenceValue (PrecisionType) override |
void | SetConvergenceWindowSize (SizeValueType) override |
const PrecisionType & | GetConvergenceValue () const override |
Private Attributes | |
double | m_CurrentGradientNorm {} |
int | m_CurrentNumberOfEvaluations {} |
double | m_CurrentParameterNorm {} |
double | m_CurrentStepSize {} |
bool | m_EstimateScalesAtEachIteration {} |
lbfgs_parameter_t | m_Parameters {} |
int | m_StatusCode {} |
using itk::LBFGS2Optimizerv4Template< TInternalComputationValueType >::ConstPointer = SmartPointer<const Self> |
Definition at line 200 of file itkLBFGS2Optimizerv4.h.
using itk::LBFGS2Optimizerv4Template< TInternalComputationValueType >::LineSearchMethodEnum = LBFGS2Optimizerv4Enums::LineSearchMethod |
Definition at line 173 of file itkLBFGS2Optimizerv4.h.
using itk::LBFGS2Optimizerv4Template< TInternalComputationValueType >::MetricType = typename Superclass::MetricType |
Definition at line 202 of file itkLBFGS2Optimizerv4.h.
using itk::LBFGS2Optimizerv4Template< TInternalComputationValueType >::ParametersType = typename Superclass::ParametersType |
Definition at line 203 of file itkLBFGS2Optimizerv4.h.
using itk::LBFGS2Optimizerv4Template< TInternalComputationValueType >::Pointer = SmartPointer<Self> |
Definition at line 199 of file itkLBFGS2Optimizerv4.h.
using itk::LBFGS2Optimizerv4Template< TInternalComputationValueType >::PrecisionType = double |
TODO: 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 correct version
Definition at line 192 of file itkLBFGS2Optimizerv4.h.
using itk::LBFGS2Optimizerv4Template< TInternalComputationValueType >::ScalesType = typename Superclass::ScalesType |
Definition at line 204 of file itkLBFGS2Optimizerv4.h.
using itk::LBFGS2Optimizerv4Template< TInternalComputationValueType >::Self = LBFGS2Optimizerv4Template |
Standard "Self" type alias.
Definition at line 197 of file itkLBFGS2Optimizerv4.h.
using itk::LBFGS2Optimizerv4Template< TInternalComputationValueType >::Superclass = GradientDescentOptimizerv4Template<TInternalComputationValueType> |
Definition at line 198 of file itkLBFGS2Optimizerv4.h.
|
protected |
|
overrideprotected |
|
inlineoverrideprivatevirtual |
Advance one step following the gradient direction. Includes transform update.
Reimplemented from itk::GradientDescentOptimizerv4Template< TInternalComputationValueType >.
Definition at line 569 of file itkLBFGS2Optimizerv4.h.
|
virtual |
Option to use ScalesEstimator for estimating scales at each iteration. The estimation overrides the scales set by SetScales(). Default is true.
|
protected |
|
staticprotected |
Function evaluation callback from libLBFGS forward to instance
|
inlineoverrideprivatevirtual |
itkGradientDecentOptimizerv4Template specific non supported methods.
Reimplemented from itk::GradientDescentOptimizerv4Template< TInternalComputationValueType >.
Definition at line 560 of file itkLBFGS2Optimizerv4.h.
|
virtual |
Get gradient norm of current iteration
|
virtual |
Get number of evaluations for current iteration
|
virtual |
Get parameter norm of current iteration
|
virtual |
Get step size of current iteration
int itk::LBFGS2Optimizerv4Template< TInternalComputationValueType >::GetDeltaConvergenceDistance | ( | ) | const |
Set/Get 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::LBFGS2Optimizerv4Template< TInternalComputationValueType >::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
.
|
virtual |
Option to use ScalesEstimator for estimating scales at each iteration. The estimation overrides the scales set by SetScales(). Default is true.
int itk::LBFGS2Optimizerv4Template< TInternalComputationValueType >::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.
LineSearchMethodEnum itk::LBFGS2Optimizerv4Template< TInternalComputationValueType >::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::LBFGS2Optimizerv4Template< TInternalComputationValueType >::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::LBFGS2Optimizerv4Template< TInternalComputationValueType >::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 should be greater than the ftol
parameter (1e-4
) and smaller than 1.0
.
PrecisionType itk::LBFGS2Optimizerv4Template< TInternalComputationValueType >::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::LBFGS2Optimizerv4Template< TInternalComputationValueType >::GetMaximumIterations | ( | ) | const |
The maximum number of iterations. The lbfgs() function terminates an optimization process with LBFGSERR_MAXIMUMITERATION
status code when the iteration count exceeds this parameter. Setting this parameter to zero continues an optimization process until a convergence or error. The default value is 0
.
int itk::LBFGS2Optimizerv4Template< TInternalComputationValueType >::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::LBFGS2Optimizerv4Template< TInternalComputationValueType >::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::LBFGS2Optimizerv4Template< TInternalComputationValueType >::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).
|
overridevirtual |
Reimplemented from itk::GradientDescentOptimizerv4Template< TInternalComputationValueType >.
|
inlineoverridevirtual |
Aliased to Set/Get MaximumIterations to match base class interface.
Reimplemented from itk::ObjectToObjectOptimizerBaseTemplate< TInternalComputationValueType >.
Definition at line 302 of file itkLBFGS2Optimizerv4.h.
PrecisionType itk::LBFGS2Optimizerv4Template< TInternalComputationValueType >::GetOrthantwiseCoefficient | ( | ) | const |
Coefficient 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::LBFGS2Optimizerv4Template< TInternalComputationValueType >::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::LBFGS2Optimizerv4Template< TInternalComputationValueType >::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::LBFGS2Optimizerv4Template< TInternalComputationValueType >::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
.
|
overridevirtual |
Get the reason for termination
Reimplemented from itk::GradientDescentOptimizerBasev4Template< TInternalComputationValueType >.
PrecisionType itk::LBFGS2Optimizerv4Template< TInternalComputationValueType >::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 |
Method for creation through the object factory.
|
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::GradientDescentOptimizerv4Template< TInternalComputationValueType >.
|
overridevirtual |
Resume optimization. This runs the optimization loop, and allows continuation of stopped optimization
Reimplemented from itk::GradientDescentOptimizerv4Template< TInternalComputationValueType >.
|
inlineoverrideprivatevirtual |
itkGradientDecentOptimizerv4Template specific non supported methods.
Reimplemented from itk::GradientDescentOptimizerv4Template< TInternalComputationValueType >.
Definition at line 555 of file itkLBFGS2Optimizerv4.h.
void itk::LBFGS2Optimizerv4Template< TInternalComputationValueType >::SetDeltaConvergenceDistance | ( | int | nPast | ) |
Set/Get 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::LBFGS2Optimizerv4Template< TInternalComputationValueType >::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
.
|
virtual |
Option to use ScalesEstimator for estimating scales at each iteration. The estimation overrides the scales set by SetScales(). Default is true.
void itk::LBFGS2Optimizerv4Template< TInternalComputationValueType >::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::LBFGS2Optimizerv4Template< TInternalComputationValueType >::SetLineSearch | ( | const LineSearchMethodEnum & | 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::LBFGS2Optimizerv4Template< TInternalComputationValueType >::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::LBFGS2Optimizerv4Template< TInternalComputationValueType >::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 should be greater than the ftol
parameter (1e-4
) and smaller than 1.0
.
void itk::LBFGS2Optimizerv4Template< TInternalComputationValueType >::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::LBFGS2Optimizerv4Template< TInternalComputationValueType >::SetMaximumIterations | ( | int | maxIterations | ) |
The maximum number of iterations. The lbfgs() function terminates an optimization process with LBFGSERR_MAXIMUMITERATION
status code when the iteration count exceeds this parameter. Setting this parameter to zero continues an optimization process until a convergence or error. The default value is 0
.
void itk::LBFGS2Optimizerv4Template< TInternalComputationValueType >::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::LBFGS2Optimizerv4Template< TInternalComputationValueType >::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).
|
inlineoverrideprivate |
itkGradientDecentOptimizerv4Template specific non supported methods.
Definition at line 550 of file itkLBFGS2Optimizerv4.h.
void itk::LBFGS2Optimizerv4Template< TInternalComputationValueType >::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).
|
inlineoverridevirtual |
Aliased to Set/Get MaximumIterations to match base class interface.
Reimplemented from itk::ObjectToObjectOptimizerBaseTemplate< TInternalComputationValueType >.
Definition at line 307 of file itkLBFGS2Optimizerv4.h.
void itk::LBFGS2Optimizerv4Template< TInternalComputationValueType >::SetOrthantwiseCoefficient | ( | PrecisionType | orthant_c | ) |
Coefficient 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::LBFGS2Optimizerv4Template< TInternalComputationValueType >::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::LBFGS2Optimizerv4Template< TInternalComputationValueType >::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.
void itk::LBFGS2Optimizerv4Template< TInternalComputationValueType >::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
.
void itk::LBFGS2Optimizerv4Template< TInternalComputationValueType >::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
.
|
overridevirtual |
Start optimization with an initial value.
Reimplemented from itk::GradientDescentOptimizerv4Template< TInternalComputationValueType >.
|
protected |
Update the progress as reported from libLBFSG and notify itkObject
|
staticprotected |
Progress callback from libLBFGS forwards it to the specific instance
|
private |
Definition at line 542 of file itkLBFGS2Optimizerv4.h.
|
private |
Definition at line 543 of file itkLBFGS2Optimizerv4.h.
|
private |
Definition at line 541 of file itkLBFGS2Optimizerv4.h.
|
private |
Definition at line 540 of file itkLBFGS2Optimizerv4.h.
|
private |
Definition at line 539 of file itkLBFGS2Optimizerv4.h.
|
private |
Definition at line 537 of file itkLBFGS2Optimizerv4.h.
|
private |
Definition at line 544 of file itkLBFGS2Optimizerv4.h.