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

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

#include <ROL_TrustRegionModel.hpp>

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

Public Member Functions

virtual ~TrustRegionModel ()
 
 TrustRegionModel (Objective< Real > &obj, BoundConstraint< Real > &bnd, const Vector< Real > &x, const Vector< Real > &g, const Ptr< Secant< Real > > &secant=nullPtr, const bool useSecantPrecond=false, const bool useSecantHessVec=false)
 
virtual void update (Objective< Real > &obj, BoundConstraint< Real > &bnd, const Vector< Real > &x, const Vector< Real > &g, const Ptr< Secant< Real > > &secant=nullPtr)
 
virtual Real value (const Vector< Real > &s, Real &tol)
 Compute value.
 
virtual void gradient (Vector< Real > &g, const Vector< Real > &s, Real &tol)
 Compute gradient.
 
virtual void hessVec (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &s, Real &tol)
 Apply Hessian approximation to vector.
 
virtual void invHessVec (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &s, Real &tol)
 Apply inverse Hessian approximation to vector.
 
virtual void precond (Vector< Real > &Pv, const Vector< Real > &v, const Vector< Real > &s, Real &tol)
 Apply preconditioner to vector.
 
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
 
virtual const Ptr< BoundConstraint< Real > > getBoundConstraint (void) const
 
virtual void dualTransform (Vector< Real > &tv, const Vector< Real > &v)
 
virtual void primalTransform (Vector< Real > &tv, const Vector< Real > &v)
 
virtual void updatePredictedReduction (Real &pred, const Vector< Real > &s)
 
virtual void updateActualReduction (Real &ared, const Vector< Real > &s)
 
- 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.
 
virtual void update (const Vector< Real > &x, bool flag=true, int iter=-1)
 Update objective function.
 
virtual Real dirDeriv (const Vector< Real > &x, const Vector< Real > &d, Real &tol)
 Compute directional derivative.
 
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.
 
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.
 
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.
 
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.
 
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.
 
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.
 
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.
 
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.
 
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.
 
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.
 
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 Member Functions

void initialize (const Vector< Real > &s)
 

Private Attributes

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

Detailed Description

template<class Real>
class ROL::TrustRegionModel< Real >

Provides the interface to evaluate trust-region model functions.

ROL::TrustRegionModel 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 67 of file ROL_TrustRegionModel.hpp.

Constructor & Destructor Documentation

◆ ~TrustRegionModel()

template<class Real >
virtual ROL::TrustRegionModel< Real >::~TrustRegionModel ( )
inlinevirtual

Definition at line 123 of file ROL_TrustRegionModel.hpp.

◆ TrustRegionModel()

template<class Real >
ROL::TrustRegionModel< Real >::TrustRegionModel ( Objective< Real > & obj,
BoundConstraint< Real > & bnd,
const Vector< Real > & x,
const Vector< Real > & g,
const Ptr< Secant< Real > > & secant = nullPtr,
const bool useSecantPrecond = false,
const bool useSecantHessVec = false )
inline

Definition at line 125 of file ROL_TrustRegionModel.hpp.

Member Function Documentation

◆ initialize()

template<class Real >
void ROL::TrustRegionModel< Real >::initialize ( const Vector< Real > & s)
inlineprivate

◆ applyHessian()

template<class Real >
void ROL::TrustRegionModel< Real >::applyHessian ( Vector< Real > & hv,
const Vector< Real > & v,
Real & tol )
inlineprotected

◆ applyInvHessian()

template<class Real >
void ROL::TrustRegionModel< Real >::applyInvHessian ( Vector< Real > & hv,
const Vector< Real > & v,
Real & tol )
inlineprotected

◆ applyPrecond()

template<class Real >
void ROL::TrustRegionModel< Real >::applyPrecond ( Vector< Real > & Pv,
const Vector< Real > & v,
Real & tol )
inlineprotected

◆ update()

template<class Real >
virtual void ROL::TrustRegionModel< Real >::update ( Objective< Real > & obj,
BoundConstraint< Real > & bnd,
const Vector< Real > & x,
const Vector< Real > & g,
const Ptr< Secant< Real > > & secant = nullPtr )
inlinevirtual

◆ value()

template<class Real >
virtual Real ROL::TrustRegionModel< 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 >.

Reimplemented in ROL::ColemanLiModel< Real >, and ROL::KelleySachsModel< Real >.

Definition at line 151 of file ROL_TrustRegionModel.hpp.

References ROL::TrustRegionModel< Real >::applyHessian(), ROL::Vector< Real >::dual(), ROL::TrustRegionModel< Real >::dual_, ROL::TrustRegionModel< Real >::g_, and ROL::TrustRegionModel< Real >::initialize().

Referenced by ROL::TrustRegion< Real >::update().

◆ gradient()

template<class Real >
virtual void ROL::TrustRegionModel< 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 >.

Reimplemented in ROL::ColemanLiModel< Real >, and ROL::KelleySachsModel< Real >.

Definition at line 159 of file ROL_TrustRegionModel.hpp.

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

◆ hessVec()

template<class Real >
virtual void ROL::TrustRegionModel< 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 >.

Reimplemented in ROL::ColemanLiModel< Real >, and ROL::KelleySachsModel< Real >.

Definition at line 164 of file ROL_TrustRegionModel.hpp.

References ROL::TrustRegionModel< Real >::applyHessian().

Referenced by ROL::CauchyPoint< Real >::cauchypoint_unc(), ROL::DogLeg< Real >::run(), ROL::DoubleDogLeg< Real >::run(), and ROL::TruncatedCG< Real >::run().

◆ invHessVec()

template<class Real >
virtual void ROL::TrustRegionModel< Real >::invHessVec ( Vector< Real > & hv,
const Vector< Real > & v,
const Vector< Real > & x,
Real & tol )
inlinevirtual

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 >.

Reimplemented in ROL::KelleySachsModel< Real >.

Definition at line 168 of file ROL_TrustRegionModel.hpp.

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

Referenced by ROL::DogLeg< Real >::run(), and ROL::DoubleDogLeg< Real >::run().

◆ precond()

template<class Real >
virtual void ROL::TrustRegionModel< Real >::precond ( Vector< Real > & Pv,
const Vector< Real > & v,
const Vector< Real > & x,
Real & tol )
inlinevirtual

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 >.

Reimplemented in ROL::KelleySachsModel< Real >.

Definition at line 172 of file ROL_TrustRegionModel.hpp.

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

Referenced by ROL::TruncatedCG< Real >::run().

◆ getGradient()

template<class Real >
virtual const Ptr< const Vector< Real > > ROL::TrustRegionModel< Real >::getGradient ( void ) const
inlinevirtual

◆ getIterate()

template<class Real >
virtual const Ptr< const Vector< Real > > ROL::TrustRegionModel< Real >::getIterate ( void ) const
inlinevirtual

◆ getObjective()

template<class Real >
virtual const Ptr< Objective< Real > > ROL::TrustRegionModel< Real >::getObjective ( void ) const
inlinevirtual

Definition at line 190 of file ROL_TrustRegionModel.hpp.

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

◆ getBoundConstraint()

template<class Real >
virtual const Ptr< BoundConstraint< Real > > ROL::TrustRegionModel< Real >::getBoundConstraint ( void ) const
inlinevirtual

◆ dualTransform()

template<class Real >
virtual void ROL::TrustRegionModel< Real >::dualTransform ( Vector< Real > & tv,
const Vector< Real > & v )
inlinevirtual

◆ primalTransform()

template<class Real >
virtual void ROL::TrustRegionModel< Real >::primalTransform ( Vector< Real > & tv,
const Vector< Real > & v )
inlinevirtual

◆ updatePredictedReduction()

template<class Real >
virtual void ROL::TrustRegionModel< Real >::updatePredictedReduction ( Real & pred,
const Vector< Real > & s )
inlinevirtual

Reimplemented in ROL::ColemanLiModel< Real >.

Definition at line 212 of file ROL_TrustRegionModel.hpp.

Referenced by ROL::TrustRegion< Real >::update().

◆ updateActualReduction()

template<class Real >
virtual void ROL::TrustRegionModel< Real >::updateActualReduction ( Real & ared,
const Vector< Real > & s )
inlinevirtual

Reimplemented in ROL::ColemanLiModel< Real >.

Definition at line 214 of file ROL_TrustRegionModel.hpp.

Referenced by ROL::TrustRegion< Real >::update().

Member Data Documentation

◆ obj_

template<class Real >
Ptr<Objective<Real> > ROL::TrustRegionModel< Real >::obj_
private

◆ bnd_

template<class Real >
Ptr<BoundConstraint<Real> > ROL::TrustRegionModel< Real >::bnd_
private

◆ x_

template<class Real >
Ptr<const Vector<Real> > ROL::TrustRegionModel< Real >::x_
private

◆ g_

template<class Real >
Ptr<const Vector<Real> > ROL::TrustRegionModel< Real >::g_
private

◆ dual_

template<class Real >
Ptr<Vector<Real> > ROL::TrustRegionModel< Real >::dual_
private

◆ secant_

template<class Real >
Ptr<Secant<Real> > ROL::TrustRegionModel< Real >::secant_
private

◆ useSecantPrecond_

template<class Real >
const bool ROL::TrustRegionModel< Real >::useSecantPrecond_
private

◆ useSecantHessVec_

template<class Real >
const bool ROL::TrustRegionModel< Real >::useSecantHessVec_
private

◆ init_

template<class Real >
bool ROL::TrustRegionModel< Real >::init_
private

Definition at line 78 of file ROL_TrustRegionModel.hpp.

Referenced by ROL::TrustRegionModel< Real >::initialize().


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