ROL
ROL_RiskNeutralConstraint.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_RISKNEUTRALCONSTRAINT_HPP
11#define ROL_RISKNEUTRALCONSTRAINT_HPP
12
13#include "ROL_Ptr.hpp"
14#include "ROL_Vector.hpp"
15#include "ROL_Constraint.hpp"
17
18namespace ROL {
19
20template<class Real>
21class RiskNeutralConstraint : public Constraint<Real> {
22private:
23 const Ptr<Constraint<Real>> con_;
24 const Ptr<SampleGenerator<Real>> xsampler_;
25 const Ptr<BatchManager<Real>> cbman_;
26
27 Ptr<Vector<Real>> conVec_;
28 Ptr<Vector<Real>> optVec_;
29
31
32 void init(const Vector<Real> &c, const Vector<Real> &x) {
33 if (!initialized_) {
34 conVec_ = c.clone();
35 optVec_ = x.dual().clone();
36 initialized_ = true;
37 }
38 }
39
40public:
42 const Ptr<SampleGenerator<Real>> &xsampler,
43 const Ptr<BatchManager<Real>> &cbman)
44 : con_(con), xsampler_(xsampler), cbman_(cbman), initialized_(false) {}
45
46 void update( const Vector<Real> &x, bool flag = true, int iter = -1 ) {
47 con_->update(x,flag,iter);
48 }
49
50 void update( const Vector<Real> &x, UpdateType type, int iter = -1 ) {
51 con_->update(x,type,iter);
52 }
53
54 void value(Vector<Real> &c, const Vector<Real> &x, Real &tol ) {
55 init(c,x);
56 conVec_->zero();
57 for ( int i = 0; i < xsampler_->numMySamples(); ++i ) {
58 con_->setParameter(xsampler_->getMyPoint(i));
59 con_->value(c,x,tol);
60 conVec_->axpy(xsampler_->getMyWeight(i),c);
61 }
62 c.zero();
63 cbman_->sumAll(*conVec_,c);
64 }
65
66 void applyJacobian(Vector<Real> &jv, const Vector<Real> &v, const Vector<Real> &x, Real &tol) {
67 init(jv,x);
68 conVec_->zero();
69 for ( int i = 0; i < xsampler_->numMySamples(); ++i ) {
70 con_->setParameter(xsampler_->getMyPoint(i));
71 con_->applyJacobian(jv,v,x,tol);
72 conVec_->axpy(xsampler_->getMyWeight(i),jv);
73 }
74 jv.zero();
75 cbman_->sumAll(*conVec_,jv);
76 }
77
78 void applyAdjointJacobian(Vector<Real> &ajv, const Vector<Real> &v, const Vector<Real> &x, Real &tol) {
79 init(v.dual(),x);
80 optVec_->zero();
81 for ( int i = 0; i < xsampler_->numMySamples(); ++i ) {
82 con_->setParameter(xsampler_->getMyPoint(i));
83 con_->applyAdjointJacobian(ajv,v,x,tol);
84 optVec_->axpy(xsampler_->getMyWeight(i),ajv);
85 }
86 ajv.zero();
87 xsampler_->sumAll(*optVec_,ajv);
88 }
89
90 void applyAdjointHessian(Vector<Real> &ahuv, const Vector<Real> &u, const Vector<Real> &v, const Vector<Real> &x, Real &tol) {
91 init(u.dual(),x);
92 optVec_->zero();
93 for ( int i = 0; i < xsampler_->numMySamples(); ++i ) {
94 con_->setParameter(xsampler_->getMyPoint(i));
95 con_->applyAdjointHessian(ahuv,u,v,x,tol);
96 optVec_->axpy(xsampler_->getMyWeight(i),ahuv);
97 }
98 ahuv.zero();
99 xsampler_->sumAll(*optVec_,ahuv);
100 }
101
102};
103
104}
105
106#endif
Defines the general constraint operator interface.
RiskNeutralConstraint(const Ptr< Constraint< Real > > &con, const Ptr< SampleGenerator< Real > > &xsampler, const Ptr< BatchManager< Real > > &cbman)
const Ptr< Constraint< Real > > con_
void value(Vector< Real > &c, const Vector< Real > &x, Real &tol)
Evaluate the constraint operator at .
void init(const Vector< Real > &c, const Vector< Real > &x)
const Ptr< BatchManager< Real > > cbman_
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 i...
void applyJacobian(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply the constraint Jacobian at , , to vector .
void update(const Vector< Real > &x, UpdateType type, int iter=-1)
Update constraint function.
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 .
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 ,...
const Ptr< SampleGenerator< Real > > xsampler_
Defines the linear algebra or vector space interface.
virtual const Vector & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis,...
virtual void zero()
Set to zero vector.
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.