|
ROL
|
Provides the interface to evaluate trust-region model functions. More...
#include <ROL_TrustRegionModel.hpp>
Public Member Functions | |
| virtual | ~TrustRegionModel () |
| TrustRegionModel (Objective< Real > &obj, BoundConstraint< Real > &bnd, const Vector< Real > &x, const Vector< Real > &g, const Ptr< Secant< Real > > &secant=nullPtr, const bool useSecantPrecond=false, const bool useSecantHessVec=false) | |
| virtual void | update (Objective< Real > &obj, BoundConstraint< Real > &bnd, const Vector< Real > &x, const Vector< Real > &g, const Ptr< Secant< Real > > &secant=nullPtr) |
| virtual Real | value (const Vector< Real > &s, Real &tol) |
| virtual void | gradient (Vector< Real > &g, const Vector< Real > &s, Real &tol) |
| virtual void | hessVec (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &s, Real &tol) |
| virtual void | invHessVec (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &s, Real &tol) |
| virtual void | precond (Vector< Real > &Pv, const Vector< Real > &v, const Vector< Real > &s, Real &tol) |
| 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 |
| virtual const Ptr< BoundConstraint< Real > > | getBoundConstraint (void) const |
| virtual void | dualTransform (Vector< Real > &tv, const Vector< Real > &v) |
| virtual void | primalTransform (Vector< Real > &tv, const Vector< Real > &v) |
| virtual void | updatePredictedReduction (Real &pred, const Vector< Real > &s) |
| virtual void | updateActualReduction (Real &ared, const Vector< Real > &s) |
| 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 Member Functions | |
| void | initialize (const Vector< Real > &s) |
Private Attributes | |
| Ptr< Objective< Real > > | obj_ |
| Ptr< BoundConstraint< Real > > | bnd_ |
| Ptr< const Vector< Real > > | x_ |
| Ptr< const Vector< Real > > | g_ |
| Ptr< Vector< Real > > | dual_ |
| Ptr< Secant< Real > > | secant_ |
| const bool | useSecantPrecond_ |
| const bool | useSecantHessVec_ |
| bool | init_ |
Provides the interface to evaluate trust-region model functions.
ROL::TrustRegionModel 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 33 of file ROL_TrustRegionModel.hpp.
|
inlinevirtual |
Definition at line 89 of file ROL_TrustRegionModel.hpp.
|
inline |
Definition at line 91 of file ROL_TrustRegionModel.hpp.
References bnd_, g_, init_, obj_, ROL::ROL::Objective< Real >::Objective(), secant_, useSecantHessVec_, useSecantPrecond_, and x_.
Referenced by ROL::ColemanLiModel< Real >::ColemanLiModel(), ROL::KelleySachsModel< Real >::KelleySachsModel(), and ROL::LinMoreModel< Real >::LinMoreModel().
|
inlineprivate |
Definition at line 46 of file ROL_TrustRegionModel.hpp.
References ROL::Vector< Real >::clone(), ROL::Vector< Real >::dual(), dual_, and init_.
Referenced by value().
|
inlineprotected |
Definition at line 57 of file ROL_TrustRegionModel.hpp.
References obj_, secant_, useSecantHessVec_, and x_.
Referenced by ROL::LinMoreModel< Real >::applyFullHessian(), gradient(), ROL::ColemanLiModel< Real >::hessVec(), ROL::KelleySachsModel< Real >::hessVec(), hessVec(), and value().
|
inlineprotected |
Definition at line 66 of file ROL_TrustRegionModel.hpp.
References obj_, secant_, useSecantHessVec_, and x_.
Referenced by ROL::KelleySachsModel< Real >::invHessVec(), and invHessVec().
|
inlineprotected |
Definition at line 75 of file ROL_TrustRegionModel.hpp.
References obj_, secant_, useSecantPrecond_, and x_.
Referenced by ROL::LinMoreModel< Real >::applyFullPrecond(), ROL::KelleySachsModel< Real >::precond(), and precond().
|
inlinevirtual |
Reimplemented in ROL::ColemanLiModel< Real >.
Definition at line 104 of file ROL_TrustRegionModel.hpp.
References bnd_, g_, obj_, ROL::ROL::Objective< Real >::Objective(), secant_, and x_.
Referenced by ROL::ColemanLiModel< Real >::update().
|
inlinevirtual |
Reimplemented in ROL::ColemanLiModel< Real >, and ROL::KelleySachsModel< Real >.
Definition at line 117 of file ROL_TrustRegionModel.hpp.
References applyHessian(), ROL::Vector< Real >::dual(), dual_, g_, and initialize().
Referenced by ROL::TrustRegion< Real >::update().
|
inlinevirtual |
Reimplemented in ROL::ColemanLiModel< Real >, and ROL::KelleySachsModel< Real >.
Definition at line 125 of file ROL_TrustRegionModel.hpp.
References applyHessian(), g_, and ROL::Vector< Real >::plus().
|
inlinevirtual |
Reimplemented in ROL::ColemanLiModel< Real >, and ROL::KelleySachsModel< Real >.
Definition at line 130 of file ROL_TrustRegionModel.hpp.
References applyHessian().
Referenced by ROL::DogLeg< Real >::run(), ROL::DoubleDogLeg< Real >::run(), and ROL::TruncatedCG< Real >::run().
|
inlinevirtual |
Reimplemented in ROL::KelleySachsModel< Real >.
Definition at line 134 of file ROL_TrustRegionModel.hpp.
References applyInvHessian().
Referenced by ROL::DogLeg< Real >::run(), and ROL::DoubleDogLeg< Real >::run().
|
inlinevirtual |
Reimplemented in ROL::KelleySachsModel< Real >.
Definition at line 138 of file ROL_TrustRegionModel.hpp.
References applyPrecond().
Referenced by ROL::TruncatedCG< Real >::run().
|
inlinevirtual |
Definition at line 148 of file ROL_TrustRegionModel.hpp.
References g_.
Referenced by ROL::ColemanLiModel< Real >::computeCauchyPoint(), ROL::ColemanLiModel< Real >::constructC(), ROL::ColemanLiModel< Real >::constructInverseD(), ROL::ColemanLiModel< Real >::gradient(), ROL::KelleySachsModel< Real >::gradient(), ROL::ColemanLiModel< Real >::hessVec(), ROL::ColemanLiModel< Real >::minimize1D(), ROL::KelleySachsModel< Real >::pruneBindingConstraints(), ROL::KelleySachsModel< Real >::pruneNonbindingConstraints(), ROL::DogLeg< Real >::run(), ROL::DoubleDogLeg< Real >::run(), ROL::LinMore< Real >::run(), ROL::TruncatedCG< Real >::run(), ROL::TrustRegion< Real >::update(), ROL::ColemanLiModel< Real >::value(), and ROL::KelleySachsModel< Real >::value().
|
inlinevirtual |
Definition at line 152 of file ROL_TrustRegionModel.hpp.
References x_.
Referenced by ROL::ColemanLiModel< Real >::computeAlpha(), ROL::ColemanLiModel< Real >::computeFullReflectiveStep(), ROL::ColemanLiModel< Real >::computeReflectiveStep(), ROL::ColemanLiModel< Real >::constructInverseD(), ROL::ColemanLiModel< Real >::getScalarBounds(), ROL::ColemanLiModel< Real >::isStrictlyFeasibleStep(), ROL::KelleySachsModel< Real >::primalTransform(), ROL::KelleySachsModel< Real >::pruneBindingConstraints(), ROL::KelleySachsModel< Real >::pruneNonbindingConstraints(), and ROL::LinMore< Real >::run().
|
inlinevirtual |
Definition at line 156 of file ROL_TrustRegionModel.hpp.
References obj_.
|
inlinevirtual |
Definition at line 160 of file ROL_TrustRegionModel.hpp.
References bnd_.
Referenced by ROL::LinMoreModel< Real >::applyFreeHessian(), ROL::LinMoreModel< Real >::applyFreePrecond(), ROL::ColemanLiModel< Real >::computeAlpha(), ROL::ColemanLiModel< Real >::computeFullReflectiveStep(), ROL::ColemanLiModel< Real >::computeReflectiveStep(), ROL::ColemanLiModel< Real >::constructC(), ROL::ColemanLiModel< Real >::constructInverseD(), ROL::LinMore< Real >::dbreakpt(), ROL::LinMore< Real >::dgpstep(), ROL::ColemanLiModel< Real >::getScalarBounds(), ROL::ColemanLiModel< Real >::isStrictlyFeasibleStep(), ROL::KelleySachsModel< Real >::primalTransform(), ROL::KelleySachsModel< Real >::pruneBindingConstraints(), ROL::KelleySachsModel< Real >::pruneNonbindingConstraints(), and ROL::LinMore< Real >::run().
|
inlinevirtual |
Reimplemented in ROL::ColemanLiModel< Real >, and ROL::KelleySachsModel< Real >.
Definition at line 170 of file ROL_TrustRegionModel.hpp.
References ROL::Vector< Real >::set().
Referenced by ROL::DogLeg< Real >::run(), ROL::DoubleDogLeg< Real >::run(), ROL::TruncatedCG< Real >::run(), and ROL::TrustRegion< Real >::update().
|
inlinevirtual |
Reimplemented in ROL::ColemanLiModel< Real >, and ROL::KelleySachsModel< Real >.
Definition at line 174 of file ROL_TrustRegionModel.hpp.
References ROL::Vector< Real >::set().
Referenced by ROL::DogLeg< Real >::run(), ROL::DoubleDogLeg< Real >::run(), and ROL::TruncatedCG< Real >::run().
|
inlinevirtual |
Reimplemented in ROL::ColemanLiModel< Real >.
Definition at line 178 of file ROL_TrustRegionModel.hpp.
Referenced by ROL::TrustRegion< Real >::update().
|
inlinevirtual |
Reimplemented in ROL::ColemanLiModel< Real >.
Definition at line 180 of file ROL_TrustRegionModel.hpp.
Referenced by ROL::TrustRegion< Real >::update().
|
private |
Definition at line 35 of file ROL_TrustRegionModel.hpp.
Referenced by applyHessian(), applyInvHessian(), applyPrecond(), getObjective(), TrustRegionModel(), and update().
|
private |
Definition at line 36 of file ROL_TrustRegionModel.hpp.
Referenced by getBoundConstraint(), TrustRegionModel(), and update().
|
private |
Definition at line 37 of file ROL_TrustRegionModel.hpp.
Referenced by applyHessian(), applyInvHessian(), applyPrecond(), getIterate(), TrustRegionModel(), and update().
|
private |
Definition at line 37 of file ROL_TrustRegionModel.hpp.
Referenced by getGradient(), gradient(), TrustRegionModel(), update(), and value().
|
private |
Definition at line 38 of file ROL_TrustRegionModel.hpp.
Referenced by initialize(), and value().
|
private |
Definition at line 39 of file ROL_TrustRegionModel.hpp.
Referenced by applyHessian(), applyInvHessian(), applyPrecond(), TrustRegionModel(), and update().
|
private |
Definition at line 41 of file ROL_TrustRegionModel.hpp.
Referenced by applyPrecond(), and TrustRegionModel().
|
private |
Definition at line 42 of file ROL_TrustRegionModel.hpp.
Referenced by applyHessian(), applyInvHessian(), and TrustRegionModel().
|
private |
Definition at line 44 of file ROL_TrustRegionModel.hpp.
Referenced by initialize(), and TrustRegionModel().