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

Implements the computation of optimization steps with composite-step trust-region methods. More...

#include <ROL_CompositeStep.hpp>

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

Public Member Functions

virtual ~CompositeStep ()
 
 CompositeStep (ROL::ParameterList &parlist)
 
void initialize (Vector< Real > &x, const Vector< Real > &g, Vector< Real > &l, const Vector< Real > &c, Objective< Real > &obj, Constraint< Real > &con, AlgorithmState< Real > &algo_state)
 Initialize step. More...
 
void compute (Vector< Real > &s, const Vector< Real > &x, const Vector< Real > &l, Objective< Real > &obj, Constraint< Real > &con, AlgorithmState< Real > &algo_state)
 Compute step. More...
 
void update (Vector< Real > &x, Vector< Real > &l, const Vector< Real > &s, Objective< Real > &obj, Constraint< Real > &con, AlgorithmState< Real > &algo_state)
 Update step, if successful. More...
 
void compute (Vector< Real > &s, const Vector< Real > &x, Objective< Real > &obj, BoundConstraint< Real > &con, AlgorithmState< Real > &algo_state)
 Compute step for bound constraints; here only to satisfy the interface requirements, does nothing, needs refactoring. More...
 
void update (Vector< Real > &x, const Vector< Real > &s, Objective< Real > &obj, BoundConstraint< Real > &con, AlgorithmState< Real > &algo_state)
 Update step, for bound constraints; here only to satisfy the interface requirements, does nothing, needs refactoring. More...
 
std::string printHeader (void) const
 Print iterate header. More...
 
std::string printName (void) const
 Print step name. More...
 
std::string print (AlgorithmState< Real > &algo_state, bool pHeader=false) const
 Print iterate status. More...
 
void computeLagrangeMultiplier (Vector< Real > &l, const Vector< Real > &x, const Vector< Real > &gf, Constraint< Real > &con)
 Compute Lagrange multipliers by solving the least-squares problem minimizing the gradient of the Lagrangian, via the augmented system formulation. More...
 
void computeQuasinormalStep (Vector< Real > &n, const Vector< Real > &c, const Vector< Real > &x, Real delta, Constraint< Real > &con)
 Compute quasi-normal step by minimizing the norm of the linearized constraint. More...
 
void solveTangentialSubproblem (Vector< Real > &t, Vector< Real > &tCP, Vector< Real > &Wg, const Vector< Real > &x, const Vector< Real > &g, const Vector< Real > &n, const Vector< Real > &l, Real delta, Objective< Real > &obj, Constraint< Real > &con)
 Solve tangential subproblem. More...
 
void accept (Vector< Real > &s, Vector< Real > &n, Vector< Real > &t, Real f_new, Vector< Real > &c_new, Vector< Real > &gf_new, Vector< Real > &l_new, Vector< Real > &g_new, const Vector< Real > &x, const Vector< Real > &l, Real f, const Vector< Real > &gf, const Vector< Real > &c, const Vector< Real > &g, Vector< Real > &tCP, Vector< Real > &Wg, Objective< Real > &obj, Constraint< Real > &con, AlgorithmState< Real > &algo_state)
 Check acceptance of subproblem solutions, adjust merit function penalty parameter, ensure global convergence. More...
 
- Public Member Functions inherited from ROL::Step< Real >
virtual ~Step ()
 
 Step (void)
 
virtual void initialize (Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, BoundConstraint< Real > &con, AlgorithmState< Real > &algo_state)
 Initialize step with bound constraint. More...
 
virtual void initialize (Vector< Real > &x, const Vector< Real > &s, const Vector< Real > &g, Objective< Real > &obj, BoundConstraint< Real > &con, AlgorithmState< Real > &algo_state)
 Initialize step with bound constraint. More...
 
virtual void initialize (Vector< Real > &x, const Vector< Real > &g, Vector< Real > &l, const Vector< Real > &c, Objective< Real > &obj, Constraint< Real > &con, BoundConstraint< Real > &bnd, AlgorithmState< Real > &algo_state)
 Initialize step with equality constraint. More...
 
virtual void compute (Vector< Real > &s, const Vector< Real > &x, const Vector< Real > &l, Objective< Real > &obj, Constraint< Real > &con, BoundConstraint< Real > &bnd, AlgorithmState< Real > &algo_state)
 Compute step (equality constraints). More...
 
virtual void update (Vector< Real > &x, Vector< Real > &l, const Vector< Real > &s, Objective< Real > &obj, Constraint< Real > &con, BoundConstraint< Real > &bnd, AlgorithmState< Real > &algo_state)
 Update step, if successful (equality constraints). More...
 
const ROL::Ptr< const
StepState< Real > > 
getStepState (void) const
 Get state for step object. More...
 
void reset (const Real searchSize=1.0)
 Get state for step object. More...
 

Private Member Functions

template<typename T >
int sgn (T val)
 
void printInfoLS (const std::vector< Real > &res) const
 
Real setTolOSS (const Real intol) const
 

Private Attributes

ROL::Ptr< Vector< Real > > xvec_
 
ROL::Ptr< Vector< Real > > gvec_
 
ROL::Ptr< Vector< Real > > cvec_
 
ROL::Ptr< Vector< Real > > lvec_
 
int flagCG_
 
int flagAC_
 
int iterCG_
 
int maxiterCG_
 
int maxiterOSS_
 
Real tolCG_
 
Real tolOSS_
 
bool tolOSSfixed_
 
Real lmhtol_
 
Real qntol_
 
Real pgtol_
 
Real projtol_
 
Real tangtol_
 
Real tntmax_
 
Real zeta_
 
Real Delta_
 
Real penalty_
 
Real eta_
 
bool useConHess_
 
Real ared_
 
Real pred_
 
Real snorm_
 
Real nnorm_
 
Real tnorm_
 
bool infoQN_
 
bool infoLM_
 
bool infoTS_
 
bool infoAC_
 
bool infoLS_
 
bool infoALL_
 
int totalIterCG_
 
int totalProj_
 
int totalNegCurv_
 
int totalRef_
 
int totalCallLS_
 
int totalIterLS_
 

Additional Inherited Members

- Protected Member Functions inherited from ROL::Step< Real >
ROL::Ptr< StepState< Real > > getState (void)
 

Detailed Description

template<class Real>
class ROL::CompositeStep< Real >

Implements the computation of optimization steps with composite-step trust-region methods.

Definition at line 63 of file ROL_CompositeStep.hpp.

Constructor & Destructor Documentation

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

Definition at line 148 of file ROL_CompositeStep.hpp.

template<class Real >
ROL::CompositeStep< Real >::CompositeStep ( ROL::ParameterList &  parlist)
inline

Member Function Documentation

template<class Real >
template<typename T >
int ROL::CompositeStep< Real >::sgn ( val)
inlineprivate
template<class Real >
void ROL::CompositeStep< Real >::printInfoLS ( const std::vector< Real > &  res) const
inlineprivate
template<class Real >
Real ROL::CompositeStep< Real >::setTolOSS ( const Real  intol) const
inlineprivate
template<class Real >
void ROL::CompositeStep< Real >::initialize ( Vector< Real > &  x,
const Vector< Real > &  g,
Vector< Real > &  l,
const Vector< Real > &  c,
Objective< Real > &  obj,
Constraint< Real > &  con,
AlgorithmState< Real > &  algo_state 
)
inlinevirtual
template<class Real >
void ROL::CompositeStep< Real >::compute ( Vector< Real > &  s,
const Vector< Real > &  x,
const Vector< Real > &  l,
Objective< Real > &  obj,
Constraint< Real > &  con,
AlgorithmState< Real > &  algo_state 
)
inlinevirtual
template<class Real >
void ROL::CompositeStep< Real >::update ( Vector< Real > &  x,
Vector< Real > &  l,
const Vector< Real > &  s,
Objective< Real > &  obj,
Constraint< Real > &  con,
AlgorithmState< Real > &  algo_state 
)
inlinevirtual
template<class Real >
void ROL::CompositeStep< Real >::compute ( Vector< Real > &  s,
const Vector< Real > &  x,
Objective< Real > &  obj,
BoundConstraint< Real > &  con,
AlgorithmState< Real > &  algo_state 
)
inlinevirtual

Compute step for bound constraints; here only to satisfy the interface requirements, does nothing, needs refactoring.

Reimplemented from ROL::Step< Real >.

Definition at line 372 of file ROL_CompositeStep.hpp.

template<class Real >
void ROL::CompositeStep< Real >::update ( Vector< Real > &  x,
const Vector< Real > &  s,
Objective< Real > &  obj,
BoundConstraint< Real > &  con,
AlgorithmState< Real > &  algo_state 
)
inlinevirtual

Update step, for bound constraints; here only to satisfy the interface requirements, does nothing, needs refactoring.

Reimplemented from ROL::Step< Real >.

Definition at line 379 of file ROL_CompositeStep.hpp.

template<class Real >
std::string ROL::CompositeStep< Real >::printHeader ( void  ) const
inlinevirtual

Print iterate header.

Reimplemented from ROL::Step< Real >.

Definition at line 385 of file ROL_CompositeStep.hpp.

Referenced by ROL::CompositeStep< Real >::print().

template<class Real >
std::string ROL::CompositeStep< Real >::printName ( void  ) const
inlinevirtual

Print step name.

Reimplemented from ROL::Step< Real >.

Definition at line 406 of file ROL_CompositeStep.hpp.

Referenced by ROL::CompositeStep< Real >::print().

template<class Real >
std::string ROL::CompositeStep< Real >::print ( AlgorithmState< Real > &  algo_state,
bool  pHeader = false 
) const
inlinevirtual
template<class Real >
void ROL::CompositeStep< Real >::computeLagrangeMultiplier ( Vector< Real > &  l,
const Vector< Real > &  x,
const Vector< Real > &  gf,
Constraint< Real > &  con 
)
inline

Compute Lagrange multipliers by solving the least-squares problem minimizing the gradient of the Lagrangian, via the augmented system formulation.

Parameters
[out]lis the updated Lagrange multiplier; a dual constraint-space vector
[in]xis the current iterate; an optimization-space vector
[in]gfis the gradient of the objective function; a dual optimization-space vector
[in]conis the equality constraint object

On return ... fill the blanks.

Definition at line 468 of file ROL_CompositeStep.hpp.

References ROL::Constraint< Real >::applyAdjointJacobian(), ROL::CompositeStep< Real >::cvec_, ROL::CompositeStep< Real >::gvec_, ROL::CompositeStep< Real >::infoLM_, ROL::CompositeStep< Real >::lmhtol_, ROL::CompositeStep< Real >::lvec_, ROL::Vector< Real >::plus(), ROL::CompositeStep< Real >::printInfoLS(), ROL::CompositeStep< Real >::setTolOSS(), ROL::Constraint< Real >::solveAugmentedSystem(), ROL::CompositeStep< Real >::totalCallLS_, ROL::CompositeStep< Real >::totalIterLS_, and ROL::CompositeStep< Real >::xvec_.

Referenced by ROL::CompositeStep< Real >::accept(), ROL::CompositeStep< Real >::initialize(), and ROL::CompositeStep< Real >::update().

template<class Real >
void ROL::CompositeStep< Real >::computeQuasinormalStep ( Vector< Real > &  n,
const Vector< Real > &  c,
const Vector< Real > &  x,
Real  delta,
Constraint< Real > &  con 
)
inline

Compute quasi-normal step by minimizing the norm of the linearized constraint.

  Compute an approximate solution of the problem

\[ \begin{array}{rl} \min_{n} & \|c'(x_k)n + c(x_k)\|^2_{\mathcal{X}} \\ \mbox{subject to} & \|n\|_{\mathcal{X}} \le \delta \end{array} \]

The approximate solution is computed using Powell's dogleg method. The dogleg path is computed using the Cauchy point and a full Newton step. The path's intersection with the trust-region constraint gives the quasi-normal step.

Parameters
[out]nis the quasi-normal step; an optimization-space vector
[in]cis the value of equality constraints; a constraint-space vector
[in]xis the current iterate; an optimization-space vector
[in]deltais the trust-region radius for the quasi-normal step
[in]conis the equality constraint object

Definition at line 535 of file ROL_CompositeStep.hpp.

References ROL::Constraint< Real >::applyAdjointJacobian(), ROL::Constraint< Real >::applyJacobian(), ROL::Vector< Real >::axpy(), ROL::CompositeStep< Real >::cvec_, ROL::Vector< Real >::dual(), ROL::CompositeStep< Real >::gvec_, ROL::CompositeStep< Real >::infoQN_, ROL::CompositeStep< Real >::lvec_, ROL::CompositeStep< Real >::printInfoLS(), ROL::CompositeStep< Real >::qntol_, ROL::Vector< Real >::scale(), ROL::Vector< Real >::set(), ROL::CompositeStep< Real >::setTolOSS(), ROL::Constraint< Real >::solveAugmentedSystem(), ROL::CompositeStep< Real >::totalCallLS_, ROL::CompositeStep< Real >::totalIterLS_, ROL::CompositeStep< Real >::xvec_, and zero.

Referenced by ROL::CompositeStep< Real >::accept(), and ROL::CompositeStep< Real >::compute().

template<class Real >
void ROL::CompositeStep< Real >::solveTangentialSubproblem ( Vector< Real > &  t,
Vector< Real > &  tCP,
Vector< Real > &  Wg,
const Vector< Real > &  x,
const Vector< Real > &  g,
const Vector< Real > &  n,
const Vector< Real > &  l,
Real  delta,
Objective< Real > &  obj,
Constraint< Real > &  con 
)
inline

Solve tangential subproblem.

Parameters
[out]tis the solution of the tangential subproblem; an optimization-space vector
[out]tCPis the Cauchy point for the tangential subproblem; an optimization-space vector
[out]Wgis the dual of the projected gradient of the Lagrangian; an optimization-space vector
[in]xis the current iterate; an optimization-space vector
[in]gis the gradient of the Lagrangian; a dual optimization-space vector
[in]nis the quasi-normal step; an optimization-space vector
[in]lis the Lagrange multiplier; a dual constraint-space vector
[in]deltais the trust-region radius for the tangential subproblem
[in]conis the equality constraint object

Definition at line 645 of file ROL_CompositeStep.hpp.

References ROL::Constraint< Real >::applyAdjointHessian(), ROL::Constraint< Real >::applyJacobian(), ROL::CompositeStep< Real >::cvec_, ROL::Vector< Real >::dot(), ROL::Vector< Real >::dual(), ROL::CompositeStep< Real >::flagCG_, ROL::CompositeStep< Real >::gvec_, ROL::Objective< Real >::hessVec(), ROL::CompositeStep< Real >::infoTS_, ROL::CompositeStep< Real >::iterCG_, ROL::CompositeStep< Real >::lvec_, ROL::CompositeStep< Real >::maxiterCG_, ROL::Vector< Real >::norm(), ROL::CompositeStep< Real >::pgtol_, ROL::Vector< Real >::plus(), ROL::CompositeStep< Real >::printInfoLS(), ROL::CompositeStep< Real >::projtol_, ROL::Vector< Real >::set(), ROL::CompositeStep< Real >::setTolOSS(), ROL::CompositeStep< Real >::sgn(), ROL::Constraint< Real >::solveAugmentedSystem(), ROL::CompositeStep< Real >::tolCG_, ROL::CompositeStep< Real >::totalCallLS_, ROL::CompositeStep< Real >::totalIterLS_, ROL::CompositeStep< Real >::useConHess_, ROL::CompositeStep< Real >::xvec_, zero, and ROL::Vector< Real >::zero().

Referenced by ROL::CompositeStep< Real >::accept(), and ROL::CompositeStep< Real >::compute().

template<class Real >
void ROL::CompositeStep< Real >::accept ( Vector< Real > &  s,
Vector< Real > &  n,
Vector< Real > &  t,
Real  f_new,
Vector< Real > &  c_new,
Vector< Real > &  gf_new,
Vector< Real > &  l_new,
Vector< Real > &  g_new,
const Vector< Real > &  x,
const Vector< Real > &  l,
Real  f,
const Vector< Real > &  gf,
const Vector< Real > &  c,
const Vector< Real > &  g,
Vector< Real > &  tCP,
Vector< Real > &  Wg,
Objective< Real > &  obj,
Constraint< Real > &  con,
AlgorithmState< Real > &  algo_state 
)
inline

Check acceptance of subproblem solutions, adjust merit function penalty parameter, ensure global convergence.

Definition at line 972 of file ROL_CompositeStep.hpp.

References ROL::Constraint< Real >::applyAdjointHessian(), ROL::Constraint< Real >::applyAdjointJacobian(), ROL::Constraint< Real >::applyJacobian(), ROL::CompositeStep< Real >::ared_, ROL::CompositeStep< Real >::computeLagrangeMultiplier(), ROL::CompositeStep< Real >::computeQuasinormalStep(), ROL::CompositeStep< Real >::cvec_, ROL::CompositeStep< Real >::Delta_, ROL::Vector< Real >::dot(), ROL::Vector< Real >::dual(), ROL::CompositeStep< Real >::eta_, ROL::CompositeStep< Real >::flagCG_, ROL::Objective< Real >::gradient(), ROL::CompositeStep< Real >::gvec_, ROL::Objective< Real >::hessVec(), ROL::CompositeStep< Real >::infoAC_, ROL::AlgorithmState< Real >::iter, ROL::CompositeStep< Real >::iterCG_, ROL::CompositeStep< Real >::lmhtol_, ROL::CompositeStep< Real >::lvec_, ROL::CompositeStep< Real >::nnorm_, ROL::Vector< Real >::norm(), ROL::CompositeStep< Real >::penalty_, ROL::CompositeStep< Real >::pgtol_, ROL::Vector< Real >::plus(), ROL::CompositeStep< Real >::pred_, ROL::CompositeStep< Real >::printInfoLS(), ROL::CompositeStep< Real >::projtol_, ROL::CompositeStep< Real >::qntol_, ROL::Vector< Real >::set(), ROL::CompositeStep< Real >::setTolOSS(), ROL::CompositeStep< Real >::snorm_, ROL::Constraint< Real >::solveAugmentedSystem(), ROL::CompositeStep< Real >::solveTangentialSubproblem(), ROL::CompositeStep< Real >::tangtol_, ROL::CompositeStep< Real >::tnorm_, ROL::CompositeStep< Real >::tntmax_, ROL::CompositeStep< Real >::tolOSSfixed_, ROL::CompositeStep< Real >::totalCallLS_, ROL::CompositeStep< Real >::totalIterCG_, ROL::CompositeStep< Real >::totalIterLS_, ROL::CompositeStep< Real >::totalNegCurv_, ROL::CompositeStep< Real >::totalProj_, ROL::CompositeStep< Real >::totalRef_, ROL::Objective< Real >::update(), ROL::Constraint< Real >::update(), ROL::CompositeStep< Real >::useConHess_, ROL::Objective< Real >::value(), ROL::Constraint< Real >::value(), ROL::CompositeStep< Real >::xvec_, zero, and ROL::CompositeStep< Real >::zeta_.

Referenced by ROL::CompositeStep< Real >::compute().

Member Data Documentation

template<class Real >
ROL::Ptr<Vector<Real> > ROL::CompositeStep< Real >::xvec_
private
template<class Real >
ROL::Ptr<Vector<Real> > ROL::CompositeStep< Real >::gvec_
private
template<class Real >
ROL::Ptr<Vector<Real> > ROL::CompositeStep< Real >::cvec_
private
template<class Real >
ROL::Ptr<Vector<Real> > ROL::CompositeStep< Real >::lvec_
private
template<class Real >
int ROL::CompositeStep< Real >::flagCG_
private
template<class Real >
int ROL::CompositeStep< Real >::flagAC_
private
template<class Real >
int ROL::CompositeStep< Real >::iterCG_
private
template<class Real >
int ROL::CompositeStep< Real >::maxiterCG_
private
template<class Real >
int ROL::CompositeStep< Real >::maxiterOSS_
private

Definition at line 79 of file ROL_CompositeStep.hpp.

template<class Real >
Real ROL::CompositeStep< Real >::tolCG_
private
template<class Real >
Real ROL::CompositeStep< Real >::tolOSS_
private
template<class Real >
bool ROL::CompositeStep< Real >::tolOSSfixed_
private
template<class Real >
Real ROL::CompositeStep< Real >::lmhtol_
private
template<class Real >
Real ROL::CompositeStep< Real >::qntol_
private
template<class Real >
Real ROL::CompositeStep< Real >::pgtol_
private
template<class Real >
Real ROL::CompositeStep< Real >::projtol_
private
template<class Real >
Real ROL::CompositeStep< Real >::tangtol_
private
template<class Real >
Real ROL::CompositeStep< Real >::tntmax_
private
template<class Real >
Real ROL::CompositeStep< Real >::zeta_
private
template<class Real >
Real ROL::CompositeStep< Real >::Delta_
private
template<class Real >
Real ROL::CompositeStep< Real >::penalty_
private
template<class Real >
Real ROL::CompositeStep< Real >::eta_
private
template<class Real >
bool ROL::CompositeStep< Real >::useConHess_
private
template<class Real >
Real ROL::CompositeStep< Real >::ared_
private
template<class Real >
Real ROL::CompositeStep< Real >::pred_
private
template<class Real >
Real ROL::CompositeStep< Real >::snorm_
private
template<class Real >
Real ROL::CompositeStep< Real >::nnorm_
private
template<class Real >
Real ROL::CompositeStep< Real >::tnorm_
private
template<class Real >
bool ROL::CompositeStep< Real >::infoQN_
private
template<class Real >
bool ROL::CompositeStep< Real >::infoLM_
private
template<class Real >
bool ROL::CompositeStep< Real >::infoTS_
private
template<class Real >
bool ROL::CompositeStep< Real >::infoAC_
private
template<class Real >
bool ROL::CompositeStep< Real >::infoLS_
private
template<class Real >
bool ROL::CompositeStep< Real >::infoALL_
private

Definition at line 111 of file ROL_CompositeStep.hpp.

Referenced by ROL::CompositeStep< Real >::CompositeStep().

template<class Real >
int ROL::CompositeStep< Real >::totalIterCG_
private
template<class Real >
int ROL::CompositeStep< Real >::totalProj_
private
template<class Real >
int ROL::CompositeStep< Real >::totalNegCurv_
private
template<class Real >
int ROL::CompositeStep< Real >::totalRef_
private
template<class Real >
int ROL::CompositeStep< Real >::totalCallLS_
private
template<class Real >
int ROL::CompositeStep< Real >::totalIterLS_
private

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