ROL
|
Provides the interface for the progressive hedging objective. More...
#include <ROL_PH_Objective.hpp>
Public Member Functions | |
PH_Objective (const Ptr< Objective< Real >> &obj, const Ptr< Vector< Real >> &x, const Real penaltyParam, ParameterList &parlist) | |
void | setData (const Ptr< const Vector< Real >> &xbar, const Ptr< const Vector< Real >> &w, const Real penaltyParam) |
void | update (const Vector< Real > &x, bool flag=true, int iter=-1) |
Update objective function. More... | |
Real | value (const Vector< Real > &x, Real &tol) |
Compute value. More... | |
void | gradient (Vector< Real > &g, const Vector< Real > &x, Real &tol) |
Compute gradient. More... | |
void | hessVec (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol) |
Apply Hessian approximation to vector. More... | |
void | setParameter (const std::vector< Real > ¶m) |
Public Member Functions inherited from ROL::Objective< Real > | |
virtual | ~Objective () |
virtual Real | dirDeriv (const Vector< Real > &x, const Vector< Real > &d, Real &tol) |
Compute directional derivative. More... | |
virtual void | invHessVec (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol) |
Apply inverse Hessian approximation to vector. More... | |
virtual void | precond (Vector< Real > &Pv, const Vector< Real > &v, const Vector< Real > &x, Real &tol) |
Apply preconditioner to vector. 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... | |
Private Attributes | |
Ptr< Objective< Real > > | obj_ |
Ptr< const Vector< Real > > | xbar_ |
Ptr< const Vector< Real > > | w_ |
Ptr< Vector< Real > > | xprimal_ |
Real | penaltyParam_ |
Additional Inherited Members | |
Protected Member Functions inherited from ROL::Objective< Real > | |
const std::vector< Real > | getParameter (void) const |
Provides the interface for the progressive hedging objective.
Definition at line 63 of file ROL_PH_Objective.hpp.
|
inline |
Definition at line 72 of file ROL_PH_Objective.hpp.
References ROL::PH_Objective< Real >::obj_, and ROL::PH_Objective< Real >::xprimal_.
|
inline |
Definition at line 114 of file ROL_PH_Objective.hpp.
References ROL::PH_Objective< Real >::penaltyParam_, ROL::PH_Objective< Real >::w_, and ROL::PH_Objective< Real >::xbar_.
|
inlinevirtual |
Update objective function.
This function updates the objective function at new iterations.
[in] | x | is the new iterate. |
[in] | flag | is true if the iterate has changed. |
[in] | iter | is the outer algorithm iterations count. |
Reimplemented from ROL::Objective< Real >.
Definition at line 122 of file ROL_PH_Objective.hpp.
References ROL::PH_Objective< Real >::obj_.
|
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 >.
Definition at line 126 of file ROL_PH_Objective.hpp.
References ROL::Vector< Real >::dot(), ROL::PH_Objective< Real >::obj_, ROL::PH_Objective< Real >::penaltyParam_, ROL::PH_Objective< Real >::w_, ROL::PH_Objective< Real >::xbar_, and ROL::PH_Objective< Real >::xprimal_.
|
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 >.
Definition at line 136 of file ROL_PH_Objective.hpp.
References ROL::PH_Objective< Real >::obj_, ROL::PH_Objective< Real >::penaltyParam_, ROL::Vector< Real >::plus(), ROL::PH_Objective< Real >::w_, ROL::PH_Objective< Real >::xbar_, and ROL::PH_Objective< Real >::xprimal_.
|
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 >.
Definition at line 144 of file ROL_PH_Objective.hpp.
References ROL::Vector< Real >::axpy(), ROL::Vector< Real >::dual(), ROL::PH_Objective< Real >::obj_, and ROL::PH_Objective< Real >::penaltyParam_.
|
inlinevirtual |
Reimplemented from ROL::Objective< Real >.
Definition at line 149 of file ROL_PH_Objective.hpp.
References ROL::PH_Objective< Real >::obj_, and ROL::Objective< Real >::setParameter().
|
private |
Definition at line 65 of file ROL_PH_Objective.hpp.
Referenced by ROL::PH_Objective< Real >::gradient(), ROL::PH_Objective< Real >::hessVec(), ROL::PH_Objective< Real >::PH_Objective(), ROL::PH_Objective< Real >::setParameter(), ROL::PH_Objective< Real >::update(), and ROL::PH_Objective< Real >::value().
|
private |
Definition at line 66 of file ROL_PH_Objective.hpp.
Referenced by ROL::PH_Objective< Real >::gradient(), ROL::PH_Objective< Real >::setData(), and ROL::PH_Objective< Real >::value().
|
private |
Definition at line 66 of file ROL_PH_Objective.hpp.
Referenced by ROL::PH_Objective< Real >::gradient(), ROL::PH_Objective< Real >::setData(), and ROL::PH_Objective< Real >::value().
|
private |
Definition at line 67 of file ROL_PH_Objective.hpp.
Referenced by ROL::PH_Objective< Real >::gradient(), ROL::PH_Objective< Real >::PH_Objective(), and ROL::PH_Objective< Real >::value().
|
private |
Definition at line 68 of file ROL_PH_Objective.hpp.
Referenced by ROL::PH_Objective< Real >::gradient(), ROL::PH_Objective< Real >::hessVec(), ROL::PH_Objective< Real >::setData(), and ROL::PH_Objective< Real >::value().