ROL
|
Provides the interface to compute optimization steps with trust regions. More...
#include <ROL_TrustRegionStep.hpp>
Public Member Functions | |
virtual | ~TrustRegionStep () |
TrustRegionStep (ROL::ParameterList &parlist) | |
Constructor. More... | |
TrustRegionStep (ROL::Ptr< Secant< Real > > &secant, ROL::ParameterList &parlist) | |
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. 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 | |
void | parseParameterList (ROL::ParameterList &parlist) |
Parse input ParameterList. More... | |
void | updateGradient (Vector< Real > &x, Objective< Real > &obj, BoundConstraint< Real > &bnd, AlgorithmState< Real > &algo_state) |
Update gradient to iteratively satisfy inexactness condition. More... | |
Real | computeCriticalityMeasure (const Vector< Real > &g, const Vector< Real > &x, BoundConstraint< Real > &bnd) |
Compute the criticality measure. More... | |
Private Attributes | |
Ptr< Vector< Real > > | xnew_ |
Container for updated iteration vector. More... | |
Ptr< Vector< Real > > | xold_ |
Container for previous iteration vector. More... | |
Ptr< Vector< Real > > | gp_ |
Container for previous gradient vector. More... | |
Ptr< TrustRegion< Real > > | trustRegion_ |
Container for trust-region solver object. More... | |
Ptr< TrustRegionModel< Real > > | model_ |
Container for trust-region model. More... | |
ETrustRegion | etr_ |
Trust-region subproblem solver type. More... | |
ETrustRegionModel | TRmodel_ |
Trust-region subproblem model type. More... | |
Real | delMax_ |
Maximum trust-region radius. More... | |
ETrustRegionFlag | TRflag_ |
Trust-region exit flag. More... | |
int | SPflag_ |
Subproblem solver termination flag. More... | |
int | SPiter_ |
Subproblem solver iteration count. More... | |
bool | bndActive_ |
Flag whether bound is activated. More... | |
Ptr< Secant< Real > > | secant_ |
Container for secant approximation. More... | |
ESecant | esec_ |
Secant type. More... | |
bool | useSecantHessVec_ |
Flag whether to use a secant Hessian. More... | |
bool | useSecantPrecond_ |
Flag whether to use a secant preconditioner. More... | |
Real | scaleEps_ |
Scaling for epsilon-active sets. More... | |
bool | useProjectedGrad_ |
Flag whether to use the projected gradient criticality measure. More... | |
Real | alpha_init_ |
Initial line-search parameter for projected methods. More... | |
int | max_fval_ |
Maximum function evaluations in line-search for projected methods. More... | |
Real | mu_ |
Post-Smoothing tolerance for projected methods. More... | |
Real | beta_ |
Post-Smoothing rate for projected methods. More... | |
Real | stepBackMax_ |
Real | stepBackScale_ |
bool | singleReflect_ |
std::vector< bool > | useInexact_ |
Flags for inexact (0) objective function, (1) gradient, (2) Hessian. More... | |
Real | scale0_ |
Scale for inexact gradient computation. More... | |
Real | scale1_ |
Scale for inexact gradient computation. More... | |
int | verbosity_ |
Print additional information to screen if > 0. More... | |
Additional Inherited Members | |
Protected Member Functions inherited from ROL::Step< Real > | |
ROL::Ptr< StepState< Real > > | getState (void) |
Provides the interface to compute optimization steps with trust regions.
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 trust-region 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\) the trial step \(s_k\) is computed by approximately solving the trust-region subproblem
\[ \min_{s} \frac{1}{2}\langle B_k s, s\rangle_{\mathcal{X}} + \langle g_k,s\rangle_{\mathcal{X}} \quad\text{s.t.}\quad \|s\|_{\mathcal{X}} \le \Delta_k \]
where \(B_k\in L(\mathcal{X},\mathcal{X})\), \(g_k\approx\nabla f(x_k)\), and \(\Delta_k > 0\). The approximate minimizer \(s_k\) must satisfy the fraction of Cauchy decrease condition
\[ -\frac{1}{2}\langle B_k s, s\rangle_{\mathcal{X}} - \langle g_k,s\rangle_{\mathcal{X}} \ge \kappa_0 \|g_k\|_{\mathcal{X}} \min\left\{\,\Delta_k,\, \frac{\|g_k\|_{\mathcal{X}}}{1+\|B_k\|_{L(\mathcal{X},\mathcal{X}})}\,\right\} \]
for some \(\kappa_0>0\) independent of \(k\). ROL's trust-region algorithm allows for both inexact objective function and gradient evaluation. The user must ensure that the inexact objective function, \(f_k\) satisfies
\[ |(f(x_k+s_k)-f_k(x_k+s_k)) - (f(x_k)-f_k(x_k))| \le \eta_1 \min\{\,-\frac{1}{2}\langle B_k s, s\rangle_{\mathcal{X}} - \langle g_k,s\rangle_{\mathcal{X}}, \,r_k\,\} \]
where \(\eta_1\) is the step acceptance threshold and \(r_k\) is a user-defined forcing sequence of positive numbers converging to zero. The inexact gradient, \(g_k\), must satisfy
\[ \|g_k-\nabla J(x_k)\|_{\mathcal{X}} \le \kappa_1\min\{\,\|g_k\|_{\mathcal{X}},\,\Delta_k\,\} \]
where \(\kappa_1 > 0\) is independent of \(k\).
For bound constrained problems, ROL employs projected Newton-type methods. For these methods, ROL requires 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 gradient approximation \(g_k\), we define the binding set as
\[ \mathcal{B}_k = \{\, \xi\in\Xi\,:\,x_k(\xi) = a(\xi) \;\text{and}\; -g_k(\xi) < 0 \,\}\cap \{\, \xi\in\Xi\,:\,x_k(\xi) = b(\xi) \;\text{and}\; -g_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, ROL prunes the variables in the binding set and runs a standard trust-region subproblem solver on the free variables. ROL then must perform a projected search to ensure the fraction of Cauchy decrease condition is satisfied.
TrustRegionStep implements a number of algorithms for both bound constrained and unconstrained optimization. These algorithms are: Cauchy Point, Dogleg, Double Dogleg, and Truncated CG. Each of these methods can be run using a secant approximation of the Hessian. These methods are chosen through the ETrustRegion enum.
Definition at line 95 of file ROL_TrustRegionStep.hpp.
|
inlinevirtual |
Definition at line 270 of file ROL_TrustRegionStep.hpp.
|
inline |
Constructor.
Standard constructor to build a TrustRegionStep object. Algorithmic specifications are passed in through a ROL::ParameterList.
[in] | parlist | is a parameter list containing algorithmic specifications |
Definition at line 279 of file ROL_TrustRegionStep.hpp.
References ROL::TrustRegionStep< Real >::esec_, ROL::TrustRegionStep< Real >::parseParameterList(), ROL::TrustRegionStep< Real >::secant_, ROL::StringToESecant(), ROL::TrustRegionStep< Real >::useSecantHessVec_, and ROL::TrustRegionStep< Real >::useSecantPrecond_.
|
inline |
Constructor.
Constructor to build a TrustRegionStep object with a user-defined secant object. Algorithmic specifications are passed in through a ROL::ParameterList.
[in] | secant | is a user-defined secant object |
[in] | parlist | is a parameter list containing algorithmic specifications |
Definition at line 312 of file ROL_TrustRegionStep.hpp.
References ROL::TrustRegionStep< Real >::parseParameterList(), ROL::TrustRegionStep< Real >::secant_, ROL::TrustRegionStep< Real >::useSecantHessVec_, and ROL::TrustRegionStep< Real >::useSecantPrecond_.
|
inlineprivate |
Parse input ParameterList.
This function sets trust region specific parameters specified in the user supplied ParameterList.
[in] | parlist | is the user-supplied ParameterList. |
Definition at line 149 of file ROL_TrustRegionStep.hpp.
References ROL::TrustRegionStep< Real >::alpha_init_, ROL::TrustRegionStep< Real >::beta_, ROL::TrustRegionStep< Real >::delMax_, ROL::TrustRegionStep< Real >::etr_, ROL::Step< Real >::getState(), ROL::TrustRegionStep< Real >::max_fval_, ROL::TrustRegionStep< Real >::mu_, ROL::TrustRegionStep< Real >::scale0_, ROL::TrustRegionStep< Real >::scale1_, ROL::TrustRegionStep< Real >::scaleEps_, ROL::TrustRegionStep< Real >::singleReflect_, ROL::TrustRegionStep< Real >::stepBackMax_, ROL::TrustRegionStep< Real >::stepBackScale_, ROL::StringToETrustRegion(), ROL::StringToETrustRegionModel(), ROL::TrustRegionStep< Real >::TRmodel_, ROL::TrustRegionStep< Real >::trustRegion_, ROL::TrustRegionStep< Real >::useInexact_, ROL::TrustRegionStep< Real >::useProjectedGrad_, and ROL::TrustRegionStep< Real >::verbosity_.
Referenced by ROL::TrustRegionStep< Real >::TrustRegionStep().
|
inlineprivate |
Update gradient to iteratively satisfy inexactness condition.
This function attempts to ensure that the inexact gradient condition,
\[ \|g_k-\nabla J(x_k)\|_{\mathcal{X}} \le \kappa_1\min\{\,\|g_k\|_{\mathcal{X}},\,\Delta_k\,\}, \]
is satisfied. This function works under the assumption that the gradient function returns a gradient approximation which satisfies the error tolerance prescribed by the tol input parameter.
[in] | x | is the current optimization variable. |
[in] | obj | is the objective function. |
[in] | bnd | is the bound constraint. |
[in,out] | algo_state | is the algorithm state. |
Definition at line 199 of file ROL_TrustRegionStep.hpp.
References ROL::TrustRegionStep< Real >::computeCriticalityMeasure(), ROL::Step< Real >::getState(), ROL::AlgorithmState< Real >::gnorm, ROL::Objective< Real >::gradient(), ROL::AlgorithmState< Real >::ngrad, ROL::TrustRegionStep< Real >::scale0_, and ROL::TrustRegionStep< Real >::useInexact_.
Referenced by ROL::TrustRegionStep< Real >::initialize(), and ROL::TrustRegionStep< Real >::update().
|
inlineprivate |
Compute the criticality measure.
This function computes either the norm of the gradient projected onto the tangent cone or the norm of \(x_k - P_{[a,b]}(x_k-g_k)\).
[in] | g | is the current gradient. |
[in] | x | is the current iterate. |
[in] | bnd | is the bound constraint. |
Definition at line 243 of file ROL_TrustRegionStep.hpp.
References ROL::BoundConstraint< Real >::computeProjectedGradient(), ROL::Vector< Real >::dual(), ROL::TrustRegionStep< Real >::gp_, ROL::BoundConstraint< Real >::isActivated(), ROL::Vector< Real >::norm(), ROL::BoundConstraint< Real >::project(), ROL::TrustRegionStep< Real >::useProjectedGrad_, and ROL::TrustRegionStep< Real >::xnew_.
Referenced by ROL::TrustRegionStep< Real >::updateGradient().
|
inlinevirtual |
Initialize step.
This function initializes the information necessary to run the trust-region algorithm.
[in] | x | is the initial guess for the optimization vector. |
[in] | obj | is the objective function. |
[in] | bnd | is the bound constraint. |
[in] | algo_state | is the algorithm state. |
Reimplemented from ROL::Step< Real >.
Definition at line 348 of file ROL_TrustRegionStep.hpp.
References ROL::TrustRegionStep< Real >::bndActive_, ROL::Vector< Real >::clone(), ROL::TrustRegionStep< Real >::delMax_, ROL::TrustRegionStep< Real >::etr_, ROL::Step< Real >::getState(), ROL::AlgorithmState< Real >::gnorm, ROL::TrustRegionStep< Real >::gp_, ROL::Objective< Real >::hessVec(), ROL::Objective< Real >::invHessVec(), ROL::BoundConstraint< Real >::isActivated(), ROL::isValidTrustRegionSubproblem(), ROL::AlgorithmState< Real >::iter, ROL::TrustRegionStep< Real >::model_, ROL::AlgorithmState< Real >::nfval, ROL::BoundConstraint< Real >::project(), ROL::BoundConstraint< Real >::projectInterior(), ROL::TrustRegionStep< Real >::secant_, ROL::TrustRegionStep< Real >::singleReflect_, ROL::AlgorithmState< Real >::snorm, ROL::TrustRegionStep< Real >::stepBackMax_, ROL::TrustRegionStep< Real >::stepBackScale_, ROL::TrustRegionStep< Real >::TRmodel_, ROL::TrustRegionStep< Real >::trustRegion_, ROL::TRUSTREGION_DOGLEG, ROL::TRUSTREGION_DOUBLEDOGLEG, ROL::TRUSTREGION_MODEL_COLEMANLI, ROL::TRUSTREGION_MODEL_KELLEYSACHS, ROL::TRUSTREGION_MODEL_LINMORE, ROL::Objective< Real >::update(), ROL::TrustRegionStep< Real >::updateGradient(), ROL::TrustRegionStep< Real >::useSecantHessVec_, ROL::TrustRegionStep< Real >::useSecantPrecond_, ROL::Objective< Real >::value(), ROL::AlgorithmState< Real >::value, ROL::TrustRegionStep< Real >::xnew_, ROL::TrustRegionStep< Real >::xold_, and zero.
Referenced by main().
|
inlinevirtual |
Compute step.
Computes a trial step, \(s_k\) by solving the trust-region subproblem. The trust-region subproblem solver is defined by the enum ETrustRegion.
[out] | s | is the computed trial step |
[in] | x | is the current iterate |
[in] | obj | is the objective function |
[in] | bnd | are the bound constraints |
[in] | algo_state | contains the current state of the algorithm |
Reimplemented from ROL::Step< Real >.
Definition at line 519 of file ROL_TrustRegionStep.hpp.
References ROL::Step< Real >::getState(), ROL::AlgorithmState< Real >::gnorm, ROL::BoundConstraint< Real >::isActivated(), ROL::TrustRegionStep< Real >::model_, ROL::TrustRegionStep< Real >::scaleEps_, ROL::TrustRegionStep< Real >::secant_, ROL::AlgorithmState< Real >::snorm, ROL::TrustRegionStep< Real >::SPflag_, ROL::TrustRegionStep< Real >::SPiter_, ROL::TrustRegionStep< Real >::TRmodel_, ROL::TrustRegionStep< Real >::trustRegion_, ROL::TRUSTREGION_MODEL_COLEMANLI, and ROL::TRUSTREGION_MODEL_KELLEYSACHS.
Referenced by main().
|
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.
[in,out] | x | is the updated iterate |
[in] | s | is the computed trial step |
[in] | obj | is the objective function |
[in] | bnd | are the bound constraints |
[in] | algo_state | contains the current state of the algorithm |
Reimplemented from ROL::Step< Real >.
Definition at line 553 of file ROL_TrustRegionStep.hpp.
References ROL::Step< Real >::getState(), ROL::TrustRegionStep< Real >::gp_, ROL::BoundConstraint< Real >::isActivated(), ROL::AlgorithmState< Real >::iter, ROL::AlgorithmState< Real >::iterateVec, ROL::TrustRegionStep< Real >::model_, ROL::AlgorithmState< Real >::nfval, ROL::AlgorithmState< Real >::ngrad, ROL::TrustRegionStep< Real >::secant_, ROL::AlgorithmState< Real >::snorm, ROL::TrustRegionStep< Real >::SPflag_, ROL::TrustRegionStep< Real >::SPiter_, ROL::TrustRegionStep< Real >::TRflag_, ROL::TrustRegionStep< Real >::trustRegion_, ROL::TRUSTREGION_FLAG_POSPREDNEG, ROL::TRUSTREGION_FLAG_SUCCESS, ROL::TrustRegionStep< Real >::updateGradient(), ROL::TrustRegionStep< Real >::useInexact_, ROL::TrustRegionStep< Real >::useSecantHessVec_, ROL::TrustRegionStep< Real >::useSecantPrecond_, ROL::AlgorithmState< Real >::value, ROL::TrustRegionStep< Real >::xnew_, and ROL::TrustRegionStep< Real >::xold_.
Referenced by main().
|
inlinevirtual |
Print iterate header.
This function produces a string containing header information.
Reimplemented from ROL::Step< Real >.
Definition at line 618 of file ROL_TrustRegionStep.hpp.
References ROL::CG_FLAG_SUCCESS, ROL::CG_FLAG_UNDEFINED, ROL::ECGFlagToString(), ROL::TrustRegionStep< Real >::etr_, ROL::ETrustRegionFlagToString(), ROL::NumberToString(), ROL::TRUSTREGION_FLAG_SUCCESS, ROL::TRUSTREGION_FLAG_UNDEFINED, ROL::TRUSTREGION_LINMORE, ROL::TRUSTREGION_TRUNCATEDCG, and ROL::TrustRegionStep< Real >::verbosity_.
Referenced by ROL::TrustRegionStep< Real >::print().
|
inlinevirtual |
Print step name.
This function produces a string containing the algorithmic step information.
Reimplemented from ROL::Step< Real >.
Definition at line 678 of file ROL_TrustRegionStep.hpp.
References ROL::TrustRegionStep< Real >::bndActive_, ROL::TrustRegionStep< Real >::esec_, ROL::ESecantToString(), ROL::TrustRegionStep< Real >::etr_, ROL::ETrustRegionModelToString(), ROL::ETrustRegionToString(), ROL::TrustRegionStep< Real >::TRmodel_, ROL::TrustRegionStep< Real >::useSecantHessVec_, and ROL::TrustRegionStep< Real >::useSecantPrecond_.
Referenced by ROL::TrustRegionStep< Real >::print().
|
inlinevirtual |
Print iterate status.
This function prints the iteration status.
[in] | algo_state | is the current state of the algorithm |
[in] | printHeader | if ste to true will print the header at each iteration |
Reimplemented from ROL::Step< Real >.
Definition at line 708 of file ROL_TrustRegionStep.hpp.
References ROL::TrustRegionStep< Real >::etr_, ROL::Step< Real >::getStepState(), ROL::AlgorithmState< Real >::gnorm, ROL::AlgorithmState< Real >::iter, ROL::AlgorithmState< Real >::nfval, ROL::AlgorithmState< Real >::ngrad, ROL::TrustRegionStep< Real >::printHeader(), ROL::TrustRegionStep< Real >::printName(), ROL::AlgorithmState< Real >::snorm, ROL::TrustRegionStep< Real >::SPflag_, ROL::TrustRegionStep< Real >::SPiter_, ROL::TrustRegionStep< Real >::TRflag_, ROL::TRUSTREGION_LINMORE, ROL::TRUSTREGION_TRUNCATEDCG, and ROL::AlgorithmState< Real >::value.
|
private |
Container for updated iteration vector.
Definition at line 99 of file ROL_TrustRegionStep.hpp.
Referenced by ROL::TrustRegionStep< Real >::computeCriticalityMeasure(), ROL::TrustRegionStep< Real >::initialize(), and ROL::TrustRegionStep< Real >::update().
|
private |
Container for previous iteration vector.
Definition at line 100 of file ROL_TrustRegionStep.hpp.
Referenced by ROL::TrustRegionStep< Real >::initialize(), and ROL::TrustRegionStep< Real >::update().
|
private |
Container for previous gradient vector.
Definition at line 101 of file ROL_TrustRegionStep.hpp.
Referenced by ROL::TrustRegionStep< Real >::computeCriticalityMeasure(), ROL::TrustRegionStep< Real >::initialize(), and ROL::TrustRegionStep< Real >::update().
|
private |
Container for trust-region solver object.
Definition at line 104 of file ROL_TrustRegionStep.hpp.
Referenced by ROL::TrustRegionStep< Real >::compute(), ROL::TrustRegionStep< Real >::initialize(), ROL::TrustRegionStep< Real >::parseParameterList(), and ROL::TrustRegionStep< Real >::update().
|
private |
Container for trust-region model.
Definition at line 105 of file ROL_TrustRegionStep.hpp.
Referenced by ROL::TrustRegionStep< Real >::compute(), ROL::TrustRegionStep< Real >::initialize(), and ROL::TrustRegionStep< Real >::update().
|
private |
Trust-region subproblem solver type.
Definition at line 106 of file ROL_TrustRegionStep.hpp.
Referenced by ROL::TrustRegionStep< Real >::initialize(), ROL::TrustRegionStep< Real >::parseParameterList(), ROL::TrustRegionStep< Real >::print(), ROL::TrustRegionStep< Real >::printHeader(), and ROL::TrustRegionStep< Real >::printName().
|
private |
Trust-region subproblem model type.
Definition at line 107 of file ROL_TrustRegionStep.hpp.
Referenced by ROL::TrustRegionStep< Real >::compute(), ROL::TrustRegionStep< Real >::initialize(), ROL::TrustRegionStep< Real >::parseParameterList(), and ROL::TrustRegionStep< Real >::printName().
|
private |
Maximum trust-region radius.
Definition at line 108 of file ROL_TrustRegionStep.hpp.
Referenced by ROL::TrustRegionStep< Real >::initialize(), and ROL::TrustRegionStep< Real >::parseParameterList().
|
private |
Trust-region exit flag.
Definition at line 109 of file ROL_TrustRegionStep.hpp.
Referenced by ROL::TrustRegionStep< Real >::print(), and ROL::TrustRegionStep< Real >::update().
|
private |
Subproblem solver termination flag.
Definition at line 110 of file ROL_TrustRegionStep.hpp.
Referenced by ROL::TrustRegionStep< Real >::compute(), ROL::TrustRegionStep< Real >::print(), and ROL::TrustRegionStep< Real >::update().
|
private |
Subproblem solver iteration count.
Definition at line 111 of file ROL_TrustRegionStep.hpp.
Referenced by ROL::TrustRegionStep< Real >::compute(), ROL::TrustRegionStep< Real >::print(), and ROL::TrustRegionStep< Real >::update().
|
private |
Flag whether bound is activated.
Definition at line 112 of file ROL_TrustRegionStep.hpp.
Referenced by ROL::TrustRegionStep< Real >::initialize(), and ROL::TrustRegionStep< Real >::printName().
|
private |
Container for secant approximation.
Definition at line 115 of file ROL_TrustRegionStep.hpp.
Referenced by ROL::TrustRegionStep< Real >::compute(), ROL::TrustRegionStep< Real >::initialize(), ROL::TrustRegionStep< Real >::TrustRegionStep(), and ROL::TrustRegionStep< Real >::update().
|
private |
Secant type.
Definition at line 116 of file ROL_TrustRegionStep.hpp.
Referenced by ROL::TrustRegionStep< Real >::printName(), and ROL::TrustRegionStep< Real >::TrustRegionStep().
|
private |
Flag whether to use a secant Hessian.
Definition at line 117 of file ROL_TrustRegionStep.hpp.
Referenced by ROL::TrustRegionStep< Real >::initialize(), ROL::TrustRegionStep< Real >::printName(), ROL::TrustRegionStep< Real >::TrustRegionStep(), and ROL::TrustRegionStep< Real >::update().
|
private |
Flag whether to use a secant preconditioner.
Definition at line 118 of file ROL_TrustRegionStep.hpp.
Referenced by ROL::TrustRegionStep< Real >::initialize(), ROL::TrustRegionStep< Real >::printName(), ROL::TrustRegionStep< Real >::TrustRegionStep(), and ROL::TrustRegionStep< Real >::update().
|
private |
Scaling for epsilon-active sets.
Definition at line 121 of file ROL_TrustRegionStep.hpp.
Referenced by ROL::TrustRegionStep< Real >::compute(), and ROL::TrustRegionStep< Real >::parseParameterList().
|
private |
Flag whether to use the projected gradient criticality measure.
Definition at line 122 of file ROL_TrustRegionStep.hpp.
Referenced by ROL::TrustRegionStep< Real >::computeCriticalityMeasure(), and ROL::TrustRegionStep< Real >::parseParameterList().
|
private |
Initial line-search parameter for projected methods.
Definition at line 125 of file ROL_TrustRegionStep.hpp.
Referenced by ROL::TrustRegionStep< Real >::parseParameterList().
|
private |
Maximum function evaluations in line-search for projected methods.
Definition at line 126 of file ROL_TrustRegionStep.hpp.
Referenced by ROL::TrustRegionStep< Real >::parseParameterList().
|
private |
Post-Smoothing tolerance for projected methods.
Definition at line 127 of file ROL_TrustRegionStep.hpp.
Referenced by ROL::TrustRegionStep< Real >::parseParameterList().
|
private |
Post-Smoothing rate for projected methods.
Definition at line 128 of file ROL_TrustRegionStep.hpp.
Referenced by ROL::TrustRegionStep< Real >::parseParameterList().
|
private |
Definition at line 131 of file ROL_TrustRegionStep.hpp.
Referenced by ROL::TrustRegionStep< Real >::initialize(), and ROL::TrustRegionStep< Real >::parseParameterList().
|
private |
Definition at line 132 of file ROL_TrustRegionStep.hpp.
Referenced by ROL::TrustRegionStep< Real >::initialize(), and ROL::TrustRegionStep< Real >::parseParameterList().
|
private |
Definition at line 133 of file ROL_TrustRegionStep.hpp.
Referenced by ROL::TrustRegionStep< Real >::initialize(), and ROL::TrustRegionStep< Real >::parseParameterList().
|
private |
Flags for inexact (0) objective function, (1) gradient, (2) Hessian.
Definition at line 136 of file ROL_TrustRegionStep.hpp.
Referenced by ROL::TrustRegionStep< Real >::parseParameterList(), ROL::TrustRegionStep< Real >::update(), and ROL::TrustRegionStep< Real >::updateGradient().
|
private |
Scale for inexact gradient computation.
Definition at line 137 of file ROL_TrustRegionStep.hpp.
Referenced by ROL::TrustRegionStep< Real >::parseParameterList(), and ROL::TrustRegionStep< Real >::updateGradient().
|
private |
Scale for inexact gradient computation.
Definition at line 138 of file ROL_TrustRegionStep.hpp.
Referenced by ROL::TrustRegionStep< Real >::parseParameterList().
|
private |
Print additional information to screen if > 0.
Definition at line 141 of file ROL_TrustRegionStep.hpp.
Referenced by ROL::TrustRegionStep< Real >::parseParameterList(), and ROL::TrustRegionStep< Real >::printHeader().