|
| | LineSearchAlgorithm (ParameterList &parlist, const Ptr< Secant< Real > > &secant=nullPtr, const Ptr< DescentDirection_U< Real > > &descent=nullPtr, const Ptr< LineSearch_U< Real > > &lineSearch=nullPtr) |
| | Constructor.
|
| 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 |
| void | writeHeader (std::ostream &os) const override |
| | Print iterate header.
|
| void | writeName (std::ostream &os) const override |
| | Print step name.
|
| void | writeOutput (std::ostream &os, bool print_header=false) const override |
| | Print iterate status.
|
| virtual | ~Algorithm () |
| | Algorithm () |
| | Constructor, given a step and a status test.
|
| 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.
|
| 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.
|
| 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.
|
| 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.
|
| virtual void | run (Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, std::ostream &outStream=std::cout)=0 |
| | 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.
|
| virtual void | writeExitStatus (std::ostream &os) const |
| Ptr< const AlgorithmState< Real > > | getState () const |
| void | reset () |
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.