ROL
|
Provides the interface for the progressive hedging probability objective. More...
#include <ROL_PH_ProbObjective.hpp>
Public Member Functions | |
PH_ProbObjective (const Ptr< Objective< Real >> &obj, ParameterList &parlist) | |
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 Member Functions | |
void | getValue (const Vector< Real > &x, Real &tol) |
void | getGradient (const Vector< Real > &x, Real &tol) |
Real | smoothHeaviside (const Real x, const int deriv=0) const |
Private Attributes | |
const Ptr< Objective< Real > > | obj_ |
Real | threshold_ |
Real | eps_ |
bool | isValueComputed_ |
Real | val_ |
bool | isGradientInitialized_ |
bool | isGradientComputed_ |
Ptr< Vector< Real > > | g_ |
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 probability objective.
Definition at line 58 of file ROL_PH_ProbObjective.hpp.
|
inline |
Definition at line 109 of file ROL_PH_ProbObjective.hpp.
References ROL::PH_ProbObjective< Real >::eps_, and ROL::PH_ProbObjective< Real >::threshold_.
|
inlineprivate |
Definition at line 71 of file ROL_PH_ProbObjective.hpp.
References ROL::PH_ProbObjective< Real >::isValueComputed_, ROL::PH_ProbObjective< Real >::obj_, and ROL::PH_ProbObjective< Real >::val_.
Referenced by ROL::PH_ProbObjective< Real >::gradient(), ROL::PH_ProbObjective< Real >::hessVec(), and ROL::PH_ProbObjective< Real >::value().
|
inlineprivate |
Definition at line 78 of file ROL_PH_ProbObjective.hpp.
References ROL::Vector< Real >::dual(), ROL::PH_ProbObjective< Real >::g_, ROL::PH_ProbObjective< Real >::isGradientComputed_, ROL::PH_ProbObjective< Real >::isGradientInitialized_, and ROL::PH_ProbObjective< Real >::obj_.
Referenced by ROL::PH_ProbObjective< Real >::gradient(), and ROL::PH_ProbObjective< Real >::hessVec().
|
inlineprivate |
Definition at line 89 of file ROL_PH_ProbObjective.hpp.
References ROL::PH_ProbObjective< Real >::eps_.
Referenced by ROL::PH_ProbObjective< Real >::gradient(), ROL::PH_ProbObjective< Real >::hessVec(), and ROL::PH_ProbObjective< Real >::value().
|
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 120 of file ROL_PH_ProbObjective.hpp.
References ROL::PH_ProbObjective< Real >::isGradientComputed_, ROL::PH_ProbObjective< Real >::isValueComputed_, and ROL::PH_ProbObjective< 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_ProbObjective.hpp.
References ROL::PH_ProbObjective< Real >::getValue(), ROL::PH_ProbObjective< Real >::smoothHeaviside(), ROL::PH_ProbObjective< Real >::threshold_, and ROL::PH_ProbObjective< Real >::val_.
|
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 135 of file ROL_PH_ProbObjective.hpp.
References ROL::PH_ProbObjective< Real >::g_, ROL::PH_ProbObjective< Real >::getGradient(), ROL::PH_ProbObjective< Real >::getValue(), ROL::Vector< Real >::scale(), ROL::Vector< Real >::set(), ROL::PH_ProbObjective< Real >::smoothHeaviside(), ROL::PH_ProbObjective< Real >::threshold_, ROL::PH_ProbObjective< Real >::val_, and ROL::Vector< Real >::zero().
|
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 147 of file ROL_PH_ProbObjective.hpp.
References ROL::Vector< Real >::axpy(), ROL::Vector< Real >::dot(), ROL::PH_ProbObjective< Real >::g_, ROL::PH_ProbObjective< Real >::getGradient(), ROL::PH_ProbObjective< Real >::getValue(), ROL::PH_ProbObjective< Real >::obj_, ROL::Vector< Real >::scale(), ROL::PH_ProbObjective< Real >::smoothHeaviside(), ROL::PH_ProbObjective< Real >::threshold_, ROL::PH_ProbObjective< Real >::val_, and ROL::Vector< Real >::zero().
|
inlinevirtual |
Reimplemented from ROL::Objective< Real >.
Definition at line 165 of file ROL_PH_ProbObjective.hpp.
References ROL::PH_ProbObjective< Real >::obj_, and ROL::Objective< Real >::setParameter().
|
private |
|
private |
Definition at line 61 of file ROL_PH_ProbObjective.hpp.
Referenced by ROL::PH_ProbObjective< Real >::gradient(), ROL::PH_ProbObjective< Real >::hessVec(), ROL::PH_ProbObjective< Real >::PH_ProbObjective(), and ROL::PH_ProbObjective< Real >::value().
|
private |
Definition at line 62 of file ROL_PH_ProbObjective.hpp.
Referenced by ROL::PH_ProbObjective< Real >::PH_ProbObjective(), and ROL::PH_ProbObjective< Real >::smoothHeaviside().
|
private |
Definition at line 64 of file ROL_PH_ProbObjective.hpp.
Referenced by ROL::PH_ProbObjective< Real >::getValue(), and ROL::PH_ProbObjective< Real >::update().
|
private |
Definition at line 65 of file ROL_PH_ProbObjective.hpp.
Referenced by ROL::PH_ProbObjective< Real >::getValue(), ROL::PH_ProbObjective< Real >::gradient(), ROL::PH_ProbObjective< Real >::hessVec(), and ROL::PH_ProbObjective< Real >::value().
|
private |
Definition at line 67 of file ROL_PH_ProbObjective.hpp.
Referenced by ROL::PH_ProbObjective< Real >::getGradient().
|
private |
Definition at line 68 of file ROL_PH_ProbObjective.hpp.
Referenced by ROL::PH_ProbObjective< Real >::getGradient(), and ROL::PH_ProbObjective< Real >::update().
|
private |
Definition at line 69 of file ROL_PH_ProbObjective.hpp.
Referenced by ROL::PH_ProbObjective< Real >::getGradient(), ROL::PH_ProbObjective< Real >::gradient(), and ROL::PH_ProbObjective< Real >::hessVec().