ROL
|
Defines the equality constraint operator interface for StdVectors. More...
#include <ROL_StdConstraint.hpp>
Public Member Functions | |
virtual | ~StdConstraint () |
void | update (const Vector< Real > &x, bool flag=true, int iter=-1) override |
Update constraint functions. x is the optimization variable, flag = true if optimization variable is changed, iter is the outer algorithm iterations count. More... | |
virtual void | update (const std::vector< Real > &x, bool flag=true, int iter=-1) |
void | update (const Vector< Real > &x, UpdateType type, int iter=-1) override |
Update constraint function. More... | |
virtual void | update (const std::vector< Real > &x, UpdateType type, int iter=-1) |
void | value (Vector< Real > &c, const Vector< Real > &x, Real &tol) override |
Evaluate the constraint operator \(c:\mathcal{X} \rightarrow \mathcal{C}\) at \(x\). More... | |
virtual void | value (std::vector< Real > &c, const std::vector< Real > &x, Real &tol)=0 |
void | applyJacobian (Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &x, Real &tol) override |
Apply the constraint Jacobian at \(x\), \(c'(x) \in L(\mathcal{X}, \mathcal{C})\), to vector \(v\). More... | |
virtual void | applyJacobian (std::vector< Real > &jv, const std::vector< Real > &v, const std::vector< Real > &x, Real &tol) |
void | applyAdjointJacobian (Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &x, Real &tol) override |
Apply the adjoint of the the constraint Jacobian at \(x\), \(c'(x)^* \in L(\mathcal{C}^*, \mathcal{X}^*)\), to vector \(v\). More... | |
virtual void | applyAdjointJacobian (std::vector< Real > &ajv, const std::vector< Real > &v, const std::vector< Real > &x, Real &tol) |
void | applyAdjointHessian (Vector< Real > &ahuv, const Vector< Real > &u, const Vector< Real > &v, const Vector< Real > &x, Real &tol) override |
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 void | applyAdjointHessian (std::vector< Real > &ahuv, const std::vector< Real > &u, const std::vector< Real > &v, const std::vector< Real > &x, Real &tol) |
std::vector< Real > | solveAugmentedSystem (Vector< Real > &v1, Vector< Real > &v2, const Vector< Real > &b1, const Vector< Real > &b2, const Vector< Real > &x, Real &tol) override |
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 or Riesz operator, and \(0 : \mathcal{C}^* \rightarrow \mathcal{C}\) is a zero operator. More... | |
virtual std::vector< Real > | solveAugmentedSystem (std::vector< Real > &v1, std::vector< Real > &v2, const std::vector< Real > &b1, const std::vector< Real > &b2, const std::vector< Real > &x, Real tol) |
void | applyPreconditioner (Vector< Real > &pv, const Vector< Real > &v, const Vector< Real > &x, const Vector< Real > &g, Real &tol) override |
Apply a constraint preconditioner at \(x\), \(P(x) \in L(\mathcal{C}, \mathcal{C}^*)\), to vector \(v\). Ideally, this preconditioner satisfies the following relationship:
\[ \left[c'(x) \circ R \circ c'(x)^* \circ P(x)\right] v = v \,, \] where R is the appropriate Riesz map in \(L(\mathcal{X}^*, \mathcal{X})\). It is used by the solveAugmentedSystem method. More... | |
virtual void | applyPreconditioner (std::vector< Real > &pv, const std::vector< Real > &v, const std::vector< Real > &x, const std::vector< Real > &g, Real &tol) |
Public Member Functions inherited from ROL::Constraint< Real > | |
virtual | ~Constraint (void) |
Constraint (void) | |
virtual void | applyAdjointJacobian (Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &x, const Vector< Real > &dualv, Real &tol) |
Apply the adjoint of the the constraint Jacobian at \(x\), \(c'(x)^* \in L(\mathcal{C}^*, \mathcal{X}^*)\), to vector \(v\). More... | |
void | activate (void) |
Turn on constraints. More... | |
void | deactivate (void) |
Turn off constraints. More... | |
bool | isActivated (void) |
Check if constraints are on. More... | |
virtual std::vector < std::vector< Real > > | checkApplyJacobian (const Vector< Real > &x, const Vector< Real > &v, const Vector< Real > &jv, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1) |
Finite-difference check for the constraint Jacobian application. More... | |
virtual std::vector < std::vector< Real > > | checkApplyJacobian (const Vector< Real > &x, const Vector< Real > &v, const Vector< Real > &jv, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1) |
Finite-difference check for the constraint Jacobian application. More... | |
virtual std::vector < std::vector< Real > > | checkApplyAdjointJacobian (const Vector< Real > &x, const Vector< Real > &v, const Vector< Real > &c, const Vector< Real > &ajv, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS) |
Finite-difference check for the application of the adjoint of constraint Jacobian. More... | |
virtual Real | checkAdjointConsistencyJacobian (const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &x, const bool printToStream=true, std::ostream &outStream=std::cout) |
virtual Real | checkAdjointConsistencyJacobian (const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &x, const Vector< Real > &dualw, const Vector< Real > &dualv, const bool printToStream=true, std::ostream &outStream=std::cout) |
virtual std::vector < std::vector< Real > > | checkApplyAdjointHessian (const Vector< Real > &x, const Vector< Real > &u, const Vector< Real > &v, const Vector< Real > &hv, const std::vector< Real > &step, const bool printToScreen=true, std::ostream &outStream=std::cout, const int order=1) |
Finite-difference check for the application of the adjoint of constraint Hessian. More... | |
virtual std::vector < std::vector< Real > > | checkApplyAdjointHessian (const Vector< Real > &x, const Vector< Real > &u, const Vector< Real > &v, const Vector< Real > &hv, const bool printToScreen=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1) |
Finite-difference check for the application of the adjoint of constraint Hessian. More... | |
virtual void | setParameter (const std::vector< Real > ¶m) |
Additional Inherited Members | |
Protected Member Functions inherited from ROL::Constraint< Real > | |
const std::vector< Real > | getParameter (void) const |
Defines the equality constraint operator interface for StdVectors.
Definition at line 25 of file ROL_StdConstraint.hpp.
|
inlinevirtual |
Definition at line 27 of file ROL_StdConstraint.hpp.
|
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::Constraint< Real >.
Definition at line 16 of file ROL_StdConstraint_Def.hpp.
References ROL::StdVector< Real, Element >::getVector(), and ROL::update().
|
inlinevirtual |
Definition at line 31 of file ROL_StdConstraint.hpp.
|
overridevirtual |
Update constraint function.
This function updates the constraint function at new iterations.
[in] | x | is the new iterate. |
[in] | type | is the type of update requested. |
[in] | iter | is the outer algorithm iterations count. |
Reimplemented from ROL::Constraint< Real >.
Definition at line 22 of file ROL_StdConstraint_Def.hpp.
References ROL::StdVector< Real, Element >::getVector(), and ROL::update().
|
inlinevirtual |
Definition at line 33 of file ROL_StdConstraint.hpp.
|
overridevirtual |
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::Constraint< Real >.
Definition at line 28 of file ROL_StdConstraint_Def.hpp.
References ROL::StdVector< Real, Element >::getVector(), and ROL::value.
|
pure virtual |
Implemented in ROL::ZOO::Constraint_Cantilever< Real >, ROL::ZOO::Constraint_CylinderHead< Real >, ROL::ZOO::Constraint_HS42b< Real >, ROL::ZOO::Constraint_HS14b< Real >, ROL::ZOO::Constraint_CantileverBeam< Real >, ROL::ZOO::Constraint_HS63b< Real >, ROL::ZOO::Constraint_HS55< Real >, ROL::ZOO::Constraint_HS42a< Real >, ROL::ZOO::Constraint_HS49< Real >, ROL::ZOO::Constraint_HS50< Real >, ROL::ZOO::Constraint_HS51< Real >, ROL::ZOO::Constraint_HS52< Real >, ROL::ZOO::Constraint_HS48< Real >, ROL::ZOO::Constraint_HS53< Real >, ROL::ZOO::Constraint_Cubic< Real >, ROL::ZOO::Constraint_HS14a< Real >, ROL::ZOO::Constraint_HS9< Real >, ROL::ZOO::Constraint_HS28< Real >, ROL::ZOO::Constraint_HS41< Real >, ROL::ZOO::Constraint_Quartic< Real >, ROL::ZOO::Constraint_HS63a< Real >, ROL::ZOO::Constraint_HS21< Real >, con2d< Real >, and con2d< Real >.
|
overridevirtual |
Apply the constraint Jacobian at \(x\), \(c'(x) \in L(\mathcal{X}, \mathcal{C})\), to vector \(v\).
@param[out] jv is the result of applying the constraint Jacobian to @b v at @b x; a constraint-space vector @param[in] v is an optimization-space vector @param[in] x is the constraint argument; an optimization-space vector @param[in,out] tol is a tolerance for inexact evaluations; currently unused On return, \form#91, where
\(v \in \mathcal{X}\), \(\mathsf{jv} \in \mathcal{C}\).
The default implementation is a finite-difference approximation.
Reimplemented from ROL::Constraint< Real >.
Definition at line 35 of file ROL_StdConstraint_Def.hpp.
References ROL::Constraint< Real >::applyJacobian(), applyJacobian(), and ROL::StdVector< Real, Element >::getVector().
|
virtual |
Reimplemented in ROL::ZOO::Constraint_Cantilever< Real >, ROL::ZOO::Constraint_CantileverBeam< Real >, ROL::ZOO::Constraint_CylinderHead< Real >, ROL::ZOO::Constraint_HS42b< Real >, ROL::ZOO::Constraint_HS14b< Real >, ROL::ZOO::Constraint_HS63b< Real >, ROL::ZOO::Constraint_HS55< Real >, ROL::ZOO::Constraint_HS50< Real >, ROL::ZOO::Constraint_HS51< Real >, ROL::ZOO::Constraint_HS52< Real >, ROL::ZOO::Constraint_HS42a< Real >, ROL::ZOO::Constraint_HS49< Real >, ROL::ZOO::Constraint_HS53< Real >, ROL::ZOO::Constraint_HS48< Real >, ROL::ZOO::Constraint_HS14a< Real >, ROL::ZOO::Constraint_HS9< Real >, ROL::ZOO::Constraint_Cubic< Real >, ROL::ZOO::Constraint_HS28< Real >, ROL::ZOO::Constraint_HS41< Real >, ROL::ZOO::Constraint_Quartic< Real >, ROL::ZOO::Constraint_HS63a< Real >, ROL::ZOO::Constraint_HS21< Real >, con2d< Real >, and con2d< Real >.
Definition at line 50 of file ROL_StdConstraint_Def.hpp.
|
overridevirtual |
Apply the adjoint of the the constraint Jacobian at \(x\), \(c'(x)^* \in L(\mathcal{C}^*, \mathcal{X}^*)\), to vector \(v\).
@param[out] ajv is the result of applying the adjoint of the constraint Jacobian to @b v at @b x; a dual optimization-space vector @param[in] v is a dual constraint-space vector @param[in] x is the constraint argument; an optimization-space vector @param[in,out] tol is a tolerance for inexact evaluations; currently unused On return, \form#95, where
\(v \in \mathcal{C}^*\), \(\mathsf{ajv} \in \mathcal{X}^*\).
The default implementation is a finite-difference approximation.
Reimplemented from ROL::Constraint< Real >.
Definition at line 58 of file ROL_StdConstraint_Def.hpp.
References ROL::Constraint< Real >::applyAdjointJacobian(), and ROL::StdVector< Real, Element >::getVector().
|
virtual |
Reimplemented in ROL::ZOO::Constraint_Cantilever< Real >, ROL::ZOO::Constraint_CantileverBeam< Real >, ROL::ZOO::Constraint_CylinderHead< Real >, ROL::ZOO::Constraint_HS42b< Real >, ROL::ZOO::Constraint_HS14b< Real >, ROL::ZOO::Constraint_HS63b< Real >, ROL::ZOO::Constraint_HS55< Real >, ROL::ZOO::Constraint_HS50< Real >, ROL::ZOO::Constraint_HS51< Real >, ROL::ZOO::Constraint_HS52< Real >, ROL::ZOO::Constraint_HS53< Real >, ROL::ZOO::Constraint_HS49< Real >, ROL::ZOO::Constraint_HS48< Real >, ROL::ZOO::Constraint_HS42a< Real >, ROL::ZOO::Constraint_HS14a< Real >, ROL::ZOO::Constraint_HS9< Real >, ROL::ZOO::Constraint_HS28< Real >, ROL::ZOO::Constraint_Cubic< Real >, ROL::ZOO::Constraint_HS41< Real >, ROL::ZOO::Constraint_Quartic< Real >, ROL::ZOO::Constraint_HS63a< Real >, ROL::ZOO::Constraint_HS21< Real >, con2d< Real >, and con2d< Real >.
Definition at line 73 of file ROL_StdConstraint_Def.hpp.
|
overridevirtual |
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 \).
@param[out] ahuv is the result of applying the derivative of the adjoint of the constraint Jacobian at @b x to vector @b u in direction @b v; a dual optimization-space vector @param[in] u is the direction vector; a dual constraint-space vector @param[in] v is an optimization-space vector @param[in] x is the constraint argument; an optimization-space vector @param[in,out] tol is a tolerance for inexact evaluations; currently unused On return, \form#100, 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::Constraint< Real >.
Definition at line 81 of file ROL_StdConstraint_Def.hpp.
References ROL::Constraint< Real >::applyAdjointHessian(), and ROL::StdVector< Real, Element >::getVector().
|
virtual |
Reimplemented in ROL::ZOO::Constraint_HS42b< Real >, ROL::ZOO::Constraint_HS14b< Real >, ROL::ZOO::Constraint_HS63b< Real >, ROL::ZOO::Constraint_HS55< Real >, ROL::ZOO::Constraint_HS50< Real >, ROL::ZOO::Constraint_HS51< Real >, ROL::ZOO::Constraint_HS52< Real >, ROL::ZOO::Constraint_HS53< Real >, ROL::ZOO::Constraint_HS49< Real >, ROL::ZOO::Constraint_HS48< Real >, ROL::ZOO::Constraint_HS42a< Real >, ROL::ZOO::Constraint_HS14a< Real >, ROL::ZOO::Constraint_HS28< Real >, ROL::ZOO::Constraint_HS41< Real >, ROL::ZOO::Constraint_HS9< Real >, ROL::ZOO::Constraint_HS63a< Real >, and ROL::ZOO::Constraint_HS21< Real >.
Definition at line 99 of file ROL_StdConstraint_Def.hpp.
|
overridevirtual |
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 or Riesz 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::Constraint< Real >.
Definition at line 108 of file ROL_StdConstraint_Def.hpp.
References ROL::StdVector< Real, Element >::getVector(), and ROL::Constraint< Real >::solveAugmentedSystem().
|
virtual |
Definition at line 128 of file ROL_StdConstraint_Def.hpp.
|
overridevirtual |
Apply a constraint preconditioner at \(x\), \(P(x) \in L(\mathcal{C}, \mathcal{C}^*)\), to vector \(v\). Ideally, this preconditioner satisfies the following relationship:
\[ \left[c'(x) \circ R \circ c'(x)^* \circ P(x)\right] v = v \,, \]
where R is the appropriate Riesz map in \(L(\mathcal{X}^*, \mathcal{X})\). It is used by the solveAugmentedSystem method.
@param[out] pv is the result of applying the constraint preconditioner to @b v at @b x; a dual constraint-space vector @param[in] v is a constraint-space vector @param[in] x is the preconditioner argument; an optimization-space vector @param[in] g is the preconditioner argument; a dual optimization-space vector, unused @param[in,out] tol is a tolerance for inexact evaluations On return, \form#114, where
\(v \in \mathcal{C}\), \(\mathsf{pv} \in \mathcal{C}^*\).
The default implementation is the Riesz map in \(L(\mathcal{C}, \mathcal{C}^*)\).
Reimplemented from ROL::Constraint< Real >.
Definition at line 139 of file ROL_StdConstraint_Def.hpp.
References ROL::Constraint< Real >::applyPreconditioner(), and ROL::StdVector< Real, Element >::getVector().
|
virtual |
Definition at line 157 of file ROL_StdConstraint_Def.hpp.