10 #ifndef ROL_CONSTRAINT_H
11 #define ROL_CONSTRAINT_H
276 const std::vector<Real> &steps,
277 const bool printToStream =
true,
278 std::ostream & outStream = std::cout,
279 const int order = 1 ) ;
290 const bool printToStream =
true,
291 std::ostream & outStream = std::cout,
293 const int order = 1 ) ;
304 const bool printToStream =
true,
305 std::ostream & outStream = std::cout,
323 const bool printToStream =
true,
324 std::ostream & outStream = std::cout) {
333 const bool printToStream =
true,
334 std::ostream & outStream = std::cout);
346 const std::vector<Real> &step,
347 const bool printToScreen =
true,
348 std::ostream & outStream = std::cout,
349 const int order = 1 ) ;
359 const bool printToScreen =
true,
360 std::ostream & outStream = std::cout,
362 const int order = 1 ) ;
375 param_.assign(param.begin(),param.end());
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 ...
virtual const Vector & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
virtual void applyJacobian(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply the constraint Jacobian at , , to vector .
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 , , to vector . Ideally, this preconditioner satisfies the follo...
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 void update(const Vector< Real > &x, UpdateType type, int iter=-1)
Update constraint function.
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.
Contains definitions of custom data types in ROL.
void deactivate(void)
Turn off constraints.
const std::vector< Real > getParameter(void) const
void activate(void)
Turn on constraints.
Defines the linear algebra or vector space interface.
virtual void value(Vector< Real > &c, const Vector< Real > &x, Real &tol)=0
Evaluate the constraint operator at .
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 to vector in direction ...
std::vector< Real > param_
virtual void setParameter(const std::vector< Real > ¶m)
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 > > 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.
#define ROL_NUM_CHECKDERIV_STEPS
Number of steps for derivative checks.
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 , , to vector .
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 where , , , , is an identity or Riesz operator...
bool isActivated(void)
Check if constraints are on.
virtual void set(const Vector &x)
Set where .
virtual ~Constraint(void)
Defines the general constraint operator interface.