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

Defines the reduced time-dependent objective function interface for simulation-based optimization. More...

#include <ROL_ReducedDynamicObjective.hpp>

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

Public Member Functions

 ReducedDynamicObjective (const Ptr< DynamicObjective< Real >> &obj, const Ptr< DynamicConstraint< Real >> &con, const Ptr< Vector< Real >> &u0, const Ptr< Vector< Real >> &zvec, const Ptr< Vector< Real >> &cvec, const std::vector< TimeStamp< Real >> &timeStamp, ROL::ParameterList &pl, const Ptr< std::ostream > &stream=nullPtr)
 
Ptr< Vector< Real > > makeDynamicVector (const Vector< Real > &x) const
 
void update (const Vector< Real > &x, bool flag=true, int iter=-1)
 Update objective function. More...
 
Real value (const Vector< Real > &x, Real &tol)
 Compute value. More...
 
void gradient (Vector< Real > &g, const Vector< Real > &x, Real &tol)
 Compute gradient. More...
 
void hessVec (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
 Apply Hessian approximation to vector. More...
 
- Public Member Functions inherited from ROL::Objective< Real >
virtual ~Objective ()
 
virtual Real dirDeriv (const Vector< Real > &x, const Vector< Real > &d, Real &tol)
 Compute directional derivative. More...
 
virtual void invHessVec (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
 Apply inverse Hessian approximation to vector. More...
 
virtual void precond (Vector< Real > &Pv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
 Apply preconditioner to vector. More...
 
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)
 

Private Types

using size_type = typename std::vector< Real >::size_type
 

Private Member Functions

PartitionedVector< Real > & partition (Vector< Real > &x) const
 
const PartitionedVector< Real > & partition (const Vector< Real > &x) const
 
void throwError (const int err, const std::string &sfunc, const std::string &func, const int line) const
 
Real solveState (const Vector< Real > &x)
 
Real updateSketch (const Vector< Real > &x, const Real tol)
 
void setTerminalCondition (Vector< Real > &l, const Vector< Real > &uold, const Vector< Real > &unew, const Vector< Real > &z, const TimeStamp< Real > &ts)
 
void computeAdjointRHS (Vector< Real > &rhs, const Vector< Real > &l, const Vector< Real > &uold, const Vector< Real > &unew, const Vector< Real > &z, const TimeStamp< Real > &ts)
 
void advanceAdjoint (Vector< Real > &l, Vector< Real > &rhs, const Vector< Real > &uold, const Vector< Real > &unew, const Vector< Real > &z, const TimeStamp< Real > &ts)
 
void solveAdjoint (const Vector< Real > &x)
 
void updateGradient (Vector< Real > &g, const Vector< Real > &l, const Vector< Real > &uold, const Vector< Real > &unew, const Vector< Real > &z, const TimeStamp< Real > &ts)
 
void advanceStateSens (Vector< Real > &wnew, const Vector< Real > &v, const Vector< Real > &wold, const Vector< Real > &uold, const Vector< Real > &unew, const Vector< Real > &z, const TimeStamp< Real > &ts)
 
void setTerminalConditionHess (Vector< Real > &p, const Vector< Real > &v, const Vector< Real > &l, const Vector< Real > &wold, const Vector< Real > &wnew, const Vector< Real > &uold, const Vector< Real > &unew, const Vector< Real > &z, const TimeStamp< Real > &ts)
 
void computeOldStateHessLag (Vector< Real > &Hv, const Vector< Real > &l, const Vector< Real > &wold, const Vector< Real > &wnew, const Vector< Real > &uold, const Vector< Real > &unew, const Vector< Real > &z, const TimeStamp< Real > &ts, const bool sumInto=false)
 
void computeOldMixedHessLag (Vector< Real > &Hv, const Vector< Real > &v, const Vector< Real > &l, const Vector< Real > &uold, const Vector< Real > &unew, const Vector< Real > &z, const TimeStamp< Real > &ts, const bool sumInto=false)
 
void computeNewStateJacobian (Vector< Real > &Hv, const Vector< Real > &p, const Vector< Real > &uold, const Vector< Real > &unew, const Vector< Real > &z, const TimeStamp< Real > &ts, const bool sumInto=false)
 
void computeNewStateHessLag (Vector< Real > &Hv, const Vector< Real > &l, const Vector< Real > &wold, const Vector< Real > &wnew, const Vector< Real > &uold, const Vector< Real > &unew, const Vector< Real > &z, const TimeStamp< Real > &ts, const bool sumInto=false)
 
void computeNewMixedHessLag (Vector< Real > &Hv, const Vector< Real > &v, const Vector< Real > &l, const Vector< Real > &uold, const Vector< Real > &unew, const Vector< Real > &z, const TimeStamp< Real > &ts, const bool sumInto=false)
 
void advanceAdjointSens (Vector< Real > &p, Vector< Real > &rhs, const Vector< Real > &uold, const Vector< Real > &unew, const Vector< Real > &z, const TimeStamp< Real > &ts)
 
void computeControlHessLag (Vector< Real > &Hv, const Vector< Real > &v, const Vector< Real > &l, const Vector< Real > &uold, const Vector< Real > &unew, const Vector< Real > &z, const TimeStamp< Real > &ts)
 
void addMixedHessLag (Vector< Real > &Hv, const Vector< Real > &l, const Vector< Real > &wold, const Vector< Real > &wnew, const Vector< Real > &uold, const Vector< Real > &unew, const Vector< Real > &z, const TimeStamp< Real > &ts)
 
void addAdjointSens (Vector< Real > &Hv, const Vector< Real > &p, const Vector< Real > &uold, const Vector< Real > &unew, const Vector< Real > &z, const TimeStamp< Real > &ts)
 

Private Attributes

const Ptr< DynamicObjective
< Real > > 
obj_
 
const Ptr< DynamicConstraint
< Real > > 
con_
 
const Ptr< Vector< Real > > u0_
 
const std::vector< TimeStamp
< Real > > 
timeStamp_
 
const size_type Nt_
 
const bool useSketch_
 
size_type rankState_
 
Ptr< Sketch< Real > > stateSketch_
 
size_type rankAdjoint_
 
Ptr< Sketch< Real > > adjointSketch_
 
size_type rankStateSens_
 
Ptr< Sketch< Real > > stateSensSketch_
 
const bool useInexact_
 
const size_type updateFactor_
 
const size_type maxRank_
 
const bool syncHessRank_
 
std::vector< Ptr< Vector< Real > > > uhist_
 
std::vector< Ptr< Vector< Real > > > lhist_
 
std::vector< Ptr< Vector< Real > > > whist_
 
std::vector< Ptr< Vector< Real > > > phist_
 
Ptr< Vector< Real > > cprimal_
 
Ptr< Vector< Real > > crhs_
 
Ptr< Vector< Real > > udual_
 
Ptr< Vector< Real > > rhs_
 
Ptr< Vector< Real > > zdual_
 
Real val_
 
bool isValueComputed_
 
bool isStateComputed_
 
bool isAdjointComputed_
 
bool useHessian_
 
bool useSymHess_
 
Ptr< std::ostream > stream_
 
const bool print_
 
const int freq_
 

Additional Inherited Members

- Protected Member Functions inherited from ROL::Objective< Real >
const std::vector< Real > getParameter (void) const
 

Detailed Description

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

Defines the reduced time-dependent objective function interface for simulation-based optimization.

This objective function implements the implicitly defined objective function given by

\[ F(z) := \sum_{n=1}^{N_t} f_n(u_{n-1}(z),u_n(z),z_n) \]

where \(f_n:\mathcal{U}\times\mathcal{U}\times\mathcal{Z}\to\mathbb{R}\), and \(u_n\in\mathcal{U}\) solves the system of equations

\[ c_n(u_{n-1},u_n,z_n) = 0,\quad n=1,\ldots,N_t \]

with \(u_0\) provided.

Disclaimer: This is currently only set up for single step time integrators and piecewise-constant-in-time controls.

Definition at line 80 of file ROL_ReducedDynamicObjective.hpp.

Member Typedef Documentation

template<typename Real >
using ROL::ReducedDynamicObjective< Real >::size_type = typename std::vector<Real>::size_type
private

Definition at line 81 of file ROL_ReducedDynamicObjective.hpp.

Constructor & Destructor Documentation

template<typename Real >
ROL::ReducedDynamicObjective< Real >::ReducedDynamicObjective ( const Ptr< DynamicObjective< Real >> &  obj,
const Ptr< DynamicConstraint< Real >> &  con,
const Ptr< Vector< Real >> &  u0,
const Ptr< Vector< Real >> &  zvec,
const Ptr< Vector< Real >> &  cvec,
const std::vector< TimeStamp< Real >> &  timeStamp,
ROL::ParameterList &  pl,
const Ptr< std::ostream > &  stream = nullPtr 
)
inline

Member Function Documentation

template<typename Real >
PartitionedVector<Real>& ROL::ReducedDynamicObjective< Real >::partition ( Vector< Real > &  x) const
inlineprivate
template<typename Real >
const PartitionedVector<Real>& ROL::ReducedDynamicObjective< Real >::partition ( const Vector< Real > &  x) const
inlineprivate

Definition at line 131 of file ROL_ReducedDynamicObjective.hpp.

template<typename Real >
void ROL::ReducedDynamicObjective< Real >::throwError ( const int  err,
const std::string &  sfunc,
const std::string &  func,
const int  line 
) const
inlineprivate
template<typename Real >
Ptr<Vector<Real> > ROL::ReducedDynamicObjective< Real >::makeDynamicVector ( const Vector< Real > &  x) const
inline
template<typename Real >
void ROL::ReducedDynamicObjective< Real >::update ( const Vector< Real > &  x,
bool  flag = true,
int  iter = -1 
)
inlinevirtual
template<typename Real >
Real ROL::ReducedDynamicObjective< Real >::value ( const Vector< Real > &  x,
Real &  tol 
)
inlinevirtual
template<typename Real >
void ROL::ReducedDynamicObjective< 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 >.

Definition at line 315 of file ROL_ReducedDynamicObjective.hpp.

References ROL::ReducedDynamicObjective< Real >::adjointSketch_, ROL::ReducedDynamicObjective< Real >::advanceAdjoint(), ROL::ReducedDynamicObjective< Real >::computeAdjointRHS(), ROL::ReducedDynamicObjective< Real >::con_, ROL::PartitionedVector< Real >::get(), ROL::ReducedDynamicObjective< Real >::isAdjointComputed_, ROL::ReducedDynamicObjective< Real >::lhist_, ROL::ReducedDynamicObjective< Real >::Nt_, ROL::ReducedDynamicObjective< Real >::obj_, ROL::ReducedDynamicObjective< Real >::partition(), ROL::ReducedDynamicObjective< Real >::rankState_, ROL::ReducedDynamicObjective< Real >::rhs_, ROL::ReducedDynamicObjective< Real >::setTerminalCondition(), ROL::ReducedDynamicObjective< Real >::solveState(), ROL::ReducedDynamicObjective< Real >::stateSketch_, ROL::ReducedDynamicObjective< Real >::stream_, ROL::ReducedDynamicObjective< Real >::throwError(), ROL::ReducedDynamicObjective< Real >::timeStamp_, ROL::ReducedDynamicObjective< Real >::u0_, ROL::ReducedDynamicObjective< Real >::uhist_, ROL::ReducedDynamicObjective< Real >::updateGradient(), ROL::ReducedDynamicObjective< Real >::updateSketch(), ROL::ReducedDynamicObjective< Real >::useInexact_, and ROL::ReducedDynamicObjective< Real >::useSketch_.

template<typename Real >
void ROL::ReducedDynamicObjective< 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 >.

Definition at line 410 of file ROL_ReducedDynamicObjective.hpp.

References ROL::ReducedDynamicObjective< Real >::addAdjointSens(), ROL::ReducedDynamicObjective< Real >::addMixedHessLag(), ROL::ReducedDynamicObjective< Real >::adjointSketch_, ROL::ReducedDynamicObjective< Real >::advanceAdjointSens(), ROL::ReducedDynamicObjective< Real >::advanceStateSens(), ROL::ReducedDynamicObjective< Real >::computeControlHessLag(), ROL::ReducedDynamicObjective< Real >::computeNewMixedHessLag(), ROL::ReducedDynamicObjective< Real >::computeNewStateHessLag(), ROL::ReducedDynamicObjective< Real >::computeNewStateJacobian(), ROL::ReducedDynamicObjective< Real >::computeOldMixedHessLag(), ROL::ReducedDynamicObjective< Real >::computeOldStateHessLag(), ROL::ReducedDynamicObjective< Real >::con_, ROL::PartitionedVector< Real >::get(), ROL::Objective< Real >::hessVec(), ROL::ReducedDynamicObjective< Real >::lhist_, ROL::ReducedDynamicObjective< Real >::Nt_, ROL::ReducedDynamicObjective< Real >::obj_, ROL::ReducedDynamicObjective< Real >::partition(), ROL::ReducedDynamicObjective< Real >::phist_, ROL::ReducedDynamicObjective< Real >::rhs_, ROL::ReducedDynamicObjective< Real >::setTerminalConditionHess(), ROL::ReducedDynamicObjective< Real >::solveAdjoint(), ROL::ReducedDynamicObjective< Real >::solveState(), ROL::ReducedDynamicObjective< Real >::stateSensSketch_, ROL::ReducedDynamicObjective< Real >::stateSketch_, ROL::ReducedDynamicObjective< Real >::throwError(), ROL::ReducedDynamicObjective< Real >::timeStamp_, ROL::ReducedDynamicObjective< Real >::u0_, ROL::ReducedDynamicObjective< Real >::uhist_, ROL::ReducedDynamicObjective< Real >::useHessian_, ROL::ReducedDynamicObjective< Real >::useSketch_, ROL::ReducedDynamicObjective< Real >::useSymHess_, and ROL::ReducedDynamicObjective< Real >::whist_.

template<typename Real >
Real ROL::ReducedDynamicObjective< Real >::solveState ( const Vector< Real > &  x)
inlineprivate
template<typename Real >
Real ROL::ReducedDynamicObjective< Real >::updateSketch ( const Vector< Real > &  x,
const Real  tol 
)
inlineprivate
template<typename Real >
void ROL::ReducedDynamicObjective< Real >::setTerminalCondition ( Vector< Real > &  l,
const Vector< Real > &  uold,
const Vector< Real > &  unew,
const Vector< Real > &  z,
const TimeStamp< Real > &  ts 
)
inlineprivate
template<typename Real >
void ROL::ReducedDynamicObjective< Real >::computeAdjointRHS ( Vector< Real > &  rhs,
const Vector< Real > &  l,
const Vector< Real > &  uold,
const Vector< Real > &  unew,
const Vector< Real > &  z,
const TimeStamp< Real > &  ts 
)
inlineprivate
template<typename Real >
void ROL::ReducedDynamicObjective< Real >::advanceAdjoint ( Vector< Real > &  l,
Vector< Real > &  rhs,
const Vector< Real > &  uold,
const Vector< Real > &  unew,
const Vector< Real > &  z,
const TimeStamp< Real > &  ts 
)
inlineprivate
template<typename Real >
void ROL::ReducedDynamicObjective< Real >::solveAdjoint ( const Vector< Real > &  x)
inlineprivate
template<typename Real >
void ROL::ReducedDynamicObjective< Real >::updateGradient ( Vector< Real > &  g,
const Vector< Real > &  l,
const Vector< Real > &  uold,
const Vector< Real > &  unew,
const Vector< Real > &  z,
const TimeStamp< Real > &  ts 
)
inlineprivate
template<typename Real >
void ROL::ReducedDynamicObjective< Real >::advanceStateSens ( Vector< Real > &  wnew,
const Vector< Real > &  v,
const Vector< Real > &  wold,
const Vector< Real > &  uold,
const Vector< Real > &  unew,
const Vector< Real > &  z,
const TimeStamp< Real > &  ts 
)
inlineprivate
template<typename Real >
void ROL::ReducedDynamicObjective< Real >::setTerminalConditionHess ( Vector< Real > &  p,
const Vector< Real > &  v,
const Vector< Real > &  l,
const Vector< Real > &  wold,
const Vector< Real > &  wnew,
const Vector< Real > &  uold,
const Vector< Real > &  unew,
const Vector< Real > &  z,
const TimeStamp< Real > &  ts 
)
inlineprivate
template<typename Real >
void ROL::ReducedDynamicObjective< Real >::computeOldStateHessLag ( Vector< Real > &  Hv,
const Vector< Real > &  l,
const Vector< Real > &  wold,
const Vector< Real > &  wnew,
const Vector< Real > &  uold,
const Vector< Real > &  unew,
const Vector< Real > &  z,
const TimeStamp< Real > &  ts,
const bool  sumInto = false 
)
inlineprivate
template<typename Real >
void ROL::ReducedDynamicObjective< Real >::computeOldMixedHessLag ( Vector< Real > &  Hv,
const Vector< Real > &  v,
const Vector< Real > &  l,
const Vector< Real > &  uold,
const Vector< Real > &  unew,
const Vector< Real > &  z,
const TimeStamp< Real > &  ts,
const bool  sumInto = false 
)
inlineprivate
template<typename Real >
void ROL::ReducedDynamicObjective< Real >::computeNewStateJacobian ( Vector< Real > &  Hv,
const Vector< Real > &  p,
const Vector< Real > &  uold,
const Vector< Real > &  unew,
const Vector< Real > &  z,
const TimeStamp< Real > &  ts,
const bool  sumInto = false 
)
inlineprivate
template<typename Real >
void ROL::ReducedDynamicObjective< Real >::computeNewStateHessLag ( Vector< Real > &  Hv,
const Vector< Real > &  l,
const Vector< Real > &  wold,
const Vector< Real > &  wnew,
const Vector< Real > &  uold,
const Vector< Real > &  unew,
const Vector< Real > &  z,
const TimeStamp< Real > &  ts,
const bool  sumInto = false 
)
inlineprivate
template<typename Real >
void ROL::ReducedDynamicObjective< Real >::computeNewMixedHessLag ( Vector< Real > &  Hv,
const Vector< Real > &  v,
const Vector< Real > &  l,
const Vector< Real > &  uold,
const Vector< Real > &  unew,
const Vector< Real > &  z,
const TimeStamp< Real > &  ts,
const bool  sumInto = false 
)
inlineprivate
template<typename Real >
void ROL::ReducedDynamicObjective< Real >::advanceAdjointSens ( Vector< Real > &  p,
Vector< Real > &  rhs,
const Vector< Real > &  uold,
const Vector< Real > &  unew,
const Vector< Real > &  z,
const TimeStamp< Real > &  ts 
)
inlineprivate
template<typename Real >
void ROL::ReducedDynamicObjective< Real >::computeControlHessLag ( Vector< Real > &  Hv,
const Vector< Real > &  v,
const Vector< Real > &  l,
const Vector< Real > &  uold,
const Vector< Real > &  unew,
const Vector< Real > &  z,
const TimeStamp< Real > &  ts 
)
inlineprivate
template<typename Real >
void ROL::ReducedDynamicObjective< Real >::addMixedHessLag ( Vector< Real > &  Hv,
const Vector< Real > &  l,
const Vector< Real > &  wold,
const Vector< Real > &  wnew,
const Vector< Real > &  uold,
const Vector< Real > &  unew,
const Vector< Real > &  z,
const TimeStamp< Real > &  ts 
)
inlineprivate
template<typename Real >
void ROL::ReducedDynamicObjective< Real >::addAdjointSens ( Vector< Real > &  Hv,
const Vector< Real > &  p,
const Vector< Real > &  uold,
const Vector< Real > &  unew,
const Vector< Real > &  z,
const TimeStamp< Real > &  ts 
)
inlineprivate

Member Data Documentation

template<typename Real >
const Ptr<DynamicObjective<Real> > ROL::ReducedDynamicObjective< Real >::obj_
private
template<typename Real >
const Ptr<DynamicConstraint<Real> > ROL::ReducedDynamicObjective< Real >::con_
private
template<typename Real >
const Ptr<Vector<Real> > ROL::ReducedDynamicObjective< Real >::u0_
private
template<typename Real >
const std::vector<TimeStamp<Real> > ROL::ReducedDynamicObjective< Real >::timeStamp_
private
template<typename Real >
const size_type ROL::ReducedDynamicObjective< Real >::Nt_
private
template<typename Real >
const bool ROL::ReducedDynamicObjective< Real >::useSketch_
private
template<typename Real >
size_type ROL::ReducedDynamicObjective< Real >::rankState_
private
template<typename Real >
Ptr<Sketch<Real> > ROL::ReducedDynamicObjective< Real >::stateSketch_
private
template<typename Real >
size_type ROL::ReducedDynamicObjective< Real >::rankAdjoint_
private
template<typename Real >
Ptr<Sketch<Real> > ROL::ReducedDynamicObjective< Real >::adjointSketch_
private
template<typename Real >
size_type ROL::ReducedDynamicObjective< Real >::rankStateSens_
private
template<typename Real >
Ptr<Sketch<Real> > ROL::ReducedDynamicObjective< Real >::stateSensSketch_
private
template<typename Real >
const bool ROL::ReducedDynamicObjective< Real >::useInexact_
private
template<typename Real >
const size_type ROL::ReducedDynamicObjective< Real >::updateFactor_
private
template<typename Real >
const size_type ROL::ReducedDynamicObjective< Real >::maxRank_
private
template<typename Real >
const bool ROL::ReducedDynamicObjective< Real >::syncHessRank_
private
template<typename Real >
std::vector<Ptr<Vector<Real> > > ROL::ReducedDynamicObjective< Real >::uhist_
private
template<typename Real >
std::vector<Ptr<Vector<Real> > > ROL::ReducedDynamicObjective< Real >::lhist_
private
template<typename Real >
std::vector<Ptr<Vector<Real> > > ROL::ReducedDynamicObjective< Real >::whist_
private
template<typename Real >
std::vector<Ptr<Vector<Real> > > ROL::ReducedDynamicObjective< Real >::phist_
private
template<typename Real >
Ptr<Vector<Real> > ROL::ReducedDynamicObjective< Real >::cprimal_
private
template<typename Real >
Ptr<Vector<Real> > ROL::ReducedDynamicObjective< Real >::crhs_
private
template<typename Real >
Ptr<Vector<Real> > ROL::ReducedDynamicObjective< Real >::udual_
private
template<typename Real >
Ptr<Vector<Real> > ROL::ReducedDynamicObjective< Real >::rhs_
private
template<typename Real >
Ptr<Vector<Real> > ROL::ReducedDynamicObjective< Real >::zdual_
private
template<typename Real >
Real ROL::ReducedDynamicObjective< Real >::val_
private
template<typename Real >
bool ROL::ReducedDynamicObjective< Real >::isValueComputed_
private
template<typename Real >
bool ROL::ReducedDynamicObjective< Real >::isStateComputed_
private
template<typename Real >
bool ROL::ReducedDynamicObjective< Real >::isAdjointComputed_
private
template<typename Real >
bool ROL::ReducedDynamicObjective< Real >::useHessian_
private
template<typename Real >
bool ROL::ReducedDynamicObjective< Real >::useSymHess_
private
template<typename Real >
Ptr<std::ostream> ROL::ReducedDynamicObjective< Real >::stream_
private
template<typename Real >
const bool ROL::ReducedDynamicObjective< Real >::print_
private
template<typename Real >
const int ROL::ReducedDynamicObjective< Real >::freq_
private

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