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) |
Compute value. More... | |
virtual void | gradient (Vector< Real > &g, const Vector< Real > &s, Real &tol) |
Compute gradient. More... | |
virtual void | hessVec (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &s, Real &tol) |
Apply Hessian approximation to vector. More... | |
virtual void | invHessVec (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &s, Real &tol) |
Apply inverse Hessian approximation to vector. More... | |
virtual void | precond (Vector< Real > &Pv, const Vector< Real > &v, const Vector< Real > &s, Real &tol) |
Apply preconditioner to vector. More... | |
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::Objective< Real > | |
virtual | ~Objective () |
Objective () | |
virtual void | update (const Vector< Real > &x, UpdateType type, int iter=-1) |
Update objective function. More... | |
virtual void | update (const Vector< Real > &x, bool flag=true, int iter=-1) |
Update objective function. More... | |
virtual Real | dirDeriv (const Vector< Real > &x, const Vector< Real > &d, Real &tol) |
Compute directional derivative. More... | |
virtual void | prox (Vector< Real > &Pv, const Vector< Real > &v, Real t, Real &tol) |
Compute the proximity operator. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
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::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.
|
inlineprivate |
Definition at line 46 of file ROL_TrustRegionModel.hpp.
References ROL::Vector< Real >::dual(), ROL::TrustRegionModel< Real >::dual_, and ROL::TrustRegionModel< Real >::init_.
Referenced by ROL::TrustRegionModel< Real >::value().
|
inlineprotected |
Definition at line 57 of file ROL_TrustRegionModel.hpp.
References ROL::TrustRegionModel< Real >::obj_, ROL::TrustRegionModel< Real >::secant_, ROL::TrustRegionModel< Real >::useSecantHessVec_, and ROL::TrustRegionModel< Real >::x_.
Referenced by ROL::LinMoreModel< Real >::applyFullHessian(), ROL::TrustRegionModel< Real >::gradient(), ROL::TrustRegionModel< Real >::hessVec(), ROL::KelleySachsModel< Real >::hessVec(), ROL::ColemanLiModel< Real >::hessVec(), and ROL::TrustRegionModel< Real >::value().
|
inlineprotected |
Definition at line 66 of file ROL_TrustRegionModel.hpp.
References ROL::TrustRegionModel< Real >::obj_, ROL::TrustRegionModel< Real >::secant_, ROL::TrustRegionModel< Real >::useSecantHessVec_, and ROL::TrustRegionModel< Real >::x_.
Referenced by ROL::TrustRegionModel< Real >::invHessVec(), and ROL::KelleySachsModel< Real >::invHessVec().
|
inlineprotected |
Definition at line 75 of file ROL_TrustRegionModel.hpp.
References ROL::TrustRegionModel< Real >::obj_, ROL::TrustRegionModel< Real >::secant_, ROL::TrustRegionModel< Real >::useSecantPrecond_, and ROL::TrustRegionModel< Real >::x_.
Referenced by ROL::LinMoreModel< Real >::applyFullPrecond(), ROL::TrustRegionModel< Real >::precond(), and ROL::KelleySachsModel< Real >::precond().
|
inlinevirtual |
Reimplemented in ROL::ColemanLiModel< Real >.
Definition at line 104 of file ROL_TrustRegionModel.hpp.
References ROL::TrustRegionModel< Real >::bnd_, ROL::TrustRegionModel< Real >::g_, ROL::TrustRegionModel< Real >::obj_, ROL::TrustRegionModel< Real >::secant_, and ROL::TrustRegionModel< Real >::x_.
Referenced by ROL::ColemanLiModel< Real >::update().
|
inlinevirtual |
Compute value.
This function returns the objective function value.
[in] | x | is the current iterate. |
[in] | tol | is a tolerance for inexact objective function computation. |
Implements ROL::Objective< Real >.
Reimplemented in ROL::ColemanLiModel< Real >, and ROL::KelleySachsModel< Real >.
Definition at line 117 of file ROL_TrustRegionModel.hpp.
References ROL::TrustRegionModel< Real >::applyHessian(), ROL::Vector< Real >::dual(), ROL::TrustRegionModel< Real >::dual_, ROL::TrustRegionModel< Real >::g_, and ROL::TrustRegionModel< Real >::initialize().
Referenced by ROL::TrustRegion< Real >::update().
|
inlinevirtual |
Compute gradient.
This function returns the objective function gradient.
[out] | g | is the gradient. |
[in] | x | is the current iterate. |
[in] | tol | is a tolerance for inexact objective function computation. |
The default implementation is a finite-difference approximation based on the function value. This requires the definition of a basis \(\{\phi_i\}\) for the optimization vectors x and the definition of a basis \(\{\psi_j\}\) for the dual optimization vectors (gradient vectors g). The bases must be related through the Riesz map, i.e., \( R \{\phi_i\} = \{\psi_j\}\), and this must be reflected in the implementation of the ROL::Vector::dual() method.
Reimplemented from ROL::Objective< Real >.
Reimplemented in ROL::ColemanLiModel< Real >, and ROL::KelleySachsModel< Real >.
Definition at line 125 of file ROL_TrustRegionModel.hpp.
References ROL::TrustRegionModel< Real >::applyHessian(), ROL::TrustRegionModel< Real >::g_, and ROL::Vector< Real >::plus().
|
inlinevirtual |
Apply Hessian approximation to vector.
This function applies the Hessian of the objective function to the vector \(v\).
[out] | hv | is the the action of the Hessian on \(v\). |
[in] | v | is the direction vector. |
[in] | x | is the current iterate. |
[in] | tol | is a tolerance for inexact objective function computation. |
Reimplemented from ROL::Objective< Real >.
Reimplemented in ROL::ColemanLiModel< Real >, and ROL::KelleySachsModel< Real >.
Definition at line 130 of file ROL_TrustRegionModel.hpp.
References ROL::TrustRegionModel< Real >::applyHessian().
Referenced by ROL::DogLeg< Real >::run(), ROL::DoubleDogLeg< Real >::run(), and ROL::TruncatedCG< Real >::run().
|
inlinevirtual |
Apply inverse Hessian approximation to vector.
This function applies the inverse Hessian of the objective function to the vector \(v\).
[out] | hv | is the action of the inverse Hessian on \(v\). |
[in] | v | is the direction vector. |
[in] | x | is the current iterate. |
[in] | tol | is a tolerance for inexact objective function computation. |
Reimplemented from ROL::Objective< Real >.
Reimplemented in ROL::KelleySachsModel< Real >.
Definition at line 134 of file ROL_TrustRegionModel.hpp.
References ROL::TrustRegionModel< Real >::applyInvHessian().
Referenced by ROL::DogLeg< Real >::run(), and ROL::DoubleDogLeg< Real >::run().
|
inlinevirtual |
Apply preconditioner to vector.
This function applies a preconditioner for the Hessian of the objective function to the vector \(v\).
[out] | Pv | is the action of the Hessian preconditioner on \(v\). |
[in] | v | is the direction vector. |
[in] | x | is the current iterate. |
[in] | tol | is a tolerance for inexact objective function computation. |
Reimplemented from ROL::Objective< Real >.
Reimplemented in ROL::KelleySachsModel< Real >.
Definition at line 138 of file ROL_TrustRegionModel.hpp.
References ROL::TrustRegionModel< Real >::applyPrecond().
Referenced by ROL::TruncatedCG< Real >::run().
|
inlinevirtual |
Definition at line 148 of file ROL_TrustRegionModel.hpp.
References ROL::TrustRegionModel< Real >::g_.
Referenced by ROL::ColemanLiModel< Real >::computeCauchyPoint(), ROL::ColemanLiModel< Real >::constructC(), ROL::ColemanLiModel< Real >::constructInverseD(), ROL::ColemanLiModel< 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::TruncatedCG< Real >::run(), ROL::LinMore< Real >::run(), ROL::TrustRegion< Real >::update(), and ROL::ColemanLiModel< Real >::value().
|
inlinevirtual |
Definition at line 152 of file ROL_TrustRegionModel.hpp.
References ROL::TrustRegionModel< Real >::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 ROL::TrustRegionModel< Real >::obj_.
|
inlinevirtual |
Definition at line 160 of file ROL_TrustRegionModel.hpp.
References ROL::TrustRegionModel< Real >::bnd_.
Referenced by ROL::LinMoreModel< Real >::applyFreeHessian(), ROL::LinMoreModel< Real >::applyFreePrecond(), ROL::ColemanLiModel< Real >::constructC(), ROL::ColemanLiModel< Real >::constructInverseD(), ROL::LinMore< Real >::dbreakpt(), ROL::LinMore< Real >::dgpstep(), ROL::ColemanLiModel< Real >::getScalarBounds(), ROL::KelleySachsModel< Real >::primalTransform(), ROL::KelleySachsModel< Real >::pruneBindingConstraints(), and ROL::KelleySachsModel< Real >::pruneNonbindingConstraints().
|
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 ROL::TrustRegionModel< Real >::applyHessian(), ROL::TrustRegionModel< Real >::applyInvHessian(), ROL::TrustRegionModel< Real >::applyPrecond(), ROL::TrustRegionModel< Real >::getObjective(), and ROL::TrustRegionModel< Real >::update().
|
private |
Definition at line 36 of file ROL_TrustRegionModel.hpp.
Referenced by ROL::TrustRegionModel< Real >::getBoundConstraint(), and ROL::TrustRegionModel< Real >::update().
|
private |
Definition at line 37 of file ROL_TrustRegionModel.hpp.
Referenced by ROL::TrustRegionModel< Real >::applyHessian(), ROL::TrustRegionModel< Real >::applyInvHessian(), ROL::TrustRegionModel< Real >::applyPrecond(), ROL::TrustRegionModel< Real >::getIterate(), and ROL::TrustRegionModel< Real >::update().
|
private |
Definition at line 37 of file ROL_TrustRegionModel.hpp.
Referenced by ROL::TrustRegionModel< Real >::getGradient(), ROL::TrustRegionModel< Real >::gradient(), ROL::TrustRegionModel< Real >::update(), and ROL::TrustRegionModel< Real >::value().
|
private |
Definition at line 38 of file ROL_TrustRegionModel.hpp.
Referenced by ROL::TrustRegionModel< Real >::initialize(), and ROL::TrustRegionModel< Real >::value().
|
private |
Definition at line 39 of file ROL_TrustRegionModel.hpp.
Referenced by ROL::TrustRegionModel< Real >::applyHessian(), ROL::TrustRegionModel< Real >::applyInvHessian(), ROL::TrustRegionModel< Real >::applyPrecond(), and ROL::TrustRegionModel< Real >::update().
|
private |
Definition at line 41 of file ROL_TrustRegionModel.hpp.
Referenced by ROL::TrustRegionModel< Real >::applyPrecond().
|
private |
Definition at line 42 of file ROL_TrustRegionModel.hpp.
Referenced by ROL::TrustRegionModel< Real >::applyHessian(), and ROL::TrustRegionModel< Real >::applyInvHessian().
|
private |
Definition at line 44 of file ROL_TrustRegionModel.hpp.
Referenced by ROL::TrustRegionModel< Real >::initialize().