ROL
|
Implements an equality constraint function that evaluates to zero on the surface of a bounded parallelpiped and is positive in the interior. More...
#include <ROL_BinaryConstraint.hpp>
Classes | |
class | BoundsCheck |
Public Member Functions | |
BinaryConstraint (const ROL::Ptr< const Vector< Real > > &lo, const ROL::Ptr< const Vector< Real > > &up, Real gamma) | |
BinaryConstraint (const BoundConstraint< Real > &bnd, Real gamma) | |
BinaryConstraint (const ROL::Ptr< const BoundConstraint< Real > > &bnd, Real gamma) | |
void | value (Vector< Real > &c, const Vector< Real > &x, Real &tol) override |
Evaluate the constraint operator \(c:\mathcal{X} \rightarrow \mathcal{C}\) at \(x\). | |
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\). | |
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\). | |
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 \). | |
void | setPenalty (Real gamma) |
![]() | |
virtual | ~Constraint (void) |
Constraint (void) | |
virtual void | update (const Vector< Real > &x, UpdateType type, int iter=-1) |
Update constraint function. | |
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. | |
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\). | |
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 | |
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: | |
void | activate (void) |
Turn on constraints. | |
void | deactivate (void) |
Turn off constraints. | |
bool | isActivated (void) |
Check if constraints are on. | |
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. | |
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. | |
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. | |
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. | |
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. | |
virtual void | setParameter (const std::vector< Real > ¶m) |
Private Attributes | |
const Ptr< const Vector< Real > > | lo_ |
const Ptr< const Vector< Real > > | up_ |
Ptr< Vector< Real > > | d_ |
Real | gamma_ |
Additional Inherited Members | |
![]() | |
const std::vector< Real > | getParameter (void) const |
Implements an equality constraint function that evaluates to zero on the surface of a bounded parallelpiped and is positive in the interior.
Definition at line 61 of file ROL_BinaryConstraint.hpp.
ROL::BinaryConstraint< Real >::BinaryConstraint | ( | const ROL::Ptr< const Vector< Real > > & | lo, |
const ROL::Ptr< const Vector< Real > > & | up, | ||
Real | gamma ) |
Definition at line 50 of file ROL_BinaryConstraint_Def.hpp.
ROL::BinaryConstraint< Real >::BinaryConstraint | ( | const BoundConstraint< Real > & | bnd, |
Real | gamma ) |
Definition at line 55 of file ROL_BinaryConstraint_Def.hpp.
ROL::BinaryConstraint< Real >::BinaryConstraint | ( | const ROL::Ptr< const BoundConstraint< Real > > & | bnd, |
Real | gamma ) |
Definition at line 59 of file ROL_BinaryConstraint_Def.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 63 of file ROL_BinaryConstraint_Def.hpp.
References ROL::Vector< Real >::applyBinary(), ROL::Vector< Real >::axpy(), ROL::Vector< Real >::scale(), and ROL::Vector< Real >::set().
|
overridevirtual |
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::Constraint< Real >.
Definition at line 74 of file ROL_BinaryConstraint_Def.hpp.
References ROL::Vector< Real >::applyBinary(), ROL::Vector< Real >::axpy(), ROL::Vector< Real >::scale(), and ROL::Vector< Real >::set().
|
overridevirtual |
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::Constraint< Real >.
Definition at line 86 of file ROL_BinaryConstraint_Def.hpp.
References applyJacobian().
|
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 \).
[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::Constraint< Real >.
Definition at line 91 of file ROL_BinaryConstraint_Def.hpp.
References ROL::Vector< Real >::applyBinary(), ROL::Vector< Real >::axpy(), ROL::Vector< Real >::scale(), and ROL::Vector< Real >::set().
void ROL::BinaryConstraint< Real >::setPenalty | ( | Real | gamma | ) |
Definition at line 104 of file ROL_BinaryConstraint_Def.hpp.
|
private |
Definition at line 63 of file ROL_BinaryConstraint.hpp.
|
private |
Definition at line 64 of file ROL_BinaryConstraint.hpp.
|
private |
Definition at line 65 of file ROL_BinaryConstraint.hpp.
|
private |
Definition at line 66 of file ROL_BinaryConstraint.hpp.