ROL
ROL_ConstraintFromObjective_Def.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Rapid Optimization Library (ROL) Package
4//
5// Copyright 2014 NTESS and the ROL contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef ROL_CONSTRAINTFROMOBJECTIVE_DEF_H
11#define ROL_CONSTRAINTFROMOBJECTIVE_DEF_H
12
13namespace ROL {
14
15template<typename Real>
17 obj_(obj), dualVector_(nullPtr), offset_(offset), isDualInitialized_(false) {}
18
19template<typename Real>
20const Ptr<Objective<Real>> ConstraintFromObjective<Real>::getObjective(void) const { return obj_; }
21
22template<typename Real>
23void ConstraintFromObjective<Real>::setParameter( const std::vector<Real> &param ) {
24 obj_->setParameter(param);
26}
27
28template<typename Real>
30 obj_->update(x,type,iter);
31}
32
33template<typename Real>
34void ConstraintFromObjective<Real>::update( const Vector<Real>& x, bool flag, int iter ) {
35 obj_->update(x,flag,iter);
36}
37
38template<typename Real>
40 setValue(c, obj_->value(x,tol) - offset_ );
41}
42
43template<typename Real>
45 if ( !isDualInitialized_ ) {
46 dualVector_ = x.dual().clone();
47 isDualInitialized_ = true;
48 }
49 obj_->gradient(*dualVector_,x,tol);
50 //setValue(jv,v.dot(dualVector_->dual()));
51 setValue(jv,v.apply(*dualVector_));
52}
53
54template<typename Real>
56 obj_->gradient(ajv,x,tol);
57 ajv.scale(getValue(v));
58}
59
60template<typename Real>
62 obj_->hessVec(ahuv,v,x,tol);
63 ahuv.scale(getValue(u));
64}
65
66template<typename Real>
68 return dynamic_cast<const SingletonVector<Real>&>(x).getValue();
69}
70
71template<typename Real>
73 dynamic_cast<SingletonVector<Real>&>(x).setValue(val);
74}
75
76} // namespace ROL
77
78#endif // ROL_CONSTRAINTFROMOBJECTIVE_H
const Ptr< Obj > obj_
void setValue(Vector< Real > &x, Real val)
void update(const Vector< Real > &x, UpdateType type, int iter=-1) override
Update constraint function.
void value(Vector< Real > &c, const Vector< Real > &x, Real &tol) override
Evaluate the constraint operator at .
ConstraintFromObjective(const Ptr< Objective< Real > > &obj, const Real offset=0)
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 ,...
void setParameter(const std::vector< Real > &param) override
const Ptr< Objective< Real > > getObjective(void) const
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 .
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 setParameter(const std::vector< Real > &param)
Provides the interface to evaluate objective functions.
Defines the linear algebra or vector space interface.
virtual Real apply(const Vector< Real > &x) const
Apply to a dual vector. This is equivalent to the call .
virtual void scale(const Real alpha)=0
Compute where .
virtual const Vector & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis,...
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.