|
ROL
|
Provides the interface to evaluate trust-region model functions. More...
#include <ROL_TrustRegionModel_U.hpp>
Public Member Functions | |
| virtual | ~TrustRegionModel_U () |
| TrustRegionModel_U (ParameterList &list, const Ptr< Secant< Real > > &secant=nullPtr, ESecantMode mode=SECANTMODE_BOTH) | |
| void | initialize (const Vector< Real > &x, const Vector< Real > &g) |
| void | validate (Objective< Real > &obj, const Vector< Real > &x, const Vector< Real > &g, ETrustRegionU etr) |
| virtual void | setData (Objective< Real > &obj, const Vector< Real > &x, const Vector< Real > &g, Real &tol) |
| void | update (const Vector< Real > &x, const Vector< Real > &s, const Vector< Real > &gold, const Vector< Real > &gnew, const Real snorm, const int iter) |
| virtual Real | value (const Vector< Real > &s, Real &tol) override |
| virtual void | gradient (Vector< Real > &g, const Vector< Real > &s, Real &tol) override |
| virtual void | hessVec (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &s, Real &tol) override |
| virtual void | invHessVec (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &s, Real &tol) override |
| virtual void | precond (Vector< Real > &Pv, const Vector< Real > &v, const Vector< Real > &s, Real &tol) override |
| virtual const Ptr< const Vector< Real > > | getGradient (void) const |
| virtual const Ptr< const Vector< Real > > | getIterate (void) const |
| virtual const Ptr< Objective< Real > > | getObjective (void) const |
| Public Member Functions inherited from ROL::ROL::Objective< Real > | |
| virtual | ~Objective () |
| Objective () | |
| virtual void | update (const Vector< Real > &x, UpdateType type, int iter=-1) |
| Update objective function. | |
| virtual void | update (const Vector< Real > &x, bool flag=true, int iter=-1) |
| Update objective function. | |
| virtual Real | value (const Vector< Real > &x, Real &tol)=0 |
| Compute value. | |
| virtual void | gradient (Vector< Real > &g, const Vector< Real > &x, Real &tol) |
| Compute gradient. | |
| virtual Real | dirDeriv (const Vector< Real > &x, const Vector< Real > &d, Real &tol) |
| Compute directional derivative. | |
| virtual void | hessVec (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol) |
| Apply Hessian approximation to vector. | |
| virtual void | invHessVec (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol) |
| Apply inverse Hessian approximation to vector. | |
| virtual void | precond (Vector< Real > &Pv, const Vector< Real > &v, const Vector< Real > &x, Real &tol) |
| Apply preconditioner to vector. | |
| virtual void | prox (Vector< Real > &Pv, const Vector< Real > &v, Real t, Real &tol) |
| Compute the proximity operator. | |
| virtual void | proxJacVec (Vector< Real > &Jv, const Vector< Real > &v, const Vector< Real > &x, Real t, Real &tol) |
| Apply the Jacobian of the proximity operator. | |
| virtual std::vector< std::vector< Real > > | checkGradient (const Vector< Real > &x, const Vector< Real > &d, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1) |
| Finite-difference gradient check. | |
| virtual std::vector< std::vector< Real > > | checkGradient (const Vector< Real > &x, const Vector< Real > &g, const Vector< Real > &d, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1) |
| Finite-difference gradient check. | |
| virtual std::vector< std::vector< Real > > | checkGradient (const Vector< Real > &x, const Vector< Real > &d, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1) |
| Finite-difference gradient check with specified step sizes. | |
| virtual std::vector< std::vector< Real > > | checkGradient (const Vector< Real > &x, const Vector< Real > &g, const Vector< Real > &d, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1) |
| Finite-difference gradient check with specified step sizes. | |
| virtual std::vector< std::vector< Real > > | checkHessVec (const Vector< Real > &x, const Vector< Real > &v, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1) |
| Finite-difference Hessian-applied-to-vector check. | |
| virtual std::vector< std::vector< Real > > | checkHessVec (const Vector< Real > &x, const Vector< Real > &hv, const Vector< Real > &v, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1) |
| Finite-difference Hessian-applied-to-vector check. | |
| virtual std::vector< std::vector< Real > > | checkHessVec (const Vector< Real > &x, const Vector< Real > &v, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1) |
| Finite-difference Hessian-applied-to-vector check with specified step sizes. | |
| virtual std::vector< std::vector< Real > > | checkHessVec (const Vector< Real > &x, const Vector< Real > &hv, const Vector< Real > &v, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1) |
| Finite-difference Hessian-applied-to-vector check with specified step sizes. | |
| virtual std::vector< Real > | checkHessSym (const Vector< Real > &x, const Vector< Real > &v, const Vector< Real > &w, const bool printToStream=true, std::ostream &outStream=std::cout) |
| Hessian symmetry check. | |
| virtual std::vector< Real > | checkHessSym (const Vector< Real > &x, const Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &w, const bool printToStream=true, std::ostream &outStream=std::cout) |
| Hessian symmetry check. | |
| virtual std::vector< std::vector< Real > > | checkProxJacVec (const Vector< Real > &x, const Vector< Real > &v, Real t=Real(1), bool printToStream=true, std::ostream &outStream=std::cout, int numSteps=ROL_NUM_CHECKDERIV_STEPS) |
| Finite-difference proximity operator Jacobian-applied-to-vector check. | |
| virtual void | setParameter (const std::vector< Real > ¶m) |
Protected Member Functions | |
| void | applyHessian (Vector< Real > &hv, const Vector< Real > &v, Real &tol) |
| void | applyInvHessian (Vector< Real > &hv, const Vector< Real > &v, Real &tol) |
| void | applyPrecond (Vector< Real > &Pv, const Vector< Real > &v, Real &tol) |
| Protected Member Functions inherited from ROL::ROL::Objective< Real > | |
| const std::vector< Real > | getParameter (void) const |
Private Attributes | |
| Ptr< Objective< Real > > | obj_ |
| Ptr< const Vector< Real > > | x_ |
| Ptr< const Vector< Real > > | g_ |
| Ptr< Vector< Real > > | dual_ |
| Real | tol_ |
| Ptr< Secant< Real > > | secant_ |
| bool | useSecantPrecond_ |
| bool | useSecantHessVec_ |
Provides the interface to evaluate trust-region model functions.
ROL::TrustRegionModel_U provides the interface to implement a number of trust-region models for unconstrained and constrained optimization. The default implementation is the standard quadratic trust region model for unconstrained optimization.
Definition at line 32 of file ROL_TrustRegionModel_U.hpp.
|
inlinevirtual |
Definition at line 79 of file ROL_TrustRegionModel_U.hpp.
|
inline |
Definition at line 81 of file ROL_TrustRegionModel_U.hpp.
References g_, obj_, secant_, ROL::SecantFactory(), ROL::SECANTMODE_BOTH, useSecantHessVec_, useSecantPrecond_, and x_.
|
inlineprotected |
Definition at line 47 of file ROL_TrustRegionModel_U.hpp.
References obj_, secant_, tol_, useSecantHessVec_, and x_.
Referenced by gradient(), hessVec(), and value().
|
inlineprotected |
Definition at line 56 of file ROL_TrustRegionModel_U.hpp.
References obj_, secant_, tol_, useSecantHessVec_, and x_.
Referenced by invHessVec().
|
inlineprotected |
Definition at line 65 of file ROL_TrustRegionModel_U.hpp.
References obj_, secant_, tol_, useSecantPrecond_, and x_.
Referenced by precond().
|
inline |
Definition at line 91 of file ROL_TrustRegionModel_U.hpp.
References ROL::Vector< Real >::clone(), and dual_.
|
inline |
Definition at line 99 of file ROL_TrustRegionModel_U.hpp.
References ROL::Vector< Real >::clone(), ROL::Objective< Real >::invHessVec(), ROL::ROL::Objective< Real >::Objective(), ROL::ROL_EPSILON(), ROL::TRUSTREGION_U_DOGLEG, ROL::TRUSTREGION_U_DOUBLEDOGLEG, and useSecantHessVec_.
|
inlinevirtual |
Definition at line 117 of file ROL_TrustRegionModel_U.hpp.
References g_, obj_, ROL::ROL::Objective< Real >::Objective(), tol_, and x_.
Referenced by ROL::TRUtils::initialRadius().
|
inline |
Definition at line 127 of file ROL_TrustRegionModel_U.hpp.
References secant_, useSecantHessVec_, and useSecantPrecond_.
|
inlineoverridevirtual |
Definition at line 138 of file ROL_TrustRegionModel_U.hpp.
References applyHessian(), dual_, and g_.
|
inlineoverridevirtual |
Definition at line 145 of file ROL_TrustRegionModel_U.hpp.
References applyHessian(), g_, and ROL::Vector< Real >::plus().
|
inlineoverridevirtual |
Definition at line 150 of file ROL_TrustRegionModel_U.hpp.
References applyHessian().
Referenced by ROL::TypeB::KelleySachsAlgorithm< Real >::applyFreeHessian(), ROL::TypeB::LinMoreAlgorithm< Real >::applyFreeHessian(), ROL::TypeB::ColemanLiAlgorithm< Real >::applyHessian(), ROL::TypeB::LinMoreAlgorithm< Real >::dcauchy(), ROL::TypeB::TrustRegionSPGAlgorithm< Real >::dcauchy(), ROL::TypeP::TrustRegionAlgorithm< Real >::dcauchy(), ROL::TypeP::TrustRegionAlgorithm< Real >::dncg(), ROL::TypeB::LinMoreAlgorithm< Real >::dprsrch(), ROL::TypeB::TrustRegionSPGAlgorithm< Real >::dpsg(), ROL::TypeB::TrustRegionSPGAlgorithm< Real >::dpsg_simple(), ROL::TypeP::TrustRegionAlgorithm< Real >::dspg(), ROL::TypeP::TrustRegionAlgorithm< Real >::dspg2(), ROL::TRUtils::initialRadius(), ROL::CauchyPoint_U< Real >::solve(), ROL::DogLeg_U< Real >::solve(), ROL::DoubleDogLeg_U< Real >::solve(), ROL::SPGTrustRegion_U< Real >::solve(), and ROL::TruncatedCG_U< Real >::solve().
|
inlineoverridevirtual |
Definition at line 154 of file ROL_TrustRegionModel_U.hpp.
References applyInvHessian().
Referenced by ROL::DogLeg_U< Real >::solve(), and ROL::DoubleDogLeg_U< Real >::solve().
|
inlineoverridevirtual |
Definition at line 158 of file ROL_TrustRegionModel_U.hpp.
References applyPrecond().
Referenced by ROL::TypeB::KelleySachsAlgorithm< Real >::applyFreePrecond(), ROL::TypeB::LinMoreAlgorithm< Real >::applyFreePrecond(), ROL::TypeB::ColemanLiAlgorithm< Real >::applyPrecond(), and ROL::TruncatedCG_U< Real >::solve().
|
inlinevirtual |
Definition at line 168 of file ROL_TrustRegionModel_U.hpp.
References g_.
Referenced by ROL::CauchyPoint_U< Real >::solve(), ROL::DogLeg_U< Real >::solve(), ROL::DoubleDogLeg_U< Real >::solve(), ROL::SPGTrustRegion_U< Real >::solve(), and ROL::TruncatedCG_U< Real >::solve().
|
inlinevirtual |
Definition at line 172 of file ROL_TrustRegionModel_U.hpp.
References x_.
|
inlinevirtual |
Definition at line 176 of file ROL_TrustRegionModel_U.hpp.
References obj_.
|
private |
Definition at line 34 of file ROL_TrustRegionModel_U.hpp.
Referenced by applyHessian(), applyInvHessian(), applyPrecond(), getObjective(), setData(), and TrustRegionModel_U().
|
private |
Definition at line 35 of file ROL_TrustRegionModel_U.hpp.
Referenced by applyHessian(), applyInvHessian(), applyPrecond(), getIterate(), setData(), and TrustRegionModel_U().
|
private |
Definition at line 35 of file ROL_TrustRegionModel_U.hpp.
Referenced by getGradient(), gradient(), setData(), TrustRegionModel_U(), and value().
|
private |
Definition at line 36 of file ROL_TrustRegionModel_U.hpp.
Referenced by initialize(), and value().
|
private |
Definition at line 37 of file ROL_TrustRegionModel_U.hpp.
Referenced by applyHessian(), applyInvHessian(), applyPrecond(), and setData().
|
private |
Definition at line 39 of file ROL_TrustRegionModel_U.hpp.
Referenced by applyHessian(), applyInvHessian(), applyPrecond(), TrustRegionModel_U(), and update().
|
private |
Definition at line 40 of file ROL_TrustRegionModel_U.hpp.
Referenced by applyPrecond(), TrustRegionModel_U(), and update().
|
private |
Definition at line 41 of file ROL_TrustRegionModel_U.hpp.
Referenced by applyHessian(), applyInvHessian(), TrustRegionModel_U(), update(), and validate().