|
ROL
|
Provides an interface to run equality constrained optimization algorithms using the Composite-Step Trust-Region Sequential Quadratic Programming (SQP) method. More...
#include <ROL_TypeE_CompositeStepAlgorithm.hpp>
Inheritance diagram for ROL::TypeE::CompositeStepAlgorithm< Real >:Public Member Functions | |
| CompositeStepAlgorithm (ParameterList &list) | |
| virtual void | run (Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, Constraint< Real > &econ, Vector< Real > &emul, const Vector< Real > &eres, std::ostream &outStream=std::cout) override |
| Run algorithm on equality constrained problems (Type-E). This general interface supports the use of dual optimization vector spaces, where the user does not define the dual() method. | |
| virtual void | writeHeader (std::ostream &os) const override |
| Print iterate header. | |
| virtual void | writeName (std::ostream &os) const override |
| Print step name. | |
| virtual void | writeOutput (std::ostream &os, const bool print_header=false) const override |
| Print iterate status. | |
Public Member Functions inherited from ROL::TypeE::Algorithm< Real > | |
| 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 equality constrained problems (Type-E). This is the primary Type-E interface. | |
| virtual void | run (Vector< Real > &x, Objective< Real > &obj, Constraint< Real > &econ, Vector< Real > &emul, std::ostream &outStream=std::cout) |
| Run algorithm on equality constrained problems (Type-E). This is the primary Type-E interface. | |
| virtual void | run (Vector< Real > &x, Objective< Real > &obj, Constraint< Real > &econ, Vector< Real > &emul, Constraint< Real > &linear_econ, Vector< Real > &linear_emul, std::ostream &outStream=std::cout) |
| Run algorithm on equality constrained problems with explicit linear equality constraints (Type-E). This is the primary Type-E with explicit linear equality constraints interface. | |
| virtual void | run (Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, Constraint< Real > &econ, Vector< Real > &emul, const Vector< Real > &eres, Constraint< Real > &linear_econ, Vector< Real > &linear_emul, const Vector< Real > &linear_eres, std::ostream &outStream=std::cout) |
| Run algorithm on equality constrained problems with explicit linear equality constraints (Type-E). 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 () |
Private Member Functions | |
| void | initialize (Vector< Real > &x, const Vector< Real > &g, Vector< Real > &l, const Vector< Real > &c, Objective< Real > &obj, Constraint< Real > &con, std::ostream &outStream=std::cout) |
| Initialize algorithm by computing a few quantities. | |
| void | computeTrial (Vector< Real > &s, const Vector< Real > &x, const Vector< Real > &l, Objective< Real > &obj, Constraint< Real > &con, std::ostream &os) |
| Compute trial step. | |
| void | updateRadius (Vector< Real > &x, Vector< Real > &l, const Vector< Real > &s, Objective< Real > &obj, Constraint< Real > &con, std::ostream &os) |
| Update trust-region radius, take step, etc. | |
| void | computeLagrangeMultiplier (Vector< Real > &l, const Vector< Real > &x, const Vector< Real > &gf, Constraint< Real > &con, std::ostream &os) |
| Compute Lagrange multipliers by solving the least-squares problem minimizing the gradient of the Lagrangian, via the augmented system formulation. | |
| void | computeQuasinormalStep (Vector< Real > &n, const Vector< Real > &c, const Vector< Real > &x, Real delta, Constraint< Real > &con, std::ostream &os) |
| Compute quasi-normal step by minimizing the norm of the linearized constraint. | |
| 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, std::ostream &os) |
| Solve tangential subproblem. | |
| 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, std::ostream &os) |
| Check acceptance of subproblem solutions, adjust merit function penalty parameter, ensure global convergence. | |
| template<typename T > | |
| int | sgn (T val) const |
| void | printInfoLS (const std::vector< Real > &res, std::ostream &os) const |
| Real | setTolOSS (const Real intol) const |
Private Attributes | |
| ParameterList | list_ |
| 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_ |
| int | verbosity_ |
| bool | printHeader_ |
Additional Inherited Members | |
Protected Member Functions inherited from ROL::TypeE::Algorithm< Real > | |
| void | initialize (const Vector< Real > &x, const Vector< Real > &g, const Vector< Real > &mul, const Vector< Real > &c) |
Protected Attributes inherited from ROL::TypeE::Algorithm< Real > | |
| const Ptr< CombinedStatusTest< Real > > | status_ |
| const Ptr< AlgorithmState< Real > > | state_ |
Provides an interface to run equality constrained optimization algorithms using the Composite-Step Trust-Region Sequential Quadratic Programming (SQP) method.
Definition at line 24 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
| ROL::TypeE::CompositeStepAlgorithm< Real >::CompositeStepAlgorithm | ( | ParameterList & | list | ) |
Definition at line 17 of file ROL_TypeE_CompositeStepAlgorithm_Def.hpp.
References ROL::TypeE::CompositeStepAlgorithm< Real >::Delta_, ROL::TypeE::CompositeStepAlgorithm< Real >::eta_, ROL::TypeE::CompositeStepAlgorithm< Real >::flagAC_, ROL::TypeE::CompositeStepAlgorithm< Real >::flagCG_, ROL::TypeE::CompositeStepAlgorithm< Real >::infoAC_, ROL::TypeE::CompositeStepAlgorithm< Real >::infoALL_, ROL::TypeE::CompositeStepAlgorithm< Real >::infoLM_, ROL::TypeE::CompositeStepAlgorithm< Real >::infoLS_, ROL::TypeE::CompositeStepAlgorithm< Real >::infoQN_, ROL::TypeE::CompositeStepAlgorithm< Real >::infoTS_, ROL::TypeE::CompositeStepAlgorithm< Real >::iterCG_, ROL::TypeE::CompositeStepAlgorithm< Real >::list_, ROL::TypeE::CompositeStepAlgorithm< Real >::lmhtol_, ROL::TypeE::CompositeStepAlgorithm< Real >::maxiterCG_, ROL::TypeE::CompositeStepAlgorithm< Real >::nnorm_, ROL::TypeE::CompositeStepAlgorithm< Real >::penalty_, ROL::TypeE::CompositeStepAlgorithm< Real >::pgtol_, ROL::TypeE::CompositeStepAlgorithm< Real >::printHeader_, ROL::TypeE::CompositeStepAlgorithm< Real >::projtol_, ROL::TypeE::CompositeStepAlgorithm< Real >::qntol_, ROL::TypeE::Algorithm< Real >::reset(), ROL::TypeE::CompositeStepAlgorithm< Real >::snorm_, ROL::TypeE::Algorithm< Real >::status_, ROL::TypeE::CompositeStepAlgorithm< Real >::tangtol_, ROL::TypeE::CompositeStepAlgorithm< Real >::tnorm_, ROL::TypeE::CompositeStepAlgorithm< Real >::tntmax_, ROL::TypeE::CompositeStepAlgorithm< Real >::tolCG_, ROL::TypeE::CompositeStepAlgorithm< Real >::tolOSS_, ROL::TypeE::CompositeStepAlgorithm< Real >::tolOSSfixed_, ROL::TypeE::CompositeStepAlgorithm< Real >::totalCallLS_, ROL::TypeE::CompositeStepAlgorithm< Real >::totalIterCG_, ROL::TypeE::CompositeStepAlgorithm< Real >::totalIterLS_, ROL::TypeE::CompositeStepAlgorithm< Real >::totalNegCurv_, ROL::TypeE::CompositeStepAlgorithm< Real >::totalProj_, ROL::TypeE::CompositeStepAlgorithm< Real >::totalRef_, ROL::TypeE::CompositeStepAlgorithm< Real >::useConHess_, ROL::TypeE::CompositeStepAlgorithm< Real >::verbosity_, zero, and ROL::TypeE::CompositeStepAlgorithm< Real >::zeta_.
|
private |
Initialize algorithm by computing a few quantities.
Definition at line 987 of file ROL_TypeE_CompositeStepAlgorithm_Def.hpp.
References ROL::Constraint< Real >::applyAdjointJacobian(), ROL::Vector< Real >::clone(), ROL::Objective< Real >::gradient(), ROL::Initial, ROL::TypeE::Algorithm< Real >::initialize(), ROL::Constraint< Real >::update(), ROL::Objective< Real >::update(), ROL::Objective< Real >::value(), and ROL::Constraint< Real >::value().
|
private |
Compute trial step.
Definition at line 82 of file ROL_TypeE_CompositeStepAlgorithm_Def.hpp.
References ROL::Constraint< Real >::applyAdjointJacobian(), ROL::Objective< Real >::gradient(), ROL::Objective< Real >::value(), and ROL::Constraint< Real >::value().
|
private |
Update trust-region radius, take step, etc.
Definition at line 916 of file ROL_TypeE_CompositeStepAlgorithm_Def.hpp.
References ROL::Accept, ROL::Constraint< Real >::applyAdjointJacobian(), ROL::Objective< Real >::gradient(), ROL::Vector< Real >::plus(), ROL::Revert, ROL::Constraint< Real >::update(), ROL::Objective< Real >::update(), ROL::Objective< Real >::value(), ROL::Constraint< Real >::value(), and zero.
|
private |
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 |
Definition at line 133 of file ROL_TypeE_CompositeStepAlgorithm_Def.hpp.
References ROL::Constraint< Real >::applyAdjointJacobian(), ROL::Vector< Real >::plus(), and ROL::Constraint< Real >::solveAugmentedSystem().
|
private |
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 183 of file ROL_TypeE_CompositeStepAlgorithm_Def.hpp.
References ROL::Constraint< Real >::applyAdjointJacobian(), ROL::Constraint< Real >::applyJacobian(), ROL::Vector< Real >::axpy(), ROL::Vector< Real >::dual(), ROL::Vector< Real >::scale(), ROL::Vector< Real >::set(), ROL::Constraint< Real >::solveAugmentedSystem(), and zero.
|
private |
Solve tangential subproblem.
@param[out] t is the solution of the tangential subproblem; an optimization-space vector @param[out] tCP is the Cauchy point for the tangential subproblem; an optimization-space vector @param[out] Wg is the dual of the projected gradient of the Lagrangian; an optimization-space vector @param[in] x is the current iterate; an optimization-space vector @param[in] g is the gradient of the Lagrangian; a dual optimization-space vector @param[in] n is the quasi-normal step; an optimization-space vector @param[in] l is the Lagrange multiplier; a dual constraint-space vector @param[in] delta is the trust-region radius for the tangential subproblem @param[in] obj is the objective function @param[in] con is the equality constraint object
Definition at line 286 of file ROL_TypeE_CompositeStepAlgorithm_Def.hpp.
References ROL::Constraint< Real >::applyAdjointHessian(), ROL::Constraint< Real >::applyJacobian(), ROL::Vector< Real >::dot(), ROL::Vector< Real >::dual(), ROL::Objective< Real >::hessVec(), ROL::Vector< Real >::norm(), ROL::Vector< Real >::plus(), ROL::Vector< Real >::set(), ROL::Constraint< Real >::solveAugmentedSystem(), zero, and ROL::Vector< Real >::zero().
|
private |
Check acceptance of subproblem solutions, adjust merit function penalty parameter, ensure global convergence.
Definition at line 620 of file ROL_TypeE_CompositeStepAlgorithm_Def.hpp.
References ROL::Vector< Real >::apply(), ROL::Constraint< Real >::applyAdjointHessian(), ROL::Constraint< Real >::applyAdjointJacobian(), ROL::Constraint< Real >::applyJacobian(), ROL::Vector< Real >::dot(), ROL::Objective< Real >::gradient(), ROL::Objective< Real >::hessVec(), ROL::Vector< Real >::norm(), ROL::Vector< Real >::plus(), ROL::Vector< Real >::set(), ROL::Constraint< Real >::solveAugmentedSystem(), ROL::Trial, ROL::Constraint< Real >::update(), ROL::Objective< Real >::update(), ROL::Objective< Real >::value(), ROL::Constraint< Real >::value(), and zero.
|
private |
Definition at line 1159 of file ROL_TypeE_CompositeStepAlgorithm_Def.hpp.
|
private |
Definition at line 1165 of file ROL_TypeE_CompositeStepAlgorithm_Def.hpp.
|
private |
Definition at line 1181 of file ROL_TypeE_CompositeStepAlgorithm_Def.hpp.
|
overridevirtual |
Run algorithm on equality constrained problems (Type-E). This general interface supports the use of dual optimization vector spaces, where the user does not define the dual() method.
Implements ROL::TypeE::Algorithm< Real >.
Definition at line 1030 of file ROL_TypeE_CompositeStepAlgorithm_Def.hpp.
References ROL::Vector< Real >::clone(), and ROL::TypeE::Algorithm< Real >::writeExitStatus().
Referenced by main().
|
overridevirtual |
Print iterate header.
Reimplemented from ROL::TypeE::Algorithm< Real >.
Definition at line 1059 of file ROL_TypeE_CompositeStepAlgorithm_Def.hpp.
|
overridevirtual |
Print step name.
Reimplemented from ROL::TypeE::Algorithm< Real >.
Definition at line 1101 of file ROL_TypeE_CompositeStepAlgorithm_Def.hpp.
|
overridevirtual |
Print iterate status.
Reimplemented from ROL::TypeE::Algorithm< Real >.
Definition at line 1110 of file ROL_TypeE_CompositeStepAlgorithm_Def.hpp.
|
private |
Definition at line 27 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::CompositeStepAlgorithm().
|
private |
Definition at line 30 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
|
private |
Definition at line 31 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
|
private |
Definition at line 32 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
|
private |
Definition at line 33 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
|
private |
Definition at line 36 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::CompositeStepAlgorithm().
|
private |
Definition at line 37 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::CompositeStepAlgorithm().
|
private |
Definition at line 38 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::CompositeStepAlgorithm().
|
private |
Definition at line 41 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::CompositeStepAlgorithm().
|
private |
Definition at line 42 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
|
private |
Definition at line 43 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::CompositeStepAlgorithm().
|
private |
Definition at line 44 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::CompositeStepAlgorithm().
|
private |
Definition at line 45 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::CompositeStepAlgorithm().
|
private |
Definition at line 48 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::CompositeStepAlgorithm().
|
private |
Definition at line 49 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::CompositeStepAlgorithm().
|
private |
Definition at line 50 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::CompositeStepAlgorithm().
|
private |
Definition at line 51 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::CompositeStepAlgorithm().
|
private |
Definition at line 52 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::CompositeStepAlgorithm().
|
private |
Definition at line 53 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::CompositeStepAlgorithm().
|
private |
Definition at line 56 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::CompositeStepAlgorithm().
|
private |
Definition at line 57 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::CompositeStepAlgorithm().
|
private |
Definition at line 58 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::CompositeStepAlgorithm().
|
private |
Definition at line 59 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::CompositeStepAlgorithm().
|
private |
Definition at line 60 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::CompositeStepAlgorithm().
|
private |
Definition at line 62 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
|
private |
Definition at line 63 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
|
private |
Definition at line 64 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::CompositeStepAlgorithm().
|
private |
Definition at line 65 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::CompositeStepAlgorithm().
|
private |
Definition at line 66 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::CompositeStepAlgorithm().
|
private |
Definition at line 69 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::CompositeStepAlgorithm().
|
private |
Definition at line 70 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::CompositeStepAlgorithm().
|
private |
Definition at line 71 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::CompositeStepAlgorithm().
|
private |
Definition at line 72 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::CompositeStepAlgorithm().
|
private |
Definition at line 73 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::CompositeStepAlgorithm().
|
private |
Definition at line 74 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::CompositeStepAlgorithm().
|
private |
Definition at line 77 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::CompositeStepAlgorithm().
|
private |
Definition at line 78 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::CompositeStepAlgorithm().
|
private |
Definition at line 79 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::CompositeStepAlgorithm().
|
private |
Definition at line 80 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::CompositeStepAlgorithm().
|
private |
Definition at line 81 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::CompositeStepAlgorithm().
|
private |
Definition at line 82 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::CompositeStepAlgorithm().
|
private |
Definition at line 85 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::CompositeStepAlgorithm().
|
private |
Definition at line 86 of file ROL_TypeE_CompositeStepAlgorithm.hpp.
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::CompositeStepAlgorithm().