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

Defines the time dependent constraint operator interface for simulation-based optimization. More...

#include <ROL_Constraint_TimeSimOpt.hpp>

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

Public Member Functions

 Constraint_TimeSimOpt ()
 
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...
 
virtual void update_1_old (const Vector< Real > &u_old, bool flag=true, int iter=-1)
 Update constraint functions with respect to Sim variable. u_old is the state variable flag = true if optimization variable is changed, iter is the outer algorithm iterations count. More...
 
virtual void update_1_new (const Vector< Real > &u_new, bool flag=true, int iter=-1)
 Update constraint functions with respect to Sim variable. u_new is the state variable flag = true if optimization variable is changed, iter is the outer algorithm iterations count. More...
 
virtual void update_2 (const Vector< Real > &z, bool flag=true, int iter=-1) override
 Update constraint functions with respect to Opt variable. z is the control variable, flag = true if optimization variable is changed, iter is the outer algorithm iterations count. More...
 
virtual void value (Vector< Real > &c, const Vector< Real > &u_old, const Vector< Real > &u_new, const Vector< Real > &z, Real &tol)=0
 Evaluate the constraint operator \(c:\mathcal{U_o}\times\mathcal{U_n}\times\mathcal{Z} \rightarrow \mathcal{C}\) at \((u,z)\). More...
 
virtual void solve (Vector< Real > &c, const Vector< Real > &u_old, Vector< Real > &u_new, const Vector< Real > &z, Real &tol)=0
 
virtual void applyJacobian_1_old (Vector< Real > &jv, const Vector< Real > &v_old, const Vector< Real > &u_old, const Vector< Real > &u_new, const Vector< Real > &z, Real &tol)=0
 
virtual void applyJacobian_1_new (Vector< Real > &jv, const Vector< Real > &v_new, const Vector< Real > &u_old, const Vector< Real > &u_new, const Vector< Real > &z, Real &tol)=0
 
virtual void applyInverseJacobian_1_new (Vector< Real > &ijv, const Vector< Real > &v_new, const Vector< Real > &u_old, const Vector< Real > &u_new, const Vector< Real > &z, Real &tol)=0
 
virtual void applyJacobian_2 (Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &u_old, const Vector< Real > &u_new, const Vector< Real > &z, Real &tol)=0
 
virtual void applyAdjointJacobian_1_old (Vector< Real > &ajv_old, const Vector< Real > &dualv, const Vector< Real > &u_old, const Vector< Real > &u_new, const Vector< Real > &z, Real &tol)=0
 
virtual void applyAdjointJacobian_1_new (Vector< Real > &ajv_new, const Vector< Real > &dualv, const Vector< Real > &u_old, const Vector< Real > &u_new, const Vector< Real > &z, Real &tol)=0
 
virtual void applyInverseAdjointJacobian_1_new (Vector< Real > &iajv, const Vector< Real > &v_new, const Vector< Real > &u_old, const Vector< Real > &u_new, const Vector< Real > &z, Real &tol)=0
 
virtual void applyAdjointJacobian_2_time (Vector< Real > &ajv, const Vector< Real > &dualv, const Vector< Real > &u_old, const Vector< Real > &u_new, const Vector< Real > &z, Real &tol)=0
 
virtual void applyAdjointHessian_11_old (Vector< Real > &ahwv_old, const Vector< Real > &w, const Vector< Real > &v_new, const Vector< Real > &u_old, const Vector< Real > &u_new, const Vector< Real > &z, Real &tol)=0
 
virtual void applyAdjointHessian_11_new (Vector< Real > &ahwv_new, const Vector< Real > &w, const Vector< Real > &v_new, const Vector< Real > &u_old, const Vector< Real > &u_new, const Vector< Real > &z, Real &tol)=0
 
virtual void update (const Vector< Real > &u, const Vector< Real > &z, bool flag=true, int iter=-1) override
 Update constraint functions. x is the optimization variable, flag = true if optimization variable is changed, iter is the outer algorithm iterations count. More...
 
virtual void value (Vector< Real > &c, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
 Evaluate the constraint operator \(c:\mathcal{U}\times\mathcal{Z} \rightarrow \mathcal{C}\) at \((u,z)\). More...
 
virtual void solve (Vector< Real > &c, Vector< Real > &u, const Vector< Real > &z, Real &tol) override
 Given \(z\), solve \(c(u,z)=0\) for \(u\). More...
 
virtual void applyJacobian_1 (Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
 Apply the partial constraint Jacobian at \((u,z)\), \(c_u(u,z) \in L(\mathcal{U}, \mathcal{C})\), to the vector \(v\). More...
 
virtual void applyJacobian_2 (Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
 Apply the partial constraint Jacobian at \((u,z)\), \(c_z(u,z) \in L(\mathcal{Z}, \mathcal{C})\), to the vector \(v\). More...
 
virtual void applyInverseJacobian_1 (Vector< Real > &ijv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) overridefinal
 Apply the inverse partial constraint Jacobian at \((u,z)\), \(c_u(u,z)^{-1} \in L(\mathcal{C}, \mathcal{U})\), to the vector \(v\). More...
 
virtual void applyAdjointJacobian_1 (Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
 Apply the adjoint of the partial constraint Jacobian at \((u,z)\), \(c_u(u,z)^* \in L(\mathcal{C}^*, \mathcal{U}^*)\), to the vector \(v\). This is the primary interface. More...
 
virtual void applyAdjointJacobian_2 (Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
 Apply the adjoint of the partial constraint Jacobian at \((u,z)\), \(c_z(u,z)^* \in L(\mathcal{C}^*, \mathcal{Z}^*)\), to vector \(v\). This is the primary interface. More...
 
virtual void applyInverseAdjointJacobian_1 (Vector< Real > &iajv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) overridefinal
 Apply the inverse of the adjoint of the partial constraint Jacobian at \((u,z)\), \(c_u(u,z)^{-*} \in L(\mathcal{U}^*, \mathcal{C}^*)\), to the vector \(v\). More...
 
virtual void applyAdjointHessian_11 (Vector< Real > &ahwv, const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
 Apply the simulation-space derivative of the adjoint of the constraint simulation-space Jacobian at \((u,z)\) to the vector \(w\) in the direction \(v\), according to \(v\mapsto c_{uu}(u,z)(v,\cdot)^*w\). More...
 
virtual void applyAdjointHessian_12 (Vector< Real > &ahwv, const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
 Apply the optimization-space derivative of the adjoint of the constraint simulation-space Jacobian at \((u,z)\) to the vector \(w\) in the direction \(v\), according to \(v\mapsto c_{uz}(u,z)(v,\cdot)^*w\). More...
 
virtual void applyAdjointHessian_21 (Vector< Real > &ahwv, const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
 Apply the simulation-space derivative of the adjoint of the constraint optimization-space Jacobian at \((u,z)\) to the vector \(w\) in the direction \(v\), according to \(v\mapsto c_{zu}(u,z)(v,\cdot)^*w\). More...
 
virtual void applyAdjointHessian_22 (Vector< Real > &ahwv, const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
 Apply the optimization-space derivative of the adjoint of the constraint optimization-space Jacobian at \((u,z)\) to the vector \(w\) in the direction \(v\), according to \(v\mapsto c_{zz}(u,z)(v,\cdot)^*w\). More...
 
virtual Real checkSolve (const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, const ROL::Vector< Real > &c, const bool printToStream=true, std::ostream &outStream=std::cout) override
 
virtual Real checkInverseJacobian_1_new (const ROL::Vector< Real > &c, const ROL::Vector< Real > &u_new, const ROL::Vector< Real > &u_old, const ROL::Vector< Real > &z, const ROL::Vector< Real > &v_new, const bool printToStream=true, std::ostream &outStream=std::cout)
 
virtual Real checkInverseAdjointJacobian_1_new (const ROL::Vector< Real > &c, const ROL::Vector< Real > &u_new, const ROL::Vector< Real > &u_old, const ROL::Vector< Real > &z, const ROL::Vector< Real > &v_new, const bool printToStream=true, std::ostream &outStream=std::cout)
 
std::vector< std::vector< Real > > checkApplyJacobian_1_new (const Vector< Real > &u_new, const Vector< Real > &u_old, const Vector< Real > &z, const Vector< Real > &v, const Vector< Real > &jv, 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 > > checkApplyJacobian_1_new (const Vector< Real > &u_new, const Vector< Real > &u_old, const Vector< Real > &z, const Vector< Real > &v, const Vector< Real > &jv, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
 
- Public Member Functions inherited from ROL::ROL::Constraint_SimOpt< Real >
 Constraint_SimOpt ()
 
virtual void update_1 (const Vector< Real > &u, bool flag=true, int iter=-1)
 Update constraint functions with respect to Sim variable. x is the optimization variable, flag = true if optimization variable is changed, iter is the outer algorithm iterations count. More...
 
virtual void setSolveParameters (ParameterList &parlist)
 Set solve parameters. More...
 
virtual void applyAdjointJacobian_1 (Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &dualv, Real &tol)
 Apply the adjoint of the partial constraint Jacobian at \((u,z)\), \(c_u(u,z)^* \in L(\mathcal{C}^*, \mathcal{U}^*)\), to the vector \(v\). This is the secondary interface, for use with dual spaces where the user does not define the dual() operation. More...
 
virtual void applyAdjointJacobian_2 (Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &dualv, Real &tol)
 Apply the adjoint of the partial constraint Jacobian at \((u,z)\), \(c_z(u,z)^* \in L(\mathcal{C}^*, \mathcal{Z}^*)\), to vector \(v\). This is the secondary interface, for use with dual spaces where the user does not define the dual() operation. More...
 
virtual std::vector< Real > solveAugmentedSystem (Vector< Real > &v1, Vector< Real > &v2, const Vector< Real > &b1, const Vector< Real > &b2, const Vector< Real > &x, Real &tol)
 Approximately solves the augmented system

\[ \begin{pmatrix} I & c'(x)^* \\ c'(x) & 0 \end{pmatrix} \begin{pmatrix} v_{1} \\ v_{2} \end{pmatrix} = \begin{pmatrix} b_{1} \\ b_{2} \end{pmatrix} \]

where \(v_{1} \in \mathcal{X}\), \(v_{2} \in \mathcal{C}^*\), \(b_{1} \in \mathcal{X}^*\), \(b_{2} \in \mathcal{C}\), \(I : \mathcal{X} \rightarrow \mathcal{X}^*\) is an identity operator, and \(0 : \mathcal{C}^* \rightarrow \mathcal{C}\) is a zero operator. More...

 
virtual void applyPreconditioner (Vector< Real > &pv, const Vector< Real > &v, const Vector< Real > &x, const Vector< Real > &g, Real &tol)
 Apply a constraint preconditioner at \(x\), \(P(x) \in L(\mathcal{C}, \mathcal{C})\), to vector \(v\). In general, this preconditioner satisfies the following relationship:

\[ c'(x) c'(x)^* P(x) v \approx v \,. \]

It is used by the solveAugmentedSystem method. More...

 
virtual void update (const Vector< Real > &x, bool flag=true, int iter=-1)
 Update constraint functions. x is the optimization variable, flag = true if optimization variable is changed, iter is the outer algorithm iterations count. More...
 
virtual void value (Vector< Real > &c, const Vector< Real > &x, Real &tol)
 Evaluate the constraint operator \(c:\mathcal{X} \rightarrow \mathcal{C}\) at \(x\). More...
 
virtual void applyJacobian (Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
 Apply the constraint Jacobian at \(x\), \(c'(x) \in L(\mathcal{X}, \mathcal{C})\), to vector \(v\). More...
 
virtual void applyAdjointJacobian (Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
 Apply the adjoint of the the constraint Jacobian at \(x\), \(c'(x)^* \in L(\mathcal{C}^*, \mathcal{X}^*)\), to vector \(v\). More...
 
virtual void applyAdjointHessian (Vector< Real > &ahwv, const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
 Apply the derivative of the adjoint of the constraint Jacobian at \(x\) to vector \(u\) in direction \(v\), according to \( v \mapsto c''(x)(v,\cdot)^*u \). More...
 
virtual Real checkAdjointConsistencyJacobian_1 (const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const bool printToStream=true, std::ostream &outStream=std::cout)
 Check the consistency of the Jacobian and its adjoint. This is the primary interface. More...
 
virtual Real checkAdjointConsistencyJacobian_1 (const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &dualw, const Vector< Real > &dualv, const bool printToStream=true, std::ostream &outStream=std::cout)
 Check the consistency of the Jacobian and its adjoint. This is the secondary interface, for use with dual spaces where the user does not define the dual() operation. More...
 
virtual Real checkAdjointConsistencyJacobian_2 (const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const bool printToStream=true, std::ostream &outStream=std::cout)
 Check the consistency of the Jacobian and its adjoint. This is the primary interface. More...
 
virtual Real checkAdjointConsistencyJacobian_2 (const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &dualw, const Vector< Real > &dualv, const bool printToStream=true, std::ostream &outStream=std::cout)
 Check the consistency of the Jacobian and its adjoint. This is the secondary interface, for use with dual spaces where the user does not define the dual() operation. More...
 
virtual Real checkInverseJacobian_1 (const Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const bool printToStream=true, std::ostream &outStream=std::cout)
 
virtual Real checkInverseAdjointJacobian_1 (const Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const bool printToStream=true, std::ostream &outStream=std::cout)
 
std::vector< std::vector< Real > > checkApplyJacobian_1 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &v, const Vector< Real > &jv, 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 > > checkApplyJacobian_1 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &v, const Vector< Real > &jv, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
 
std::vector< std::vector< Real > > checkApplyJacobian_2 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &v, const Vector< Real > &jv, 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 > > checkApplyJacobian_2 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &v, const Vector< Real > &jv, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
 
std::vector< std::vector< Real > > checkApplyAdjointHessian_11 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &p, const Vector< Real > &v, const Vector< Real > &hv, 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 > > checkApplyAdjointHessian_11 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &p, const Vector< Real > &v, const Vector< Real > &hv, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
 
std::vector< std::vector< Real > > checkApplyAdjointHessian_21 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &p, const Vector< Real > &v, const Vector< Real > &hv, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
 \( u\in U \), \( z\in Z \), \( p\in C^\ast \), \( v \in U \), \( hv \in U^\ast \) More...
 
std::vector< std::vector< Real > > checkApplyAdjointHessian_21 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &p, const Vector< Real > &v, const Vector< Real > &hv, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
 \( u\in U \), \( z\in Z \), \( p\in C^\ast \), \( v \in U \), \( hv \in U^\ast \) More...
 
std::vector< std::vector< Real > > checkApplyAdjointHessian_12 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &p, const Vector< Real > &v, const Vector< Real > &hv, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
 \( u\in U \), \( z\in Z \), \( p\in C^\ast \), \( v \in U \), \( hv \in U^\ast \) More...
 
std::vector< std::vector< Real > > checkApplyAdjointHessian_12 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &p, const Vector< Real > &v, const Vector< Real > &hv, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
 
std::vector< std::vector< Real > > checkApplyAdjointHessian_22 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &p, const Vector< Real > &v, const Vector< Real > &hv, 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 > > checkApplyAdjointHessian_22 (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &p, const Vector< Real > &v, const Vector< Real > &hv, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
 
- Public Member Functions inherited from ROL::Constraint< Real >
virtual ~Constraint (void)
 
 Constraint (void)
 
virtual void applyAdjointJacobian (Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &x, const Vector< Real > &dualv, Real &tol)
 Apply the adjoint of the the constraint Jacobian at \(x\), \(c'(x)^* \in L(\mathcal{C}^*, \mathcal{X}^*)\), to vector \(v\). More...
 
void activate (void)
 Turn on constraints. More...
 
void deactivate (void)
 Turn off constraints. More...
 
bool isActivated (void)
 Check if constraints are on. More...
 
virtual std::vector
< std::vector< Real > > 
checkApplyJacobian (const Vector< Real > &x, const Vector< Real > &v, const Vector< Real > &jv, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
 Finite-difference check for the constraint Jacobian application. More...
 
virtual std::vector
< std::vector< Real > > 
checkApplyJacobian (const Vector< Real > &x, const Vector< Real > &v, const Vector< Real > &jv, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
 Finite-difference check for the constraint Jacobian application. More...
 
virtual std::vector
< std::vector< Real > > 
checkApplyAdjointJacobian (const Vector< Real > &x, const Vector< Real > &v, const Vector< Real > &c, const Vector< Real > &ajv, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS)
 Finite-difference check for the application of the adjoint of constraint Jacobian. More...
 
virtual Real checkAdjointConsistencyJacobian (const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &x, const bool printToStream=true, std::ostream &outStream=std::cout)
 
virtual Real checkAdjointConsistencyJacobian (const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &x, const Vector< Real > &dualw, const Vector< Real > &dualv, const bool printToStream=true, std::ostream &outStream=std::cout)
 
virtual std::vector
< std::vector< Real > > 
checkApplyAdjointHessian (const Vector< Real > &x, const Vector< Real > &u, const Vector< Real > &v, const Vector< Real > &hv, const std::vector< Real > &step, const bool printToScreen=true, std::ostream &outStream=std::cout, const int order=1)
 Finite-difference check for the application of the adjoint of constraint Hessian. More...
 
virtual std::vector
< std::vector< Real > > 
checkApplyAdjointHessian (const Vector< Real > &x, const Vector< Real > &u, const Vector< Real > &v, const Vector< Real > &hv, const bool printToScreen=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
 Finite-difference check for the application of the adjoint of constraint Hessian. More...
 
virtual void setParameter (const std::vector< Real > &param)
 

Protected Member Functions

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

Private Member Functions

Vector< Real > & getNewVector (Vector< Real > &x) const
 
const Vector< Real > & getNewVector (const Vector< Real > &x) const
 
Vector< Real > & getOldVector (Vector< Real > &x) const
 
const Vector< Real > & getOldVector (const Vector< Real > &x) const
 

Private Attributes

VectorWorkspace< Real > workspace_
 

Additional Inherited Members

- Protected Attributes inherited from ROL::ROL::Constraint_SimOpt< Real >
Real atol_
 
Real rtol_
 
Real stol_
 
Real factor_
 
Real decr_
 
int maxit_
 
bool print_
 
bool zero_
 
int solverType_
 
bool firstSolve_
 

Detailed Description

template<class Real>
class ROL::Constraint_TimeSimOpt< Real >

Defines the time dependent constraint operator interface for simulation-based optimization.

This constraint interface inherits from ROL_Constraint_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. The space \(\mathcal{U_n}\) is the ``new'' variables that can be determined by satisfying constraints in the form

\[ c(u_o,u_n,z) = 0 \,. \]

where \(u_0 \in \mathcal{U_o},\; u_n\in\mathcal{U_n},\) and \(z\in\mathcal{Z}\). In this way this constraint defines a sequence of state variables. The basic operator interface, to be implemented by the user, requires:

The user may also overload:

Definition at line 100 of file ROL_Constraint_TimeSimOpt.hpp.

Constructor & Destructor Documentation

template<class Real >
ROL::Constraint_TimeSimOpt< Real >::Constraint_TimeSimOpt ( )
inline

Definition at line 136 of file ROL_Constraint_TimeSimOpt.hpp.

Member Function Documentation

template<class Real >
Vector<Real>& ROL::Constraint_TimeSimOpt< Real >::getNewVector ( Vector< Real > &  x) const
inlineprivate
template<class Real >
const Vector<Real>& ROL::Constraint_TimeSimOpt< Real >::getNewVector ( const Vector< Real > &  x) const
inlineprivate
template<class Real >
Vector<Real>& ROL::Constraint_TimeSimOpt< Real >::getOldVector ( Vector< Real > &  x) const
inlineprivate
template<class Real >
const Vector<Real>& ROL::Constraint_TimeSimOpt< Real >::getOldVector ( const Vector< Real > &  x) const
inlineprivate
template<class Real >
VectorWorkspace<Real>& ROL::Constraint_TimeSimOpt< Real >::getVectorWorkspace ( ) const
inlineprotected
template<class Real >
virtual void ROL::Constraint_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 149 of file ROL_Constraint_TimeSimOpt.hpp.

References ROL::Constraint_TimeSimOpt< Real >::update_1_new(), ROL::Constraint_TimeSimOpt< Real >::update_1_old(), and ROL::Constraint_TimeSimOpt< Real >::update_2().

Referenced by ROL::Constraint_TimeSimOpt< Real >::checkApplyJacobian_1_new(), ROL::Constraint_TimeSimOpt< Real >::checkInverseAdjointJacobian_1_new(), ROL::Constraint_TimeSimOpt< Real >::checkInverseJacobian_1_new(), ROL::Constraint_TimeSimOpt< Real >::checkSolve(), and ROL::Constraint_TimeSimOpt< Real >::update().

template<class Real >
virtual void ROL::Constraint_TimeSimOpt< Real >::update_1_old ( const Vector< Real > &  u_old,
bool  flag = true,
int  iter = -1 
)
inlinevirtual

Update constraint functions with respect to Sim variable. u_old is the state variable flag = true if optimization variable is changed, iter is the outer algorithm iterations count.

Definition at line 163 of file ROL_Constraint_TimeSimOpt.hpp.

Referenced by ROL::Constraint_TimeSimOpt< Real >::update().

template<class Real >
virtual void ROL::Constraint_TimeSimOpt< Real >::update_1_new ( const Vector< Real > &  u_new,
bool  flag = true,
int  iter = -1 
)
inlinevirtual

Update constraint functions with respect to Sim variable. u_new is the state variable flag = true if optimization variable is changed, iter is the outer algorithm iterations count.

Definition at line 170 of file ROL_Constraint_TimeSimOpt.hpp.

Referenced by ROL::Constraint_TimeSimOpt< Real >::update().

template<class Real >
virtual void ROL::Constraint_TimeSimOpt< Real >::update_2 ( const Vector< Real > &  z,
bool  flag = true,
int  iter = -1 
)
inlineoverridevirtual

Update constraint functions with respect to Opt variable. z is the control variable, flag = true if optimization variable is changed, iter is the outer algorithm iterations count.

Reimplemented from ROL::ROL::Constraint_SimOpt< Real >.

Definition at line 177 of file ROL_Constraint_TimeSimOpt.hpp.

Referenced by ROL::Constraint_TimeSimOpt< Real >::update().

template<class Real >
virtual void ROL::Constraint_TimeSimOpt< Real >::value ( Vector< Real > &  c,
const Vector< Real > &  u_old,
const Vector< Real > &  u_new,
const Vector< Real > &  z,
Real &  tol 
)
pure virtual

Evaluate the constraint operator \(c:\mathcal{U_o}\times\mathcal{U_n}\times\mathcal{Z} \rightarrow \mathcal{C}\) at \((u,z)\).

  @param[out]      c      is the result of evaluating the constraint operator at @b  \form#196; a constraint-space vector
  @param[in]       u_old  is the constraint argument; a simulation-space vector from the previous interval
  @param[in]       u_new  is the constraint argument; a simulation-space vector from the current interval
  @param[in]       z      is the constraint argument; an optimization-space vector
  @param[in,out]   tol    is a tolerance for inexact evaluations; currently unused

  On return, \form#197,
  where \form#73, \form#240, 

\(\mathsf{u_n} \in \mathcal{U_n}\), and $ \(\mathsf{z} \in\mathcal{Z}\).


Referenced by ROL::Constraint_TimeSimOpt< Real >::checkApplyJacobian_1_new(), ROL::Constraint_TimeSimOpt< Real >::checkSolve(), and ROL::Constraint_TimeSimOpt< Real >::value().

template<class Real >
virtual void ROL::Constraint_TimeSimOpt< Real >::solve ( Vector< Real > &  c,
const Vector< Real > &  u_old,
Vector< Real > &  u_new,
const Vector< Real > &  z,
Real &  tol 
)
pure virtual
template<class Real >
virtual void ROL::Constraint_TimeSimOpt< Real >::applyJacobian_1_old ( Vector< Real > &  jv,
const Vector< Real > &  v_old,
const Vector< Real > &  u_old,
const Vector< Real > &  u_new,
const Vector< Real > &  z,
Real &  tol 
)
pure virtual
template<class Real >
virtual void ROL::Constraint_TimeSimOpt< Real >::applyJacobian_1_new ( Vector< Real > &  jv,
const Vector< Real > &  v_new,
const Vector< Real > &  u_old,
const Vector< Real > &  u_new,
const Vector< Real > &  z,
Real &  tol 
)
pure virtual
template<class Real >
virtual void ROL::Constraint_TimeSimOpt< Real >::applyInverseJacobian_1_new ( Vector< Real > &  ijv,
const Vector< Real > &  v_new,
const Vector< Real > &  u_old,
const Vector< Real > &  u_new,
const Vector< Real > &  z,
Real &  tol 
)
pure virtual
template<class Real >
virtual void ROL::Constraint_TimeSimOpt< Real >::applyJacobian_2 ( Vector< Real > &  jv,
const Vector< Real > &  v,
const Vector< Real > &  u_old,
const Vector< Real > &  u_new,
const Vector< Real > &  z,
Real &  tol 
)
pure virtual
template<class Real >
virtual void ROL::Constraint_TimeSimOpt< Real >::applyAdjointJacobian_1_old ( Vector< Real > &  ajv_old,
const Vector< Real > &  dualv,
const Vector< Real > &  u_old,
const Vector< Real > &  u_new,
const Vector< Real > &  z,
Real &  tol 
)
pure virtual
template<class Real >
virtual void ROL::Constraint_TimeSimOpt< Real >::applyAdjointJacobian_1_new ( Vector< Real > &  ajv_new,
const Vector< Real > &  dualv,
const Vector< Real > &  u_old,
const Vector< Real > &  u_new,
const Vector< Real > &  z,
Real &  tol 
)
pure virtual
template<class Real >
virtual void ROL::Constraint_TimeSimOpt< Real >::applyInverseAdjointJacobian_1_new ( Vector< Real > &  iajv,
const Vector< Real > &  v_new,
const Vector< Real > &  u_old,
const Vector< Real > &  u_new,
const Vector< Real > &  z,
Real &  tol 
)
pure virtual
template<class Real >
virtual void ROL::Constraint_TimeSimOpt< Real >::applyAdjointJacobian_2_time ( Vector< Real > &  ajv,
const Vector< Real > &  dualv,
const Vector< Real > &  u_old,
const Vector< Real > &  u_new,
const Vector< Real > &  z,
Real &  tol 
)
pure virtual
template<class Real >
virtual void ROL::Constraint_TimeSimOpt< Real >::applyAdjointHessian_11_old ( Vector< Real > &  ahwv_old,
const Vector< Real > &  w,
const Vector< Real > &  v_new,
const Vector< Real > &  u_old,
const Vector< Real > &  u_new,
const Vector< Real > &  z,
Real &  tol 
)
pure virtual
template<class Real >
virtual void ROL::Constraint_TimeSimOpt< Real >::applyAdjointHessian_11_new ( Vector< Real > &  ahwv_new,
const Vector< Real > &  w,
const Vector< Real > &  v_new,
const Vector< Real > &  u_old,
const Vector< Real > &  u_new,
const Vector< Real > &  z,
Real &  tol 
)
pure virtual
template<class Real >
virtual void ROL::Constraint_TimeSimOpt< Real >::update ( const Vector< Real > &  u,
const Vector< Real > &  z,
bool  flag = true,
int  iter = -1 
)
inlineoverridevirtual

Update constraint functions. x is the optimization variable, flag = true if optimization variable is changed, iter is the outer algorithm iterations count.

Reimplemented from ROL::ROL::Constraint_SimOpt< Real >.

Definition at line 273 of file ROL_Constraint_TimeSimOpt.hpp.

References ROL::Constraint_TimeSimOpt< Real >::getNewVector(), ROL::Constraint_TimeSimOpt< Real >::getOldVector(), and ROL::Constraint_TimeSimOpt< Real >::update().

template<class Real >
virtual void ROL::Constraint_TimeSimOpt< Real >::value ( Vector< Real > &  c,
const Vector< Real > &  u,
const Vector< Real > &  z,
Real &  tol 
)
inlineoverridevirtual

Evaluate the constraint operator \(c:\mathcal{U}\times\mathcal{Z} \rightarrow \mathcal{C}\) at \((u,z)\).

Parameters
[out]cis the result of evaluating the constraint operator at \((u,z)\); a constraint-space vector
[in]uis the constraint argument; a simulation-space vector
[in]zis the constraint argument; an optimization-space vector
[in,out]tolis a tolerance for inexact evaluations; currently unused

On return, \(\mathsf{c} = c(u,z)\), where \(\mathsf{c} \in \mathcal{C}\), \(\mathsf{u} \in \mathcal{U}\), and $ \(\mathsf{z} \in\mathcal{Z}\).


Implements ROL::ROL::Constraint_SimOpt< Real >.

Definition at line 282 of file ROL_Constraint_TimeSimOpt.hpp.

References ROL::Constraint_TimeSimOpt< Real >::getNewVector(), ROL::Constraint_TimeSimOpt< Real >::getOldVector(), and ROL::Constraint_TimeSimOpt< Real >::value().

template<class Real >
virtual void ROL::Constraint_TimeSimOpt< Real >::solve ( Vector< Real > &  c,
Vector< Real > &  u,
const Vector< Real > &  z,
Real &  tol 
)
inlineoverridevirtual

Given \(z\), solve \(c(u,z)=0\) for \(u\).

Parameters
[out]cis the result of evaluating the constraint operator at \((u,z)\); a constraint-space vector
[in,out]uis the solution vector; a simulation-space vector
[in]zis the constraint argument; an optimization-space vector
[in,out]tolis a tolerance for inexact evaluations; currently unused

The defualt implementation is Newton's method globalized with a backtracking line search.


Reimplemented from ROL::ROL::Constraint_SimOpt< Real >.

Definition at line 294 of file ROL_Constraint_TimeSimOpt.hpp.

References ROL::Constraint_TimeSimOpt< Real >::getNewVector(), ROL::Constraint_TimeSimOpt< Real >::getOldVector(), and ROL::Constraint_TimeSimOpt< Real >::solve().

template<class Real >
virtual void ROL::Constraint_TimeSimOpt< Real >::applyJacobian_1 ( Vector< Real > &  jv,
const Vector< Real > &  v,
const Vector< Real > &  u,
const Vector< Real > &  z,
Real &  tol 
)
inlineoverridevirtual

Apply the partial constraint Jacobian at \((u,z)\), \(c_u(u,z) \in L(\mathcal{U}, \mathcal{C})\), to the vector \(v\).

  @param[out]      jv  is the result of applying the constraint Jacobian to @b v at @b  \form#196; a constraint-space vector
  @param[in]       v   is a simulation-space vector
  @param[in]       u   is the constraint argument; an simulation-space vector
  @param[in]       z   is the constraint argument; an optimization-space vector
  @param[in,out]   tol is a tolerance for inexact evaluations; currently unused

  On return, \form#201, where

\(v \in \mathcal{U}\), \(\mathsf{jv} \in \mathcal{C}\).


Reimplemented from ROL::ROL::Constraint_SimOpt< Real >.

Definition at line 305 of file ROL_Constraint_TimeSimOpt.hpp.

References ROL::Constraint_TimeSimOpt< Real >::applyJacobian_1_new(), ROL::Constraint_TimeSimOpt< Real >::applyJacobian_1_old(), ROL::Vector< Real >::axpy(), ROL::Constraint_TimeSimOpt< Real >::getNewVector(), ROL::Constraint_TimeSimOpt< Real >::getOldVector(), and ROL::Constraint_TimeSimOpt< Real >::workspace_.

template<class Real >
virtual void ROL::Constraint_TimeSimOpt< Real >::applyJacobian_2 ( Vector< Real > &  jv,
const Vector< Real > &  v,
const Vector< Real > &  u,
const Vector< Real > &  z,
Real &  tol 
)
inlineoverridevirtual

Apply the partial constraint Jacobian at \((u,z)\), \(c_z(u,z) \in L(\mathcal{Z}, \mathcal{C})\), to the vector \(v\).

  @param[out]      jv  is the result of applying the constraint Jacobian to @b v at @b  \form#196; a constraint-space vector
  @param[in]       v   is an optimization-space vector
  @param[in]       u   is the constraint argument; a simulation-space vector
  @param[in]       z   is the constraint argument; an optimization-space vector
  @param[in,out]   tol is a tolerance for inexact evaluations; currently unused

  On return, \form#204, where

\(v \in \mathcal{Z}\), \(\mathsf{jv} \in \mathcal{C}\).


Reimplemented from ROL::ROL::Constraint_SimOpt< Real >.

Definition at line 332 of file ROL_Constraint_TimeSimOpt.hpp.

References ROL::Constraint_TimeSimOpt< Real >::applyJacobian_2(), ROL::Constraint_TimeSimOpt< Real >::getNewVector(), and ROL::Constraint_TimeSimOpt< Real >::getOldVector().

template<class Real >
virtual void ROL::Constraint_TimeSimOpt< Real >::applyInverseJacobian_1 ( Vector< Real > &  ijv,
const Vector< Real > &  v,
const Vector< Real > &  u,
const Vector< Real > &  z,
Real &  tol 
)
inlinefinaloverridevirtual

Apply the inverse partial constraint Jacobian at \((u,z)\), \(c_u(u,z)^{-1} \in L(\mathcal{C}, \mathcal{U})\), to the vector \(v\).

  @param[out]      ijv is the result of applying the inverse constraint Jacobian to @b v at @b  \form#196; a simulation-space vector
  @param[in]       v   is a constraint-space vector
  @param[in]       u   is the constraint argument; a simulation-space vector
  @param[in]       z   is the constraint argument; an optimization-space vector
  @param[in,out]   tol is a tolerance for inexact evaluations; currently unused

  On return, \form#207, where

\(v \in \mathcal{C}\), \(\mathsf{ijv} \in \mathcal{U}\).


Reimplemented from ROL::ROL::Constraint_SimOpt< Real >.

Definition at line 347 of file ROL_Constraint_TimeSimOpt.hpp.

template<class Real >
virtual void ROL::Constraint_TimeSimOpt< Real >::applyAdjointJacobian_1 ( Vector< Real > &  ajv,
const Vector< Real > &  v,
const Vector< Real > &  u,
const Vector< Real > &  z,
Real &  tol 
)
inlineoverridevirtual

Apply the adjoint of the partial constraint Jacobian at \((u,z)\), \(c_u(u,z)^* \in L(\mathcal{C}^*, \mathcal{U}^*)\), to the vector \(v\). This is the primary interface.

  @param[out]      ajv    is the result of applying the adjoint of the constraint Jacobian to @b v at @b (u,z); a dual simulation-space vector
  @param[in]       v      is a dual constraint-space vector
  @param[in]       u      is the constraint argument; a simulation-space vector
  @param[in]       z      is the constraint argument; an optimization-space vector
  @param[in,out]   tol    is a tolerance for inexact evaluations; currently unused

  On return, \form#210, where

\(v \in \mathcal{C}^*\), \(\mathsf{ajv} \in \mathcal{U}^*\).


Reimplemented from ROL::ROL::Constraint_SimOpt< Real >.

Definition at line 356 of file ROL_Constraint_TimeSimOpt.hpp.

References ROL::Constraint_TimeSimOpt< Real >::applyAdjointJacobian_1_new(), ROL::Constraint_TimeSimOpt< Real >::applyAdjointJacobian_1_old(), ROL::Constraint_TimeSimOpt< Real >::getNewVector(), and ROL::Constraint_TimeSimOpt< Real >::getOldVector().

template<class Real >
virtual void ROL::Constraint_TimeSimOpt< Real >::applyAdjointJacobian_2 ( Vector< Real > &  ajv,
const Vector< Real > &  v,
const Vector< Real > &  u,
const Vector< Real > &  z,
Real &  tol 
)
inlineoverridevirtual

Apply the adjoint of the partial constraint Jacobian at \((u,z)\), \(c_z(u,z)^* \in L(\mathcal{C}^*, \mathcal{Z}^*)\), to vector \(v\). This is the primary interface.

  @param[out]      ajv    is the result of applying the adjoint of the constraint Jacobian to @b v at @b  \form#196; a dual optimization-space vector
  @param[in]       v      is a dual constraint-space vector
  @param[in]       u      is the constraint argument; a simulation-space vector
  @param[in]       z      is the constraint argument; an optimization-space vector
  @param[in,out]   tol    is a tolerance for inexact evaluations; currently unused

  On return, \form#213, where

\(v \in \mathcal{C}^*\), \(\mathsf{ajv} \in \mathcal{Z}^*\).


Reimplemented from ROL::ROL::Constraint_SimOpt< Real >.

Definition at line 370 of file ROL_Constraint_TimeSimOpt.hpp.

References ROL::Constraint_TimeSimOpt< Real >::applyAdjointJacobian_2_time(), ROL::Constraint_TimeSimOpt< Real >::getNewVector(), and ROL::Constraint_TimeSimOpt< Real >::getOldVector().

template<class Real >
virtual void ROL::Constraint_TimeSimOpt< Real >::applyInverseAdjointJacobian_1 ( Vector< Real > &  iajv,
const Vector< Real > &  v,
const Vector< Real > &  u,
const Vector< Real > &  z,
Real &  tol 
)
inlinefinaloverridevirtual

Apply the inverse of the adjoint of the partial constraint Jacobian at \((u,z)\), \(c_u(u,z)^{-*} \in L(\mathcal{U}^*, \mathcal{C}^*)\), to the vector \(v\).

  @param[out]      iajv is the result of applying the inverse adjoint of the constraint Jacobian to @b v at @b (u,z); a dual constraint-space vector
  @param[in]       v   is a dual simulation-space vector
  @param[in]       u   is the constraint argument; a simulation-space vector
  @param[in]       z   is the constraint argument; an optimization-space vector
  @param[in,out]   tol is a tolerance for inexact evaluations; currently unused

  On return, \form#216, where

\(v \in \mathcal{U}^*\), \(\mathsf{iajv} \in \mathcal{C}^*\).


Reimplemented from ROL::ROL::Constraint_SimOpt< Real >.

Definition at line 381 of file ROL_Constraint_TimeSimOpt.hpp.

template<class Real >
virtual void ROL::Constraint_TimeSimOpt< Real >::applyAdjointHessian_11 ( Vector< Real > &  ahwv,
const Vector< Real > &  w,
const Vector< Real > &  v,
const Vector< Real > &  u,
const Vector< Real > &  z,
Real &  tol 
)
inlineoverridevirtual

Apply the simulation-space derivative of the adjoint of the constraint simulation-space Jacobian at \((u,z)\) to the vector \(w\) in the direction \(v\), according to \(v\mapsto c_{uu}(u,z)(v,\cdot)^*w\).

  @param[out]      ahwv is the result of applying the simulation-space derivative of the adjoint of the constraint simulation-space Jacobian at @b  \form#196 to the vector @b  \form#219 in direction @b  \form#219; a dual simulation-space vector
  @param[in]       w    is the direction vector; a dual constraint-space vector
  @param[in]       v    is a simulation-space vector
  @param[in]       u    is the constraint argument; a simulation-space vector
  @param[in]       z    is the constraint argument; an optimization-space vector
  @param[in,out]   tol  is a tolerance for inexact evaluations; currently unused

  On return, \form#221, where

\(w \in \mathcal{C}^*\), \(v \in \mathcal{U}\), and \(\mathsf{ahwv} \in \mathcal{U}^*\).


Reimplemented from ROL::ROL::Constraint_SimOpt< Real >.

Definition at line 407 of file ROL_Constraint_TimeSimOpt.hpp.

References ROL::Constraint_TimeSimOpt< Real >::applyAdjointHessian_11_new(), ROL::Constraint_TimeSimOpt< Real >::applyAdjointHessian_11_old(), ROL::Constraint_TimeSimOpt< Real >::getNewVector(), and ROL::Constraint_TimeSimOpt< Real >::getOldVector().

template<class Real >
virtual void ROL::Constraint_TimeSimOpt< Real >::applyAdjointHessian_12 ( Vector< Real > &  ahwv,
const Vector< Real > &  w,
const Vector< Real > &  v,
const Vector< Real > &  u,
const Vector< Real > &  z,
Real &  tol 
)
inlineoverridevirtual

Apply the optimization-space derivative of the adjoint of the constraint simulation-space Jacobian at \((u,z)\) to the vector \(w\) in the direction \(v\), according to \(v\mapsto c_{uz}(u,z)(v,\cdot)^*w\).

  @param[out]      ahwv is the result of applying the optimization-space derivative of the adjoint of the constraint simulation-space Jacobian at @b  \form#196 to the vector @b  \form#219 in direction @b  \form#219; a dual optimization-space vector
  @param[in]       w    is the direction vector; a dual constraint-space vector
  @param[in]       v    is a simulation-space vector
  @param[in]       u    is the constraint argument; a simulation-space vector
  @param[in]       z    is the constraint argument; an optimization-space vector
  @param[in,out]   tol  is a tolerance for inexact evaluations; currently unused

  On return, \form#225, where

\(w \in \mathcal{C}^*\), \(v \in \mathcal{U}\), and \(\mathsf{ahwv} \in \mathcal{Z}^*\).


Reimplemented from ROL::ROL::Constraint_SimOpt< Real >.

Definition at line 445 of file ROL_Constraint_TimeSimOpt.hpp.

References ROL::Vector< Real >::zero().

template<class Real >
virtual void ROL::Constraint_TimeSimOpt< Real >::applyAdjointHessian_21 ( Vector< Real > &  ahwv,
const Vector< Real > &  w,
const Vector< Real > &  v,
const Vector< Real > &  u,
const Vector< Real > &  z,
Real &  tol 
)
inlineoverridevirtual

Apply the simulation-space derivative of the adjoint of the constraint optimization-space Jacobian at \((u,z)\) to the vector \(w\) in the direction \(v\), according to \(v\mapsto c_{zu}(u,z)(v,\cdot)^*w\).

  @param[out]      ahwv is the result of applying the simulation-space derivative of the adjoint of the constraint optimization-space Jacobian at @b  \form#196 to the vector @b  \form#219 in direction @b  \form#219; a dual simulation-space vector 
  @param[in]       w    is the direction vector; a dual constraint-space vector 
  @param[in]       v    is a optimization-space vector 
  @param[in]       u    is the constraint argument; a simulation-space vector 
  @param[in]       z    is the constraint argument; an optimization-space vector 
  @param[in,out]   tol  is a tolerance for inexact evaluations; currently unused 

  On return, \form#228, where 

\(w \in \mathcal{C}^*\), \(v \in \mathcal{Z}\), and \(\mathsf{ahwv} \in \mathcal{U}^*\).


Reimplemented from ROL::ROL::Constraint_SimOpt< Real >.

Definition at line 472 of file ROL_Constraint_TimeSimOpt.hpp.

References ROL::Vector< Real >::zero().

template<class Real >
virtual void ROL::Constraint_TimeSimOpt< Real >::applyAdjointHessian_22 ( Vector< Real > &  ahwv,
const Vector< Real > &  w,
const Vector< Real > &  v,
const Vector< Real > &  u,
const Vector< Real > &  z,
Real &  tol 
)
inlineoverridevirtual

Apply the optimization-space derivative of the adjoint of the constraint optimization-space Jacobian at \((u,z)\) to the vector \(w\) in the direction \(v\), according to \(v\mapsto c_{zz}(u,z)(v,\cdot)^*w\).

  @param[out]      ahwv is the result of applying the optimization-space derivative of the adjoint of the constraint optimization-space Jacobian at @b  \form#196 to the vector @b  \form#219 in direction @b  \form#219; a dual optimization-space vector
  @param[in]       w    is the direction vector; a dual constraint-space vector
  @param[in]       v    is a optimization-space vector
  @param[in]       u    is the constraint argument; a simulation-space vector
  @param[in]       z    is the constraint argument; an optimization-space vector
  @param[in,out]   tol  is a tolerance for inexact evaluations; currently unused

  On return, \form#230, where

\(w \in \mathcal{C}^*\), \(v \in \mathcal{Z}\), and \(\mathsf{ahwv} \in \mathcal{Z}^*\).


Reimplemented from ROL::ROL::Constraint_SimOpt< Real >.

Definition at line 498 of file ROL_Constraint_TimeSimOpt.hpp.

References ROL::Vector< Real >::zero().

template<class Real >
virtual Real ROL::Constraint_TimeSimOpt< Real >::checkSolve ( const ROL::Vector< Real > &  u,
const ROL::Vector< Real > &  z,
const ROL::Vector< Real > &  c,
const bool  printToStream = true,
std::ostream &  outStream = std::cout 
)
inlineoverridevirtual
template<class Real >
virtual Real ROL::Constraint_TimeSimOpt< Real >::checkInverseJacobian_1_new ( const ROL::Vector< Real > &  c,
const ROL::Vector< Real > &  u_new,
const ROL::Vector< Real > &  u_old,
const ROL::Vector< Real > &  z,
const ROL::Vector< Real > &  v_new,
const bool  printToStream = true,
std::ostream &  outStream = std::cout 
)
inlinevirtual
template<class Real >
virtual Real ROL::Constraint_TimeSimOpt< Real >::checkInverseAdjointJacobian_1_new ( const ROL::Vector< Real > &  c,
const ROL::Vector< Real > &  u_new,
const ROL::Vector< Real > &  u_old,
const ROL::Vector< Real > &  z,
const ROL::Vector< Real > &  v_new,
const bool  printToStream = true,
std::ostream &  outStream = std::cout 
)
inlinevirtual
template<class Real >
std::vector<std::vector<Real> > ROL::Constraint_TimeSimOpt< Real >::checkApplyJacobian_1_new ( const Vector< Real > &  u_new,
const Vector< Real > &  u_old,
const Vector< Real > &  z,
const Vector< Real > &  v,
const Vector< Real > &  jv,
const bool  printToStream = true,
std::ostream &  outStream = std::cout,
const int  numSteps = ROL_NUM_CHECKDERIV_STEPS,
const int  order = 1 
)
inline

Definition at line 600 of file ROL_Constraint_TimeSimOpt.hpp.

template<class Real >
std::vector<std::vector<Real> > ROL::Constraint_TimeSimOpt< Real >::checkApplyJacobian_1_new ( const Vector< Real > &  u_new,
const Vector< Real > &  u_old,
const Vector< Real > &  z,
const Vector< Real > &  v,
const Vector< Real > &  jv,
const std::vector< Real > &  steps,
const bool  printToStream = true,
std::ostream &  outStream = std::cout,
const int  order = 1 
)
inline

Member Data Documentation

template<class Real >
VectorWorkspace<Real> ROL::Constraint_TimeSimOpt< Real >::workspace_
mutableprivate

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