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

Provides an interface to run unconstrained line search algorithms. More...

#include <ROL_TypeU_LineSearchAlgorithm.hpp>

+ Inheritance diagram for ROL::TypeU::LineSearchAlgorithm< Real >:

Public Member Functions

 LineSearchAlgorithm (ParameterList &parlist, const Ptr< Secant< Real >> &secant=nullPtr, const Ptr< DescentDirection_U< Real >> &descent=nullPtr, const Ptr< LineSearch_U< Real >> &lineSearch=nullPtr)
 Constructor. More...
 
void initialize (const Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, std::ostream &outStream=std::cout)
 
void run (Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, std::ostream &outStream=std::cout) override
 Run algorithm on unconstrained problems (Type-U). This general interface supports the use of dual optimization vector spaces, where the user does not define the dual() method. More...
 
void writeHeader (std::ostream &os) const override
 Print iterate header. More...
 
void writeName (std::ostream &os) const override
 Print step name. More...
 
void writeOutput (std::ostream &os, const bool print_header=false) const override
 Print iterate status. More...
 
- Public Member Functions inherited from ROL::TypeU::Algorithm< Real >
virtual ~Algorithm ()
 
 Algorithm ()
 Constructor, given a step and a status test. More...
 
void setStatusTest (const Ptr< StatusTest< Real >> &status, bool combineStatus=false)
 
virtual void run (Problem< Real > &problem, std::ostream &outStream=std::cout)
 Run algorithm on unconstrained problems (Type-U). This is the primary Type-U interface. More...
 
virtual void run (Vector< Real > &x, Objective< Real > &obj, std::ostream &outStream=std::cout)
 Run algorithm on unconstrained problems (Type-U). This is the primary Type-U interface. More...
 
virtual void run (Vector< Real > &x, Objective< Real > &obj, Constraint< Real > &linear_con, Vector< Real > &linear_mul, std::ostream &outStream=std::cout)
 Run algorithm on unconstrained problems with explicit linear equality constraints (Type-U). This is the primary Type-U with explicit linear equality constraints interface. More...
 
virtual void run (Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, Constraint< Real > &linear_con, Vector< Real > &linear_mul, const Vector< Real > &linear_c, std::ostream &outStream=std::cout)
 Run algorithm on unconstrained problems with explicit linear equality constraints (Type-U). This general interface supports the use of dual optimization vector spaces, where the user does not define the dual() method. More...
 
virtual void writeExitStatus (std::ostream &os) const
 
Ptr< const AlgorithmState< Real > > getState () const
 
void reset ()
 

Private Attributes

Ptr< DescentDirection_U< Real > > desc_
 Unglobalized step object. More...
 
Ptr< LineSearch_U< Real > > lineSearch_
 Line-search object. More...
 
EDescentU edesc_
 enum determines type of descent direction More...
 
ELineSearchU els_
 enum determines type of line search More...
 
ECurvatureConditionU 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 printHeader_
 
int ls_nfval_
 
int ls_ngrad_
 
int SPflag_
 
int SPiter_
 
std::string lineSearchName_
 
std::string descentName_
 

Additional Inherited Members

- Protected Member Functions inherited from ROL::TypeU::Algorithm< Real >
void initialize (const Vector< Real > &x, const Vector< Real > &g)
 
- Protected Attributes inherited from ROL::TypeU::Algorithm< Real >
const Ptr< CombinedStatusTest
< Real > > 
status_
 
const Ptr< AlgorithmState< Real > > state_
 

Detailed Description

template<class Real>
class ROL::TypeU::LineSearchAlgorithm< Real >

Provides an interface to run unconstrained line search algorithms.

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 optimization problems,

\[ \min_x\quad f(x). \]

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}}\).

LineSearchStep implements a number of algorithms for 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 101 of file ROL_TypeU_LineSearchAlgorithm.hpp.

Constructor & Destructor Documentation

template<typename Real >
ROL::TypeU::LineSearchAlgorithm< Real >::LineSearchAlgorithm ( ParameterList &  parlist,
const Ptr< Secant< Real >> &  secant = nullPtr,
const Ptr< DescentDirection_U< Real >> &  descent = nullPtr,
const Ptr< LineSearch_U< Real >> &  lineSearch = nullPtr 
)

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 descent direction object
[in]lineSearchis a user-defined line search object

Definition at line 54 of file ROL_TypeU_LineSearchAlgorithm_Def.hpp.

References ROL::TypeU::LineSearchAlgorithm< Real >::acceptLastAlpha_, ROL::TypeU::LineSearchAlgorithm< Real >::desc_, ROL::TypeU::LineSearchAlgorithm< Real >::descentName_, ROL::TypeU::LineSearchAlgorithm< Real >::econd_, ROL::TypeU::LineSearchAlgorithm< Real >::edesc_, ROL::TypeU::LineSearchAlgorithm< Real >::els_, ROL::TypeU::LineSearchAlgorithm< Real >::lineSearch_, ROL::TypeU::LineSearchAlgorithm< Real >::lineSearchName_, ROL::TypeU::LineSearchAlgorithm< Real >::printHeader_, ROL::TypeU::Algorithm< Real >::status_, ROL::StringToECurvatureConditionU(), ROL::StringToEDescentU(), ROL::StringToELineSearchU(), and ROL::TypeU::LineSearchAlgorithm< Real >::verbosity_.

Member Function Documentation

template<typename Real >
void ROL::TypeU::LineSearchAlgorithm< Real >::initialize ( const Vector< Real > &  x,
const Vector< Real > &  g,
Objective< Real > &  obj,
std::ostream &  outStream = std::cout 
)
template<typename Real >
void ROL::TypeU::LineSearchAlgorithm< Real >::run ( Vector< Real > &  x,
const Vector< Real > &  g,
Objective< Real > &  obj,
std::ostream &  outStream = std::cout 
)
overridevirtual

Run algorithm on unconstrained problems (Type-U). This general interface supports the use of dual optimization vector spaces, where the user does not define the dual() method.

Implements ROL::TypeU::Algorithm< Real >.

Definition at line 113 of file ROL_TypeU_LineSearchAlgorithm_Def.hpp.

References ROL::Accept, ROL::Vector< Real >::clone(), ROL::Objective< Real >::gradient(), ROL::Vector< Real >::plus(), ROL::Objective< Real >::update(), ROL::Objective< Real >::value(), and ROL::TypeU::Algorithm< Real >::writeExitStatus().

template<typename Real >
void ROL::TypeU::LineSearchAlgorithm< Real >::writeHeader ( std::ostream &  os) const
overridevirtual

Print iterate header.

Reimplemented from ROL::TypeU::Algorithm< Real >.

Definition at line 172 of file ROL_TypeU_LineSearchAlgorithm_Def.hpp.

References ROL::DESCENT_U_NEWTONKRYLOV.

template<typename Real >
void ROL::TypeU::LineSearchAlgorithm< Real >::writeName ( std::ostream &  os) const
overridevirtual

Print step name.

Reimplemented from ROL::TypeU::Algorithm< Real >.

Definition at line 213 of file ROL_TypeU_LineSearchAlgorithm_Def.hpp.

References ROL::ECurvatureConditionUToString().

template<typename Real >
void ROL::TypeU::LineSearchAlgorithm< Real >::writeOutput ( std::ostream &  os,
const bool  write_header = false 
) const
overridevirtual

Print iterate status.

Reimplemented from ROL::TypeU::Algorithm< Real >.

Definition at line 223 of file ROL_TypeU_LineSearchAlgorithm_Def.hpp.

References ROL::DESCENT_U_NEWTONKRYLOV.

Member Data Documentation

template<class Real >
Ptr<DescentDirection_U<Real> > ROL::TypeU::LineSearchAlgorithm< Real >::desc_
private

Unglobalized step object.

Definition at line 104 of file ROL_TypeU_LineSearchAlgorithm.hpp.

Referenced by ROL::TypeU::LineSearchAlgorithm< Real >::LineSearchAlgorithm().

template<class Real >
Ptr<LineSearch_U<Real> > ROL::TypeU::LineSearchAlgorithm< Real >::lineSearch_
private

Line-search object.

Definition at line 105 of file ROL_TypeU_LineSearchAlgorithm.hpp.

Referenced by ROL::TypeU::LineSearchAlgorithm< Real >::LineSearchAlgorithm().

template<class Real >
EDescentU ROL::TypeU::LineSearchAlgorithm< Real >::edesc_
private

enum determines type of descent direction

Definition at line 107 of file ROL_TypeU_LineSearchAlgorithm.hpp.

Referenced by ROL::TypeU::LineSearchAlgorithm< Real >::LineSearchAlgorithm().

template<class Real >
ELineSearchU ROL::TypeU::LineSearchAlgorithm< Real >::els_
private

enum determines type of line search

Definition at line 108 of file ROL_TypeU_LineSearchAlgorithm.hpp.

Referenced by ROL::TypeU::LineSearchAlgorithm< Real >::LineSearchAlgorithm().

template<class Real >
ECurvatureConditionU ROL::TypeU::LineSearchAlgorithm< Real >::econd_
private

enum determines type of curvature condition

Definition at line 109 of file ROL_TypeU_LineSearchAlgorithm.hpp.

Referenced by ROL::TypeU::LineSearchAlgorithm< Real >::LineSearchAlgorithm().

template<class Real >
bool ROL::TypeU::LineSearchAlgorithm< Real >::acceptLastAlpha_
private

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

Definition at line 111 of file ROL_TypeU_LineSearchAlgorithm.hpp.

Referenced by ROL::TypeU::LineSearchAlgorithm< Real >::LineSearchAlgorithm().

template<class Real >
bool ROL::TypeU::LineSearchAlgorithm< Real >::usePreviousAlpha_
private

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

Definition at line 113 of file ROL_TypeU_LineSearchAlgorithm.hpp.

template<class Real >
int ROL::TypeU::LineSearchAlgorithm< Real >::verbosity_
private
template<class Real >
bool ROL::TypeU::LineSearchAlgorithm< Real >::printHeader_
private
template<class Real >
int ROL::TypeU::LineSearchAlgorithm< Real >::ls_nfval_
private

Definition at line 117 of file ROL_TypeU_LineSearchAlgorithm.hpp.

template<class Real >
int ROL::TypeU::LineSearchAlgorithm< Real >::ls_ngrad_
private

Definition at line 117 of file ROL_TypeU_LineSearchAlgorithm.hpp.

template<class Real >
int ROL::TypeU::LineSearchAlgorithm< Real >::SPflag_
private

Definition at line 118 of file ROL_TypeU_LineSearchAlgorithm.hpp.

template<class Real >
int ROL::TypeU::LineSearchAlgorithm< Real >::SPiter_
private

Definition at line 118 of file ROL_TypeU_LineSearchAlgorithm.hpp.

template<class Real >
std::string ROL::TypeU::LineSearchAlgorithm< Real >::lineSearchName_
private
template<class Real >
std::string ROL::TypeU::LineSearchAlgorithm< Real >::descentName_
private

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