ROL
|
Defines the general constraint operator interface. More...
#include <ROL_Constraint.hpp>
Public Member Functions | |
virtual | ~Constraint (void) |
Constraint (void) | |
virtual void | update (const Vector< Real > &x, UpdateType type, int iter=-1) |
Update constraint function. 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)=0 |
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 | 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... | |
virtual void | applyAdjointHessian (Vector< Real > &ahuv, const Vector< Real > &u, 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 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 or Riesz 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\). 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... | |
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) |
Protected Member Functions | |
const std::vector< Real > | getParameter (void) const |
Private Attributes | |
bool | activated_ |
std::vector< Real > | param_ |
Defines the general constraint operator interface.
ROL's constraint interface is designed for Fréchet differentiable operators \(c:\mathcal{X} \rightarrow \mathcal{C}\), where \(\mathcal{X}\) and \(\mathcal{C}\) are Banach spaces. The constraints are of the form
\[ c(x) = 0 \,. \]
The basic operator interface, to be implemented by the user, requires:
It is strongly recommended that the user additionally overload:
The user may also overload:
Definition at line 52 of file ROL_Constraint.hpp.
|
inlinevirtual |
Definition at line 57 of file ROL_Constraint.hpp.
|
inline |
Definition at line 59 of file ROL_Constraint.hpp.
|
inlinevirtual |
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 in ROL::ROL::Constraint_SimOpt< Real >, ROL::Constraint_SimOpt< Real >, ROL::StochasticConstraint< Real >, ROL::AffineTransformConstraint< Real >, ROL::RiskNeutralConstraint< Real >, ROL::ChainRuleConstraint< Real >, ROL::Constraint_Partitioned< Real >, ROL::ConstraintFromObjective< Real >, ROL::AlmostSureConstraint< Real >, ROL::SimulatedConstraint< Real >, ROL::ElasticLinearConstraint< Real >, ROL::LinearConstraint< Real >, ROL::ROL::SimConstraint< Real >, ROL::SlacklessConstraint< Real >, ROL::SimConstraint< Real >, ROL::StdConstraint< Real >, ROL::StdConstraint< RealT >, ROL::MeanValueConstraint< Real >, and ROL::RiskLessConstraint< Real >.
Definition at line 68 of file ROL_Constraint.hpp.
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::accept(), ROL::CompositeStep< Real >::accept(), ROL::TypeE::CompositeStepAlgorithm< Real >::initialize(), ROL::CompositeStep< Real >::initialize(), ROL::CompositeStep< Real >::update(), ROL::MoreauYosidaPenaltyStep< Real >::update(), ROL::TypeE::CompositeStepAlgorithm< Real >::updateRadius(), ROL::TypeG::MoreauYosidaAlgorithm< Real >::updateState(), ROL::TypeG::InteriorPointAlgorithm< Real >::updateState(), and ROL::MoreauYosidaPenaltyStep< Real >::updateState().
|
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::ROL::Constraint_SimOpt< Real >, ROL::Constraint_SimOpt< Real >, ROL::Reduced_Constraint_SimOpt< Real >, ROL::PrimalDualInteriorPointResidual< Real >, ROL::PrimalDualInteriorPointResidual< Real >, ROL::StochasticConstraint< Real >, ROL::ChainRuleConstraint< Real >, ROL::Constraint_DynamicState< Real >, ROL::AffineTransformConstraint< Real >, ROL::Constraint_Partitioned< Real >, ROL::RiskNeutralConstraint< Real >, ROL::ConstraintFromObjective< Real >, ROL::ElasticLinearConstraint< Real >, ROL::AlmostSureConstraint< Real >, ROL::LinearConstraint< Real >, ROL::SimulatedConstraint< Real >, ROL::SlacklessConstraint< Real >, ROL::ROL::SimConstraint< Real >, ROL::SimConstraint< Real >, ROL::StdConstraint< Real >, ROL::StdConstraint< RealT >, ROL::MeanValueConstraint< Real >, and ROL::RiskLessConstraint< Real >.
Definition at line 79 of file ROL_Constraint.hpp.
|
pure virtual |
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}\).
Implemented in ROL::ROL::Constraint_SimOpt< Real >, ROL::Constraint_SimOpt< Real >, Normalization_Constraint< Real >, Normalization_Constraint< Real >, ROL::Reduced_Constraint_SimOpt< Real >, ROL::ZOO::EqualityConstraint_SimpleEqConstrained< Real, XPrim, XDual, CPrim, CDual >, ROL::PrimalDualInteriorPointResidual< Real >, ROL::PrimalDualInteriorPointResidual< Real >, ROL::ZOO::Constraint_ParaboloidCircle< Real, XPrim, XDual, CPrim, CDual >, ROL::ZOO::InequalityConstraint_HS32< Real >, ROL::ZOO::Constraint_HS39b< Real >, ROL::ZOO::InequalityConstraint_HS29< Real >, ROL::ZOO::EqualityConstraint_HS32< Real >, ROL::InteriorPoint::PrimalDualResidual< Real >, ROL::InteriorPoint::PrimalDualResidual< Real >, ROL::ZOO::Constraint_HS24< Real >, ROL::BinaryConstraint< Real >, ROL::ZOO::Constraint_HS39a< Real >, ROL::StochasticConstraint< Real >, ROL::ChainRuleConstraint< Real >, ROL::ScalarLinearConstraint< Real >, ROL::RiskNeutralConstraint< Real >, ROL::AffineTransformConstraint< Real >, ROL::Constraint_Partitioned< Real >, ROL::ReducedLinearConstraint< Real >, ROL::ConstraintFromObjective< Real >, ROL::AlmostSureConstraint< Real >, ROL::SimulatedConstraint< Real >, ROL::ElasticLinearConstraint< Real >, ROL::LinearConstraint< Real >, ROL::Constraint_DynamicState< Real >, ROL::StdConstraint< Real >, ROL::StdConstraint< RealT >, ROL::BoundToConstraint< Real >, ROL::SlacklessConstraint< Real >, ROL::ROL::SimConstraint< Real >, ROL::SimConstraint< Real >, ROL::LowerBoundToConstraint< Real >, ROL::UpperBoundToConstraint< Real >, ROL::MeanValueConstraint< Real >, and ROL::RiskLessConstraint< Real >.
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::accept(), ROL::CompositeStep< Real >::accept(), ROL::CompositeStep< Real >::compute(), ROL::TypeE::CompositeStepAlgorithm< Real >::computeTrial(), ROL::TypeE::CompositeStepAlgorithm< Real >::initialize(), ROL::InteriorPointStep< Real >::initialize(), ROL::CompositeStep< Real >::initialize(), testRandomInputs(), ROL::CompositeStep< Real >::update(), ROL::TypeE::CompositeStepAlgorithm< Real >::updateRadius(), ROL::TypeG::MoreauYosidaAlgorithm< Real >::updateState(), ROL::TypeG::InteriorPointAlgorithm< Real >::updateState(), and ROL::MoreauYosidaPenaltyStep< Real >::updateState().
|
virtual |
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 in ROL::ROL::Constraint_SimOpt< Real >, ROL::Constraint_SimOpt< Real >, Normalization_Constraint< Real >, ROL::PrimalDualInteriorPointResidual< Real >, Normalization_Constraint< Real >, ROL::Reduced_Constraint_SimOpt< Real >, ROL::ZOO::EqualityConstraint_SimpleEqConstrained< Real, XPrim, XDual, CPrim, CDual >, ROL::PrimalDualInteriorPointResidual< Real >, ROL::ZOO::Constraint_ParaboloidCircle< Real, XPrim, XDual, CPrim, CDual >, ROL::ZOO::InequalityConstraint_HS32< Real >, ROL::InteriorPoint::PrimalDualResidual< Real >, ROL::InteriorPoint::PrimalDualResidual< Real >, ROL::ZOO::Constraint_HS39b< Real >, ROL::ZOO::InequalityConstraint_HS29< Real >, ROL::ZOO::EqualityConstraint_HS32< Real >, ROL::ZOO::Constraint_HS24< Real >, ROL::BinaryConstraint< Real >, ROL::ZOO::Constraint_HS39a< Real >, ROL::StochasticConstraint< Real >, ROL::SimulatedConstraint< Real >, ROL::ChainRuleConstraint< Real >, ROL::RiskNeutralConstraint< Real >, ROL::AlmostSureConstraint< Real >, ROL::ScalarLinearConstraint< Real >, ROL::AffineTransformConstraint< Real >, ROL::Constraint_Partitioned< Real >, ROL::ReducedLinearConstraint< Real >, ROL::ConstraintFromObjective< Real >, ROL::Constraint_DynamicState< Real >, ROL::StdConstraint< Real >, ROL::StdConstraint< RealT >, ROL::ElasticLinearConstraint< Real >, ROL::LinearConstraint< Real >, ROL::BoundToConstraint< Real >, ROL::SlacklessConstraint< Real >, ROL::ROL::SimConstraint< Real >, ROL::SimConstraint< Real >, ROL::LowerBoundToConstraint< Real >, ROL::UpperBoundToConstraint< Real >, ROL::MeanValueConstraint< Real >, and ROL::RiskLessConstraint< Real >.
Definition at line 19 of file ROL_ConstraintDef.hpp.
References ROL::Vector< Real >::axpy(), ROL::Vector< Real >::clone(), ROL::Vector< Real >::norm(), ROL::Vector< Real >::scale(), ROL::Temp, ROL::update(), ROL::value, and ROL::Vector< Real >::zero().
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::accept(), ROL::CompositeStep< Real >::accept(), ROL::StdConstraint< Real >::applyJacobian(), ROL::TypeE::CompositeStepAlgorithm< Real >::computeQuasinormalStep(), ROL::CompositeStep< Real >::computeQuasinormalStep(), ROL::TypeE::CompositeStepAlgorithm< Real >::solveTangentialSubproblem(), ROL::CompositeStep< Real >::solveTangentialSubproblem(), and testRandomInputs().
|
virtual |
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 in ROL::ROL::Constraint_SimOpt< Real >, ROL::Constraint_SimOpt< Real >, Normalization_Constraint< Real >, Normalization_Constraint< Real >, ROL::ZOO::EqualityConstraint_SimpleEqConstrained< Real, XPrim, XDual, CPrim, CDual >, ROL::Reduced_Constraint_SimOpt< Real >, ROL::ZOO::Constraint_ParaboloidCircle< Real, XPrim, XDual, CPrim, CDual >, ROL::ZOO::InequalityConstraint_HS32< Real >, ROL::ZOO::Constraint_HS39b< Real >, ROL::ZOO::InequalityConstraint_HS29< Real >, ROL::SimulatedConstraint< Real >, ROL::ZOO::Constraint_HS24< Real >, ROL::ZOO::EqualityConstraint_HS32< Real >, ROL::ZOO::Constraint_HS39a< Real >, ROL::BinaryConstraint< Real >, ROL::ChainRuleConstraint< Real >, ROL::RiskNeutralConstraint< Real >, ROL::StochasticConstraint< Real >, ROL::AlmostSureConstraint< Real >, ROL::ScalarLinearConstraint< Real >, ROL::AffineTransformConstraint< Real >, ROL::Constraint_Partitioned< Real >, ROL::ReducedLinearConstraint< Real >, ROL::ConstraintFromObjective< Real >, ROL::StdConstraint< Real >, ROL::StdConstraint< RealT >, ROL::Constraint_DynamicState< Real >, ROL::ElasticLinearConstraint< Real >, ROL::LinearConstraint< Real >, ROL::ROL::SimConstraint< Real >, ROL::BoundToConstraint< Real >, ROL::SimConstraint< Real >, ROL::LowerBoundToConstraint< Real >, ROL::UpperBoundToConstraint< Real >, ROL::MeanValueConstraint< Real >, and ROL::RiskLessConstraint< Real >.
Definition at line 52 of file ROL_ConstraintDef.hpp.
References ROL::Vector< Real >::dual().
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::accept(), ROL::CompositeStep< Real >::accept(), ROL::StdConstraint< Real >::applyAdjointJacobian(), ROL::CompositeStep< Real >::compute(), ROL::TypeE::CompositeStepAlgorithm< Real >::computeLagrangeMultiplier(), ROL::CompositeStep< Real >::computeLagrangeMultiplier(), ROL::TypeE::CompositeStepAlgorithm< Real >::computeQuasinormalStep(), ROL::CompositeStep< Real >::computeQuasinormalStep(), ROL::TypeE::CompositeStepAlgorithm< Real >::computeTrial(), ROL::TypeE::AugmentedLagrangianAlgorithm< Real >::initialize(), ROL::TypeG::StabilizedLCLAlgorithm< Real >::initialize(), ROL::TypeE::StabilizedLCLAlgorithm< Real >::initialize(), ROL::TypeG::AugmentedLagrangianAlgorithm< Real >::initialize(), ROL::TypeE::CompositeStepAlgorithm< Real >::initialize(), ROL::CompositeStep< Real >::initialize(), ROL::AugmentedLagrangianStep< Real >::initialize(), ROL::TypeE::StabilizedLCLAlgorithm< Real >::run(), ROL::TypeG::StabilizedLCLAlgorithm< Real >::run(), ROL::CompositeStep< Real >::update(), ROL::TypeE::CompositeStepAlgorithm< Real >::updateRadius(), ROL::TypeG::MoreauYosidaAlgorithm< Real >::updateState(), ROL::TypeG::InteriorPointAlgorithm< Real >::updateState(), and ROL::MoreauYosidaPenaltyStep< Real >::updateState().
|
virtual |
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] dualv is a vector used for temporary variables; a constraint-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 in ROL::ElasticLinearConstraint< Real >, ROL::LinearConstraint< Real >, and ROL::SlacklessConstraint< Real >.
Definition at line 64 of file ROL_ConstraintDef.hpp.
References ROL::Vector< Real >::axpy(), ROL::Vector< Real >::basis(), ROL::Vector< Real >::clone(), ROL::Vector< Real >::dimension(), ROL::Vector< Real >::norm(), ROL::Temp, ROL::update(), ROL::value, and ROL::Vector< Real >::zero().
|
virtual |
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 in ROL::ROL::Constraint_SimOpt< Real >, ROL::Constraint_SimOpt< Real >, Normalization_Constraint< Real >, ROL::ZOO::EqualityConstraint_SimpleEqConstrained< Real, XPrim, XDual, CPrim, CDual >, Normalization_Constraint< Real >, ROL::Reduced_Constraint_SimOpt< Real >, ROL::ZOO::Constraint_ParaboloidCircle< Real, XPrim, XDual, CPrim, CDual >, ROL::ZOO::InequalityConstraint_HS32< Real >, ROL::ZOO::Constraint_HS39b< Real >, ROL::SimulatedConstraint< Real >, ROL::ZOO::InequalityConstraint_HS29< Real >, ROL::ZOO::Constraint_HS24< Real >, ROL::ZOO::EqualityConstraint_HS32< Real >, ROL::ZOO::Constraint_HS39a< Real >, ROL::AlmostSureConstraint< Real >, ROL::RiskNeutralConstraint< Real >, ROL::ChainRuleConstraint< Real >, ROL::BinaryConstraint< Real >, ROL::StochasticConstraint< Real >, ROL::ScalarLinearConstraint< Real >, ROL::Constraint_Partitioned< Real >, ROL::AffineTransformConstraint< Real >, ROL::StdConstraint< Real >, ROL::StdConstraint< RealT >, ROL::ReducedLinearConstraint< Real >, ROL::ConstraintFromObjective< Real >, ROL::Constraint_DynamicState< Real >, ROL::ElasticLinearConstraint< Real >, ROL::LinearConstraint< Real >, ROL::ROL::SimConstraint< Real >, ROL::BoundToConstraint< Real >, ROL::SlacklessConstraint< Real >, ROL::SimConstraint< Real >, ROL::LowerBoundToConstraint< Real >, ROL::UpperBoundToConstraint< Real >, ROL::MeanValueConstraint< Real >, and ROL::RiskLessConstraint< Real >.
Definition at line 132 of file ROL_ConstraintDef.hpp.
References ROL::Vector< Real >::axpy(), ROL::Vector< Real >::clone(), ROL::Vector< Real >::norm(), ROL::Vector< Real >::scale(), ROL::Temp, ROL::update(), and ROL::Vector< Real >::zero().
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::accept(), ROL::CompositeStep< Real >::accept(), ROL::StdConstraint< Real >::applyAdjointHessian(), ROL::ZOO::Constraint_ParaboloidCircle< Real, XPrim, XDual, CPrim, CDual >::applyAdjointHessian(), ROL::Reduced_Constraint_SimOpt< Real >::applyAdjointHessian(), ROL::TypeE::CompositeStepAlgorithm< Real >::solveTangentialSubproblem(), and ROL::CompositeStep< Real >::solveTangentialSubproblem().
|
virtual |
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 in ROL::ROL::Constraint_SimOpt< Real >, ROL::Constraint_SimOpt< Real >, Normalization_Constraint< Real >, ROL::ScalarLinearConstraint< Real >, ROL::StdConstraint< Real >, and ROL::StdConstraint< RealT >.
Definition at line 161 of file ROL_ConstraintDef.hpp.
References applyJacobian(), ROL::Vector< Real >::clone(), ROL::Vector< Real >::dot(), ROL::Vector< Real >::dual(), ROL::Vector< Real >::plus(), ROL::Objective_SerialSimOpt< Real >::zero, zero, and ROL::Vector< Real >::zero().
Referenced by ROL::TypeE::CompositeStepAlgorithm< Real >::accept(), ROL::CompositeStep< Real >::accept(), ROL::TypeE::CompositeStepAlgorithm< Real >::computeLagrangeMultiplier(), ROL::CompositeStep< Real >::computeLagrangeMultiplier(), ROL::TypeE::CompositeStepAlgorithm< Real >::computeQuasinormalStep(), ROL::CompositeStep< Real >::computeQuasinormalStep(), ROL::StdConstraint< Real >::solveAugmentedSystem(), Normalization_Constraint< Real >::solveAugmentedSystem(), ROL::Constraint_SimOpt< Real >::solveAugmentedSystem(), ROL::TypeE::CompositeStepAlgorithm< Real >::solveTangentialSubproblem(), and ROL::CompositeStep< Real >::solveTangentialSubproblem().
|
inlinevirtual |
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 in ROL::ROL::Constraint_SimOpt< Real >, ROL::Constraint_SimOpt< Real >, ROL::SimulatedConstraint< Real >, ROL::AlmostSureConstraint< Real >, ROL::StdConstraint< Real >, ROL::StdConstraint< RealT >, ROL::Constraint_Partitioned< Real >, ROL::Constraint_DynamicState< Real >, ROL::ROL::SimConstraint< Real >, and ROL::SimConstraint< Real >.
Definition at line 249 of file ROL_Constraint.hpp.
Referenced by ROL::StdConstraint< Real >::applyPreconditioner(), and ROL::Constraint_SimOpt< Real >::applyPreconditioner().
|
inline |
Turn on constraints.
Definition at line 259 of file ROL_Constraint.hpp.
|
inline |
Turn off constraints.
Definition at line 263 of file ROL_Constraint.hpp.
|
inline |
Check if constraints are on.
Definition at line 267 of file ROL_Constraint.hpp.
|
virtual |
Finite-difference check for the constraint Jacobian application.
Details here.
Definition at line 349 of file ROL_ConstraintDef.hpp.
References applyJacobian(), ROL::Vector< Real >::clone(), ROL::Finite_Difference_Arrays::shifts, ROL::Temp, ROL::update(), ROL::value, and ROL::Finite_Difference_Arrays::weights.
Referenced by main().
|
virtual |
Finite-difference check for the constraint Jacobian application.
Details here.
Definition at line 330 of file ROL_ConstraintDef.hpp.
|
virtual |
Finite-difference check for the application of the adjoint of constraint Jacobian.
Details here. (This function should be deprecated)
Definition at line 454 of file ROL_ConstraintDef.hpp.
References ROL::Vector< Real >::basis(), ROL::Vector< Real >::clone(), ROL::Vector< Real >::dimension(), ROL::Temp, ROL::update(), and ROL::value.
Referenced by main().
|
inlinevirtual |
Definition at line 320 of file ROL_Constraint.hpp.
Referenced by ROL::Constraint< RealT >::checkAdjointConsistencyJacobian(), and main().
|
virtual |
Definition at line 550 of file ROL_ConstraintDef.hpp.
References ROL::Vector< Real >::apply(), applyJacobian(), ROL::Vector< Real >::clone(), ROL::Temp, and ROL::update().
|
virtual |
Finite-difference check for the application of the adjoint of constraint Hessian.
Details here.
Definition at line 604 of file ROL_ConstraintDef.hpp.
References ROL::Vector< Real >::clone(), ROL::Finite_Difference_Arrays::shifts, ROL::Temp, ROL::update(), and ROL::Finite_Difference_Arrays::weights.
Referenced by main().
|
virtual |
Finite-difference check for the application of the adjoint of constraint Hessian.
Details here.
Definition at line 586 of file ROL_ConstraintDef.hpp.
|
inlineprotected |
Definition at line 369 of file ROL_Constraint.hpp.
Referenced by Constraint_BurgersControl< Real >::applyAdjointHessian_11(), Constraint_BurgersControl< Real >::applyAdjointHessian_12(), Constraint_BurgersControl< Real >::applyAdjointHessian_21(), Constraint_BurgersControl< Real >::applyAdjointHessian_22(), Constraint_BurgersControl< Real >::applyAdjointJacobian_1(), DiffusionConstraint< Real >::applyAdjointJacobian_2(), Constraint_BurgersControl< Real >::applyAdjointJacobian_2(), Constraint_BurgersControl< Real >::applyInverseAdjointJacobian_1(), DiffusionConstraint< Real >::applyInverseJacobian_1(), Constraint_BurgersControl< Real >::applyInverseJacobian_1(), Constraint_BurgersControl< Real >::applyJacobian_1(), DiffusionConstraint< Real >::applyJacobian_1(), DiffusionConstraint< Real >::applyJacobian_2(), Constraint_BurgersControl< Real >::applyJacobian_2(), Constraint_BurgersControl< Real >::compute_pde_jacobian(), Constraint_BurgersControl< Real >::compute_residual(), DiffusionConstraint< Real >::solve(), ROL::CompositeConstraint_SimOpt< Real >::solveConRed(), DiffusionConstraint< Real >::value(), and Constraint_BurgersControl< Real >::value().
|
inlinevirtual |
Reimplemented in ROL::Reduced_Constraint_SimOpt< Real >, ROL::CompositeConstraint_SimOpt< Real >, ROL::Constraint_Partitioned< Real >, ROL::ConstraintFromObjective< Real >, ROL::ROL::SimConstraint< Real >, ROL::SlacklessConstraint< Real >, ROL::SimConstraint< Real >, and ROL::RiskLessConstraint< Real >.
Definition at line 374 of file ROL_Constraint.hpp.
Referenced by ROL::RiskLessConstraint< Real >::setParameter(), ROL::SlacklessConstraint< Real >::setParameter(), ROL::ConstraintFromObjective< Real >::setParameter(), ROL::Constraint_Partitioned< Real >::setParameter(), ROL::CompositeConstraint_SimOpt< Real >::setParameter(), and ROL::Reduced_Constraint_SimOpt< Real >::setParameter().
|
private |
Definition at line 54 of file ROL_Constraint.hpp.
Referenced by ROL::Constraint< RealT >::activate(), ROL::Constraint< RealT >::deactivate(), and ROL::Constraint< RealT >::isActivated().
|
private |
Definition at line 366 of file ROL_Constraint.hpp.
Referenced by ROL::Constraint< RealT >::getParameter(), and ROL::Constraint< RealT >::setParameter().