ROL
TrustRegionModel_U Class Reference

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

#include <ROL_TrustRegionModel_U.hpp>

+ Inheritance diagram for TrustRegionModel_U:

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)
 
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
 
virtual void gradient (Vector< Real > &g, const Vector< Real > &s, Real &tol) override
 
virtual void hessVec (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &s, Real &tol) override
 
virtual void invHessVec (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &s, Real &tol) override
 
virtual void precond (Vector< Real > &Pv, const Vector< Real > &v, const Vector< Real > &s, Real &tol) override
 
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.
 
virtual void update (const Vector< Real > &x, bool flag=true, int iter=-1)
 Update objective function.
 
virtual Real value (const Vector< Real > &x, Real &tol)=0
 Compute value.
 
virtual void gradient (Vector< Real > &g, const Vector< Real > &x, Real &tol)
 Compute gradient.
 
virtual Real dirDeriv (const Vector< Real > &x, const Vector< Real > &d, Real &tol)
 Compute directional derivative.
 
virtual void hessVec (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
 Apply Hessian approximation to vector.
 
virtual void invHessVec (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
 Apply inverse Hessian approximation to vector.
 
virtual void precond (Vector< Real > &Pv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
 Apply preconditioner to vector.
 
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 Attributes

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

Detailed Description

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.


Constructor & Destructor Documentation

◆ ~TrustRegionModel_U()

virtual TrustRegionModel_U::~TrustRegionModel_U ( )
inlinevirtual

Definition at line 112 of file ROL_TrustRegionModel_U.hpp.

◆ TrustRegionModel_U()

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

Definition at line 114 of file ROL_TrustRegionModel_U.hpp.

References mode.

Member Function Documentation

◆ applyHessian()

void TrustRegionModel_U::applyHessian ( Vector< Real > & hv,
const Vector< Real > & v,
Real & tol )
inlineprotected

◆ applyInvHessian()

void TrustRegionModel_U::applyInvHessian ( Vector< Real > & hv,
const Vector< Real > & v,
Real & tol )
inlineprotected

Definition at line 89 of file ROL_TrustRegionModel_U.hpp.

References obj_.

Referenced by TrustRegionModel_U::invHessVec().

◆ applyPrecond()

void TrustRegionModel_U::applyPrecond ( Vector< Real > & Pv,
const Vector< Real > & v,
Real & tol )
inlineprotected

Definition at line 98 of file ROL_TrustRegionModel_U.hpp.

References obj_.

Referenced by TrustRegionModel_U::precond().

◆ initialize()

void TrustRegionModel_U::initialize ( const Vector< Real > & x,
const Vector< Real > & g )
inline

Definition at line 126 of file ROL_TrustRegionModel_U.hpp.

◆ validate()

void TrustRegionModel_U::validate ( Objective< Real > & obj,
const Vector< Real > & x,
const Vector< Real > & g,
ETrustRegionU etr )
inline

Definition at line 134 of file ROL_TrustRegionModel_U.hpp.

◆ setData()

virtual void TrustRegionModel_U::setData ( Objective< Real > & obj,
const Vector< Real > & x,
const Vector< Real > & g )
inlinevirtual

Definition at line 152 of file ROL_TrustRegionModel_U.hpp.

References obj_.

Referenced by ROL::TRUtils::initialRadius().

◆ update()

void TrustRegionModel_U::update ( const Vector< Real > & x,
const Vector< Real > & s,
const Vector< Real > & gold,
const Vector< Real > & gnew,
const Real snorm,
const int iter )
inline

Definition at line 160 of file ROL_TrustRegionModel_U.hpp.

References iter, and snorm.

◆ value()

virtual Real TrustRegionModel_U::value ( const Vector< Real > & s,
Real & tol )
inlineoverridevirtual

Definition at line 172 of file ROL_TrustRegionModel_U.hpp.

References TrustRegionModel_U::applyHessian().

◆ gradient()

virtual void TrustRegionModel_U::gradient ( Vector< Real > & g,
const Vector< Real > & s,
Real & tol )
inlineoverridevirtual

Definition at line 180 of file ROL_TrustRegionModel_U.hpp.

References TrustRegionModel_U::applyHessian().

◆ hessVec()

◆ invHessVec()

virtual void TrustRegionModel_U::invHessVec ( Vector< Real > & hv,
const Vector< Real > & v,
const Vector< Real > & s,
Real & tol )
inlineoverridevirtual

◆ precond()

virtual void TrustRegionModel_U::precond ( Vector< Real > & Pv,
const Vector< Real > & v,
const Vector< Real > & s,
Real & tol )
inlineoverridevirtual

◆ getGradient()

virtual const Ptr< const Vector< Real > > TrustRegionModel_U::getGradient ( void ) const
inlinevirtual

◆ getIterate()

virtual const Ptr< const Vector< Real > > TrustRegionModel_U::getIterate ( void ) const
inlinevirtual

Definition at line 207 of file ROL_TrustRegionModel_U.hpp.

◆ getObjective()

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

Definition at line 211 of file ROL_TrustRegionModel_U.hpp.

References obj_.

Member Data Documentation

◆ obj_

Ptr<Objective<Real> > TrustRegionModel_U::obj_
private

Definition at line 68 of file ROL_TrustRegionModel_U.hpp.

◆ x_

Ptr<const Vector<Real> > TrustRegionModel_U::x_
private

Definition at line 69 of file ROL_TrustRegionModel_U.hpp.

◆ g_

Ptr<const Vector<Real> > TrustRegionModel_U::g_
private

Definition at line 69 of file ROL_TrustRegionModel_U.hpp.

◆ dual_

Ptr<Vector<Real> > TrustRegionModel_U::dual_
private

Definition at line 70 of file ROL_TrustRegionModel_U.hpp.

◆ secant_

Ptr<Secant<Real> > TrustRegionModel_U::secant_
private

Definition at line 72 of file ROL_TrustRegionModel_U.hpp.

◆ useSecantPrecond_

bool TrustRegionModel_U::useSecantPrecond_
private

Definition at line 73 of file ROL_TrustRegionModel_U.hpp.

◆ useSecantHessVec_

bool TrustRegionModel_U::useSecantHessVec_
private

Definition at line 74 of file ROL_TrustRegionModel_U.hpp.


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