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

Provides an interface to run trust-region methods for unconstrained optimization algorithms. More...

#include <ROL_TypeU_TrustRegionAlgorithm.hpp>

+ Inheritance diagram for ROL::TypeU::TrustRegionAlgorithm< Real >:

Public Member Functions

 TrustRegionAlgorithm (ParameterList &parlist, const Ptr< Secant< Real >> &secant=nullPtr)
 
void run (Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, std::ostream &outStream=std::cout) override
 Run algorithm on unconstrained problems (Type-U). This general interface supports the use of dual optimization vector spaces, where the user does not define the dual() method. More...
 
void writeHeader (std::ostream &os) const override
 Print iterate header. More...
 
void writeName (std::ostream &os) const override
 Print step name. More...
 
void writeOutput (std::ostream &os, const bool print_header=false) const override
 Print iterate status. More...
 
- Public Member Functions inherited from ROL::TypeU::Algorithm< Real >
virtual ~Algorithm ()
 
 Algorithm ()
 Constructor, given a step and a status test. More...
 
void setStatusTest (const Ptr< StatusTest< Real >> &status, bool combineStatus=false)
 
virtual void run (Problem< Real > &problem, std::ostream &outStream=std::cout)
 Run algorithm on unconstrained problems (Type-U). This is the primary Type-U interface. More...
 
virtual void run (Vector< Real > &x, Objective< Real > &obj, std::ostream &outStream=std::cout)
 Run algorithm on unconstrained problems (Type-U). This is the primary Type-U interface. More...
 
virtual void run (Vector< Real > &x, Objective< Real > &obj, Constraint< Real > &linear_con, Vector< Real > &linear_mul, std::ostream &outStream=std::cout)
 Run algorithm on unconstrained problems with explicit linear equality constraints (Type-U). This is the primary Type-U with explicit linear equality constraints interface. More...
 
virtual void run (Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, Constraint< Real > &linear_con, Vector< Real > &linear_mul, const Vector< Real > &linear_c, std::ostream &outStream=std::cout)
 Run algorithm on unconstrained problems with explicit linear equality constraints (Type-U). This general interface supports the use of dual optimization vector spaces, where the user does not define the dual() method. More...
 
virtual void writeExitStatus (std::ostream &os) const
 
Ptr< const AlgorithmState< Real > > getState () const
 
void reset ()
 

Private Member Functions

void initialize (const Vector< Real > &x, const Vector< Real > &g, Vector< Real > &Bg, Objective< Real > &obj, std::ostream &outStream=std::cout)
 
Real computeValue (const Vector< Real > &x, Objective< Real > &obj, Real pRed)
 
void computeGradient (const Vector< Real > &x, Objective< Real > &obj, bool accept)
 Compute gradient to iteratively satisfy inexactness condition. More...
 

Private Attributes

Ptr< TrustRegion_U< Real > > solver_
 Container for trust-region solver object. More...
 
Ptr< TrustRegionModel_U< Real > > model_
 Container for trust-region model. More...
 
ETrustRegionU etr_
 Trust-region subproblem solver type. More...
 
Real delMax_
 Maximum trust-region radius. More...
 
Real eta0_
 Step acceptance threshold. More...
 
Real eta1_
 Radius decrease threshold. More...
 
Real eta2_
 Radius increase threshold. More...
 
Real gamma0_
 Radius decrease rate (negative rho). More...
 
Real gamma1_
 Radius decrease rate (positive rho). More...
 
Real gamma2_
 Radius increase rate. More...
 
Real TRsafe_
 Safeguard size for numerically evaluating ratio. More...
 
Real eps_
 Safeguard for numerically evaluating ratio. More...
 
TRUtils::ETRFlag TRflag_
 Trust-region exit flag. More...
 
int SPflag_
 Subproblem solver termination flag. More...
 
int SPiter_
 Subproblem solver iteration count. More...
 
bool useNM_
 
int NMstorage_
 
ESecant esec_
 Secant type. More...
 
bool useSecantPrecond_
 
bool useSecantHessVec_
 
std::vector< bool > useInexact_
 Flags for inexact (0) objective function, (1) gradient, (2) Hessian. More...
 
Real scale0_
 Scale for inexact gradient computation. More...
 
Real scale1_
 Scale for inexact gradient computation. More...
 
Real scale_
 
Real omega_
 
Real force_
 
Real forceFactor_
 
int updateIter_
 
Real gtol_
 
int verbosity_
 Print additional information to screen if > 0. More...
 
bool printHeader_
 Print header at every iteration. More...
 

Additional Inherited Members

- Protected Member Functions inherited from ROL::TypeU::Algorithm< Real >
void initialize (const Vector< Real > &x, const Vector< Real > &g)
 
- Protected Attributes inherited from ROL::TypeU::Algorithm< Real >
const Ptr< CombinedStatusTest
< Real > > 
status_
 
const Ptr< AlgorithmState< Real > > state_
 

Detailed Description

template<typename Real>
class ROL::TypeU::TrustRegionAlgorithm< Real >

Provides an interface to run trust-region methods for unconstrained optimization algorithms.

Definition at line 62 of file ROL_TypeU_TrustRegionAlgorithm.hpp.

Constructor & Destructor Documentation

template<typename Real >
ROL::TypeU::TrustRegionAlgorithm< Real >::TrustRegionAlgorithm ( ParameterList &  parlist,
const Ptr< Secant< Real >> &  secant = nullPtr 
)

Definition at line 53 of file ROL_TypeU_TrustRegionAlgorithm_Def.hpp.

References ROL::TypeU::TrustRegionAlgorithm< Real >::delMax_, ROL::TypeU::TrustRegionAlgorithm< Real >::eps_, ROL::TypeU::TrustRegionAlgorithm< Real >::esec_, ROL::TypeU::TrustRegionAlgorithm< Real >::eta0_, ROL::TypeU::TrustRegionAlgorithm< Real >::eta1_, ROL::TypeU::TrustRegionAlgorithm< Real >::eta2_, ROL::TypeU::TrustRegionAlgorithm< Real >::etr_, ROL::TypeU::TrustRegionAlgorithm< Real >::force_, ROL::TypeU::TrustRegionAlgorithm< Real >::forceFactor_, ROL::TypeU::TrustRegionAlgorithm< Real >::gamma0_, ROL::TypeU::TrustRegionAlgorithm< Real >::gamma1_, ROL::TypeU::TrustRegionAlgorithm< Real >::gamma2_, ROL::TypeU::TrustRegionAlgorithm< Real >::model_, ROL::TypeU::TrustRegionAlgorithm< Real >::NMstorage_, ROL::TypeU::TrustRegionAlgorithm< Real >::omega_, ROL::TypeU::TrustRegionAlgorithm< Real >::printHeader_, ROL::TypeU::TrustRegionAlgorithm< Real >::scale0_, ROL::TypeU::TrustRegionAlgorithm< Real >::scale1_, ROL::TypeU::TrustRegionAlgorithm< Real >::scale_, ROL::TypeU::TrustRegionAlgorithm< Real >::solver_, ROL::TypeU::Algorithm< Real >::state_, ROL::TypeU::Algorithm< Real >::status_, ROL::StringToESecant(), ROL::StringToETrustRegionU(), ROL::TypeU::TrustRegionAlgorithm< Real >::TRsafe_, ROL::TypeU::TrustRegionAlgorithm< Real >::updateIter_, ROL::TypeU::TrustRegionAlgorithm< Real >::useInexact_, ROL::TypeU::TrustRegionAlgorithm< Real >::useNM_, ROL::TypeU::TrustRegionAlgorithm< Real >::useSecantHessVec_, ROL::TypeU::TrustRegionAlgorithm< Real >::useSecantPrecond_, and ROL::TypeU::TrustRegionAlgorithm< Real >::verbosity_.

Member Function Documentation

template<typename Real >
void ROL::TypeU::TrustRegionAlgorithm< Real >::run ( Vector< Real > &  x,
const Vector< Real > &  g,
Objective< Real > &  obj,
std::ostream &  outStream = std::cout 
)
overridevirtual

Run algorithm on unconstrained problems (Type-U). This general interface supports the use of dual optimization vector spaces, where the user does not define the dual() method.

Implements ROL::TypeU::Algorithm< Real >.

Definition at line 189 of file ROL_TypeU_TrustRegionAlgorithm_Def.hpp.

References ROL::Accept, ROL::Vector< Real >::clone(), ROL::TRUtils::NPOSPREDNEG, ROL::Vector< Real >::plus(), ROL::TRUtils::POSPREDNEG, ROL::Revert, ROL::Vector< Real >::set(), ROL::TRUtils::SUCCESS, ROL::TRUtils::TRNAN, ROL::Objective< Real >::update(), ROL::TypeU::Algorithm< Real >::writeExitStatus(), and zero.

Referenced by main().

template<typename Real >
void ROL::TypeU::TrustRegionAlgorithm< Real >::writeHeader ( std::ostream &  os) const
overridevirtual
template<typename Real >
void ROL::TypeU::TrustRegionAlgorithm< Real >::writeName ( std::ostream &  os) const
overridevirtual

Print step name.

Reimplemented from ROL::TypeU::Algorithm< Real >.

Definition at line 330 of file ROL_TypeU_TrustRegionAlgorithm_Def.hpp.

References ROL::ESecantToString(), and ROL::ETrustRegionUToString().

template<typename Real >
void ROL::TypeU::TrustRegionAlgorithm< Real >::writeOutput ( std::ostream &  os,
const bool  write_header = false 
) const
overridevirtual

Print iterate status.

Reimplemented from ROL::TypeU::Algorithm< Real >.

Definition at line 351 of file ROL_TypeU_TrustRegionAlgorithm_Def.hpp.

References ROL::TRUSTREGION_U_SPG, and ROL::TRUSTREGION_U_TRUNCATEDCG.

template<typename Real >
void ROL::TypeU::TrustRegionAlgorithm< Real >::initialize ( const Vector< Real > &  x,
const Vector< Real > &  g,
Vector< Real > &  Bg,
Objective< Real > &  obj,
std::ostream &  outStream = std::cout 
)
private
template<typename Real >
Real ROL::TypeU::TrustRegionAlgorithm< Real >::computeValue ( const Vector< Real > &  x,
Objective< Real > &  obj,
Real  pRed 
)
private
template<typename Real >
void ROL::TypeU::TrustRegionAlgorithm< Real >::computeGradient ( const Vector< Real > &  x,
Objective< Real > &  obj,
bool  accept 
)
private

Compute gradient to iteratively satisfy inexactness condition.

This function attempts to ensure that the inexact gradient condition,

\[ \|g_k-\nabla J(x_k)\|_{\mathcal{X}} \le \kappa_1\min\{\,\|g_k\|_{\mathcal{X}},\,\Delta_k\,\}, \]

is satisfied. This function works under the assumption that the gradient function returns a gradient approximation which satisfies the error tolerance prescribed by the tol input parameter.

Parameters
[in]xis the current optimization variable.
[in]objis the objective function.

Definition at line 165 of file ROL_TypeU_TrustRegionAlgorithm_Def.hpp.

References ROL::Objective< Real >::gradient().

Member Data Documentation

template<typename Real>
Ptr<TrustRegion_U<Real> > ROL::TypeU::TrustRegionAlgorithm< Real >::solver_
private

Container for trust-region solver object.

Definition at line 65 of file ROL_TypeU_TrustRegionAlgorithm.hpp.

Referenced by ROL::TypeU::TrustRegionAlgorithm< Real >::TrustRegionAlgorithm().

template<typename Real>
Ptr<TrustRegionModel_U<Real> > ROL::TypeU::TrustRegionAlgorithm< Real >::model_
private

Container for trust-region model.

Definition at line 66 of file ROL_TypeU_TrustRegionAlgorithm.hpp.

Referenced by ROL::TypeU::TrustRegionAlgorithm< Real >::TrustRegionAlgorithm().

template<typename Real>
ETrustRegionU ROL::TypeU::TrustRegionAlgorithm< Real >::etr_
private

Trust-region subproblem solver type.

Definition at line 67 of file ROL_TypeU_TrustRegionAlgorithm.hpp.

Referenced by ROL::TypeU::TrustRegionAlgorithm< Real >::TrustRegionAlgorithm().

template<typename Real>
Real ROL::TypeU::TrustRegionAlgorithm< Real >::delMax_
private

Maximum trust-region radius.

Definition at line 68 of file ROL_TypeU_TrustRegionAlgorithm.hpp.

Referenced by ROL::TypeU::TrustRegionAlgorithm< Real >::TrustRegionAlgorithm().

template<typename Real>
Real ROL::TypeU::TrustRegionAlgorithm< Real >::eta0_
private

Step acceptance threshold.

Definition at line 69 of file ROL_TypeU_TrustRegionAlgorithm.hpp.

Referenced by ROL::TypeU::TrustRegionAlgorithm< Real >::TrustRegionAlgorithm().

template<typename Real>
Real ROL::TypeU::TrustRegionAlgorithm< Real >::eta1_
private

Radius decrease threshold.

Definition at line 70 of file ROL_TypeU_TrustRegionAlgorithm.hpp.

Referenced by ROL::TypeU::TrustRegionAlgorithm< Real >::TrustRegionAlgorithm().

template<typename Real>
Real ROL::TypeU::TrustRegionAlgorithm< Real >::eta2_
private

Radius increase threshold.

Definition at line 71 of file ROL_TypeU_TrustRegionAlgorithm.hpp.

Referenced by ROL::TypeU::TrustRegionAlgorithm< Real >::TrustRegionAlgorithm().

template<typename Real>
Real ROL::TypeU::TrustRegionAlgorithm< Real >::gamma0_
private

Radius decrease rate (negative rho).

Definition at line 72 of file ROL_TypeU_TrustRegionAlgorithm.hpp.

Referenced by ROL::TypeU::TrustRegionAlgorithm< Real >::TrustRegionAlgorithm().

template<typename Real>
Real ROL::TypeU::TrustRegionAlgorithm< Real >::gamma1_
private

Radius decrease rate (positive rho).

Definition at line 73 of file ROL_TypeU_TrustRegionAlgorithm.hpp.

Referenced by ROL::TypeU::TrustRegionAlgorithm< Real >::TrustRegionAlgorithm().

template<typename Real>
Real ROL::TypeU::TrustRegionAlgorithm< Real >::gamma2_
private

Radius increase rate.

Definition at line 74 of file ROL_TypeU_TrustRegionAlgorithm.hpp.

Referenced by ROL::TypeU::TrustRegionAlgorithm< Real >::TrustRegionAlgorithm().

template<typename Real>
Real ROL::TypeU::TrustRegionAlgorithm< Real >::TRsafe_
private

Safeguard size for numerically evaluating ratio.

Definition at line 75 of file ROL_TypeU_TrustRegionAlgorithm.hpp.

Referenced by ROL::TypeU::TrustRegionAlgorithm< Real >::TrustRegionAlgorithm().

template<typename Real>
Real ROL::TypeU::TrustRegionAlgorithm< Real >::eps_
private

Safeguard for numerically evaluating ratio.

Definition at line 76 of file ROL_TypeU_TrustRegionAlgorithm.hpp.

Referenced by ROL::TypeU::TrustRegionAlgorithm< Real >::TrustRegionAlgorithm().

template<typename Real>
TRUtils::ETRFlag ROL::TypeU::TrustRegionAlgorithm< Real >::TRflag_
private

Trust-region exit flag.

Definition at line 77 of file ROL_TypeU_TrustRegionAlgorithm.hpp.

template<typename Real>
int ROL::TypeU::TrustRegionAlgorithm< Real >::SPflag_
private

Subproblem solver termination flag.

Definition at line 78 of file ROL_TypeU_TrustRegionAlgorithm.hpp.

template<typename Real>
int ROL::TypeU::TrustRegionAlgorithm< Real >::SPiter_
private

Subproblem solver iteration count.

Definition at line 79 of file ROL_TypeU_TrustRegionAlgorithm.hpp.

template<typename Real>
bool ROL::TypeU::TrustRegionAlgorithm< Real >::useNM_
private
template<typename Real>
int ROL::TypeU::TrustRegionAlgorithm< Real >::NMstorage_
private
template<typename Real>
ESecant ROL::TypeU::TrustRegionAlgorithm< Real >::esec_
private
template<typename Real>
bool ROL::TypeU::TrustRegionAlgorithm< Real >::useSecantPrecond_
private
template<typename Real>
bool ROL::TypeU::TrustRegionAlgorithm< Real >::useSecantHessVec_
private
template<typename Real>
std::vector<bool> ROL::TypeU::TrustRegionAlgorithm< Real >::useInexact_
private

Flags for inexact (0) objective function, (1) gradient, (2) Hessian.

Definition at line 91 of file ROL_TypeU_TrustRegionAlgorithm.hpp.

Referenced by ROL::TypeU::TrustRegionAlgorithm< Real >::TrustRegionAlgorithm().

template<typename Real>
Real ROL::TypeU::TrustRegionAlgorithm< Real >::scale0_
private

Scale for inexact gradient computation.

Definition at line 92 of file ROL_TypeU_TrustRegionAlgorithm.hpp.

Referenced by ROL::TypeU::TrustRegionAlgorithm< Real >::TrustRegionAlgorithm().

template<typename Real>
Real ROL::TypeU::TrustRegionAlgorithm< Real >::scale1_
private

Scale for inexact gradient computation.

Definition at line 93 of file ROL_TypeU_TrustRegionAlgorithm.hpp.

Referenced by ROL::TypeU::TrustRegionAlgorithm< Real >::TrustRegionAlgorithm().

template<typename Real>
Real ROL::TypeU::TrustRegionAlgorithm< Real >::scale_
private
template<typename Real>
Real ROL::TypeU::TrustRegionAlgorithm< Real >::omega_
private
template<typename Real>
Real ROL::TypeU::TrustRegionAlgorithm< Real >::force_
private
template<typename Real>
Real ROL::TypeU::TrustRegionAlgorithm< Real >::forceFactor_
private
template<typename Real>
int ROL::TypeU::TrustRegionAlgorithm< Real >::updateIter_
private
template<typename Real>
Real ROL::TypeU::TrustRegionAlgorithm< Real >::gtol_
private

Definition at line 96 of file ROL_TypeU_TrustRegionAlgorithm.hpp.

template<typename Real>
int ROL::TypeU::TrustRegionAlgorithm< Real >::verbosity_
private

Print additional information to screen if > 0.

Definition at line 99 of file ROL_TypeU_TrustRegionAlgorithm.hpp.

Referenced by ROL::TypeU::TrustRegionAlgorithm< Real >::TrustRegionAlgorithm().

template<typename Real>
bool ROL::TypeU::TrustRegionAlgorithm< Real >::printHeader_
private

Print header at every iteration.

Definition at line 100 of file ROL_TypeU_TrustRegionAlgorithm.hpp.

Referenced by ROL::TypeU::TrustRegionAlgorithm< Real >::TrustRegionAlgorithm().


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