ROL
|
#include <ROL_ScaledObjective.hpp>
Public Member Functions | |
ScaledObjective (const Ptr< Objective< Real >> &obj, Real scale) | |
void | update (const Vector< Real > &x, UpdateType type, int iter=-1) override |
Update objective function. More... | |
void | setParameter (const std::vector< Real > ¶m) override |
Real | value (const Vector< Real > &x, Real &tol) override |
Compute value. More... | |
void | gradient (Vector< Real > &g, const Vector< Real > &x, Real &tol) override |
Compute gradient. More... | |
void | hessVec (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol) override |
Apply Hessian approximation to vector. More... | |
void | invHessVec (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol) override |
Apply inverse Hessian approximation to vector. More... | |
void | precond (Vector< Real > &Pv, const Vector< Real > &v, const Vector< Real > &x, Real &tol) override |
Apply preconditioner to vector. More... | |
Public Member Functions inherited from ROL::Objective< Real > | |
virtual | ~Objective () |
Objective () | |
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... | |
Private Attributes | |
const Ptr< Objective< Real > > | obj_ |
const Real | scale_ |
Additional Inherited Members | |
Protected Member Functions inherited from ROL::Objective< Real > | |
const std::vector< Real > | getParameter (void) const |
Definition at line 18 of file ROL_ScaledObjective.hpp.
|
inline |
Definition at line 24 of file ROL_ScaledObjective.hpp.
|
overridevirtual |
Update objective function.
This function updates the objective function at new iterations.
[in] | x | is the new iterate. |
[in] | type | is the type of update requested. |
[in] | iter | is the outer algorithm iterations count. |
Reimplemented from ROL::Objective< Real >.
Definition at line 16 of file ROL_ScaledObjective_Def.hpp.
References obj_.
|
overridevirtual |
Reimplemented from ROL::Objective< Real >.
Definition at line 21 of file ROL_ScaledObjective_Def.hpp.
References obj_, and ROL::Objective< Real >::setParameter().
|
overridevirtual |
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 >.
Definition at line 27 of file ROL_ScaledObjective_Def.hpp.
References obj_.
|
overridevirtual |
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 >.
Definition at line 32 of file ROL_ScaledObjective_Def.hpp.
References obj_, and ROL::Vector< Real >::scale().
|
overridevirtual |
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 >.
Definition at line 38 of file ROL_ScaledObjective_Def.hpp.
References obj_, and ROL::Vector< Real >::scale().
|
overridevirtual |
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 >.
Definition at line 44 of file ROL_ScaledObjective_Def.hpp.
References obj_, and ROL::Vector< Real >::scale().
|
overridevirtual |
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 >.
Definition at line 50 of file ROL_ScaledObjective_Def.hpp.
References obj_, and ROL::Vector< Real >::scale().
|
private |
Definition at line 20 of file ROL_ScaledObjective.hpp.
|
private |
Definition at line 21 of file ROL_ScaledObjective.hpp.