44 #ifndef ROL_CONSTRAINTFROMOBJECTIVE_DEF_H
45 #define ROL_CONSTRAINTFROMOBJECTIVE_DEF_H
49 template<
typename Real>
51 obj_(obj), dualVector_(nullPtr), offset_(offset), isDualInitialized_(false) {}
53 template<
typename Real>
56 template<
typename Real>
58 obj_->setParameter(param);
62 template<
typename Real>
64 obj_->update(x,type,iter);
67 template<
typename Real>
69 obj_->update(x,flag,iter);
72 template<
typename Real>
74 setValue(c,
obj_->value(x,tol) - offset_ );
77 template<
typename Real>
79 if ( !isDualInitialized_ ) {
81 isDualInitialized_ =
true;
83 obj_->gradient(*dualVector_,x,tol);
85 setValue(jv,v.
apply(*dualVector_));
88 template<
typename Real>
90 obj_->gradient(ajv,x,tol);
91 ajv.
scale(getValue(v));
94 template<
typename Real>
96 obj_->hessVec(ahuv,v,x,tol);
97 ahuv.
scale(getValue(u));
100 template<
typename Real>
105 template<
typename Real>
112 #endif // ROL_CONSTRAINTFROMOBJECTIVE_H
Provides the interface to evaluate objective functions.
virtual const Vector & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
void applyJacobian(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &x, Real &tol) override
Apply the constraint Jacobian at , , to vector .
virtual void scale(const Real alpha)=0
Compute where .
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
virtual Real apply(const Vector< Real > &x) const
Apply to a dual vector. This is equivalent to the call .
void setValue(Vector< Real > &x, Real val)
const Ptr< Objective< Real > > getObjective(void) const
void value(Vector< Real > &c, const Vector< Real > &x, Real &tol) override
Evaluate the constraint operator at .
Defines the linear algebra or vector space interface.
Real getValue(const Vector< Real > &x)
void update(const Vector< Real > &x, UpdateType type, int iter=-1) override
Update constraint function.
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 to vector in direction ...
virtual void setParameter(const std::vector< Real > ¶m)
ConstraintFromObjective(const Ptr< Objective< Real >> &obj, const Real offset=0)
void setParameter(const std::vector< Real > ¶m) override
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 , , to vector .