ROL
|
Defines a composite equality constraint operator interface for simulation-based optimization. More...
#include <ROL_CompositeConstraint_SimOpt.hpp>
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, const bool storage=true, const bool isConRedParametrized=false) | |
void | update (const Vector< Real > &u, const Vector< Real > &z, 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... | |
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... | |
void | update_2 (const Vector< Real > &z, bool flag=true, int iter=-1) |
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 | value (Vector< Real > &c, const Vector< Real > &u, const Vector< Real > &z, Real &tol) |
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) |
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) |
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) |
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) |
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) |
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) |
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) |
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) |
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) |
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) |
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) |
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 > ¶m) |
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 | 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 (const Vector< Real > &z, Real &tol) |
void | applySens (Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &z, Real &tol) |
void | applyAdjointSens (Vector< Real > &ajv, const Vector< Real > &v, 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< SimController< Real > > | stateStore_ |
bool | updateFlag_ |
int | updateIter_ |
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_ |
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 75 of file ROL_CompositeConstraint_SimOpt.hpp.
|
inline |
Definition at line 138 of file ROL_CompositeConstraint_SimOpt.hpp.
References ROL::Vector< Real >::clone(), ROL::Vector< Real >::dual(), ROL::CompositeConstraint_SimOpt< Real >::dualRed_, ROL::CompositeConstraint_SimOpt< Real >::dualZ1_, ROL::CompositeConstraint_SimOpt< Real >::dualZ_, ROL::CompositeConstraint_SimOpt< Real >::primRed_, ROL::CompositeConstraint_SimOpt< Real >::primZ_, ROL::CompositeConstraint_SimOpt< Real >::stateStore_, and ROL::CompositeConstraint_SimOpt< Real >::Sz_.
|
inlineprivate |
Definition at line 95 of file ROL_CompositeConstraint_SimOpt.hpp.
References ROL::CompositeConstraint_SimOpt< Real >::conRed_, ROL::Constraint< Real >::getParameter(), ROL::CompositeConstraint_SimOpt< Real >::primRed_, ROL::CompositeConstraint_SimOpt< Real >::stateStore_, ROL::CompositeConstraint_SimOpt< Real >::storage_, ROL::CompositeConstraint_SimOpt< Real >::Sz_, ROL::CompositeConstraint_SimOpt< Real >::updateFlag_, and ROL::CompositeConstraint_SimOpt< Real >::updateIter_.
Referenced by ROL::CompositeConstraint_SimOpt< Real >::applyAdjointHessian_11(), ROL::CompositeConstraint_SimOpt< Real >::applyAdjointHessian_12(), ROL::CompositeConstraint_SimOpt< Real >::applyAdjointJacobian_1(), ROL::CompositeConstraint_SimOpt< Real >::applyAdjointJacobian_2(), ROL::CompositeConstraint_SimOpt< Real >::applyAdjointSens(), ROL::CompositeConstraint_SimOpt< Real >::applyInverseAdjointJacobian_1(), ROL::CompositeConstraint_SimOpt< Real >::applyInverseJacobian_1(), ROL::CompositeConstraint_SimOpt< Real >::applyJacobian_1(), ROL::CompositeConstraint_SimOpt< Real >::applySens(), ROL::CompositeConstraint_SimOpt< Real >::solve(), ROL::CompositeConstraint_SimOpt< Real >::update_2(), and ROL::CompositeConstraint_SimOpt< Real >::value().
|
inlineprivate |
Definition at line 119 of file ROL_CompositeConstraint_SimOpt.hpp.
References ROL::CompositeConstraint_SimOpt< Real >::conRed_, ROL::CompositeConstraint_SimOpt< Real >::primRed_, ROL::Vector< Real >::scale(), ROL::CompositeConstraint_SimOpt< Real >::solveConRed(), and ROL::CompositeConstraint_SimOpt< Real >::Sz_.
Referenced by ROL::CompositeConstraint_SimOpt< Real >::applyAdjointHessian_21(), ROL::CompositeConstraint_SimOpt< Real >::applyAdjointHessian_22(), and ROL::CompositeConstraint_SimOpt< Real >::applyJacobian_2().
|
inlineprivate |
Definition at line 128 of file ROL_CompositeConstraint_SimOpt.hpp.
References ROL::CompositeConstraint_SimOpt< Real >::conRed_, ROL::CompositeConstraint_SimOpt< Real >::dualRed_, ROL::Vector< Real >::scale(), ROL::CompositeConstraint_SimOpt< Real >::solveConRed(), and ROL::CompositeConstraint_SimOpt< Real >::Sz_.
Referenced by ROL::CompositeConstraint_SimOpt< Real >::applyAdjointHessian_12(), ROL::CompositeConstraint_SimOpt< Real >::applyAdjointHessian_22(), and ROL::CompositeConstraint_SimOpt< Real >::applyAdjointJacobian_2().
|
inlinevirtual |
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 155 of file ROL_CompositeConstraint_SimOpt.hpp.
References ROL::CompositeConstraint_SimOpt< Real >::update_1(), and ROL::CompositeConstraint_SimOpt< Real >::update_2().
|
inlinevirtual |
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 161 of file ROL_CompositeConstraint_SimOpt.hpp.
References ROL::CompositeConstraint_SimOpt< Real >::conVal_, and ROL::CompositeConstraint_SimOpt< Real >::Sz_.
Referenced by ROL::CompositeConstraint_SimOpt< Real >::update().
|
inlinevirtual |
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 167 of file ROL_CompositeConstraint_SimOpt.hpp.
References ROL::CompositeConstraint_SimOpt< Real >::solveConRed(), ROL::CompositeConstraint_SimOpt< Real >::stateStore_, ROL::CompositeConstraint_SimOpt< Real >::updateFlag_, and ROL::CompositeConstraint_SimOpt< Real >::updateIter_.
Referenced by ROL::CompositeConstraint_SimOpt< Real >::update().
|
inlinevirtual |
Evaluate the constraint operator \(c:\mathcal{U}\times\mathcal{Z} \rightarrow \mathcal{C}\) at \((u,z)\).
[out] | c | is the result of evaluating the constraint operator at \((u,z)\); a constraint-space vector |
[in] | u | is the constraint argument; a simulation-space vector |
[in] | z | is the constraint argument; an optimization-space vector |
[in,out] | tol | is 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 177 of file ROL_CompositeConstraint_SimOpt.hpp.
References ROL::CompositeConstraint_SimOpt< Real >::conVal_, ROL::CompositeConstraint_SimOpt< Real >::solveConRed(), and ROL::CompositeConstraint_SimOpt< Real >::Sz_.
|
inlinevirtual |
Given \(z\), solve \(c(u,z)=0\) for \(u\).
[out] | c | is the result of evaluating the constraint operator at \((u,z)\); a constraint-space vector |
[in,out] | u | is the solution vector; a simulation-space vector |
[in] | z | is the constraint argument; an optimization-space vector |
[in,out] | tol | is 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 182 of file ROL_CompositeConstraint_SimOpt.hpp.
References ROL::CompositeConstraint_SimOpt< Real >::conVal_, ROL::CompositeConstraint_SimOpt< Real >::solveConRed(), and ROL::CompositeConstraint_SimOpt< Real >::Sz_.
|
inlinevirtual |
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 187 of file ROL_CompositeConstraint_SimOpt.hpp.
References ROL::CompositeConstraint_SimOpt< Real >::conVal_, ROL::CompositeConstraint_SimOpt< Real >::solveConRed(), and ROL::CompositeConstraint_SimOpt< Real >::Sz_.
|
inlinevirtual |
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 193 of file ROL_CompositeConstraint_SimOpt.hpp.
References ROL::CompositeConstraint_SimOpt< Real >::applySens(), ROL::CompositeConstraint_SimOpt< Real >::conVal_, ROL::CompositeConstraint_SimOpt< Real >::primZ_, and ROL::CompositeConstraint_SimOpt< Real >::Sz_.
|
inlinevirtual |
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 199 of file ROL_CompositeConstraint_SimOpt.hpp.
References ROL::CompositeConstraint_SimOpt< Real >::conVal_, ROL::CompositeConstraint_SimOpt< Real >::solveConRed(), and ROL::CompositeConstraint_SimOpt< Real >::Sz_.
|
inlinevirtual |
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 205 of file ROL_CompositeConstraint_SimOpt.hpp.
References ROL::CompositeConstraint_SimOpt< Real >::conVal_, ROL::CompositeConstraint_SimOpt< Real >::solveConRed(), and ROL::CompositeConstraint_SimOpt< Real >::Sz_.
|
inlinevirtual |
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 211 of file ROL_CompositeConstraint_SimOpt.hpp.
References ROL::CompositeConstraint_SimOpt< Real >::applyAdjointSens(), ROL::CompositeConstraint_SimOpt< Real >::conVal_, ROL::CompositeConstraint_SimOpt< Real >::dualZ_, ROL::CompositeConstraint_SimOpt< Real >::solveConRed(), and ROL::CompositeConstraint_SimOpt< Real >::Sz_.
|
inlinevirtual |
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 218 of file ROL_CompositeConstraint_SimOpt.hpp.
References ROL::CompositeConstraint_SimOpt< Real >::conVal_, ROL::CompositeConstraint_SimOpt< Real >::solveConRed(), and ROL::CompositeConstraint_SimOpt< Real >::Sz_.
|
inlinevirtual |
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 224 of file ROL_CompositeConstraint_SimOpt.hpp.
References ROL::CompositeConstraint_SimOpt< Real >::conVal_, and ROL::CompositeConstraint_SimOpt< Real >::solveConRed().
|
inlinevirtual |
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 230 of file ROL_CompositeConstraint_SimOpt.hpp.
References ROL::CompositeConstraint_SimOpt< Real >::applyAdjointSens(), ROL::CompositeConstraint_SimOpt< Real >::conVal_, ROL::CompositeConstraint_SimOpt< Real >::dualZ_, ROL::CompositeConstraint_SimOpt< Real >::solveConRed(), and ROL::CompositeConstraint_SimOpt< Real >::Sz_.
|
inlinevirtual |
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 237 of file ROL_CompositeConstraint_SimOpt.hpp.
References ROL::CompositeConstraint_SimOpt< Real >::applySens(), ROL::CompositeConstraint_SimOpt< Real >::conVal_, ROL::CompositeConstraint_SimOpt< Real >::primZ_, and ROL::CompositeConstraint_SimOpt< Real >::Sz_.
|
inlinevirtual |
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 243 of file ROL_CompositeConstraint_SimOpt.hpp.
References ROL::CompositeConstraint_SimOpt< Real >::applyAdjointSens(), ROL::CompositeConstraint_SimOpt< Real >::applySens(), ROL::Vector< Real >::axpy(), ROL::CompositeConstraint_SimOpt< Real >::conRed_, ROL::CompositeConstraint_SimOpt< Real >::conVal_, ROL::CompositeConstraint_SimOpt< Real >::dualRed_, ROL::CompositeConstraint_SimOpt< Real >::dualZ1_, ROL::CompositeConstraint_SimOpt< Real >::dualZ_, ROL::Vector< Real >::plus(), ROL::CompositeConstraint_SimOpt< Real >::primZ_, ROL::CompositeConstraint_SimOpt< Real >::Sz_, and ROL::Vector< Real >::zero().
|
inlinevirtual |
Reimplemented from ROL::Constraint< Real >.
Definition at line 269 of file ROL_CompositeConstraint_SimOpt.hpp.
References ROL::CompositeConstraint_SimOpt< Real >::conRed_, ROL::CompositeConstraint_SimOpt< Real >::conVal_, ROL::CompositeConstraint_SimOpt< Real >::isConRedParametrized_, and ROL::Constraint< Real >::setParameter().
|
private |
Definition at line 78 of file ROL_CompositeConstraint_SimOpt.hpp.
Referenced by ROL::CompositeConstraint_SimOpt< Real >::applyAdjointHessian_11(), ROL::CompositeConstraint_SimOpt< Real >::applyAdjointHessian_12(), ROL::CompositeConstraint_SimOpt< Real >::applyAdjointHessian_21(), ROL::CompositeConstraint_SimOpt< Real >::applyAdjointHessian_22(), ROL::CompositeConstraint_SimOpt< Real >::applyAdjointJacobian_1(), ROL::CompositeConstraint_SimOpt< Real >::applyAdjointJacobian_2(), ROL::CompositeConstraint_SimOpt< Real >::applyInverseAdjointJacobian_1(), ROL::CompositeConstraint_SimOpt< Real >::applyInverseJacobian_1(), ROL::CompositeConstraint_SimOpt< Real >::applyJacobian_1(), ROL::CompositeConstraint_SimOpt< Real >::applyJacobian_2(), ROL::CompositeConstraint_SimOpt< Real >::setParameter(), ROL::CompositeConstraint_SimOpt< Real >::solve(), ROL::CompositeConstraint_SimOpt< Real >::update_1(), and ROL::CompositeConstraint_SimOpt< Real >::value().
|
private |
Definition at line 79 of file ROL_CompositeConstraint_SimOpt.hpp.
Referenced by ROL::CompositeConstraint_SimOpt< Real >::applyAdjointHessian_22(), ROL::CompositeConstraint_SimOpt< Real >::applyAdjointSens(), ROL::CompositeConstraint_SimOpt< Real >::applySens(), ROL::CompositeConstraint_SimOpt< Real >::setParameter(), and ROL::CompositeConstraint_SimOpt< Real >::solveConRed().
|
private |
Definition at line 81 of file ROL_CompositeConstraint_SimOpt.hpp.
Referenced by ROL::CompositeConstraint_SimOpt< Real >::applyAdjointHessian_12(), ROL::CompositeConstraint_SimOpt< Real >::applyAdjointHessian_21(), ROL::CompositeConstraint_SimOpt< Real >::applyAdjointHessian_22(), ROL::CompositeConstraint_SimOpt< Real >::applyAdjointJacobian_1(), ROL::CompositeConstraint_SimOpt< Real >::applyAdjointJacobian_2(), ROL::CompositeConstraint_SimOpt< Real >::applyAdjointSens(), ROL::CompositeConstraint_SimOpt< Real >::applyInverseAdjointJacobian_1(), ROL::CompositeConstraint_SimOpt< Real >::applyInverseJacobian_1(), ROL::CompositeConstraint_SimOpt< Real >::applyJacobian_1(), ROL::CompositeConstraint_SimOpt< Real >::applyJacobian_2(), ROL::CompositeConstraint_SimOpt< Real >::applySens(), ROL::CompositeConstraint_SimOpt< Real >::CompositeConstraint_SimOpt(), ROL::CompositeConstraint_SimOpt< Real >::solve(), ROL::CompositeConstraint_SimOpt< Real >::solveConRed(), ROL::CompositeConstraint_SimOpt< Real >::update_1(), and ROL::CompositeConstraint_SimOpt< Real >::value().
|
private |
Definition at line 82 of file ROL_CompositeConstraint_SimOpt.hpp.
Referenced by ROL::CompositeConstraint_SimOpt< Real >::applySens(), ROL::CompositeConstraint_SimOpt< Real >::CompositeConstraint_SimOpt(), and ROL::CompositeConstraint_SimOpt< Real >::solveConRed().
|
private |
|
private |
Definition at line 84 of file ROL_CompositeConstraint_SimOpt.hpp.
Referenced by ROL::CompositeConstraint_SimOpt< Real >::applyAdjointHessian_21(), ROL::CompositeConstraint_SimOpt< Real >::applyAdjointHessian_22(), ROL::CompositeConstraint_SimOpt< Real >::applyJacobian_2(), and ROL::CompositeConstraint_SimOpt< Real >::CompositeConstraint_SimOpt().
|
private |
Definition at line 85 of file ROL_CompositeConstraint_SimOpt.hpp.
Referenced by ROL::CompositeConstraint_SimOpt< Real >::applyAdjointHessian_12(), ROL::CompositeConstraint_SimOpt< Real >::applyAdjointHessian_22(), ROL::CompositeConstraint_SimOpt< Real >::applyAdjointJacobian_2(), and ROL::CompositeConstraint_SimOpt< Real >::CompositeConstraint_SimOpt().
|
private |
Definition at line 86 of file ROL_CompositeConstraint_SimOpt.hpp.
Referenced by ROL::CompositeConstraint_SimOpt< Real >::applyAdjointHessian_22(), and ROL::CompositeConstraint_SimOpt< Real >::CompositeConstraint_SimOpt().
|
private |
Definition at line 88 of file ROL_CompositeConstraint_SimOpt.hpp.
Referenced by ROL::CompositeConstraint_SimOpt< Real >::CompositeConstraint_SimOpt(), ROL::CompositeConstraint_SimOpt< Real >::solveConRed(), and ROL::CompositeConstraint_SimOpt< Real >::update_2().
|
private |
Definition at line 90 of file ROL_CompositeConstraint_SimOpt.hpp.
Referenced by ROL::CompositeConstraint_SimOpt< Real >::solveConRed(), and ROL::CompositeConstraint_SimOpt< Real >::update_2().
|
private |
Definition at line 91 of file ROL_CompositeConstraint_SimOpt.hpp.
Referenced by ROL::CompositeConstraint_SimOpt< Real >::solveConRed(), and ROL::CompositeConstraint_SimOpt< Real >::update_2().
|
private |
Definition at line 93 of file ROL_CompositeConstraint_SimOpt.hpp.
Referenced by ROL::CompositeConstraint_SimOpt< Real >::solveConRed().
|
private |
Definition at line 93 of file ROL_CompositeConstraint_SimOpt.hpp.
Referenced by ROL::CompositeConstraint_SimOpt< Real >::setParameter().