ROL
Public Member Functions | Protected Member Functions | Private Attributes | List of all members
ROL::Constraint< Real > Class Template Referenceabstract

Defines the general constraint operator interface. More...

#include <ROL_Constraint.hpp>

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

Public Member Functions

virtual ~Constraint (void)
 
 Constraint (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 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 > &param)
 

Protected Member Functions

const std::vector< Real > getParameter (void) const
 

Private Attributes

bool activated_
 
std::vector< Real > param_
 

Detailed Description

template<class Real>
class ROL::Constraint< Real >

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 85 of file ROL_Constraint.hpp.

Constructor & Destructor Documentation

template<class Real>
virtual ROL::Constraint< Real >::~Constraint ( void  )
inlinevirtual

Definition at line 90 of file ROL_Constraint.hpp.

template<class Real>
ROL::Constraint< Real >::Constraint ( void  )
inline

Definition at line 92 of file ROL_Constraint.hpp.

Member Function Documentation

template<class Real>
virtual void ROL::Constraint< Real >::update ( const Vector< Real > &  x,
bool  flag = true,
int  iter = -1 
)
inlinevirtual
template<class Real>
virtual void ROL::Constraint< Real >::value ( Vector< Real > &  c,
const Vector< Real > &  x,
Real &  tol 
)
pure virtual

Evaluate the constraint operator \(c:\mathcal{X} \rightarrow \mathcal{C}\) at \(x\).

Parameters
[out]cis the result of evaluating the constraint operator at x; a constraint-space vector
[in]xis the constraint argument; an optimization-space vector
[in,out]tolis 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::BinaryConstraint< Real >, ROL::ZOO::InequalityConstraint_HS29< Real >, ROL::InteriorPoint::PrimalDualResidual< Real >, ROL::InteriorPoint::PrimalDualResidual< Real >, ROL::ZOO::EqualityConstraint_HS32< Real >, ROL::ZOO::Constraint_HS24< Real >, ROL::Constraint_Partitioned< Real >, ROL::ZOO::Constraint_HS39a< Real >, ROL::ConstraintFromObjective< Real >, ROL::AffineTransformConstraint< Real >, ROL::StochasticConstraint< Real >, ROL::ScalarLinearConstraint< Real >, ROL::LinearConstraint< Real >, ROL::RiskNeutralConstraint< Real >, ROL::BoundToConstraint< Real >, ROL::StdConstraint< Real >, ROL::AlmostSureConstraint< Real >, ROL::LowerBoundToConstraint< Real >, ROL::UpperBoundToConstraint< Real >, ROL::SimulatedConstraint< Real >, ROL::Constraint_DynamicState< Real >, ROL::ROL::Constraint_State< Real >, ROL::Constraint_State< Real >, and ROL::RiskLessConstraint< Real >.

Referenced by ROL::CompositeStep< Real >::accept(), ROL::CompositeStep< Real >::compute(), ROL::InteriorPointStep< Real >::initialize(), ROL::CompositeStep< Real >::initialize(), ROL::CompositeStep< Real >::update(), and ROL::MoreauYosidaPenaltyStep< Real >::updateState().

template<class Real >
void ROL::Constraint< Real >::applyJacobian ( Vector< Real > &  jv,
const Vector< Real > &  v,
const Vector< Real > &  x,
Real &  tol 
)
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#77, 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::BinaryConstraint< Real >, ROL::ZOO::InequalityConstraint_HS29< Real >, ROL::Constraint_Partitioned< Real >, ROL::ZOO::EqualityConstraint_HS32< Real >, ROL::ZOO::Constraint_HS24< Real >, ROL::ZOO::Constraint_HS39a< Real >, ROL::ConstraintFromObjective< Real >, ROL::AffineTransformConstraint< Real >, ROL::SimulatedConstraint< Real >, ROL::StochasticConstraint< Real >, ROL::ScalarLinearConstraint< Real >, ROL::RiskNeutralConstraint< Real >, ROL::LinearConstraint< Real >, ROL::AlmostSureConstraint< Real >, ROL::StdConstraint< Real >, ROL::BoundToConstraint< Real >, ROL::LowerBoundToConstraint< Real >, ROL::UpperBoundToConstraint< Real >, ROL::Constraint_DynamicState< Real >, ROL::ROL::Constraint_State< Real >, ROL::RiskLessConstraint< Real >, and ROL::Constraint_State< Real >.

Definition at line 54 of file ROL_ConstraintDef.hpp.

References ROL::Vector< Real >::axpy(), ROL::Vector< Real >::clone(), ROL::Vector< Real >::norm(), ROL::Vector< Real >::scale(), ROL::update(), ROL::value, and ROL::Vector< Real >::zero().

Referenced by ROL::CompositeStep< Real >::accept(), ROL::StdConstraint< Real >::applyJacobian(), ROL::CompositeStep< Real >::computeQuasinormalStep(), and ROL::CompositeStep< Real >::solveTangentialSubproblem().

template<class Real >
void ROL::Constraint< Real >::applyAdjointJacobian ( Vector< Real > &  ajv,
const Vector< Real > &  v,
const Vector< Real > &  x,
Real &  tol 
)
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#81, 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::BinaryConstraint< Real >, ROL::ZOO::Constraint_HS39b< Real >, ROL::ZOO::InequalityConstraint_HS29< Real >, ROL::Constraint_Partitioned< Real >, ROL::ZOO::Constraint_HS24< Real >, ROL::SimulatedConstraint< Real >, ROL::ZOO::EqualityConstraint_HS32< Real >, ROL::ConstraintFromObjective< Real >, ROL::ZOO::Constraint_HS39a< Real >, ROL::AffineTransformConstraint< Real >, ROL::StdConstraint< Real >, ROL::RiskNeutralConstraint< Real >, ROL::AlmostSureConstraint< Real >, ROL::StochasticConstraint< Real >, ROL::ScalarLinearConstraint< Real >, ROL::LinearConstraint< Real >, ROL::BoundToConstraint< Real >, ROL::UpperBoundToConstraint< Real >, ROL::LowerBoundToConstraint< Real >, ROL::Constraint_DynamicState< Real >, ROL::RiskLessConstraint< Real >, ROL::ROL::Constraint_State< Real >, and ROL::Constraint_State< Real >.

Definition at line 87 of file ROL_ConstraintDef.hpp.

References ROL::Vector< Real >::dual().

Referenced by ROL::CompositeStep< Real >::accept(), ROL::StdConstraint< Real >::applyAdjointJacobian(), ROL::CompositeStep< Real >::compute(), ROL::CompositeStep< Real >::computeLagrangeMultiplier(), ROL::CompositeStep< Real >::computeQuasinormalStep(), ROL::CompositeStep< Real >::initialize(), ROL::AugmentedLagrangianStep< Real >::initialize(), ROL::CompositeStep< Real >::update(), and ROL::MoreauYosidaPenaltyStep< Real >::updateState().

template<class Real >
void ROL::Constraint< Real >::applyAdjointJacobian ( Vector< Real > &  ajv,
const Vector< Real > &  v,
const Vector< Real > &  x,
const Vector< Real > &  dualv,
Real &  tol 
)
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#81, where

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

The default implementation is a finite-difference approximation.


Definition at line 99 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 >::dual(), ROL::Vector< Real >::norm(), ROL::update(), ROL::value, and ROL::Vector< Real >::zero().

template<class Real >
void ROL::Constraint< Real >::applyAdjointHessian ( Vector< Real > &  ahuv,
const Vector< Real > &  u,
const Vector< Real > &  v,
const Vector< Real > &  x,
Real &  tol 
)
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#86, 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::BinaryConstraint< Real >, ROL::SimulatedConstraint< Real >, ROL::Constraint_Partitioned< Real >, ROL::ZOO::InequalityConstraint_HS29< Real >, ROL::ZOO::Constraint_HS24< Real >, ROL::ZOO::EqualityConstraint_HS32< Real >, ROL::ZOO::Constraint_HS39a< Real >, ROL::StdConstraint< Real >, ROL::AlmostSureConstraint< Real >, ROL::RiskNeutralConstraint< Real >, ROL::AffineTransformConstraint< Real >, ROL::ScalarLinearConstraint< Real >, ROL::StochasticConstraint< Real >, ROL::LinearConstraint< Real >, ROL::BoundToConstraint< Real >, ROL::UpperBoundToConstraint< Real >, ROL::LowerBoundToConstraint< Real >, ROL::Constraint_DynamicState< Real >, ROL::RiskLessConstraint< Real >, ROL::ROL::Constraint_State< Real >, and ROL::Constraint_State< Real >.

Definition at line 166 of file ROL_ConstraintDef.hpp.

References ROL::Vector< Real >::axpy(), ROL::Vector< Real >::clone(), ROL::Vector< Real >::norm(), ROL::Vector< Real >::scale(), ROL::update(), and ROL::Vector< Real >::zero().

Referenced by ROL::CompositeStep< Real >::accept(), ROL::StdConstraint< Real >::applyAdjointHessian(), ROL::ZOO::Constraint_ParaboloidCircle< Real, XPrim, XDual, CPrim, CDual >::applyAdjointHessian(), ROL::Reduced_Constraint_SimOpt< Real >::applyAdjointHessian(), and ROL::CompositeStep< Real >::solveTangentialSubproblem().

template<class Real >
std::vector< Real > ROL::Constraint< Real >::solveAugmentedSystem ( Vector< Real > &  v1,
Vector< Real > &  v2,
const Vector< Real > &  b1,
const Vector< Real > &  b2,
const Vector< Real > &  x,
Real &  tol 
)
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.

Parameters
[out]v1is the optimization-space component of the result
[out]v2is the dual constraint-space component of the result
[in]b1is the dual optimization-space component of the right-hand side
[in]b2is the constraint-space component of the right-hand side
[in]xis the constraint argument; an optimization-space vector
[in,out]tolis 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::StdConstraint< Real >, and ROL::ScalarLinearConstraint< Real >.

Definition at line 195 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::CompositeStep< Real >::accept(), ROL::CompositeStep< Real >::computeLagrangeMultiplier(), ROL::CompositeStep< Real >::computeQuasinormalStep(), ROL::StdConstraint< Real >::solveAugmentedSystem(), Normalization_Constraint< Real >::solveAugmentedSystem(), ROL::Constraint_SimOpt< Real >::solveAugmentedSystem(), and ROL::CompositeStep< Real >::solveTangentialSubproblem().

template<class Real>
virtual void ROL::Constraint< Real >::applyPreconditioner ( Vector< Real > &  pv,
const Vector< Real > &  v,
const Vector< Real > &  x,
const Vector< Real > &  g,
Real &  tol 
)
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#100, 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::Constraint_Partitioned< Real >, ROL::StdConstraint< Real >, ROL::AlmostSureConstraint< Real >, ROL::Constraint_DynamicState< Real >, ROL::ROL::Constraint_State< Real >, and ROL::Constraint_State< Real >.

Definition at line 269 of file ROL_Constraint.hpp.

References ROL::Vector< Real >::dual(), and ROL::Vector< Real >::set().

Referenced by ROL::StdConstraint< Real >::applyPreconditioner(), and ROL::Constraint_SimOpt< Real >::applyPreconditioner().

template<class Real>
void ROL::Constraint< Real >::activate ( void  )
inline

Turn on constraints.

Definition at line 279 of file ROL_Constraint.hpp.

References ROL::Constraint< Real >::activated_.

template<class Real>
void ROL::Constraint< Real >::deactivate ( void  )
inline

Turn off constraints.

Definition at line 283 of file ROL_Constraint.hpp.

References ROL::Constraint< Real >::activated_.

template<class Real>
bool ROL::Constraint< Real >::isActivated ( void  )
inline

Check if constraints are on.

Definition at line 287 of file ROL_Constraint.hpp.

References ROL::Constraint< Real >::activated_.

template<class Real >
std::vector< std::vector< Real > > ROL::Constraint< 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 
)
virtual

Finite-difference check for the constraint Jacobian application.

Details here.

Definition at line 383 of file ROL_ConstraintDef.hpp.

References applyJacobian(), ROL::Vector< Real >::clone(), ROL::Finite_Difference_Arrays::shifts, ROL::update(), ROL::value, and ROL::Finite_Difference_Arrays::weights.

Referenced by main().

template<class Real >
std::vector< std::vector< Real > > ROL::Constraint< 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 
)
virtual

Finite-difference check for the constraint Jacobian application.

Details here.

Definition at line 364 of file ROL_ConstraintDef.hpp.

template<class Real >
std::vector< std::vector< Real > > ROL::Constraint< 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 
)
virtual

Finite-difference check for the application of the adjoint of constraint Jacobian.

Details here. (This function should be deprecated)

Definition at line 488 of file ROL_ConstraintDef.hpp.

References ROL::Vector< Real >::basis(), ROL::Vector< Real >::clone(), ROL::Vector< Real >::dimension(), ROL::Vector< Real >::dual(), ROL::update(), and ROL::value.

Referenced by main().

template<class Real>
virtual Real ROL::Constraint< Real >::checkAdjointConsistencyJacobian ( const Vector< Real > &  w,
const Vector< Real > &  v,
const Vector< Real > &  x,
const bool  printToStream = true,
std::ostream &  outStream = std::cout 
)
inlinevirtual

Definition at line 340 of file ROL_Constraint.hpp.

References ROL::Vector< Real >::dual().

Referenced by main().

template<class Real >
Real ROL::Constraint< 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
template<class Real >
std::vector< std::vector< Real > > ROL::Constraint< 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 
)
virtual

Finite-difference check for the application of the adjoint of constraint Hessian.

Details here.

Definition at line 635 of file ROL_ConstraintDef.hpp.

References ROL::Vector< Real >::clone(), ROL::Finite_Difference_Arrays::shifts, ROL::update(), and ROL::Finite_Difference_Arrays::weights.

Referenced by main().

template<class Real >
std::vector< std::vector< Real > > ROL::Constraint< 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 
)
virtual

Finite-difference check for the application of the adjoint of constraint Hessian.

Details here.

Definition at line 617 of file ROL_ConstraintDef.hpp.

template<class Real>
const std::vector<Real> ROL::Constraint< Real >::getParameter ( void  ) const
inlineprotected
template<class Real>
virtual void ROL::Constraint< Real >::setParameter ( const std::vector< Real > &  param)
inlinevirtual

Member Data Documentation

template<class Real>
bool ROL::Constraint< Real >::activated_
private
template<class Real>
std::vector<Real> ROL::Constraint< Real >::param_
private

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