ROL
|
Provides the interface to evaluate projected trust-region model functions from the Kelley-Sachs bound constrained trust-region algorithm. More...
#include <ROL_KelleySachsModel.hpp>
Public Member Functions | |
KelleySachsModel (Objective< Real > &obj, BoundConstraint< Real > &bnd, const Vector< Real > &x, const Vector< Real > &g, const Real eps) | |
KelleySachsModel (Objective< Real > &obj, BoundConstraint< Real > &bnd, const Vector< Real > &x, const Vector< Real > &g, const Real eps, const Teuchos::RCP< Secant< Real > > &secant, const bool useSecantPrecond, const bool useSecantHessVec) | |
Real | value (const Vector< Real > &s, Real &tol) |
Compute value. More... | |
void | gradient (Vector< Real > &g, const Vector< Real > &s, Real &tol) |
Compute gradient. More... | |
void | hessVec (Vector< Real > &Hv, const Vector< Real > &v, const Vector< Real > &s, Real &tol) |
Apply Hessian approximation to vector. More... | |
void | invHessVec (Vector< Real > &Hv, const Vector< Real > &v, const Vector< Real > &s, Real &tol) |
Apply inverse Hessian approximation to vector. More... | |
void | precond (Vector< Real > &Mv, const Vector< Real > &v, const Vector< Real > &x, Real &tol) |
Apply preconditioner to vector. More... | |
void | dualTransform (Vector< Real > &tv, const Vector< Real > &v) |
void | primalTransform (Vector< Real > &tv, const Vector< Real > &v) |
const Teuchos::RCP< BoundConstraint< Real > > | getBoundConstraint (void) const |
![]() | |
virtual | ~TrustRegionModel () |
TrustRegionModel (Objective< Real > &obj, const Vector< Real > &x, const Vector< Real > &g, const bool init=true) | |
TrustRegionModel (Objective< Real > &obj, const Vector< Real > &x, const Vector< Real > &g, const Teuchos::RCP< Secant< Real > > &secant, const bool useSecantPrecond, const bool useSecantHessVec) | |
virtual void | updatePredictedReduction (Real &pred, const Vector< Real > &s) |
virtual void | updateActualReduction (Real &ared, const Vector< Real > &s) |
virtual const Teuchos::RCP< const Vector< Real > > | getGradient (void) const |
virtual const Teuchos::RCP< const Vector< Real > > | getIterate (void) const |
virtual const Teuchos::RCP< Objective< Real > > | getObjective (void) const |
virtual void | update (const Vector< Real > &x, bool flag=true, int iter=-1) |
Update objective function. More... | |
![]() | |
virtual | ~Objective () |
virtual Real | dirDeriv (const Vector< Real > &x, const Vector< Real > &d, Real &tol) |
Compute directional derivative. 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 void | setParameter (const std::vector< Real > ¶m) |
Private Attributes | |
Teuchos::RCP< BoundConstraint< Real > > | bnd_ |
Teuchos::RCP< Secant< Real > > | secant_ |
Teuchos::RCP< Vector< Real > > | dual_ |
Teuchos::RCP< Vector< Real > > | prim_ |
const bool | useSecantPrecond_ |
const bool | useSecantHessVec_ |
Real | eps_ |
Additional Inherited Members | |
![]() | |
const std::vector< Real > | getParameter (void) const |
Provides the interface to evaluate projected trust-region model functions from the Kelley-Sachs bound constrained trust-region algorithm.
Definition at line 61 of file ROL_KelleySachsModel.hpp.
|
inline |
Definition at line 74 of file ROL_KelleySachsModel.hpp.
References ROL::KelleySachsModel< Real >::bnd_, ROL::Vector< Real >::clone(), ROL::KelleySachsModel< Real >::dual_, and ROL::KelleySachsModel< Real >::prim_.
|
inline |
Definition at line 83 of file ROL_KelleySachsModel.hpp.
References ROL::KelleySachsModel< Real >::bnd_, ROL::Vector< Real >::clone(), ROL::KelleySachsModel< Real >::dual_, and ROL::KelleySachsModel< Real >::prim_.
|
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. |
Reimplemented from ROL::TrustRegionModel< Real >.
Definition at line 94 of file ROL_KelleySachsModel.hpp.
References ROL::KelleySachsModel< Real >::bnd_, ROL::Vector< Real >::dual(), ROL::KelleySachsModel< Real >::dual_, ROL::KelleySachsModel< Real >::eps_, ROL::TrustRegionModel< Real >::getGradient(), ROL::TrustRegionModel< Real >::getIterate(), ROL::KelleySachsModel< Real >::hessVec(), and ROL::KelleySachsModel< Real >::prim_.
|
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::TrustRegionModel< Real >.
Definition at line 107 of file ROL_KelleySachsModel.hpp.
References ROL::KelleySachsModel< Real >::bnd_, ROL::KelleySachsModel< Real >::eps_, ROL::TrustRegionModel< Real >::getGradient(), ROL::TrustRegionModel< Real >::getIterate(), ROL::KelleySachsModel< Real >::hessVec(), ROL::Vector< Real >::plus(), and ROL::KelleySachsModel< Real >::prim_.
|
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::TrustRegionModel< Real >.
Definition at line 119 of file ROL_KelleySachsModel.hpp.
References ROL::KelleySachsModel< Real >::bnd_, ROL::KelleySachsModel< Real >::dual_, ROL::KelleySachsModel< Real >::eps_, ROL::TrustRegionModel< Real >::getGradient(), ROL::TrustRegionModel< Real >::getIterate(), ROL::TrustRegionModel< Real >::getObjective(), ROL::Vector< Real >::plus(), ROL::KelleySachsModel< Real >::prim_, ROL::KelleySachsModel< Real >::secant_, and ROL::KelleySachsModel< Real >::useSecantHessVec_.
Referenced by ROL::KelleySachsModel< Real >::gradient(), and ROL::KelleySachsModel< Real >::value().
|
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::TrustRegionModel< Real >.
Definition at line 145 of file ROL_KelleySachsModel.hpp.
References ROL::KelleySachsModel< Real >::bnd_, ROL::KelleySachsModel< Real >::dual_, ROL::KelleySachsModel< Real >::eps_, ROL::TrustRegionModel< Real >::getGradient(), ROL::TrustRegionModel< Real >::getIterate(), ROL::TrustRegionModel< Real >::getObjective(), ROL::Vector< Real >::plus(), ROL::KelleySachsModel< Real >::prim_, ROL::KelleySachsModel< Real >::secant_, and ROL::KelleySachsModel< Real >::useSecantHessVec_.
|
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::TrustRegionModel< Real >.
Definition at line 171 of file ROL_KelleySachsModel.hpp.
References ROL::KelleySachsModel< Real >::bnd_, ROL::KelleySachsModel< Real >::dual_, ROL::KelleySachsModel< Real >::eps_, ROL::TrustRegionModel< Real >::getGradient(), ROL::TrustRegionModel< Real >::getIterate(), ROL::TrustRegionModel< Real >::getObjective(), ROL::Vector< Real >::plus(), ROL::KelleySachsModel< Real >::prim_, ROL::KelleySachsModel< Real >::secant_, and ROL::KelleySachsModel< Real >::useSecantPrecond_.
|
inlinevirtual |
Reimplemented from ROL::TrustRegionModel< Real >.
Definition at line 197 of file ROL_KelleySachsModel.hpp.
References ROL::KelleySachsModel< Real >::bnd_, ROL::KelleySachsModel< Real >::eps_, ROL::TrustRegionModel< Real >::getGradient(), ROL::TrustRegionModel< Real >::getIterate(), and ROL::Vector< Real >::set().
|
inlinevirtual |
Reimplemented from ROL::TrustRegionModel< Real >.
Definition at line 205 of file ROL_KelleySachsModel.hpp.
References ROL::Vector< Real >::axpy(), ROL::KelleySachsModel< Real >::bnd_, ROL::TrustRegionModel< Real >::getIterate(), ROL::Vector< Real >::plus(), and ROL::Vector< Real >::set().
|
inlinevirtual |
Reimplemented from ROL::TrustRegionModel< Real >.
Definition at line 214 of file ROL_KelleySachsModel.hpp.
References ROL::KelleySachsModel< Real >::bnd_.
|
private |
Definition at line 63 of file ROL_KelleySachsModel.hpp.
Referenced by ROL::KelleySachsModel< Real >::dualTransform(), ROL::KelleySachsModel< Real >::getBoundConstraint(), ROL::KelleySachsModel< Real >::gradient(), ROL::KelleySachsModel< Real >::hessVec(), ROL::KelleySachsModel< Real >::invHessVec(), ROL::KelleySachsModel< Real >::KelleySachsModel(), ROL::KelleySachsModel< Real >::precond(), ROL::KelleySachsModel< Real >::primalTransform(), and ROL::KelleySachsModel< Real >::value().
|
private |
Definition at line 64 of file ROL_KelleySachsModel.hpp.
Referenced by ROL::KelleySachsModel< Real >::hessVec(), ROL::KelleySachsModel< Real >::invHessVec(), and ROL::KelleySachsModel< Real >::precond().
|
private |
|
private |
Definition at line 65 of file ROL_KelleySachsModel.hpp.
Referenced by ROL::KelleySachsModel< Real >::gradient(), ROL::KelleySachsModel< Real >::hessVec(), ROL::KelleySachsModel< Real >::invHessVec(), ROL::KelleySachsModel< Real >::KelleySachsModel(), ROL::KelleySachsModel< Real >::precond(), and ROL::KelleySachsModel< Real >::value().
|
private |
Definition at line 67 of file ROL_KelleySachsModel.hpp.
Referenced by ROL::KelleySachsModel< Real >::precond().
|
private |
Definition at line 68 of file ROL_KelleySachsModel.hpp.
Referenced by ROL::KelleySachsModel< Real >::hessVec(), and ROL::KelleySachsModel< Real >::invHessVec().
|
private |
Definition at line 70 of file ROL_KelleySachsModel.hpp.
Referenced by ROL::KelleySachsModel< Real >::dualTransform(), ROL::KelleySachsModel< Real >::gradient(), ROL::KelleySachsModel< Real >::hessVec(), ROL::KelleySachsModel< Real >::invHessVec(), ROL::KelleySachsModel< Real >::precond(), and ROL::KelleySachsModel< Real >::value().