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

Defines the time-dependent objective function interface for simulation-based optimization. Computes time-local contributions of value, gradient, Hessian-vector product etc to a larger composite objective defined over the simulation time. In contrast to other objective classes Objective_TimeSimOpt has a default implementation of value which returns zero, as time-dependent simulation based optimization problems may have an objective value which depends only on the final state of the system. More...

#include <ROL_Objective_TimeSimOpt.hpp>

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

Public Member Functions

virtual void update (const Vector< Real > &u_old, const Vector< Real > &u_new, const Vector< Real > &z, bool flag=true, int iter=-1)
 Update constraint functions. u_old Is the state from the end of the previous time step. u_new Is the state from the end of this time step. z Is the control variable flag = true if optimization variable is changed, iter is the outer algorithm iterations count. More...
 
- Public Member Functions inherited from ROL::Objective_SimOpt< Real >
virtual void update (const Vector< Real > &u, const Vector< Real > &z, bool flag=true, int iter=-1)
 Update objective function. u is an iterate, z is an iterate, flag = true if the iterate has changed, iter is the outer algorithm iterations count. More...
 
void update (const Vector< Real > &x, bool flag=true, int iter=-1)
 Update objective function. More...
 
virtual Real value (const Vector< Real > &u, const Vector< Real > &z, Real &tol)=0
 Compute value. More...
 
Real value (const Vector< Real > &x, Real &tol)
 Compute value. More...
 
virtual void gradient_1 (Vector< Real > &g, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
 Compute gradient with respect to first component. More...
 
virtual void gradient_2 (Vector< Real > &g, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
 Compute gradient with respect to second component. More...
 
void gradient (Vector< Real > &g, const Vector< Real > &x, Real &tol)
 Compute gradient. More...
 
virtual void hessVec_11 (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
 Apply Hessian approximation to vector. More...
 
virtual void hessVec_12 (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
 
virtual void hessVec_21 (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
 
virtual void hessVec_22 (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
 
void hessVec (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
 Apply Hessian approximation to vector. More...
 
std::vector< std::vector< Real > > checkGradient_1 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &d, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
 
std::vector< std::vector< Real > > checkGradient_1 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &g, const Vector< Real > &d, const bool printToStream, std::ostream &outStream, const int numSteps, const int order)
 
std::vector< std::vector< Real > > checkGradient_1 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &g, const Vector< Real > &d, const std::vector< Real > &steps, const bool printToStream, std::ostream &outStream, const int order)
 
std::vector< std::vector< Real > > checkGradient_2 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &d, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
 
std::vector< std::vector< Real > > checkGradient_2 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &g, const Vector< Real > &d, const bool printToStream, std::ostream &outStream, const int numSteps, const int order)
 
std::vector< std::vector< Real > > checkGradient_2 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &g, const Vector< Real > &d, const std::vector< Real > &steps, const bool printToStream, std::ostream &outStream, const int order)
 
std::vector< std::vector< Real > > checkHessVec_11 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &v, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
 
std::vector< std::vector< Real > > checkHessVec_11 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &v, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
 
std::vector< std::vector< Real > > checkHessVec_11 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &hv, const Vector< Real > &v, const bool printToStream, std::ostream &outStream, const int numSteps, const int order)
 
std::vector< std::vector< Real > > checkHessVec_11 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &hv, const Vector< Real > &v, const std::vector< Real > &steps, const bool printToStream, std::ostream &outStream, const int order)
 
std::vector< std::vector< Real > > checkHessVec_12 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &v, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
 
std::vector< std::vector< Real > > checkHessVec_12 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &v, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
 
std::vector< std::vector< Real > > checkHessVec_12 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &hv, const Vector< Real > &v, const bool printToStream, std::ostream &outStream, const int numSteps, const int order)
 
std::vector< std::vector< Real > > checkHessVec_12 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &hv, const Vector< Real > &v, const std::vector< Real > &steps, const bool printToStream, std::ostream &outStream, const int order)
 
std::vector< std::vector< Real > > checkHessVec_21 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &v, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
 
std::vector< std::vector< Real > > checkHessVec_21 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &v, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
 
std::vector< std::vector< Real > > checkHessVec_21 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &hv, const Vector< Real > &v, const bool printToStream, std::ostream &outStream, const int numSteps, const int order)
 
std::vector< std::vector< Real > > checkHessVec_21 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &hv, const Vector< Real > &v, const std::vector< Real > &steps, const bool printToStream, std::ostream &outStream, const int order)
 
std::vector< std::vector< Real > > checkHessVec_22 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &v, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
 
std::vector< std::vector< Real > > checkHessVec_22 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &v, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
 
std::vector< std::vector< Real > > checkHessVec_22 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &hv, const Vector< Real > &v, const bool printToStream, std::ostream &outStream, const int numSteps, const int order)
 
std::vector< std::vector< Real > > checkHessVec_22 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &hv, const Vector< Real > &v, const std::vector< Real > &steps, const bool printToStream, std::ostream &outStream, const int order)
 
- 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)
 

Protected Member Functions

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

Private Member Functions

template<int I>
Vector< Real > & getVector (Vector< Real > &x) const
 
template<int I>
const Vector< Real > & getVector (const Vector< Real > &x) const
 

Private Attributes

VectorWorkspace< Real > workspace_
 

Detailed Description

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

Defines the time-dependent objective function interface for simulation-based optimization. Computes time-local contributions of value, gradient, Hessian-vector product etc to a larger composite objective defined over the simulation time. In contrast to other objective classes Objective_TimeSimOpt has a default implementation of value which returns zero, as time-dependent simulation based optimization problems may have an objective value which depends only on the final state of the system.

      This objective interface inherits from ROL_Objective_SimOpt. Though the interface
      takes two simulation space vectors from spaces

\(\mathcal{U_o}\times\mathcal{U_n}\). The space \(\mathcal{U_o}\) is ``old'' information that accounts for the initial condition on the time interval.

Comments: It may be worthwhile to provide an implementation of OptimizationProblem for TimeSimOpt which recognizes this common use case and avoids unnecessary objective calls.

NOTE: As written, this interface is step agnostic, which needs to be changed if a final-time cost is desire or if the time-distributed cost is to be weighted non-uniformly


Definition at line 80 of file ROL_Objective_TimeSimOpt.hpp.

Member Function Documentation

template<typename Real >
template<int I>
Vector<Real>& ROL::Objective_TimeSimOpt< Real >::getVector ( Vector< Real > &  x) const
inlineprivate

Definition at line 85 of file ROL_Objective_TimeSimOpt.hpp.

template<typename Real >
template<int I>
const Vector<Real>& ROL::Objective_TimeSimOpt< Real >::getVector ( const Vector< Real > &  x) const
inlineprivate

Definition at line 90 of file ROL_Objective_TimeSimOpt.hpp.

template<typename Real >
VectorWorkspace<Real>& ROL::Objective_TimeSimOpt< Real >::getVectorWorkspace ( ) const
inlineprotected
template<typename Real >
virtual void ROL::Objective_TimeSimOpt< Real >::update ( const Vector< Real > &  u_old,
const Vector< Real > &  u_new,
const Vector< Real > &  z,
bool  flag = true,
int  iter = -1 
)
inlinevirtual

Update constraint functions. u_old Is the state from the end of the previous time step. u_new Is the state from the end of this time step. z Is the control variable flag = true if optimization variable is changed, iter is the outer algorithm iterations count.

Definition at line 113 of file ROL_Objective_TimeSimOpt.hpp.

References ROL::update_1_new(), ROL::update_1_old(), and ROL::update_2().

Member Data Documentation

template<typename Real >
VectorWorkspace<Real> ROL::Objective_TimeSimOpt< Real >::workspace_
mutableprivate

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