ROL
|
Defines the equality constraint operator interface for simulation-based optimization. More...
#include <ROL_EqualityConstraint_SimOpt.hpp>
Public Member Functions | |
virtual 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... | |
virtual void | value (Vector< Real > &c, const Vector< Real > &u, const Vector< Real > &z, Real &tol)=0 |
Evaluate the constraint operator \(c:\mathcal{U}\times\mathcal{Z} \rightarrow \mathcal{C}\) at \((u,z)\). More... | |
virtual void | solve (Vector< Real > &u, const Vector< Real > &z, Real &tol) |
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) |
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) |
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) |
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) |
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_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, 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... | |
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 void | applyInverseAdjointJacobian_1 (Vector< Real > &iajv, 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... | |
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) |
Apply the adjoint of the partial constraint Hessian at \((u,z)\), \(c_{uu}(u,z)^* \in L(L(\mathcal{C}^*, \mathcal{U}^*), \mathcal{U}^*)\), to vector \(v\) in direction \(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) |
Apply the adjoint of the partial constraint Hessian at \((u,z)\), \(c_{uz}(u,z)^* \in L(L(\mathcal{C}^*, \mathcal{U}^*), \mathcal{Z}^*)\), to vector \(v\) in direction \(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) |
Apply the adjoint of the partial constraint Hessian at \((u,z)\), \(c_{zu}(u,z)^* \in L(L(\mathcal{C}^*, \mathcal{Z}^*), \mathcal{U}^*)\), to vector \(v\) in direction \(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) |
Apply the adjoint of the partial constraint Hessian at \((u,z)\), \(c_{zz}(u,z)^* \in L(L(\mathcal{C}^*, \mathcal{Z}^*), \mathcal{Z}^*)\), to vector \(v\) in direction \(w\). 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... | |
EqualityConstraint_SimOpt (void) | |
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 bool | isFeasible (const Vector< Real > &v) |
Check if the vector, v, is feasible. 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 ROL::Vector< Real > &u, const ROL::Vector< Real > &z, const ROL::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) |
![]() | |
virtual | ~EqualityConstraint () |
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... | |
EqualityConstraint (void) | |
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... | |
Defines the equality constraint operator interface for simulation-based optimization.
This equality constraint interface inherits from ROL_EqualityConstraint, 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,z) = 0 \,. \]
The basic operator interface, to be implemented by the user, requires:
The user may also overload:
Definition at line 95 of file ROL_EqualityConstraint_SimOpt.hpp.
|
inline |
Definition at line 683 of file ROL_EqualityConstraint_SimOpt.hpp.
|
inlinevirtual |
Update constraint functions. x is the optimization variable, flag = true if optimization variable is changed, iter is the outer algorithm iterations count.
Reimplemented in ROL::ParametrizedEqualityConstraint_SimOpt< Real >.
Definition at line 102 of file ROL_EqualityConstraint_SimOpt.hpp.
Referenced by ROL::EqualityConstraint_SimOpt< Real >::applyAdjointHessian_11(), ROL::EqualityConstraint_SimOpt< Real >::applyAdjointHessian_12(), ROL::EqualityConstraint_SimOpt< Real >::applyAdjointHessian_21(), ROL::EqualityConstraint_SimOpt< Real >::applyAdjointHessian_22(), ROL::EqualityConstraint_SimOpt< Real >::applyAdjointJacobian_1(), ROL::EqualityConstraint_SimOpt< Real >::applyAdjointJacobian_2(), ROL::EqualityConstraint_SimOpt< Real >::applyJacobian_1(), ROL::EqualityConstraint_SimOpt< Real >::applyJacobian_2(), and ROL::EqualityConstraint_SimOpt< Real >::update().
|
pure virtual |
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}\).
Implemented in EqualityConstraint_BurgersControl< Real >, EqualityConstraint_BurgersControl< Real >, EqualityConstraint_BurgersControl< Real >, and ROL::ParametrizedEqualityConstraint_SimOpt< Real >.
Referenced by ROL::EqualityConstraint_SimOpt< Real >::applyAdjointJacobian_1(), ROL::EqualityConstraint_SimOpt< Real >::applyAdjointJacobian_2(), ROL::EqualityConstraint_SimOpt< Real >::applyJacobian_1(), ROL::EqualityConstraint_SimOpt< Real >::applyJacobian_2(), ROL::EqualityConstraint_SimOpt< Real >::checkSolve(), and ROL::EqualityConstraint_SimOpt< Real >::value().
|
inlinevirtual |
Given \(z\), solve \(c(u,z)=0\) for \(u\).
[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 |
Reimplemented in EqualityConstraint_BurgersControl< Real >, EqualityConstraint_BurgersControl< Real >, EqualityConstraint_BurgersControl< Real >, and ROL::ParametrizedEqualityConstraint_SimOpt< Real >.
Definition at line 131 of file ROL_EqualityConstraint_SimOpt.hpp.
Referenced by ROL::EqualityConstraint_SimOpt< Real >::checkSolve().
|
inlinevirtual |
Apply the partial constraint Jacobian at \((u,z)\), \(c_u(u,z) \in L(\mathcal{U}, \mathcal{C})\), to the vector \(v\).
[out] | jv | is the result of applying the constraint Jacobian to v at \((u,z)\); a constraint-space vector |
[in] | v | is a simulation-space vector |
[in] | u | is the constraint argument; an 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{jv} = c_u(u,z)v\), where \(v \in \mathcal{U}\), \(\mathsf{jv} \in \mathcal{C}\).
Reimplemented in EqualityConstraint_BurgersControl< Real >, EqualityConstraint_BurgersControl< Real >, and EqualityConstraint_BurgersControl< Real >.
Definition at line 153 of file ROL_EqualityConstraint_SimOpt.hpp.
References ROL::Vector< Real >::axpy(), ROL::Vector< Real >::clone(), ROL::Vector< Real >::norm(), ROL::ROL_EPSILON, ROL::Vector< Real >::scale(), ROL::EqualityConstraint_SimOpt< Real >::update(), and ROL::EqualityConstraint_SimOpt< Real >::value().
Referenced by ROL::EqualityConstraint_SimOpt< Real >::applyJacobian(), ROL::EqualityConstraint_SimOpt< Real >::checkAdjointConsistencyJacobian_1(), and ROL::EqualityConstraint_SimOpt< Real >::checkInverseJacobian_1().
|
inlinevirtual |
Apply the partial constraint Jacobian at \((u,z)\), \(c_z(u,z) \in L(\mathcal{Z}, \mathcal{C})\), to the vector \(v\).
[out] | jv | is the result of applying the constraint Jacobian to v at \((u,z)\); a constraint-space vector |
[in] | v | is an optimization-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{jv} = c_z(u,z)v\), where \(v \in \mathcal{Z}\), \(\mathsf{jv} \in \mathcal{C}\).
Reimplemented in EqualityConstraint_BurgersControl< Real >, EqualityConstraint_BurgersControl< Real >, and EqualityConstraint_BurgersControl< Real >.
Definition at line 196 of file ROL_EqualityConstraint_SimOpt.hpp.
References ROL::Vector< Real >::axpy(), ROL::Vector< Real >::clone(), ROL::Vector< Real >::norm(), ROL::ROL_EPSILON, ROL::Vector< Real >::scale(), ROL::EqualityConstraint_SimOpt< Real >::update(), and ROL::EqualityConstraint_SimOpt< Real >::value().
Referenced by ROL::EqualityConstraint_SimOpt< Real >::applyJacobian(), and ROL::EqualityConstraint_SimOpt< Real >::checkAdjointConsistencyJacobian_2().
|
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\).
[out] | ijv | is the result of applying the inverse constraint Jacobian to v at \((u,z)\); a simulation-space vector |
[in] | v | is 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{ijv} = c_u(u,z)^{-1}v\), where \(v \in \mathcal{C}\), \(\mathsf{ijv} \in \mathcal{U}\).
Reimplemented in EqualityConstraint_BurgersControl< Real >, EqualityConstraint_BurgersControl< Real >, and EqualityConstraint_BurgersControl< Real >.
Definition at line 238 of file ROL_EqualityConstraint_SimOpt.hpp.
Referenced by ROL::EqualityConstraint_SimOpt< Real >::applyPreconditioner(), and ROL::EqualityConstraint_SimOpt< Real >::checkInverseJacobian_1().
|
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.
[out] | ajv | is the result of applying the adjoint of the constraint Jacobian to v at (u,z); a dual simulation-space vector |
[in] | v | is a dual 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{ajv} = c_u(u,z)^*v\), where \(v \in \mathcal{C}^*\), \(\mathsf{ajv} \in \mathcal{U}^*\).
Reimplemented in EqualityConstraint_BurgersControl< Real >, EqualityConstraint_BurgersControl< Real >, and EqualityConstraint_BurgersControl< Real >.
Definition at line 262 of file ROL_EqualityConstraint_SimOpt.hpp.
References ROL::Vector< Real >::dual().
Referenced by ROL::EqualityConstraint_SimOpt< Real >::applyAdjointHessian_11(), ROL::EqualityConstraint_SimOpt< Real >::applyAdjointHessian_21(), ROL::EqualityConstraint_SimOpt< Real >::applyAdjointJacobian(), ROL::EqualityConstraint_SimOpt< Real >::checkAdjointConsistencyJacobian_1(), and ROL::EqualityConstraint_SimOpt< Real >::checkInverseAdjointJacobian_1().
|
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 secondary interface, for use with dual spaces where the user does not define the dual() operation.
[out] | ajv | is the result of applying the adjoint of the constraint Jacobian to v at (u,z); a dual simulation-space vector |
[in] | v | is a dual constraint-space vector |
[in] | u | is the constraint argument; a simulation-space vector |
[in] | z | is the constraint argument; an optimization-space vector |
[in] | dualv | is a vector used for temporary variables; a constraint-space vector |
[in,out] | tol | is a tolerance for inexact evaluations; currently unused |
On return, \(\mathsf{ajv} = c_u(u,z)^*v\), where \(v \in \mathcal{C}^*\), \(\mathsf{ajv} \in \mathcal{U}^*\).
Definition at line 288 of file ROL_EqualityConstraint_SimOpt.hpp.
References ROL::Vector< Real >::axpy(), ROL::Vector< Real >::basis(), ROL::Vector< Real >::clone(), ROL::Vector< Real >::dimension(), ROL::Vector< Real >::dual(), ROL::Vector< Real >::norm(), ROL::ROL_EPSILON, ROL::EqualityConstraint_SimOpt< Real >::update(), ROL::EqualityConstraint_SimOpt< Real >::value(), and ROL::Vector< Real >::zero().
|
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.
[out] | ajv | is the result of applying the adjoint of the constraint Jacobian to v at \((u,z)\); a dual optimization-space vector |
[in] | v | is a dual 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{ajv} = c_z(u,z)^*v\), where \(v \in \mathcal{C}^*\), \(\mathsf{ajv} \in \mathcal{Z}^*\).
Reimplemented in EqualityConstraint_BurgersControl< Real >, EqualityConstraint_BurgersControl< Real >, and EqualityConstraint_BurgersControl< Real >.
Definition at line 333 of file ROL_EqualityConstraint_SimOpt.hpp.
References ROL::Vector< Real >::dual().
Referenced by ROL::EqualityConstraint_SimOpt< Real >::applyAdjointHessian_12(), ROL::EqualityConstraint_SimOpt< Real >::applyAdjointHessian_22(), ROL::EqualityConstraint_SimOpt< Real >::applyAdjointJacobian(), and ROL::EqualityConstraint_SimOpt< Real >::checkAdjointConsistencyJacobian_2().
|
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 secondary interface, for use with dual spaces where the user does not define the dual() operation.
[out] | ajv | is the result of applying the adjoint of the constraint Jacobian to v at \((u,z)\); a dual optimization-space vector |
[in] | v | is a dual constraint-space vector |
[in] | u | is the constraint argument; a simulation-space vector |
[in] | z | is the constraint argument; an optimization-space vector |
[in] | dualv | is a vector used for temporary variables; a constraint-space vector |
[in,out] | tol | is a tolerance for inexact evaluations; currently unused |
On return, \(\mathsf{ajv} = c_z(u,z)^*v\), where \(v \in \mathcal{C}^*\), \(\mathsf{ajv} \in \mathcal{Z}^*\).
Definition at line 359 of file ROL_EqualityConstraint_SimOpt.hpp.
References ROL::Vector< Real >::axpy(), ROL::Vector< Real >::basis(), ROL::Vector< Real >::clone(), ROL::Vector< Real >::dimension(), ROL::Vector< Real >::dual(), ROL::Vector< Real >::norm(), ROL::ROL_EPSILON, ROL::EqualityConstraint_SimOpt< Real >::update(), ROL::EqualityConstraint_SimOpt< Real >::value(), and ROL::Vector< Real >::zero().
|
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\).
[out] | iajv | is the result of applying the inverse adjoint of the constraint Jacobian to v at (u,z); a dual constraint-space vector |
[in] | v | is a dual simulation-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{iajv} = c_u(u,z)^{-*}v\), where \(v \in \mathcal{U}^*\), \(\mathsf{iajv} \in \mathcal{C}^*\).
Reimplemented in EqualityConstraint_BurgersControl< Real >, EqualityConstraint_BurgersControl< Real >, and EqualityConstraint_BurgersControl< Real >.
Definition at line 403 of file ROL_EqualityConstraint_SimOpt.hpp.
Referenced by ROL::EqualityConstraint_SimOpt< Real >::applyPreconditioner(), and ROL::EqualityConstraint_SimOpt< Real >::checkInverseAdjointJacobian_1().
|
inlinevirtual |
Apply the adjoint of the partial constraint Hessian at \((u,z)\), \(c_{uu}(u,z)^* \in L(L(\mathcal{C}^*, \mathcal{U}^*), \mathcal{U}^*)\), to vector \(v\) in direction \(w\).
[out] | ahwv | is the result of applying the adjoint of the constraint Hessian to v at \((u,z)\) in direction w; a dual simulation-space vector |
[in] | w | is the direction vector; a dual constraint-space vector |
[in] | v | is a dual simulation-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{ahwv} = c_{uu}(u,z)^*(w,v) \), where \(w \in \mathcal{C}^*\), \(v \in \mathcal{U}^*\), and \(\mathsf{ahwv} \in \mathcal{U}^*\).
Reimplemented in EqualityConstraint_BurgersControl< Real >, EqualityConstraint_BurgersControl< Real >, and EqualityConstraint_BurgersControl< Real >.
Definition at line 428 of file ROL_EqualityConstraint_SimOpt.hpp.
References ROL::EqualityConstraint_SimOpt< Real >::applyAdjointJacobian_1(), ROL::Vector< Real >::axpy(), ROL::Vector< Real >::clone(), ROL::Vector< Real >::norm(), ROL::ROL_EPSILON, ROL::Vector< Real >::scale(), and ROL::EqualityConstraint_SimOpt< Real >::update().
Referenced by ROL::EqualityConstraint_SimOpt< Real >::applyAdjointHessian().
|
inlinevirtual |
Apply the adjoint of the partial constraint Hessian at \((u,z)\), \(c_{uz}(u,z)^* \in L(L(\mathcal{C}^*, \mathcal{U}^*), \mathcal{Z}^*)\), to vector \(v\) in direction \(w\).
[out] | ahwv | is the result of applying the adjoint of the constraint Hessian to v at \((u,z)\) in direction w; a dual optimization-space vector |
[in] | w | is the direction vector; a dual constraint-space vector |
[in] | v | is a dual simulation-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{ahwv} = c_{uz}(u,z)^*(w,v) \), where \(w \in \mathcal{C}^*\), \(v \in \mathcal{U}^*\), and \(\mathsf{ahwv} \in \mathcal{Z}^*\).
Reimplemented in EqualityConstraint_BurgersControl< Real >, EqualityConstraint_BurgersControl< Real >, and EqualityConstraint_BurgersControl< Real >.
Definition at line 472 of file ROL_EqualityConstraint_SimOpt.hpp.
References ROL::EqualityConstraint_SimOpt< Real >::applyAdjointJacobian_2(), ROL::Vector< Real >::axpy(), ROL::Vector< Real >::clone(), ROL::Vector< Real >::norm(), ROL::ROL_EPSILON, ROL::Vector< Real >::scale(), and ROL::EqualityConstraint_SimOpt< Real >::update().
Referenced by ROL::EqualityConstraint_SimOpt< Real >::applyAdjointHessian().
|
inlinevirtual |
Apply the adjoint of the partial constraint Hessian at \((u,z)\), \(c_{zu}(u,z)^* \in L(L(\mathcal{C}^*, \mathcal{Z}^*), \mathcal{U}^*)\), to vector \(v\) in direction \(w\).
[out] | ahwv | is the result of applying the adjoint of the constraint Hessian to v at \((u,z)\) in direction w; a dual simulation-space vector |
[in] | w | is the direction vector; a dual constraint-space vector |
[in] | v | is a dual optimization-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{ahwv} = c_{zu}(u,z)^*(w,v) \), where \(w \in \mathcal{C}^*\), \(v \in \mathcal{Z}^*\), and \(\mathsf{ahwv} \in \mathcal{U}^*\).
Reimplemented in EqualityConstraint_BurgersControl< Real >, EqualityConstraint_BurgersControl< Real >, and EqualityConstraint_BurgersControl< Real >.
Definition at line 516 of file ROL_EqualityConstraint_SimOpt.hpp.
References ROL::EqualityConstraint_SimOpt< Real >::applyAdjointJacobian_1(), ROL::Vector< Real >::axpy(), ROL::Vector< Real >::clone(), ROL::Vector< Real >::norm(), ROL::ROL_EPSILON, ROL::Vector< Real >::scale(), and ROL::EqualityConstraint_SimOpt< Real >::update().
Referenced by ROL::EqualityConstraint_SimOpt< Real >::applyAdjointHessian().
|
inlinevirtual |
Apply the adjoint of the partial constraint Hessian at \((u,z)\), \(c_{zz}(u,z)^* \in L(L(\mathcal{C}^*, \mathcal{Z}^*), \mathcal{Z}^*)\), to vector \(v\) in direction \(w\).
[out] | ahwv | is the result of applying the adjoint of the constraint Hessian to v at \((u,z)\) in direction w; a dual optimization-space vector |
[in] | w | is the direction vector; a dual constraint-space vector |
[in] | v | is a dual optimization-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{ahwv} = c_{zz}(u,z)^*(w,v) \), where \(w \in \mathcal{C}^*\), \(v \in \mathcal{Z}^*\), and \(\mathsf{ahwv} \in \mathcal{Z}^*\).
Reimplemented in EqualityConstraint_BurgersControl< Real >, EqualityConstraint_BurgersControl< Real >, and EqualityConstraint_BurgersControl< Real >.
Definition at line 559 of file ROL_EqualityConstraint_SimOpt.hpp.
References ROL::EqualityConstraint_SimOpt< Real >::applyAdjointJacobian_2(), ROL::Vector< Real >::axpy(), ROL::Vector< Real >::clone(), ROL::Vector< Real >::norm(), ROL::ROL_EPSILON, ROL::Vector< Real >::scale(), and ROL::EqualityConstraint_SimOpt< Real >::update().
Referenced by ROL::EqualityConstraint_SimOpt< Real >::applyAdjointHessian().
|
inlinevirtual |
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.
[out] | v1 | is the optimization-space component of the result |
[out] | v2 | is the dual constraint-space component of the result |
[in] | b1 | is the dual optimization-space component of the right-hand side |
[in] | b2 | is the constraint-space component of the right-hand side |
[in] | x | is the constraint argument; an optimization-space vector |
[in,out] | tol | is the nominal relative residual tolerance |
On return, \( [\mathsf{v1} \,\, \mathsf{v2}] \) approximately solves the augmented system, where the size of the residual is governed by special stopping conditions.
The default implementation is the preconditioned generalized minimal residual (GMRES) method, which enables the use of nonsymmetric preconditioners.
Reimplemented from ROL::EqualityConstraint< Real >.
Definition at line 624 of file ROL_EqualityConstraint_SimOpt.hpp.
References ROL::EqualityConstraint< Real >::solveAugmentedSystem().
|
inlinevirtual |
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.
[out] | pv | is the result of applying the constraint preconditioner to v at x; a constraint-space vector |
[in] | v | is a constraint-space vector |
[in] | x | is the preconditioner argument; an optimization-space vector |
[in,out] | tol | is a tolerance for inexact evaluations |
On return, \(\mathsf{pv} = P(x)v\), where \(v \in \mathcal{C}\), \(\mathsf{pv} \in \mathcal{C}\).
The default implementation is a null-op.
Reimplemented from ROL::EqualityConstraint< Real >.
Definition at line 652 of file ROL_EqualityConstraint_SimOpt.hpp.
References ROL::EqualityConstraint_SimOpt< Real >::applyInverseAdjointJacobian_1(), ROL::EqualityConstraint_SimOpt< Real >::applyInverseJacobian_1(), ROL::EqualityConstraint< Real >::applyPreconditioner(), ROL::Vector_SimOpt< Real >::get_1(), ROL::Vector_SimOpt< Real >::get_2(), and ROL::Vector< Real >::set().
|
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::EqualityConstraint< Real >.
Definition at line 690 of file ROL_EqualityConstraint_SimOpt.hpp.
References ROL::Vector_SimOpt< Real >::get_1(), ROL::Vector_SimOpt< Real >::get_2(), and ROL::EqualityConstraint_SimOpt< Real >::update().
|
inlinevirtual |
Check if the vector, v, is feasible.
Reimplemented from ROL::EqualityConstraint< Real >.
Definition at line 698 of file ROL_EqualityConstraint_SimOpt.hpp.
|
inlinevirtual |
Evaluate the constraint operator \(c:\mathcal{X} \rightarrow \mathcal{C}\) at \(x\).
[out] | c | is the result of evaluating the constraint operator at x; a constraint-space vector |
[in] | x | is the constraint argument; an optimization-space vector |
[in,out] | tol | is a tolerance for inexact evaluations; currently unused |
On return, \(\mathsf{c} = c(x)\), where \(\mathsf{c} \in \mathcal{C}\), \(\mathsf{x} \in \mathcal{X}\).
Implements ROL::EqualityConstraint< Real >.
Definition at line 700 of file ROL_EqualityConstraint_SimOpt.hpp.
References ROL::Vector_SimOpt< Real >::get_1(), ROL::Vector_SimOpt< Real >::get_2(), and ROL::EqualityConstraint_SimOpt< Real >::value().
|
inlinevirtual |
Apply the constraint Jacobian at \(x\), \(c'(x) \in L(\mathcal{X}, \mathcal{C})\), to vector \(v\).
[out] | jv | is the result of applying the constraint Jacobian to v at x; a constraint-space vector |
[in] | v | is an optimization-space vector |
[in] | x | is the constraint argument; an optimization-space vector |
[in,out] | tol | is a tolerance for inexact evaluations; currently unused |
On return, \(\mathsf{jv} = c'(x)v\), where \(v \in \mathcal{X}\), \(\mathsf{jv} \in \mathcal{C}\).
The default implementation is a finite-difference approximation.
Reimplemented from ROL::EqualityConstraint< Real >.
Definition at line 709 of file ROL_EqualityConstraint_SimOpt.hpp.
References ROL::EqualityConstraint_SimOpt< Real >::applyJacobian_1(), ROL::EqualityConstraint_SimOpt< Real >::applyJacobian_2(), ROL::Vector< Real >::clone(), ROL::Vector_SimOpt< Real >::get_1(), ROL::Vector_SimOpt< Real >::get_2(), and ROL::Vector< Real >::plus().
|
inlinevirtual |
Apply the adjoint of the the constraint Jacobian at \(x\), \(c'(x)^* \in L(\mathcal{C}^*, \mathcal{X}^*)\), to vector \(v\).
[out] | ajv | is the result of applying the adjoint of the constraint Jacobian to v at x; a dual optimization-space vector |
[in] | v | is a dual constraint-space vector |
[in] | x | is the constraint argument; an optimization-space vector |
[in,out] | tol | is a tolerance for inexact evaluations; currently unused |
On return, \(\mathsf{ajv} = c'(x)^*v\), where \(v \in \mathcal{C}^*\), \(\mathsf{ajv} \in \mathcal{X}^*\).
The default implementation is a finite-difference approximation.
Reimplemented from ROL::EqualityConstraint< Real >.
Definition at line 724 of file ROL_EqualityConstraint_SimOpt.hpp.
References ROL::EqualityConstraint_SimOpt< Real >::applyAdjointJacobian_1(), ROL::EqualityConstraint_SimOpt< Real >::applyAdjointJacobian_2(), ROL::Vector_SimOpt< Real >::get_1(), ROL::Vector_SimOpt< Real >::get_2(), ROL::Vector_SimOpt< Real >::set_1(), and ROL::Vector_SimOpt< Real >::set_2().
|
inlinevirtual |
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 \).
[out] | ahuv | is the result of applying the derivative of the adjoint of the constraint Jacobian at x to vector u in direction v; a dual optimization-space vector |
[in] | u | is the direction vector; a dual constraint-space vector |
[in] | v | is an optimization-space vector |
[in] | x | is the constraint argument; an optimization-space vector |
[in,out] | tol | is a tolerance for inexact evaluations; currently unused |
On return, \( \mathsf{ahuv} = c''(x)(v,\cdot)^*u \), where \(u \in \mathcal{C}^*\), \(v \in \mathcal{X}\), and \(\mathsf{ahuv} \in \mathcal{X}^*\).
The default implementation is a finite-difference approximation based on the adjoint Jacobian.
Reimplemented from ROL::EqualityConstraint< Real >.
Definition at line 741 of file ROL_EqualityConstraint_SimOpt.hpp.
References ROL::EqualityConstraint_SimOpt< Real >::applyAdjointHessian_11(), ROL::EqualityConstraint_SimOpt< Real >::applyAdjointHessian_12(), ROL::EqualityConstraint_SimOpt< Real >::applyAdjointHessian_21(), ROL::EqualityConstraint_SimOpt< Real >::applyAdjointHessian_22(), ROL::Vector_SimOpt< Real >::get_1(), ROL::Vector_SimOpt< Real >::get_2(), ROL::Vector_SimOpt< Real >::set_1(), and ROL::Vector_SimOpt< Real >::set_2().
|
inlinevirtual |
Definition at line 770 of file ROL_EqualityConstraint_SimOpt.hpp.
References ROL::Vector< Real >::clone(), ROL::ROL_EPSILON, ROL::EqualityConstraint_SimOpt< Real >::solve(), and ROL::EqualityConstraint_SimOpt< Real >::value().
Referenced by main().
|
inlinevirtual |
Check the consistency of the Jacobian and its adjoint. This is the primary interface.
[out] | w | is a dual constraint-space vector |
[in] | v | is a simulation-space vector |
[in] | u | is the constraint argument; a simulation-space vector |
[in] | z | is the constraint argument; an optimization-space vector |
[in] | printToStream | is is a flag that turns on/off output |
[in] | outStream | is the output stream |
Definition at line 806 of file ROL_EqualityConstraint_SimOpt.hpp.
References ROL::Vector< Real >::dual().
Referenced by main().
|
inlinevirtual |
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.
[out] | w | is a dual constraint-space vector |
[in] | v | is a simulation-space vector |
[in] | u | is the constraint argument; a simulation-space vector |
[in] | z | is the constraint argument; an optimization-space vector |
[in] | dualw | is a constraint-space vector |
[in] | dualv | is a dual simulation-space vector |
[in] | printToStream | is is a flag that turns on/off output |
[in] | outStream | is the output stream |
Definition at line 831 of file ROL_EqualityConstraint_SimOpt.hpp.
References ROL::EqualityConstraint_SimOpt< Real >::applyAdjointJacobian_1(), ROL::EqualityConstraint_SimOpt< Real >::applyJacobian_1(), ROL::Vector< Real >::clone(), ROL::Vector< Real >::dot(), ROL::ROL_EPSILON, and ROL::ROL_UNDERFLOW.
|
inlinevirtual |
Check the consistency of the Jacobian and its adjoint. This is the primary interface.
[out] | w | is a dual constraint-space vector |
[in] | v | is an optimization-space vector |
[in] | u | is the constraint argument; a simulation-space vector |
[in] | z | is the constraint argument; an optimization-space vector |
[in] | printToStream | is is a flag that turns on/off output |
[in] | outStream | is the output stream |
Definition at line 872 of file ROL_EqualityConstraint_SimOpt.hpp.
References ROL::Vector< Real >::dual().
Referenced by main().
|
inlinevirtual |
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.
[out] | w | is a dual constraint-space vector |
[in] | v | is an optimization-space vector |
[in] | u | is the constraint argument; a simulation-space vector |
[in] | z | is the constraint argument; an optimization-space vector |
[in] | dualw | is a constraint-space vector |
[in] | dualv | is a dual optimization-space vector |
[in] | printToStream | is is a flag that turns on/off output |
[in] | outStream | is the output stream |
Definition at line 896 of file ROL_EqualityConstraint_SimOpt.hpp.
References ROL::EqualityConstraint_SimOpt< Real >::applyAdjointJacobian_2(), ROL::EqualityConstraint_SimOpt< Real >::applyJacobian_2(), ROL::Vector< Real >::clone(), ROL::Vector< Real >::dot(), ROL::ROL_EPSILON, and ROL::ROL_UNDERFLOW.
|
inlinevirtual |
Definition at line 924 of file ROL_EqualityConstraint_SimOpt.hpp.
References ROL::EqualityConstraint_SimOpt< Real >::applyInverseJacobian_1(), ROL::EqualityConstraint_SimOpt< Real >::applyJacobian_1(), ROL::Vector< Real >::clone(), ROL::Vector< Real >::norm(), ROL::ROL_EPSILON, and ROL::ROL_UNDERFLOW.
Referenced by main().
|
inlinevirtual |
Definition at line 952 of file ROL_EqualityConstraint_SimOpt.hpp.
References ROL::EqualityConstraint_SimOpt< Real >::applyAdjointJacobian_1(), ROL::EqualityConstraint_SimOpt< Real >::applyInverseAdjointJacobian_1(), ROL::Vector< Real >::clone(), ROL::Vector< Real >::norm(), ROL::ROL_EPSILON, and ROL::ROL_UNDERFLOW.
Referenced by main().