ROL
Public Member Functions | Private Member Functions | Private Attributes | List of all members
ROL::LineSearchStep< Real > Class Template Reference

Provides the interface to compute optimization steps with line search. More...

#include <ROL_LineSearchStep.hpp>

+ Inheritance diagram for ROL::LineSearchStep< Real >:

Public Member Functions

 LineSearchStep (ROL::ParameterList &parlist, const ROL::Ptr< LineSearch< Real > > &lineSearch=ROL::nullPtr, const ROL::Ptr< Secant< Real > > &secant=ROL::nullPtr, const ROL::Ptr< Krylov< Real > > &krylov=ROL::nullPtr, const ROL::Ptr< NonlinearCG< Real > > &nlcg=ROL::nullPtr)
 Constructor. More...
 
void initialize (Vector< Real > &x, const Vector< Real > &s, const Vector< Real > &g, Objective< Real > &obj, BoundConstraint< Real > &bnd, AlgorithmState< Real > &algo_state)
 Initialize step with bound constraint. More...
 
void compute (Vector< Real > &s, const Vector< Real > &x, Objective< Real > &obj, BoundConstraint< Real > &bnd, AlgorithmState< Real > &algo_state)
 Compute step. More...
 
void update (Vector< Real > &x, const Vector< Real > &s, Objective< Real > &obj, BoundConstraint< Real > &bnd, AlgorithmState< Real > &algo_state)
 Update step, if successful. More...
 
std::string printHeader (void) const
 Print iterate header. More...
 
std::string printName (void) const
 Print step name. More...
 
std::string print (AlgorithmState< Real > &algo_state, bool print_header=false) const
 Print iterate status. More...
 
- Public Member Functions inherited from ROL::Step< Real >
virtual ~Step ()
 
 Step (void)
 
virtual void initialize (Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, BoundConstraint< Real > &con, AlgorithmState< Real > &algo_state)
 Initialize step with bound constraint. More...
 
virtual void initialize (Vector< Real > &x, const Vector< Real > &g, Vector< Real > &l, const Vector< Real > &c, Objective< Real > &obj, Constraint< Real > &con, AlgorithmState< Real > &algo_state)
 Initialize step with equality constraint. More...
 
virtual void initialize (Vector< Real > &x, const Vector< Real > &g, Vector< Real > &l, const Vector< Real > &c, Objective< Real > &obj, Constraint< Real > &con, BoundConstraint< Real > &bnd, AlgorithmState< Real > &algo_state)
 Initialize step with equality constraint. More...
 
virtual void compute (Vector< Real > &s, const Vector< Real > &x, const Vector< Real > &l, Objective< Real > &obj, Constraint< Real > &con, AlgorithmState< Real > &algo_state)
 Compute step (equality constraints). More...
 
virtual void update (Vector< Real > &x, Vector< Real > &l, const Vector< Real > &s, Objective< Real > &obj, Constraint< Real > &con, AlgorithmState< Real > &algo_state)
 Update step, if successful (equality constraints). More...
 
virtual void compute (Vector< Real > &s, const Vector< Real > &x, const Vector< Real > &l, Objective< Real > &obj, Constraint< Real > &con, BoundConstraint< Real > &bnd, AlgorithmState< Real > &algo_state)
 Compute step (equality constraints). More...
 
virtual void update (Vector< Real > &x, Vector< Real > &l, const Vector< Real > &s, Objective< Real > &obj, Constraint< Real > &con, BoundConstraint< Real > &bnd, AlgorithmState< Real > &algo_state)
 Update step, if successful (equality constraints). More...
 
const ROL::Ptr< const
StepState< Real > > 
getStepState (void) const
 Get state for step object. More...
 
void reset (const Real searchSize=1.0)
 Get state for step object. More...
 

Private Member Functions

Real GradDotStep (const Vector< Real > &g, const Vector< Real > &s, const Vector< Real > &x, BoundConstraint< Real > &bnd, Real eps=0)
 

Private Attributes

ROL::Ptr< Step< Real > > desc_
 Unglobalized step object. More...
 
ROL::Ptr< Secant< Real > > secant_
 Secant object (used for quasi-Newton) More...
 
ROL::Ptr< Krylov< Real > > krylov_
 Krylov solver object (used for inexact Newton) More...
 
ROL::Ptr< NonlinearCG< Real > > nlcg_
 Nonlinear CG object (used for nonlinear CG) More...
 
ROL::Ptr< LineSearch< Real > > lineSearch_
 Line-search object. More...
 
ROL::Ptr< Vector< Real > > d_
 
ELineSearch els_
 enum determines type of line search More...
 
ECurvatureCondition econd_
 enum determines type of curvature condition More...
 
bool acceptLastAlpha_
 For backwards compatibility. When max function evaluations are reached take last step. More...
 
bool usePreviousAlpha_
 If true, use the previously accepted step length (if any) as the new initial step length. More...
 
int verbosity_
 
bool computeObj_
 
Real fval_
 
ROL::ParameterList parlist_
 
std::string lineSearchName_
 

Additional Inherited Members

- Protected Member Functions inherited from ROL::Step< Real >
ROL::Ptr< StepState< Real > > getState (void)
 

Detailed Description

template<class Real>
class ROL::LineSearchStep< Real >

Provides the interface to compute optimization steps with line search.

Suppose \(\mathcal{X}\) is a Hilbert space of functions mapping \(\Xi\) to \(\mathbb{R}\). For example, \(\Xi\subset\mathbb{R}^n\) and \(\mathcal{X}=L^2(\Xi)\) or \(\Xi = \{1,\ldots,n\}\) and \(\mathcal{X}=\mathbb{R}^n\). We assume \(f:\mathcal{X}\to\mathbb{R}\) is twice-continuously Fréchet differentiable and \(a,\,b\in\mathcal{X}\) with \(a\le b\) almost everywhere in \(\Xi\). Note that these line-search algorithms will also work with secant approximations of the Hessian. This step applies to unconstrained and bound constrained optimization problems,

\[ \min_x\quad f(x) \qquad\text{and}\qquad \min_x\quad f(x)\quad\text{s.t.}\quad a\le x\le b, \]

respectively.

For unconstrained problems, given the \(k\)-th iterate \(x_k\) and a descent direction \(s_k\), the line search approximately minimizes the 1D objective function \(\phi_k(t) = f(x_k + t s_k)\). The approximate minimizer \(t\) must satisfy sufficient decrease and curvature conditions into guarantee global convergence. The sufficient decrease condition (often called the Armijo condition) is

\[ \phi_k(t) \le \phi_k(0) + c_1 t \phi_k'(0) \qquad\iff\qquad f(x_k+ts_k) \le f(x_k) + c_1 t \langle \nabla f(x_k),s_k\rangle_{\mathcal{X}} \]

where \(0 < c_1 < 1\). The curvature conditions implemented in ROL include:

Name Condition Parameters
Wolfe \(\phi_k'(t) \ge c_2\phi_k'(0)\) \(c_1<c_2<1\)
Strong Wolfe \(\left|\phi_k'(t)\right| \le c_2 \left|\phi_k'(0)\right|\) \(c_1<c_2<1\)
Generalized Wolfe \(c_2\phi_k'(0)\le \phi_k'(t)\le-c_3\phi_k'(0)\) \(0<c_3<1\)
Approximate Wolfe \(c_2\phi_k'(0)\le \phi_k'(t)\le(2c_1-1)\phi_k'(0)\) \(c_1<c_2<1\)
Goldstein \(\phi_k(0)+(1-c_1)t\phi_k'(0)\le \phi_k(t)\) \(0<c_1<\frac{1}{2}\)

Note that \(\phi_k'(t) = \langle \nabla f(x_k+ts_k),s_k\rangle_{\mathcal{X}}\).

For bound constrained problems, the line-search step performs a projected search. That is, the line search approximately minimizes the 1D objective function \(\phi_k(t) = f(P_{[a,b]}(x_k+ts_k))\) where \(P_{[a,b]}\) denotes the projection onto the upper and lower bounds. Such line-search algorithms correspond to projected gradient and projected Newton-type algorithms.

For projected methods, we require the notion of an active set of an iterate \(x_k\),

\[ \mathcal{A}_k = \{\, \xi\in\Xi\,:\,x_k(\xi) = a(\xi)\,\}\cap \{\, \xi\in\Xi\,:\,x_k(\xi) = b(\xi)\,\}. \]

Given \(\mathcal{A}_k\) and a search direction \(s_k\), we define the binding set as

\[ \mathcal{B}_k = \{\, \xi\in\Xi\,:\,x_k(\xi) = a(\xi) \;\text{and}\; s_k(\xi) < 0 \,\}\cap \{\, \xi\in\Xi\,:\,x_k(\xi) = b(\xi) \;\text{and}\; s_k(\xi) > 0 \,\}. \]

The binding set contains the values of \(\xi\in\Xi\) such that if \(x_k(\xi)\) is on a bound, then \((x_k+s_k)(\xi)\) will violate bound. Using these definitions, we can redefine the sufficient decrease and curvature conditions from the unconstrained case to the case of bound constraints.

LineSearchStep implements a number of algorithms for both bound constrained and unconstrained optimization. These algorithms are: Steepest descent; Nonlinear CG; Quasi-Newton methods; Inexact Newton methods; Newton's method. These methods are chosen through the EDescent enum.

Definition at line 138 of file ROL_LineSearchStep.hpp.

Constructor & Destructor Documentation

template<class Real>
ROL::LineSearchStep< Real >::LineSearchStep ( ROL::ParameterList &  parlist,
const ROL::Ptr< LineSearch< Real > > &  lineSearch = ROL::nullPtr,
const ROL::Ptr< Secant< Real > > &  secant = ROL::nullPtr,
const ROL::Ptr< Krylov< Real > > &  krylov = ROL::nullPtr,
const ROL::Ptr< NonlinearCG< Real > > &  nlcg = ROL::nullPtr 
)
inline

Constructor.

Standard constructor to build a LineSearchStep object. Algorithmic specifications are passed in through a ROL::ParameterList. The line-search type, secant type, Krylov type, or nonlinear CG type can be set using user-defined objects.

Parameters
[in]parlistis a parameter list containing algorithmic specifications
[in]lineSearchis a user-defined line search object
[in]secantis a user-defined secant object
[in]krylovis a user-defined Krylov object
[in]nlcgis a user-defined Nonlinear CG object

Definition at line 205 of file ROL_LineSearchStep.hpp.

References ROL::LineSearchStep< Real >::acceptLastAlpha_, ROL::LineSearchStep< Real >::computeObj_, ROL::LineSearchStep< Real >::econd_, ROL::LineSearchStep< Real >::els_, ROL::LineSearchStep< Real >::lineSearch_, ROL::LineSearchStep< Real >::lineSearchName_, ROL::StringToECurvatureCondition(), ROL::StringToELineSearch(), and ROL::LineSearchStep< Real >::verbosity_.

Member Function Documentation

template<class Real>
Real ROL::LineSearchStep< Real >::GradDotStep ( const Vector< Real > &  g,
const Vector< Real > &  s,
const Vector< Real > &  x,
BoundConstraint< Real > &  bnd,
Real  eps = 0 
)
inlineprivate
template<class Real>
void ROL::LineSearchStep< Real >::initialize ( Vector< Real > &  x,
const Vector< Real > &  s,
const Vector< Real > &  g,
Objective< Real > &  obj,
BoundConstraint< Real > &  con,
AlgorithmState< Real > &  algo_state 
)
inlinevirtual
template<class Real>
void ROL::LineSearchStep< Real >::compute ( Vector< Real > &  s,
const Vector< Real > &  x,
Objective< Real > &  obj,
BoundConstraint< Real > &  bnd,
AlgorithmState< Real > &  algo_state 
)
inlinevirtual

Compute step.

Computes a trial step, \(s_k\) as defined by the enum EDescent. Once the trial step is determined, this function determines an approximate minimizer of the 1D function \(\phi_k(t) = f(x_k+ts_k)\). This approximate minimizer must satisfy sufficient decrease and curvature conditions.

Parameters
[out]sis the computed trial step
[in]xis the current iterate
[in]objis the objective function
[in]bndare the bound constraints
[in]algo_statecontains the current state of the algorithm

Reimplemented from ROL::Step< Real >.

Definition at line 318 of file ROL_LineSearchStep.hpp.

References ROL::LineSearchStep< Real >::acceptLastAlpha_, ROL::Vector< Real >::axpy(), ROL::LineSearchStep< Real >::desc_, ROL::LineSearchStep< Real >::fval_, ROL::Step< Real >::getState(), ROL::AlgorithmState< Real >::gnorm, ROL::LineSearchStep< Real >::GradDotStep(), ROL::BoundConstraint< Real >::isActivated(), ROL::LineSearchStep< Real >::lineSearch_, ROL::Vector< Real >::plus(), ROL::BoundConstraint< Real >::project(), ROL::Vector< Real >::scale(), ROL::Vector< Real >::set(), ROL::AlgorithmState< Real >::value, and zero.

Referenced by main().

template<class Real>
void ROL::LineSearchStep< Real >::update ( Vector< Real > &  x,
const Vector< Real > &  s,
Objective< Real > &  obj,
BoundConstraint< Real > &  bnd,
AlgorithmState< Real > &  algo_state 
)
inlinevirtual

Update step, if successful.

Given a trial step, \(s_k\), this function updates \(x_{k+1}=x_k+s_k\). This function also updates the secant approximation.

Parameters
[in,out]xis the updated iterate
[in]sis the computed trial step
[in]objis the objective function
[in]conare the bound constraints
[in]algo_statecontains the current state of the algorithm

Reimplemented from ROL::Step< Real >.

Definition at line 367 of file ROL_LineSearchStep.hpp.

References ROL::LineSearchStep< Real >::computeObj_, ROL::LineSearchStep< Real >::desc_, ROL::LineSearchStep< Real >::fval_, ROL::Step< Real >::getState(), ROL::AlgorithmState< Real >::nfval, ROL::AlgorithmState< Real >::ngrad, and ROL::AlgorithmState< Real >::value.

Referenced by main().

template<class Real>
std::string ROL::LineSearchStep< Real >::printHeader ( void  ) const
inlinevirtual

Print iterate header.

This function produces a string containing header information.

Reimplemented from ROL::Step< Real >.

Definition at line 386 of file ROL_LineSearchStep.hpp.

References ROL::LineSearchStep< Real >::desc_.

Referenced by ROL::LineSearchStep< Real >::print().

template<class Real>
std::string ROL::LineSearchStep< Real >::printName ( void  ) const
inlinevirtual

Print step name.

This function produces a string containing the algorithmic step information.

Reimplemented from ROL::Step< Real >.

Definition at line 401 of file ROL_LineSearchStep.hpp.

References ROL::LineSearchStep< Real >::desc_, ROL::LineSearchStep< Real >::econd_, ROL::ECurvatureConditionToString(), and ROL::LineSearchStep< Real >::lineSearchName_.

Referenced by ROL::LineSearchStep< Real >::print().

template<class Real>
std::string ROL::LineSearchStep< Real >::print ( AlgorithmState< Real > &  algo_state,
bool  print_header = false 
) const
inlinevirtual

Print iterate status.

This function prints the iteration status.

Parameters
[in]algo_stateis the current state of the algorithm
[in]printHeaderif set to true will print the header at each iteration

Reimplemented from ROL::Step< Real >.

Definition at line 417 of file ROL_LineSearchStep.hpp.

References ROL::LineSearchStep< Real >::desc_, ROL::Step< Real >::getStepState(), ROL::AlgorithmState< Real >::iter, ROL::LineSearchStep< Real >::printHeader(), and ROL::LineSearchStep< Real >::printName().

Member Data Documentation

template<class Real>
ROL::Ptr<Step<Real> > ROL::LineSearchStep< Real >::desc_
private
template<class Real>
ROL::Ptr<Secant<Real> > ROL::LineSearchStep< Real >::secant_
private

Secant object (used for quasi-Newton)

Definition at line 142 of file ROL_LineSearchStep.hpp.

Referenced by ROL::LineSearchStep< Real >::initialize().

template<class Real>
ROL::Ptr<Krylov<Real> > ROL::LineSearchStep< Real >::krylov_
private

Krylov solver object (used for inexact Newton)

Definition at line 143 of file ROL_LineSearchStep.hpp.

Referenced by ROL::LineSearchStep< Real >::initialize().

template<class Real>
ROL::Ptr<NonlinearCG<Real> > ROL::LineSearchStep< Real >::nlcg_
private

Nonlinear CG object (used for nonlinear CG)

Definition at line 144 of file ROL_LineSearchStep.hpp.

Referenced by ROL::LineSearchStep< Real >::initialize().

template<class Real>
ROL::Ptr<LineSearch<Real> > ROL::LineSearchStep< Real >::lineSearch_
private
template<class Real>
ROL::Ptr<Vector<Real> > ROL::LineSearchStep< Real >::d_
private
template<class Real>
ELineSearch ROL::LineSearchStep< Real >::els_
private

enum determines type of line search

Definition at line 149 of file ROL_LineSearchStep.hpp.

Referenced by ROL::LineSearchStep< Real >::LineSearchStep().

template<class Real>
ECurvatureCondition ROL::LineSearchStep< Real >::econd_
private

enum determines type of curvature condition

Definition at line 150 of file ROL_LineSearchStep.hpp.

Referenced by ROL::LineSearchStep< Real >::LineSearchStep(), and ROL::LineSearchStep< Real >::printName().

template<class Real>
bool ROL::LineSearchStep< Real >::acceptLastAlpha_
private

For backwards compatibility. When max function evaluations are reached take last step.

Definition at line 152 of file ROL_LineSearchStep.hpp.

Referenced by ROL::LineSearchStep< Real >::compute(), and ROL::LineSearchStep< Real >::LineSearchStep().

template<class Real>
bool ROL::LineSearchStep< Real >::usePreviousAlpha_
private

If true, use the previously accepted step length (if any) as the new initial step length.

Definition at line 154 of file ROL_LineSearchStep.hpp.

template<class Real>
int ROL::LineSearchStep< Real >::verbosity_
private

Definition at line 156 of file ROL_LineSearchStep.hpp.

Referenced by ROL::LineSearchStep< Real >::LineSearchStep().

template<class Real>
bool ROL::LineSearchStep< Real >::computeObj_
private
template<class Real>
Real ROL::LineSearchStep< Real >::fval_
private
template<class Real>
ROL::ParameterList ROL::LineSearchStep< Real >::parlist_
private

Definition at line 160 of file ROL_LineSearchStep.hpp.

Referenced by ROL::LineSearchStep< Real >::initialize().

template<class Real>
std::string ROL::LineSearchStep< Real >::lineSearchName_
private

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