44 #ifndef ROL_QUADRATICPENALTY_SIMOPT_H 45 #define ROL_QUADRATICPENALTY_SIMOPT_H 51 #include "Teuchos_RCP.hpp" 94 const Teuchos::RCP<EqualityConstraint_SimOpt<Real> >
con_;
118 if ( !isConstraintComputed_ ) {
120 con_->value(*conValue_,u,z,tol); ncval_++;
121 isConstraintComputed_ =
true;
128 const Real penaltyParameter,
132 const bool useScaling =
false,
133 const int HessianApprox = 0 )
134 : con_(con), penaltyParameter_(penaltyParameter), ncval_(0),
135 useScaling_(useScaling), HessianApprox_(HessianApprox), isConstraintComputed_(false) {
137 dualSimVector_ = simVec.
dual().clone();
138 dualOptVector_ = optVec.
dual().clone();
139 primalConVector_ = conVec.
clone();
140 conValue_ = conVec.
clone();
141 multiplier_ = multiplier.
clone();
142 primalMultiplierVector_ = multiplier.
clone();
146 con_->update(u,z,flag,iter);
154 Real cval = multiplier_->dot(conValue_->dual());
156 Real pval = conValue_->dot(*conValue_);
158 const Real half(0.5);
161 val = cval/penaltyParameter_ + half*pval;
164 val = cval + half*penaltyParameter_*pval;
173 primalMultiplierVector_->set(conValue_->dual());
175 primalMultiplierVector_->axpy(static_cast<Real>(1)/penaltyParameter_,*multiplier_);
178 primalMultiplierVector_->scale(penaltyParameter_);
179 primalMultiplierVector_->plus(*multiplier_);
181 con_->applyAdjointJacobian_1(g,*primalMultiplierVector_,u,z,tol);
188 primalMultiplierVector_->set(conValue_->dual());
190 primalMultiplierVector_->axpy(static_cast<Real>(1)/penaltyParameter_,*multiplier_);
193 primalMultiplierVector_->scale(penaltyParameter_);
194 primalMultiplierVector_->plus(*multiplier_);
196 con_->applyAdjointJacobian_2(g,*primalMultiplierVector_,u,z,tol);
202 if (HessianApprox_ < 2) {
203 con_->applyJacobian_1(*primalConVector_,v,u,z,tol);
204 con_->applyAdjointJacobian_1(hv,primalConVector_->
dual(),u,z,tol);
206 hv.
scale(penaltyParameter_);
209 if (HessianApprox_ == 0) {
213 primalMultiplierVector_->set(conValue_->dual());
215 primalMultiplierVector_->axpy(static_cast<Real>(1)/penaltyParameter_,*multiplier_);
218 primalMultiplierVector_->scale(penaltyParameter_);
219 primalMultiplierVector_->plus(*multiplier_);
221 con_->applyAdjointHessian_11(*dualSimVector_,*primalMultiplierVector_,v,u,z,tol);
222 hv.
plus(*dualSimVector_);
233 if (HessianApprox_ < 2) {
234 con_->applyJacobian_2(*primalConVector_,v,u,z,tol);
235 con_->applyAdjointJacobian_1(hv,primalConVector_->
dual(),u,z,tol);
237 hv.
scale(penaltyParameter_);
240 if (HessianApprox_ == 0) {
244 primalMultiplierVector_->set(conValue_->dual());
246 primalMultiplierVector_->axpy(static_cast<Real>(1)/penaltyParameter_,*multiplier_);
249 primalMultiplierVector_->scale(penaltyParameter_);
250 primalMultiplierVector_->plus(*multiplier_);
252 con_->applyAdjointHessian_21(*dualSimVector_,*primalMultiplierVector_,v,u,z,tol);
253 hv.
plus(*dualSimVector_);
264 if (HessianApprox_ < 2) {
265 con_->applyJacobian_1(*primalConVector_,v,u,z,tol);
266 con_->applyAdjointJacobian_2(hv,primalConVector_->
dual(),u,z,tol);
268 hv.
scale(penaltyParameter_);
271 if (HessianApprox_ == 0) {
275 primalMultiplierVector_->set(conValue_->dual());
277 primalMultiplierVector_->axpy(static_cast<Real>(1)/penaltyParameter_,*multiplier_);
280 primalMultiplierVector_->scale(penaltyParameter_);
281 primalMultiplierVector_->plus(*multiplier_);
283 con_->applyAdjointHessian_12(*dualOptVector_,*primalMultiplierVector_,v,u,z,tol);
284 hv.
plus(*dualOptVector_);
295 if (HessianApprox_ < 2) {
296 con_->applyJacobian_2(*primalConVector_,v,u,z,tol);
297 con_->applyAdjointJacobian_2(hv,primalConVector_->
dual(),u,z,tol);
299 hv.
scale(penaltyParameter_);
302 if (HessianApprox_ == 0) {
306 primalMultiplierVector_->set(conValue_->dual());
308 primalMultiplierVector_->axpy(static_cast<Real>(1)/penaltyParameter_,*multiplier_);
311 primalMultiplierVector_->scale(penaltyParameter_);
312 primalMultiplierVector_->plus(*multiplier_);
314 con_->applyAdjointHessian_22(*dualOptVector_,*primalMultiplierVector_,v,u,z,tol);
315 hv.
plus(*dualOptVector_);
325 Real tol = std::sqrt(ROL_EPSILON<Real>());
339 multiplier_->set(multiplier);
340 penaltyParameter_ = penaltyParameter;
Provides the interface to evaluate simulation-based objective functions.
virtual void scale(const Real alpha)=0
Compute where .
Teuchos::RCP< Vector< Real > > dualOptVector_
virtual void hessVec_21(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
virtual void plus(const Vector &x)=0
Compute , where .
QuadraticPenalty_SimOpt(const Teuchos::RCP< EqualityConstraint_SimOpt< Real > > &con, const Vector< Real > &multiplier, const Real penaltyParameter, const Vector< Real > &simVec, const Vector< Real > &optVec, const Vector< Real > &conVec, const bool useScaling=false, const int HessianApprox=0)
Teuchos::RCP< Vector< Real > > primalConVector_
Contains definitions of custom data types in ROL.
Teuchos::RCP< Vector< Real > > primalMultiplierVector_
virtual void hessVec_22(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
bool isConstraintComputed_
Provides the interface to evaluate the quadratic SimOpt constraint penalty.
virtual Teuchos::RCP< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
virtual void zero()
Set to zero vector.
Defines the linear algebra or vector space interface.
const Teuchos::RCP< EqualityConstraint_SimOpt< Real > > con_
Defines the equality constraint operator interface for simulation-based optimization.
virtual void hessVec_11(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Apply Hessian approximation to vector.
virtual int getNumberConstraintEvaluations(void) const
virtual const Vector & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
virtual Real value(const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Compute value.
void evaluateConstraint(const Vector< Real > &u, const Vector< Real > &z, Real &tol)
virtual void hessVec_12(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Teuchos::RCP< Vector< Real > > conValue_
virtual void getConstraintVec(Vector< Real > &c, const Vector< Real > &u, const Vector< Real > &z)
virtual void reset(const Vector< Real > &multiplier, const Real penaltyParameter)
virtual void set(const Vector &x)
Set where .
virtual void update(const Vector< Real > &u, const Vector< Real > &z, bool flag=true, int iter=-1)
Update objective function. u is an iterate, z is an iterate, flag = true if the iterate has changed...
Teuchos::RCP< Vector< Real > > multiplier_
virtual void gradient_1(Vector< Real > &g, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Compute gradient with respect to first component.
Teuchos::RCP< Vector< Real > > dualSimVector_
virtual void gradient_2(Vector< Real > &g, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Compute gradient with respect to second component.