RVLUmin::LBFGSBT< Scalar > Class Template Reference

Limited memory Broyden-Fletcher-Goldfarb-Shanno (LBFGS) quasi-Newton optimization with geometric backtracking line search globalization. More...

#include <LBFGSBT.hh>

Inheritance diagram for RVLUmin::LBFGSBT< Scalar >:

RVLAlg::Algorithm List of all members.

Public Member Functions

 LBFGSBT (Functional< Scalar > const &f, Vector< Scalar > &x, Scalar _ihs, int _mud, int _maxsamp, bool _disp, Scalar _sl1, Scalar _eta1, Scalar _eta2, Scalar _gamma1, Scalar _gamma2, Scalar _maxfrac, Scalar _minsteptol, int _maxits, Scalar _agradtol, Scalar _rgradtol, ostream &_str=cout)
 constructor, first form:
 LBFGSBT (Functional< Scalar > const &f, Vector< Scalar > &x, Table const &t, ostream &_str=cout)
 constructor, second form - scalar arguments passed in Table.
int getCount ()
void run ()
 ~LBFGSBT ()
FunctionalEvaluation< Scalar
> const & 
getFunctionalEvaluation () const
 supplied to provide access to any intermediate data that a subclass of Functional may make available.

Detailed Description

template<typename Scalar>
class RVLUmin::LBFGSBT< Scalar >

Limited memory Broyden-Fletcher-Goldfarb-Shanno (LBFGS) quasi-Newton optimization with geometric backtracking line search globalization.

For details of the LBFGS algorithm, see the paper

"Updating Quasi-Newton Matrices with Limited Storage" by Jorge Nocedal, Math. of Computation, Vol. 35, no. 151, p.p. 773--782.

Approximates solution to unconstrained optimization problem

\[ \min_{x} f(x) \]

via generalized secant update

\[ x_+ = x_c - \alpha B \nabla f(x_c). \]

$B$ is a limited-memory secant approximation to the inverse Hessian of $f$, parametrized by maximum rank. $\alpha$ is a step length.

Note that this problem setting makes sense only for functions $f$ taking real values. Accordingly, the implementation includes a compile-time check that the Scalar template parameter designates a real field (i.e. double or float).

This implementation globalizes convergence to a local min via a geometric backtracking linesearch method, impemented in RVLUmin::BacktrackingLineSearchAlg. Line search methods are described in Nocedal and Wright, Numerical Optimization, Springer 1999. For additional documentation, see RVLUmin::BacktrackingLineSearchAlg docs.

Structure and function: Data members of three types combine to implement LBFGS:

Usage: initialize solution vector (constructor argument x) externally; construct LBFGS object; call LBFGSBT::run(). On return, x has been updated to estimated solution.

Parameters for LBFGS iteration:

NOTES:

Parameters for backtracking line search (for detailed description and notes, see docs for RVLUmin::BacktrackingLineSearchAlg):

Typical use case: see functional test source.

Definition at line 141 of file LBFGSBT.hh.


Constructor & Destructor Documentation

template<typename Scalar>
RVLUmin::LBFGSBT< Scalar >::LBFGSBT ( Functional< Scalar > const &  f,
Vector< Scalar > &  x,
Scalar  _ihs,
int  _mud,
int  _maxsamp,
bool  _disp,
Scalar  _sl1,
Scalar  _eta1,
Scalar  _eta2,
Scalar  _gamma1,
Scalar  _gamma2,
Scalar  _maxfrac,
Scalar  _minsteptol,
int  _maxits,
Scalar  _agradtol,
Scalar  _rgradtol,
ostream &  _str = cout 
)

constructor, first form:

parameters:

Parameters:
f - function to be minimized (RVL::Functional)
x - solution RVL::Vector - initial guess on call, estimated solution on return
_ihs - inverse Hessian scale - overall scale factor, so initial Hessian is this Scalar multiple of identity operator
_maxits - max number of LBFGS iterations
_mud - max stored BFGS updates - stored inverse Hessian approximation has this rank (at most)
_maxsamp - max number of steps permitted in each line search
_disp - verbosity flag - false = no output, true = function value, gradient norm at each iteration, report of line search
_sl1 - first line search step
_eta1 - lower G-A parameter
_eta2 - upper G-A parameter
_gamma1 - line search backtrack factor
_gamma2 - line search extrapolation factor ("internal doubling")
_maxfrac - fraction of max step to boundary permitted
_agradtol - stopping tolerance for gradient norm, absolute
_rgradtol - stopping tolerance for gradient norm, relative to initial gradient norm
_str - verbose output unit

Definition at line 198 of file LBFGSBT.hh.

template<typename Scalar>
RVLUmin::LBFGSBT< Scalar >::LBFGSBT ( Functional< Scalar > const &  f,
Vector< Scalar > &  x,
Table const &  t,
ostream &  _str = cout 
)

constructor, second form - scalar arguments passed in Table.

Definition at line 236 of file LBFGSBT.hh.

template<typename Scalar>
RVLUmin::LBFGSBT< Scalar >::~LBFGSBT (  ) 

Definition at line 278 of file LBFGSBT.hh.


Member Function Documentation

template<typename Scalar>
int RVLUmin::LBFGSBT< Scalar >::getCount (  ) 

Definition at line 260 of file LBFGSBT.hh.

template<typename Scalar>
void RVLUmin::LBFGSBT< Scalar >::run (  )  [virtual]

Implements RVLAlg::Algorithm.

Definition at line 262 of file LBFGSBT.hh.

References RVLAlg::LoopAlg::run().

template<typename Scalar>
FunctionalEvaluation<Scalar> const& RVLUmin::LBFGSBT< Scalar >::getFunctionalEvaluation (  )  const

supplied to provide access to any intermediate data that a subclass of Functional may make available.

Can extract const reference to current Functional, as constructed by FunctionalEvaluation, via the getFcnl() method. A cast will be required to extract any further subclass attributes.

Definition at line 285 of file LBFGSBT.hh.


The documentation for this class was generated from the following file:
Generated on 5 Jan 2017 for RvlUmin by  doxygen 1.4.7