ROL
|
Implements the computation of optimization steps with composite-step trust-region methods. More...
#include <ROL_CompositeStep.hpp>
Public Member Functions | |
virtual | ~CompositeStep () |
CompositeStep (ROL::ParameterList &parlist) | |
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. More... | |
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. More... | |
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. More... | |
void | compute (Vector< Real > &s, const Vector< Real > &x, Objective< Real > &obj, BoundConstraint< Real > &con, AlgorithmState< Real > &algo_state) |
Compute step for bound constraints; here only to satisfy the interface requirements, does nothing, needs refactoring. More... | |
void | update (Vector< Real > &x, const Vector< Real > &s, Objective< Real > &obj, BoundConstraint< Real > &con, AlgorithmState< Real > &algo_state) |
Update step, for bound constraints; here only to satisfy the interface requirements, does nothing, needs refactoring. 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 pHeader=false) const |
Print iterate status. More... | |
void | computeLagrangeMultiplier (Vector< Real > &l, const Vector< Real > &x, const Vector< Real > &gf, Constraint< Real > &con) |
Compute Lagrange multipliers by solving the least-squares problem minimizing the gradient of the Lagrangian, via the augmented system formulation. More... | |
void | computeQuasinormalStep (Vector< Real > &n, const Vector< Real > &c, const Vector< Real > &x, Real delta, Constraint< Real > &con) |
Compute quasi-normal step by minimizing the norm of the linearized constraint. More... | |
void | solveTangentialSubproblem (Vector< Real > &t, Vector< Real > &tCP, Vector< Real > &Wg, const Vector< Real > &x, const Vector< Real > &g, const Vector< Real > &n, const Vector< Real > &l, Real delta, Objective< Real > &obj, Constraint< Real > &con) |
Solve tangential subproblem. More... | |
void | accept (Vector< Real > &s, Vector< Real > &n, Vector< Real > &t, Real f_new, Vector< Real > &c_new, Vector< Real > &gf_new, Vector< Real > &l_new, Vector< Real > &g_new, const Vector< Real > &x, const Vector< Real > &l, Real f, const Vector< Real > &gf, const Vector< Real > &c, const Vector< Real > &g, Vector< Real > &tCP, Vector< Real > &Wg, Objective< Real > &obj, Constraint< Real > &con, AlgorithmState< Real > &algo_state) |
Check acceptance of subproblem solutions, adjust merit function penalty parameter, ensure global convergence. 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 > &s, 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, 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, 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 | |
template<typename T > | |
int | sgn (T val) |
void | printInfoLS (const std::vector< Real > &res) const |
Real | setTolOSS (const Real intol) const |
Private Attributes | |
ROL::Ptr< Vector< Real > > | xvec_ |
ROL::Ptr< Vector< Real > > | gvec_ |
ROL::Ptr< Vector< Real > > | cvec_ |
ROL::Ptr< Vector< Real > > | lvec_ |
int | flagCG_ |
int | flagAC_ |
int | iterCG_ |
int | maxiterCG_ |
int | maxiterOSS_ |
Real | tolCG_ |
Real | tolOSS_ |
bool | tolOSSfixed_ |
Real | lmhtol_ |
Real | qntol_ |
Real | pgtol_ |
Real | projtol_ |
Real | tangtol_ |
Real | tntmax_ |
Real | zeta_ |
Real | Delta_ |
Real | penalty_ |
Real | eta_ |
bool | useConHess_ |
Real | ared_ |
Real | pred_ |
Real | snorm_ |
Real | nnorm_ |
Real | tnorm_ |
bool | infoQN_ |
bool | infoLM_ |
bool | infoTS_ |
bool | infoAC_ |
bool | infoLS_ |
bool | infoALL_ |
int | totalIterCG_ |
int | totalProj_ |
int | totalNegCurv_ |
int | totalRef_ |
int | totalCallLS_ |
int | totalIterLS_ |
Additional Inherited Members | |
Protected Member Functions inherited from ROL::Step< Real > | |
ROL::Ptr< StepState< Real > > | getState (void) |
Implements the computation of optimization steps with composite-step trust-region methods.
Definition at line 29 of file ROL_CompositeStep.hpp.
|
inlinevirtual |
Definition at line 114 of file ROL_CompositeStep.hpp.
|
inline |
Definition at line 116 of file ROL_CompositeStep.hpp.
References ROL::CompositeStep< Real >::Delta_, ROL::CompositeStep< Real >::eta_, ROL::CompositeStep< Real >::flagAC_, ROL::CompositeStep< Real >::flagCG_, ROL::CompositeStep< Real >::infoAC_, ROL::CompositeStep< Real >::infoALL_, ROL::CompositeStep< Real >::infoLM_, ROL::CompositeStep< Real >::infoLS_, ROL::CompositeStep< Real >::infoQN_, ROL::CompositeStep< Real >::infoTS_, ROL::CompositeStep< Real >::iterCG_, ROL::CompositeStep< Real >::lmhtol_, ROL::CompositeStep< Real >::maxiterCG_, ROL::CompositeStep< Real >::nnorm_, ROL::CompositeStep< Real >::penalty_, ROL::CompositeStep< Real >::pgtol_, ROL::CompositeStep< Real >::projtol_, ROL::CompositeStep< Real >::qntol_, ROL::CompositeStep< Real >::snorm_, ROL::CompositeStep< Real >::tangtol_, ROL::CompositeStep< Real >::tnorm_, ROL::CompositeStep< Real >::tntmax_, ROL::CompositeStep< Real >::tolCG_, ROL::CompositeStep< Real >::tolOSS_, ROL::CompositeStep< Real >::tolOSSfixed_, ROL::CompositeStep< Real >::totalCallLS_, ROL::CompositeStep< Real >::totalIterCG_, ROL::CompositeStep< Real >::totalIterLS_, ROL::CompositeStep< Real >::totalNegCurv_, ROL::CompositeStep< Real >::totalProj_, ROL::CompositeStep< Real >::totalRef_, ROL::CompositeStep< Real >::useConHess_, and ROL::CompositeStep< Real >::zeta_.
|
inlineprivate |
Definition at line 87 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::solveTangentialSubproblem().
|
inlineprivate |
Definition at line 91 of file ROL_CompositeStep.hpp.
References ROL::CompositeStep< Real >::infoLS_.
Referenced by ROL::CompositeStep< Real >::accept(), ROL::CompositeStep< Real >::computeLagrangeMultiplier(), ROL::CompositeStep< Real >::computeQuasinormalStep(), and ROL::CompositeStep< Real >::solveTangentialSubproblem().
|
inlineprivate |
Definition at line 105 of file ROL_CompositeStep.hpp.
References ROL::CompositeStep< Real >::tolOSS_, and ROL::CompositeStep< Real >::tolOSSfixed_.
Referenced by ROL::CompositeStep< Real >::accept(), ROL::CompositeStep< Real >::computeLagrangeMultiplier(), ROL::CompositeStep< Real >::computeQuasinormalStep(), and ROL::CompositeStep< Real >::solveTangentialSubproblem().
|
inlinevirtual |
Initialize step.
Reimplemented from ROL::Step< Real >.
Definition at line 175 of file ROL_CompositeStep.hpp.
References ROL::Constraint< Real >::applyAdjointJacobian(), ROL::Vector< Real >::clone(), ROL::AlgorithmState< Real >::cnorm, ROL::CompositeStep< Real >::computeLagrangeMultiplier(), ROL::CompositeStep< Real >::cvec_, ROL::Step< Real >::getState(), ROL::AlgorithmState< Real >::gnorm, ROL::Objective< Real >::gradient(), ROL::CompositeStep< Real >::gvec_, ROL::AlgorithmState< Real >::iter, ROL::CompositeStep< Real >::lvec_, ROL::AlgorithmState< Real >::ncval, ROL::AlgorithmState< Real >::nfval, ROL::AlgorithmState< Real >::ngrad, ROL::Objective< Real >::update(), ROL::Constraint< Real >::update(), ROL::Objective< Real >::value(), ROL::Constraint< Real >::value(), ROL::AlgorithmState< Real >::value, and ROL::CompositeStep< Real >::xvec_.
|
inlinevirtual |
Compute step.
Reimplemented from ROL::Step< Real >.
Definition at line 218 of file ROL_CompositeStep.hpp.
References ROL::CompositeStep< Real >::accept(), ROL::Constraint< Real >::applyAdjointJacobian(), ROL::CompositeStep< Real >::computeQuasinormalStep(), ROL::CompositeStep< Real >::cvec_, ROL::CompositeStep< Real >::Delta_, ROL::Objective< Real >::gradient(), ROL::CompositeStep< Real >::gvec_, ROL::CompositeStep< Real >::iterCG_, ROL::CompositeStep< Real >::lvec_, ROL::AlgorithmState< Real >::nfval, ROL::AlgorithmState< Real >::ngrad, ROL::CompositeStep< Real >::solveTangentialSubproblem(), ROL::CompositeStep< Real >::totalIterCG_, ROL::Objective< Real >::value(), ROL::Constraint< Real >::value(), ROL::CompositeStep< Real >::xvec_, and ROL::CompositeStep< Real >::zeta_.
|
inlinevirtual |
Update step, if successful.
Reimplemented from ROL::Step< Real >.
Definition at line 266 of file ROL_CompositeStep.hpp.
References ROL::Constraint< Real >::applyAdjointJacobian(), ROL::CompositeStep< Real >::ared_, ROL::AlgorithmState< Real >::cnorm, ROL::CompositeStep< Real >::computeLagrangeMultiplier(), ROL::CompositeStep< Real >::cvec_, ROL::CompositeStep< Real >::Delta_, ROL::CompositeStep< Real >::eta_, ROL::CompositeStep< Real >::flagAC_, ROL::Step< Real >::getState(), ROL::AlgorithmState< Real >::gnorm, ROL::Objective< Real >::gradient(), ROL::CompositeStep< Real >::gvec_, ROL::AlgorithmState< Real >::iter, ROL::AlgorithmState< Real >::nfval, ROL::AlgorithmState< Real >::ngrad, ROL::CompositeStep< Real >::nnorm_, ROL::Vector< Real >::plus(), ROL::CompositeStep< Real >::pred_, ROL::AlgorithmState< Real >::snorm, ROL::CompositeStep< Real >::snorm_, ROL::CompositeStep< Real >::tnorm_, ROL::Objective< Real >::update(), ROL::Constraint< Real >::update(), ROL::Objective< Real >::value(), ROL::Constraint< Real >::value(), ROL::AlgorithmState< Real >::value, and zero.
|
inlinevirtual |
Compute step for bound constraints; here only to satisfy the interface requirements, does nothing, needs refactoring.
Reimplemented from ROL::Step< Real >.
Definition at line 338 of file ROL_CompositeStep.hpp.
|
inlinevirtual |
Update step, for bound constraints; here only to satisfy the interface requirements, does nothing, needs refactoring.
Reimplemented from ROL::Step< Real >.
Definition at line 345 of file ROL_CompositeStep.hpp.
|
inlinevirtual |
Print iterate header.
Reimplemented from ROL::Step< Real >.
Definition at line 351 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::print().
|
inlinevirtual |
Print step name.
Reimplemented from ROL::Step< Real >.
Definition at line 372 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::print().
|
inlinevirtual |
Print iterate status.
Reimplemented from ROL::Step< Real >.
Definition at line 381 of file ROL_CompositeStep.hpp.
References ROL::AlgorithmState< Real >::cnorm, ROL::CompositeStep< Real >::Delta_, ROL::CompositeStep< Real >::flagAC_, ROL::CompositeStep< Real >::flagCG_, ROL::AlgorithmState< Real >::gnorm, ROL::AlgorithmState< Real >::iter, ROL::CompositeStep< Real >::iterCG_, ROL::AlgorithmState< Real >::nfval, ROL::AlgorithmState< Real >::ngrad, ROL::CompositeStep< Real >::nnorm_, ROL::CompositeStep< Real >::printHeader(), ROL::CompositeStep< Real >::printName(), ROL::AlgorithmState< Real >::snorm, ROL::CompositeStep< Real >::tnorm_, ROL::CompositeStep< Real >::totalCallLS_, ROL::CompositeStep< Real >::totalIterLS_, and ROL::AlgorithmState< Real >::value.
|
inline |
Compute Lagrange multipliers by solving the least-squares problem minimizing the gradient of the Lagrangian, via the augmented system formulation.
[out] | l | is the updated Lagrange multiplier; a dual constraint-space vector |
[in] | x | is the current iterate; an optimization-space vector |
[in] | gf | is the gradient of the objective function; a dual optimization-space vector |
[in] | con | is the equality constraint object |
On return ... fill the blanks.
Definition at line 434 of file ROL_CompositeStep.hpp.
References ROL::Constraint< Real >::applyAdjointJacobian(), ROL::CompositeStep< Real >::cvec_, ROL::CompositeStep< Real >::gvec_, ROL::CompositeStep< Real >::infoLM_, ROL::CompositeStep< Real >::lmhtol_, ROL::CompositeStep< Real >::lvec_, ROL::Vector< Real >::plus(), ROL::CompositeStep< Real >::printInfoLS(), ROL::CompositeStep< Real >::setTolOSS(), ROL::Constraint< Real >::solveAugmentedSystem(), ROL::CompositeStep< Real >::totalCallLS_, ROL::CompositeStep< Real >::totalIterLS_, and ROL::CompositeStep< Real >::xvec_.
Referenced by ROL::CompositeStep< Real >::accept(), ROL::CompositeStep< Real >::initialize(), and ROL::CompositeStep< Real >::update().
|
inline |
Compute quasi-normal step by minimizing the norm of the linearized constraint.
Compute an approximate solution of the problem
\[ \begin{array}{rl} \min_{n} & \|c'(x_k)n + c(x_k)\|^2_{\mathcal{X}} \\ \mbox{subject to} & \|n\|_{\mathcal{X}} \le \delta \end{array} \]
The approximate solution is computed using Powell's dogleg method. The dogleg path is computed using the Cauchy point and a full Newton step. The path's intersection with the trust-region constraint gives the quasi-normal step.
[out] | n | is the quasi-normal step; an optimization-space vector |
[in] | c | is the value of equality constraints; a constraint-space vector |
[in] | x | is the current iterate; an optimization-space vector |
[in] | delta | is the trust-region radius for the quasi-normal step |
[in] | con | is the equality constraint object |
Definition at line 501 of file ROL_CompositeStep.hpp.
References ROL::Constraint< Real >::applyAdjointJacobian(), ROL::Constraint< Real >::applyJacobian(), ROL::Vector< Real >::axpy(), ROL::CompositeStep< Real >::cvec_, ROL::Vector< Real >::dual(), ROL::CompositeStep< Real >::gvec_, ROL::CompositeStep< Real >::infoQN_, ROL::CompositeStep< Real >::lvec_, ROL::CompositeStep< Real >::printInfoLS(), ROL::CompositeStep< Real >::qntol_, ROL::Vector< Real >::scale(), ROL::Vector< Real >::set(), ROL::CompositeStep< Real >::setTolOSS(), ROL::Constraint< Real >::solveAugmentedSystem(), ROL::CompositeStep< Real >::totalCallLS_, ROL::CompositeStep< Real >::totalIterLS_, ROL::CompositeStep< Real >::xvec_, and zero.
Referenced by ROL::CompositeStep< Real >::accept(), and ROL::CompositeStep< Real >::compute().
|
inline |
Solve tangential subproblem.
[out] | t | is the solution of the tangential subproblem; an optimization-space vector |
[out] | tCP | is the Cauchy point for the tangential subproblem; an optimization-space vector |
[out] | Wg | is the dual of the projected gradient of the Lagrangian; an optimization-space vector |
[in] | x | is the current iterate; an optimization-space vector |
[in] | g | is the gradient of the Lagrangian; a dual optimization-space vector |
[in] | n | is the quasi-normal step; an optimization-space vector |
[in] | l | is the Lagrange multiplier; a dual constraint-space vector |
[in] | delta | is the trust-region radius for the tangential subproblem |
[in] | con | is the equality constraint object |
Definition at line 611 of file ROL_CompositeStep.hpp.
References ROL::Constraint< Real >::applyAdjointHessian(), ROL::Constraint< Real >::applyJacobian(), ROL::CompositeStep< Real >::cvec_, ROL::Vector< Real >::dot(), ROL::Vector< Real >::dual(), ROL::CompositeStep< Real >::flagCG_, ROL::CompositeStep< Real >::gvec_, ROL::Objective< Real >::hessVec(), ROL::CompositeStep< Real >::infoTS_, ROL::CompositeStep< Real >::iterCG_, ROL::CompositeStep< Real >::lvec_, ROL::CompositeStep< Real >::maxiterCG_, ROL::Vector< Real >::norm(), ROL::CompositeStep< Real >::pgtol_, ROL::Vector< Real >::plus(), ROL::CompositeStep< Real >::printInfoLS(), ROL::CompositeStep< Real >::projtol_, ROL::Vector< Real >::set(), ROL::CompositeStep< Real >::setTolOSS(), ROL::CompositeStep< Real >::sgn(), ROL::Constraint< Real >::solveAugmentedSystem(), ROL::CompositeStep< Real >::tolCG_, ROL::CompositeStep< Real >::totalCallLS_, ROL::CompositeStep< Real >::totalIterLS_, ROL::CompositeStep< Real >::useConHess_, ROL::CompositeStep< Real >::xvec_, zero, and ROL::Vector< Real >::zero().
Referenced by ROL::CompositeStep< Real >::accept(), and ROL::CompositeStep< Real >::compute().
|
inline |
Check acceptance of subproblem solutions, adjust merit function penalty parameter, ensure global convergence.
Definition at line 938 of file ROL_CompositeStep.hpp.
References ROL::Constraint< Real >::applyAdjointHessian(), ROL::Constraint< Real >::applyAdjointJacobian(), ROL::Constraint< Real >::applyJacobian(), ROL::CompositeStep< Real >::ared_, ROL::CompositeStep< Real >::computeLagrangeMultiplier(), ROL::CompositeStep< Real >::computeQuasinormalStep(), ROL::CompositeStep< Real >::cvec_, ROL::CompositeStep< Real >::Delta_, ROL::Vector< Real >::dot(), ROL::Vector< Real >::dual(), ROL::CompositeStep< Real >::eta_, ROL::CompositeStep< Real >::flagCG_, ROL::Objective< Real >::gradient(), ROL::CompositeStep< Real >::gvec_, ROL::Objective< Real >::hessVec(), ROL::CompositeStep< Real >::infoAC_, ROL::AlgorithmState< Real >::iter, ROL::CompositeStep< Real >::iterCG_, ROL::CompositeStep< Real >::lmhtol_, ROL::CompositeStep< Real >::lvec_, ROL::CompositeStep< Real >::nnorm_, ROL::Vector< Real >::norm(), ROL::CompositeStep< Real >::penalty_, ROL::CompositeStep< Real >::pgtol_, ROL::Vector< Real >::plus(), ROL::CompositeStep< Real >::pred_, ROL::CompositeStep< Real >::printInfoLS(), ROL::CompositeStep< Real >::projtol_, ROL::CompositeStep< Real >::qntol_, ROL::Vector< Real >::set(), ROL::CompositeStep< Real >::setTolOSS(), ROL::CompositeStep< Real >::snorm_, ROL::Constraint< Real >::solveAugmentedSystem(), ROL::CompositeStep< Real >::solveTangentialSubproblem(), ROL::CompositeStep< Real >::tangtol_, ROL::CompositeStep< Real >::tnorm_, ROL::CompositeStep< Real >::tntmax_, ROL::CompositeStep< Real >::tolOSSfixed_, ROL::CompositeStep< Real >::totalCallLS_, ROL::CompositeStep< Real >::totalIterCG_, ROL::CompositeStep< Real >::totalIterLS_, ROL::CompositeStep< Real >::totalNegCurv_, ROL::CompositeStep< Real >::totalProj_, ROL::CompositeStep< Real >::totalRef_, ROL::Objective< Real >::update(), ROL::Constraint< Real >::update(), ROL::CompositeStep< Real >::useConHess_, ROL::Objective< Real >::value(), ROL::Constraint< Real >::value(), ROL::CompositeStep< Real >::xvec_, zero, and ROL::CompositeStep< Real >::zeta_.
Referenced by ROL::CompositeStep< Real >::compute().
|
private |
Definition at line 33 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::accept(), ROL::CompositeStep< Real >::compute(), ROL::CompositeStep< Real >::computeLagrangeMultiplier(), ROL::CompositeStep< Real >::computeQuasinormalStep(), ROL::CompositeStep< Real >::initialize(), and ROL::CompositeStep< Real >::solveTangentialSubproblem().
|
private |
Definition at line 34 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::accept(), ROL::CompositeStep< Real >::compute(), ROL::CompositeStep< Real >::computeLagrangeMultiplier(), ROL::CompositeStep< Real >::computeQuasinormalStep(), ROL::CompositeStep< Real >::initialize(), ROL::CompositeStep< Real >::solveTangentialSubproblem(), and ROL::CompositeStep< Real >::update().
|
private |
Definition at line 35 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::accept(), ROL::CompositeStep< Real >::compute(), ROL::CompositeStep< Real >::computeLagrangeMultiplier(), ROL::CompositeStep< Real >::computeQuasinormalStep(), ROL::CompositeStep< Real >::initialize(), ROL::CompositeStep< Real >::solveTangentialSubproblem(), and ROL::CompositeStep< Real >::update().
|
private |
Definition at line 36 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::accept(), ROL::CompositeStep< Real >::compute(), ROL::CompositeStep< Real >::computeLagrangeMultiplier(), ROL::CompositeStep< Real >::computeQuasinormalStep(), ROL::CompositeStep< Real >::initialize(), and ROL::CompositeStep< Real >::solveTangentialSubproblem().
|
private |
Definition at line 39 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::accept(), ROL::CompositeStep< Real >::CompositeStep(), ROL::CompositeStep< Real >::print(), and ROL::CompositeStep< Real >::solveTangentialSubproblem().
|
private |
Definition at line 40 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::CompositeStep(), ROL::CompositeStep< Real >::print(), and ROL::CompositeStep< Real >::update().
|
private |
|
private |
Definition at line 44 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::CompositeStep(), and ROL::CompositeStep< Real >::solveTangentialSubproblem().
|
private |
Definition at line 45 of file ROL_CompositeStep.hpp.
|
private |
Definition at line 46 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::CompositeStep(), and ROL::CompositeStep< Real >::solveTangentialSubproblem().
|
private |
Definition at line 47 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::CompositeStep(), and ROL::CompositeStep< Real >::setTolOSS().
|
private |
Definition at line 48 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::accept(), ROL::CompositeStep< Real >::CompositeStep(), and ROL::CompositeStep< Real >::setTolOSS().
|
private |
Definition at line 51 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::accept(), ROL::CompositeStep< Real >::CompositeStep(), and ROL::CompositeStep< Real >::computeLagrangeMultiplier().
|
private |
Definition at line 52 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::accept(), ROL::CompositeStep< Real >::CompositeStep(), and ROL::CompositeStep< Real >::computeQuasinormalStep().
|
private |
Definition at line 53 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::accept(), ROL::CompositeStep< Real >::CompositeStep(), and ROL::CompositeStep< Real >::solveTangentialSubproblem().
|
private |
Definition at line 54 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::accept(), ROL::CompositeStep< Real >::CompositeStep(), and ROL::CompositeStep< Real >::solveTangentialSubproblem().
|
private |
Definition at line 55 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::accept(), and ROL::CompositeStep< Real >::CompositeStep().
|
private |
Definition at line 56 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::accept(), and ROL::CompositeStep< Real >::CompositeStep().
|
private |
Definition at line 59 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::accept(), ROL::CompositeStep< Real >::CompositeStep(), and ROL::CompositeStep< Real >::compute().
|
private |
Definition at line 60 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::accept(), ROL::CompositeStep< Real >::CompositeStep(), ROL::CompositeStep< Real >::compute(), ROL::CompositeStep< Real >::print(), and ROL::CompositeStep< Real >::update().
|
private |
Definition at line 61 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::accept(), and ROL::CompositeStep< Real >::CompositeStep().
|
private |
Definition at line 62 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::accept(), ROL::CompositeStep< Real >::CompositeStep(), and ROL::CompositeStep< Real >::update().
|
private |
Definition at line 63 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::accept(), ROL::CompositeStep< Real >::CompositeStep(), and ROL::CompositeStep< Real >::solveTangentialSubproblem().
|
private |
Definition at line 65 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::accept(), and ROL::CompositeStep< Real >::update().
|
private |
Definition at line 66 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::accept(), and ROL::CompositeStep< Real >::update().
|
private |
Definition at line 67 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::accept(), ROL::CompositeStep< Real >::CompositeStep(), and ROL::CompositeStep< Real >::update().
|
private |
Definition at line 68 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::accept(), ROL::CompositeStep< Real >::CompositeStep(), ROL::CompositeStep< Real >::print(), and ROL::CompositeStep< Real >::update().
|
private |
Definition at line 69 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::accept(), ROL::CompositeStep< Real >::CompositeStep(), ROL::CompositeStep< Real >::print(), and ROL::CompositeStep< Real >::update().
|
private |
Definition at line 72 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::CompositeStep(), and ROL::CompositeStep< Real >::computeQuasinormalStep().
|
private |
Definition at line 73 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::CompositeStep(), and ROL::CompositeStep< Real >::computeLagrangeMultiplier().
|
private |
Definition at line 74 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::CompositeStep(), and ROL::CompositeStep< Real >::solveTangentialSubproblem().
|
private |
Definition at line 75 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::accept(), and ROL::CompositeStep< Real >::CompositeStep().
|
private |
Definition at line 76 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::CompositeStep(), and ROL::CompositeStep< Real >::printInfoLS().
|
private |
Definition at line 77 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::CompositeStep().
|
private |
Definition at line 80 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::accept(), ROL::CompositeStep< Real >::CompositeStep(), and ROL::CompositeStep< Real >::compute().
|
private |
Definition at line 81 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::accept(), and ROL::CompositeStep< Real >::CompositeStep().
|
private |
Definition at line 82 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::accept(), and ROL::CompositeStep< Real >::CompositeStep().
|
private |
Definition at line 83 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::accept(), and ROL::CompositeStep< Real >::CompositeStep().
|
private |
Definition at line 84 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::accept(), ROL::CompositeStep< Real >::CompositeStep(), ROL::CompositeStep< Real >::computeLagrangeMultiplier(), ROL::CompositeStep< Real >::computeQuasinormalStep(), ROL::CompositeStep< Real >::print(), and ROL::CompositeStep< Real >::solveTangentialSubproblem().
|
private |
Definition at line 85 of file ROL_CompositeStep.hpp.
Referenced by ROL::CompositeStep< Real >::accept(), ROL::CompositeStep< Real >::CompositeStep(), ROL::CompositeStep< Real >::computeLagrangeMultiplier(), ROL::CompositeStep< Real >::computeQuasinormalStep(), ROL::CompositeStep< Real >::print(), and ROL::CompositeStep< Real >::solveTangentialSubproblem().