RVLUmin::BacktrackingLineSearchAlgBase< Scalar > Class Template Reference

Does a backtracking line search starting from a prescribed step, passed as argument firststep to the constructor. More...

#include <lnsrchBT.hh>

Inheritance diagram for RVLUmin::BacktrackingLineSearchAlgBase< Scalar >:

RVLUmin::LineSearchAlgBase< Scalar > RVLAlg::Algorithm RVLAlg::Terminator List of all members.

Public Member Functions

 BacktrackingLineSearchAlgBase (LineSearchAlg< Scalar > const &lsalg, FunctionalEvaluation< Scalar > &fx, Vector< Scalar > const &dx, Scalar _step, Scalar _eta1, Scalar _eta2, Scalar _gamma1, Scalar _gamma2, bool _DispFlag, Scalar _fudge, int _maxsteps, ostream &_str)
 Constructor arguments:.
 BacktrackingLineSearchAlgBase (const BacktrackingLineSearchAlgBase< Scalar > &ls)
virtual ~BacktrackingLineSearchAlgBase ()
ostream & write (ostream &str) const
bool query ()
Scalar getStep () const
void run ()
 Checks for sufficient decrease, and while not found, shrink step.

Protected Member Functions

LineSearchAlgBase< Scalar > * clone () const

Detailed Description

template<class Scalar>
class RVLUmin::BacktrackingLineSearchAlgBase< Scalar >

Does a backtracking line search starting from a prescribed step, passed as argument firststep to the constructor.

The initial step must be prescribed externally, as no scale information is available at the outset. Starts with backtrack loop: at each unsuccessful iteration, step <- step*gamma, with gamma < 1. Stops when sufficient decrease condition is satisfied. If this happens on first step, goes into internal doubling loop: at each doubling iteration, step <- step*stretch, with 1 < stretch < gamma (last inequality to avoid cycling). Internal doubling stops if function value increases.

Definition at line 61 of file lnsrchBT.hh.


Constructor & Destructor Documentation

template<class Scalar>
RVLUmin::BacktrackingLineSearchAlgBase< Scalar >::BacktrackingLineSearchAlgBase ( LineSearchAlg< Scalar > const &  lsalg,
FunctionalEvaluation< Scalar > &  fx,
Vector< Scalar > const &  dx,
Scalar  _step,
Scalar  _eta1,
Scalar  _eta2,
Scalar  _gamma1,
Scalar  _gamma2,
bool  _DispFlag,
Scalar  _fudge,
int  _maxsteps,
ostream &  _str 
)

Constructor arguments:.

Parameters:
lsalg const reference to line search algorithm
fx current FunctionalEvaluation on call, updated on return
dx const reference to search direction
_step initial step on call, last successful step on return
_eta1 lower G-A parameter > 0
_eta2 upper G-A parameter > eta1
_gamma1 trust region reduction factor < 1
_gamma2 trust region expansion factor > 1, gamma1*gamma2 < 1
_DispFlag verbosity flag
_fudge fraction of max step to permit
_maxsteps max number of line search steps
_str verbose output unit

Definition at line 106 of file lnsrchBT.hh.

template<class Scalar>
RVLUmin::BacktrackingLineSearchAlgBase< Scalar >::BacktrackingLineSearchAlgBase ( const BacktrackingLineSearchAlgBase< Scalar > &  ls  ) 

Definition at line 132 of file lnsrchBT.hh.

template<class Scalar>
virtual RVLUmin::BacktrackingLineSearchAlgBase< Scalar >::~BacktrackingLineSearchAlgBase (  )  [virtual]

Definition at line 145 of file lnsrchBT.hh.


Member Function Documentation

template<class Scalar>
LineSearchAlgBase<Scalar>* RVLUmin::BacktrackingLineSearchAlgBase< Scalar >::clone (  )  const [protected, virtual]

Implements RVLUmin::LineSearchAlgBase< Scalar >.

Definition at line 84 of file lnsrchBT.hh.

template<class Scalar>
ostream& RVLUmin::BacktrackingLineSearchAlgBase< Scalar >::write ( ostream &  str  )  const [virtual]

Reimplemented from RVLUmin::LineSearchAlgBase< Scalar >.

Definition at line 147 of file lnsrchBT.hh.

template<class Scalar>
bool RVLUmin::BacktrackingLineSearchAlgBase< Scalar >::query (  )  [virtual]

Implements RVLAlg::Terminator.

Definition at line 152 of file lnsrchBT.hh.

template<class Scalar>
Scalar RVLUmin::BacktrackingLineSearchAlgBase< Scalar >::getStep (  )  const [virtual]

Implements RVLUmin::LineSearchAlgBase< Scalar >.

Definition at line 154 of file lnsrchBT.hh.

template<class Scalar>
void RVLUmin::BacktrackingLineSearchAlgBase< Scalar >::run (  )  [virtual]

Checks for sufficient decrease, and while not found, shrink step.

Implements RVLAlg::Algorithm.

Definition at line 160 of file lnsrchBT.hh.

References RVL::Vector< Scalar >::copy(), RVLUmin::LineSearchAlgBase< Scalar >::dx, RVLUmin::LineSearchAlgBase< Scalar >::fx, RVLUmin::LineSearchAlgBase< Scalar >::getBasePoint(), RVLUmin::LineSearchAlgBase< Scalar >::getMinStep(), RVLUmin::LineSearchAlgBase< Scalar >::getSearchDirection(), RVL::Vector< Scalar >::getSpace(), and RVL::Vector< Scalar >::linComb().


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