ROL
Public Member Functions | Protected Member Functions | Private Attributes | List of all members
ROL::TrustRegionModel_U< Real > Class Template Reference

Provides the interface to evaluate trust-region model functions. More...

#include <ROL_TrustRegionModel_U.hpp>

+ Inheritance diagram for ROL::TrustRegionModel_U< Real >:

Public Member Functions

virtual ~TrustRegionModel_U ()
 
 TrustRegionModel_U (ParameterList &list, const Ptr< Secant< Real >> &secant=nullPtr, ESecantMode mode=SECANTMODE_BOTH)
 
void initialize (const Vector< Real > &x, const Vector< Real > &g)
 
void validate (Objective< Real > &obj, const Vector< Real > &x, const Vector< Real > &g, ETrustRegionU etr)
 
virtual void setData (Objective< Real > &obj, const Vector< Real > &x, const Vector< Real > &g, Real &tol)
 
void update (const Vector< Real > &x, const Vector< Real > &s, const Vector< Real > &gold, const Vector< Real > &gnew, const Real snorm, const int iter)
 
virtual Real value (const Vector< Real > &s, Real &tol) override
 Compute value. More...
 
virtual void gradient (Vector< Real > &g, const Vector< Real > &s, Real &tol) override
 Compute gradient. More...
 
virtual void hessVec (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &s, Real &tol) override
 Apply Hessian approximation to vector. More...
 
virtual void invHessVec (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &s, Real &tol) override
 Apply inverse Hessian approximation to vector. More...
 
virtual void precond (Vector< Real > &Pv, const Vector< Real > &v, const Vector< Real > &s, Real &tol) override
 Apply preconditioner to vector. More...
 
virtual const Ptr< const
Vector< Real > > 
getGradient (void) const
 
virtual const Ptr< const
Vector< Real > > 
getIterate (void) const
 
virtual const Ptr< Objective
< Real > > 
getObjective (void) const
 
- Public Member Functions inherited from ROL::Objective< Real >
virtual ~Objective ()
 
 Objective ()
 
virtual void update (const Vector< Real > &x, UpdateType type, int iter=-1)
 Update objective function. More...
 
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)
 
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 > &param)
 

Protected Member Functions

void applyHessian (Vector< Real > &hv, const Vector< Real > &v, Real &tol)
 
void applyInvHessian (Vector< Real > &hv, const Vector< Real > &v, Real &tol)
 
void applyPrecond (Vector< Real > &Pv, const Vector< Real > &v, Real &tol)
 
- Protected Member Functions inherited from ROL::Objective< Real >
const std::vector< Real > getParameter (void) const
 

Private Attributes

Ptr< Objective< Real > > obj_
 
Ptr< const Vector< Real > > x_
 
Ptr< const Vector< Real > > g_
 
Ptr< Vector< Real > > dual_
 
Real tol_
 
Ptr< Secant< Real > > secant_
 
bool useSecantPrecond_
 
bool useSecantHessVec_
 

Detailed Description

template<typename Real>
class ROL::TrustRegionModel_U< Real >

Provides the interface to evaluate trust-region model functions.

ROL::TrustRegionModel_U provides the interface to implement a number of trust-region models for unconstrained and constrained optimization. The default implementation is the standard quadratic trust region model for unconstrained optimization.


Definition at line 66 of file ROL_TrustRegionModel_U.hpp.

Constructor & Destructor Documentation

template<typename Real>
virtual ROL::TrustRegionModel_U< Real >::~TrustRegionModel_U ( )
inlinevirtual

Definition at line 113 of file ROL_TrustRegionModel_U.hpp.

template<typename Real>
ROL::TrustRegionModel_U< Real >::TrustRegionModel_U ( ParameterList &  list,
const Ptr< Secant< Real >> &  secant = nullPtr,
ESecantMode  mode = SECANTMODE_BOTH 
)
inline

Member Function Documentation

template<typename Real>
void ROL::TrustRegionModel_U< Real >::applyHessian ( Vector< Real > &  hv,
const Vector< Real > &  v,
Real &  tol 
)
inlineprotected
template<typename Real>
void ROL::TrustRegionModel_U< Real >::applyInvHessian ( Vector< Real > &  hv,
const Vector< Real > &  v,
Real &  tol 
)
inlineprotected
template<typename Real>
void ROL::TrustRegionModel_U< Real >::applyPrecond ( Vector< Real > &  Pv,
const Vector< Real > &  v,
Real &  tol 
)
inlineprotected
template<typename Real>
void ROL::TrustRegionModel_U< Real >::initialize ( const Vector< Real > &  x,
const Vector< Real > &  g 
)
inline
template<typename Real>
void ROL::TrustRegionModel_U< Real >::validate ( Objective< Real > &  obj,
const Vector< Real > &  x,
const Vector< Real > &  g,
ETrustRegionU  etr 
)
inline
template<typename Real>
virtual void ROL::TrustRegionModel_U< Real >::setData ( Objective< Real > &  obj,
const Vector< Real > &  x,
const Vector< Real > &  g,
Real &  tol 
)
inlinevirtual
template<typename Real>
void ROL::TrustRegionModel_U< Real >::update ( const Vector< Real > &  x,
const Vector< Real > &  s,
const Vector< Real > &  gold,
const Vector< Real > &  gnew,
const Real  snorm,
const int  iter 
)
inline
template<typename Real>
virtual Real ROL::TrustRegionModel_U< Real >::value ( const Vector< Real > &  x,
Real &  tol 
)
inlineoverridevirtual

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 172 of file ROL_TrustRegionModel_U.hpp.

References ROL::TrustRegionModel_U< Real >::applyHessian(), ROL::TrustRegionModel_U< Real >::dual_, and ROL::TrustRegionModel_U< Real >::g_.

template<typename Real>
virtual void ROL::TrustRegionModel_U< Real >::gradient ( Vector< Real > &  g,
const Vector< Real > &  x,
Real &  tol 
)
inlineoverridevirtual

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 179 of file ROL_TrustRegionModel_U.hpp.

References ROL::TrustRegionModel_U< Real >::applyHessian(), ROL::TrustRegionModel_U< Real >::g_, and ROL::Vector< Real >::plus().

template<typename Real>
virtual void ROL::TrustRegionModel_U< Real >::hessVec ( Vector< Real > &  hv,
const Vector< Real > &  v,
const Vector< Real > &  x,
Real &  tol 
)
inlineoverridevirtual
template<typename Real>
virtual void ROL::TrustRegionModel_U< Real >::invHessVec ( Vector< Real > &  hv,
const Vector< Real > &  v,
const Vector< Real > &  x,
Real &  tol 
)
inlineoverridevirtual

Apply inverse Hessian approximation to vector.

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

Parameters
[out]hvis the action of the inverse 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 188 of file ROL_TrustRegionModel_U.hpp.

References ROL::TrustRegionModel_U< Real >::applyInvHessian().

Referenced by ROL::DoubleDogLeg_U< Real >::solve(), and ROL::DogLeg_U< Real >::solve().

template<typename Real>
virtual void ROL::TrustRegionModel_U< Real >::precond ( Vector< Real > &  Pv,
const Vector< Real > &  v,
const Vector< Real > &  x,
Real &  tol 
)
inlineoverridevirtual

Apply preconditioner to vector.

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

Parameters
[out]Pvis the action of the Hessian preconditioner 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 192 of file ROL_TrustRegionModel_U.hpp.

References ROL::TrustRegionModel_U< Real >::applyPrecond().

Referenced by ROL::TypeB::KelleySachsAlgorithm< Real >::applyFreePrecond(), ROL::TypeB::LinMoreAlgorithm< Real >::applyFreePrecond(), ROL::TypeB::ColemanLiAlgorithm< Real >::applyPrecond(), and ROL::TruncatedCG_U< Real >::solve().

template<typename Real>
virtual const Ptr<const Vector<Real> > ROL::TrustRegionModel_U< Real >::getGradient ( void  ) const
inlinevirtual
template<typename Real>
virtual const Ptr<const Vector<Real> > ROL::TrustRegionModel_U< Real >::getIterate ( void  ) const
inlinevirtual

Definition at line 206 of file ROL_TrustRegionModel_U.hpp.

References ROL::TrustRegionModel_U< Real >::x_.

template<typename Real>
virtual const Ptr<Objective<Real> > ROL::TrustRegionModel_U< Real >::getObjective ( void  ) const
inlinevirtual

Definition at line 210 of file ROL_TrustRegionModel_U.hpp.

References ROL::TrustRegionModel_U< Real >::obj_.

Member Data Documentation

template<typename Real>
Ptr<Objective<Real> > ROL::TrustRegionModel_U< Real >::obj_
private
template<typename Real>
Ptr<const Vector<Real> > ROL::TrustRegionModel_U< Real >::x_
private
template<typename Real>
Ptr<const Vector<Real> > ROL::TrustRegionModel_U< Real >::g_
private
template<typename Real>
Ptr<Vector<Real> > ROL::TrustRegionModel_U< Real >::dual_
private
template<typename Real>
Real ROL::TrustRegionModel_U< Real >::tol_
private
template<typename Real>
Ptr<Secant<Real> > ROL::TrustRegionModel_U< Real >::secant_
private
template<typename Real>
bool ROL::TrustRegionModel_U< Real >::useSecantPrecond_
private
template<typename Real>
bool ROL::TrustRegionModel_U< Real >::useSecantHessVec_
private

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