ROL
Public Member Functions | List of all members
Objective_Sacado< Real > Class Template Reference
+ Inheritance diagram for Objective_Sacado< Real >:

Public Member Functions

 Objective_Sacado (void)
 
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...
 
- Public Member Functions inherited from ROL::Objective< Real >
virtual ~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 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)
 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)
 Finite-difference gradient check. 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)
 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)
 Finite-difference Hessian-applied-to-vector check. 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...
 

Detailed Description

template<class Real>
class Objective_Sacado< Real >

Objective Class Definition begin

Definition at line 205 of file sacado/example_01.cpp.

Constructor & Destructor Documentation

template<class Real>
Objective_Sacado< Real >::Objective_Sacado ( void  )
inline

Definition at line 209 of file sacado/example_01.cpp.

Member Function Documentation

template<class Real>
Real Objective_Sacado< Real >::value ( const Vector< Real > &  x,
Real &  tol 
)
inlinevirtual

Compute value.

This function returns the objective function value.

Parameters
[in]xis the current iterate.
[in]tolis a tolerance for inexact objective function computation.

Implements ROL::Objective< Real >.

Definition at line 212 of file sacado/example_01.cpp.

References ROL::StdVector< Real, Element >::getVector(), and Zakharov().

template<class Real>
void Objective_Sacado< Real >::gradient ( Vector< Real > &  g,
const Vector< Real > &  x,
Real &  tol 
)
inlinevirtual

Compute gradient.

This function returns the objective function gradient.

Parameters
[out]gis the gradient.
[in]xis the current iterate.
[in]tolis 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 223 of file sacado/example_01.cpp.

References objgrad().

template<class Real>
void Objective_Sacado< Real >::hessVec ( Vector< Real > &  hv,
const Vector< Real > &  v,
const Vector< Real > &  x,
Real &  tol 
)
inlinevirtual

Apply Hessian approximation to vector.

This function applies the Hessian of the objective function to the vector \(v\).

Parameters
[out]hvis the the action of the Hessian on \(v\).
[in]vis the direction vector.
[in]xis the current iterate.
[in]tolis a tolerance for inexact objective function computation.

Reimplemented from ROL::Objective< Real >.

Definition at line 235 of file sacado/example_01.cpp.

References applyHessian().


The documentation for this class was generated from the following file: