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

Defines a composite equality constraint operator interface for simulation-based optimization. More...

#include <ROL_CompositeConstraint_SimOpt.hpp>

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

Public Member Functions

 CompositeConstraint_SimOpt (const ROL::Ptr< Constraint_SimOpt< Real >> &conVal, const ROL::Ptr< Constraint_SimOpt< Real >> &conRed, const Vector< Real > &cVal, const Vector< Real > &cRed, const Vector< Real > &u, const Vector< Real > &Sz, const Vector< Real > &z, bool storage=true, bool isConRedParametrized=false)
 
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...
 
void update_1 (const Vector< Real > &u, bool flag=true, int iter=-1) override
 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...
 
void update_2 (const Vector< Real > &z, bool flag=true, int iter=-1) override
 Update constraint functions with respect to Opt variable. x is the optimization variable, flag = true if optimization variable is changed, iter is the outer algorithm iterations count. More...
 
void update (const Vector< Real > &u, const Vector< Real > &z, UpdateType type, int iter=-1) override
 
void update_1 (const Vector< Real > &u, UpdateType type, int iter=-1) override
 
void update_2 (const Vector< Real > &z, UpdateType type, int iter=-1) override
 
void solve_update (const Vector< Real > &u, const Vector< Real > &z, UpdateType type, int iter=-1) override
 Update SimOpt constraint during solve (disconnected from optimization updates). More...
 
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...
 
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...
 
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...
 
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...
 
void applyInverseJacobian_1 (Vector< Real > &ijv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
 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...
 
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...
 
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...
 
void applyInverseAdjointJacobian_1 (Vector< Real > &ijv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
 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...
 
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...
 
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...
 
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...
 
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...
 
void setParameter (const std::vector< Real > &param) override
 
- Public Member Functions inherited from ROL::ROL::Constraint_SimOpt< Real >
 Constraint_SimOpt ()
 
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 update (const Vector< Real > &x, UpdateType type, int iter=-1)
 Update constraint function. 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 checkSolve (const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &c, const bool printToStream=true, std::ostream &outStream=std::cout)
 
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...
 

Private Member Functions

void solveConRed (Vector< Real > &Sz, const Vector< Real > &z, Real &tol)
 
void applySens (Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &Sz, const Vector< Real > &z, Real &tol)
 
void applyAdjointSens (Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &Sz, const Vector< Real > &z, Real &tol)
 

Private Attributes

const ROL::Ptr
< Constraint_SimOpt< Real > > 
conVal_
 
const ROL::Ptr
< Constraint_SimOpt< Real > > 
conRed_
 
ROL::Ptr< Vector< Real > > Sz_
 
ROL::Ptr< Vector< Real > > primRed_
 
ROL::Ptr< Vector< Real > > dualRed_
 
ROL::Ptr< Vector< Real > > primZ_
 
ROL::Ptr< Vector< Real > > dualZ_
 
ROL::Ptr< Vector< Real > > dualZ1_
 
ROL::Ptr< Vector< Real > > primU_
 
ROL::Ptr< VectorController
< Real > > 
stateStore_
 
bool updateFlag_
 
bool newUpdate_
 
int updateIter_
 
UpdateType updateType_
 
const bool storage_
 
const bool isConRedParametrized_
 

Additional Inherited Members

- Protected Member Functions inherited from ROL::Constraint< Real >
const std::vector< Real > getParameter (void) const
 
- 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<typename Real>
class ROL::CompositeConstraint_SimOpt< Real >

Defines a composite equality constraint operator interface for simulation-based optimization.

This equality constraint interface inherits from ROL_Constraint_SimOpt, for the use case when \(\mathcal{X}=\mathcal{U}\times\mathcal{Z}\) where \(\mathcal{U}\) and \(\mathcal{Z}\) are Banach spaces. \(\mathcal{U}\) denotes the "simulation space" and \(\mathcal{Z}\) denotes the "optimization space" (of designs, controls, parameters). The simulation-based constraints are of the form

\[ c(u,S(z)) = 0 \]

where \(S(z)\) solves the reducible constraint

\[ c_0(S(z),z) = 0. \]


Definition at line 74 of file ROL_CompositeConstraint_SimOpt.hpp.

Constructor & Destructor Documentation

template<typename Real >
ROL::CompositeConstraint_SimOpt< Real >::CompositeConstraint_SimOpt ( const ROL::Ptr< Constraint_SimOpt< Real >> &  conVal,
const ROL::Ptr< Constraint_SimOpt< Real >> &  conRed,
const Vector< Real > &  cVal,
const Vector< Real > &  cRed,
const Vector< Real > &  u,
const Vector< Real > &  Sz,
const Vector< Real > &  z,
bool  storage = true,
bool  isConRedParametrized = false 
)

Member Function Documentation

template<typename Real >
void ROL::CompositeConstraint_SimOpt< Real >::update ( const Vector< Real > &  u,
const Vector< Real > &  z,
bool  flag = true,
int  iter = -1 
)
overridevirtual

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 74 of file ROL_CompositeConstraint_SimOpt_Def.hpp.

References ROL::update_2().

template<typename Real >
void ROL::CompositeConstraint_SimOpt< Real >::update_1 ( const Vector< Real > &  u,
bool  flag = true,
int  iter = -1 
)
overridevirtual

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.

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

Definition at line 82 of file ROL_CompositeConstraint_SimOpt_Def.hpp.

template<typename Real >
void ROL::CompositeConstraint_SimOpt< Real >::update_2 ( const Vector< Real > &  z,
bool  flag = true,
int  iter = -1 
)
overridevirtual

Update constraint functions with respect to Opt variable. 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 89 of file ROL_CompositeConstraint_SimOpt_Def.hpp.

template<typename Real >
void ROL::CompositeConstraint_SimOpt< Real >::update ( const Vector< Real > &  u,
const Vector< Real > &  z,
UpdateType  type,
int  iter = -1 
)
overridevirtual

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

Definition at line 98 of file ROL_CompositeConstraint_SimOpt_Def.hpp.

References ROL::update_2().

template<typename Real >
void ROL::CompositeConstraint_SimOpt< Real >::update_1 ( const Vector< Real > &  u,
UpdateType  type,
int  iter = -1 
)
overridevirtual

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

Definition at line 106 of file ROL_CompositeConstraint_SimOpt_Def.hpp.

template<typename Real >
void ROL::CompositeConstraint_SimOpt< Real >::update_2 ( const Vector< Real > &  z,
UpdateType  type,
int  iter = -1 
)
overridevirtual

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

Definition at line 113 of file ROL_CompositeConstraint_SimOpt_Def.hpp.

template<typename Real >
void ROL::CompositeConstraint_SimOpt< Real >::solve_update ( const Vector< Real > &  u,
const Vector< Real > &  z,
UpdateType  type,
int  iter = -1 
)
overridevirtual

Update SimOpt constraint during solve (disconnected from optimization updates).

Parameters
[in]xis the optimization variable
[in]typeis the update type
[in]iteris the solver iteration count

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

Definition at line 122 of file ROL_CompositeConstraint_SimOpt_Def.hpp.

template<typename Real >
void ROL::CompositeConstraint_SimOpt< Real >::value ( Vector< Real > &  c,
const Vector< Real > &  u,
const Vector< Real > &  z,
Real &  tol 
)
overridevirtual

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 131 of file ROL_CompositeConstraint_SimOpt_Def.hpp.

template<typename Real >
void ROL::CompositeConstraint_SimOpt< Real >::solve ( Vector< Real > &  c,
Vector< Real > &  u,
const Vector< Real > &  z,
Real &  tol 
)
overridevirtual

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 140 of file ROL_CompositeConstraint_SimOpt_Def.hpp.

template<typename Real >
void ROL::CompositeConstraint_SimOpt< Real >::applyJacobian_1 ( Vector< Real > &  jv,
const Vector< Real > &  v,
const Vector< Real > &  u,
const Vector< Real > &  z,
Real &  tol 
)
overridevirtual

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#218; 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#224, where

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


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

Definition at line 149 of file ROL_CompositeConstraint_SimOpt_Def.hpp.

template<typename Real >
void ROL::CompositeConstraint_SimOpt< Real >::applyJacobian_2 ( Vector< Real > &  jv,
const Vector< Real > &  v,
const Vector< Real > &  u,
const Vector< Real > &  z,
Real &  tol 
)
overridevirtual

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#218; 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#227, where

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


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

Definition at line 159 of file ROL_CompositeConstraint_SimOpt_Def.hpp.

template<typename Real >
void ROL::CompositeConstraint_SimOpt< Real >::applyInverseJacobian_1 ( Vector< Real > &  ijv,
const Vector< Real > &  v,
const Vector< Real > &  u,
const Vector< Real > &  z,
Real &  tol 
)
overridevirtual

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#218; 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#230, where

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


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

Definition at line 170 of file ROL_CompositeConstraint_SimOpt_Def.hpp.

template<typename Real >
void ROL::CompositeConstraint_SimOpt< Real >::applyAdjointJacobian_1 ( Vector< Real > &  ajv,
const Vector< Real > &  v,
const Vector< Real > &  u,
const Vector< Real > &  z,
Real &  tol 
)
overridevirtual

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#233, where

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


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

Definition at line 180 of file ROL_CompositeConstraint_SimOpt_Def.hpp.

template<typename Real >
void ROL::CompositeConstraint_SimOpt< Real >::applyAdjointJacobian_2 ( Vector< Real > &  ajv,
const Vector< Real > &  v,
const Vector< Real > &  u,
const Vector< Real > &  z,
Real &  tol 
)
overridevirtual

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#218; 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#236, where

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


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

Definition at line 190 of file ROL_CompositeConstraint_SimOpt_Def.hpp.

template<typename Real >
void ROL::CompositeConstraint_SimOpt< Real >::applyInverseAdjointJacobian_1 ( Vector< Real > &  iajv,
const Vector< Real > &  v,
const Vector< Real > &  u,
const Vector< Real > &  z,
Real &  tol 
)
overridevirtual

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#239, where

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


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

Definition at line 201 of file ROL_CompositeConstraint_SimOpt_Def.hpp.

template<typename Real >
void ROL::CompositeConstraint_SimOpt< Real >::applyAdjointHessian_11 ( Vector< Real > &  ahwv,
const Vector< Real > &  w,
const Vector< Real > &  v,
const Vector< Real > &  u,
const Vector< Real > &  z,
Real &  tol 
)
overridevirtual

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#218 to the vector @b  \form#242 in direction @b  \form#242; 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#244, 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 211 of file ROL_CompositeConstraint_SimOpt_Def.hpp.

template<typename Real >
void ROL::CompositeConstraint_SimOpt< Real >::applyAdjointHessian_12 ( Vector< Real > &  ahwv,
const Vector< Real > &  w,
const Vector< Real > &  v,
const Vector< Real > &  u,
const Vector< Real > &  z,
Real &  tol 
)
overridevirtual

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#218 to the vector @b  \form#242 in direction @b  \form#242; 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#248, 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 222 of file ROL_CompositeConstraint_SimOpt_Def.hpp.

template<typename Real >
void ROL::CompositeConstraint_SimOpt< Real >::applyAdjointHessian_21 ( Vector< Real > &  ahwv,
const Vector< Real > &  w,
const Vector< Real > &  v,
const Vector< Real > &  u,
const Vector< Real > &  z,
Real &  tol 
)
overridevirtual

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#218 to the vector @b  \form#242 in direction @b  \form#242; 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#251, 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 234 of file ROL_CompositeConstraint_SimOpt_Def.hpp.

template<typename Real >
void ROL::CompositeConstraint_SimOpt< Real >::applyAdjointHessian_22 ( Vector< Real > &  ahwv,
const Vector< Real > &  w,
const Vector< Real > &  v,
const Vector< Real > &  u,
const Vector< Real > &  z,
Real &  tol 
)
overridevirtual

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#218 to the vector @b  \form#242 in direction @b  \form#242; 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#253, 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 246 of file ROL_CompositeConstraint_SimOpt_Def.hpp.

References ROL::Vector< Real >::axpy(), ROL::Vector< Real >::plus(), and ROL::Vector< Real >::zero().

template<typename Real >
void ROL::CompositeConstraint_SimOpt< Real >::setParameter ( const std::vector< Real > &  param)
overridevirtual
template<typename Real >
void ROL::CompositeConstraint_SimOpt< Real >::solveConRed ( Vector< Real > &  Sz,
const Vector< Real > &  z,
Real &  tol 
)
private
template<typename Real >
void ROL::CompositeConstraint_SimOpt< Real >::applySens ( Vector< Real > &  jv,
const Vector< Real > &  v,
const Vector< Real > &  Sz,
const Vector< Real > &  z,
Real &  tol 
)
private
template<typename Real >
void ROL::CompositeConstraint_SimOpt< Real >::applyAdjointSens ( Vector< Real > &  ajv,
const Vector< Real > &  v,
const Vector< Real > &  Sz,
const Vector< Real > &  z,
Real &  tol 
)
private

Member Data Documentation

template<typename Real>
const ROL::Ptr<Constraint_SimOpt<Real> > ROL::CompositeConstraint_SimOpt< Real >::conVal_
private

Definition at line 77 of file ROL_CompositeConstraint_SimOpt.hpp.

template<typename Real>
const ROL::Ptr<Constraint_SimOpt<Real> > ROL::CompositeConstraint_SimOpt< Real >::conRed_
private

Definition at line 77 of file ROL_CompositeConstraint_SimOpt.hpp.

template<typename Real>
ROL::Ptr<Vector<Real> > ROL::CompositeConstraint_SimOpt< Real >::Sz_
private
template<typename Real>
ROL::Ptr<Vector<Real> > ROL::CompositeConstraint_SimOpt< Real >::primRed_
private
template<typename Real>
ROL::Ptr<Vector<Real> > ROL::CompositeConstraint_SimOpt< Real >::dualRed_
private
template<typename Real>
ROL::Ptr<Vector<Real> > ROL::CompositeConstraint_SimOpt< Real >::primZ_
private
template<typename Real>
ROL::Ptr<Vector<Real> > ROL::CompositeConstraint_SimOpt< Real >::dualZ_
private
template<typename Real>
ROL::Ptr<Vector<Real> > ROL::CompositeConstraint_SimOpt< Real >::dualZ1_
private
template<typename Real>
ROL::Ptr<Vector<Real> > ROL::CompositeConstraint_SimOpt< Real >::primU_
private
template<typename Real>
ROL::Ptr<VectorController<Real> > ROL::CompositeConstraint_SimOpt< Real >::stateStore_
private
template<typename Real>
bool ROL::CompositeConstraint_SimOpt< Real >::updateFlag_
private

Definition at line 83 of file ROL_CompositeConstraint_SimOpt.hpp.

template<typename Real>
bool ROL::CompositeConstraint_SimOpt< Real >::newUpdate_
private

Definition at line 83 of file ROL_CompositeConstraint_SimOpt.hpp.

template<typename Real>
int ROL::CompositeConstraint_SimOpt< Real >::updateIter_
private

Definition at line 84 of file ROL_CompositeConstraint_SimOpt.hpp.

template<typename Real>
UpdateType ROL::CompositeConstraint_SimOpt< Real >::updateType_
private

Definition at line 85 of file ROL_CompositeConstraint_SimOpt.hpp.

template<typename Real>
const bool ROL::CompositeConstraint_SimOpt< Real >::storage_
private

Definition at line 87 of file ROL_CompositeConstraint_SimOpt.hpp.

template<typename Real>
const bool ROL::CompositeConstraint_SimOpt< Real >::isConRedParametrized_
private

Definition at line 87 of file ROL_CompositeConstraint_SimOpt.hpp.


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