ROL
ROL_RiskNeutralConstraint.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ************************************************************************
3 //
4 // Rapid Optimization Library (ROL) Package
5 // Copyright (2014) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact lead developers:
38 // Drew Kouri (dpkouri@sandia.gov) and
39 // Denis Ridzal (dridzal@sandia.gov)
40 //
41 // ************************************************************************
42 // @HEADER
43 
44 #ifndef ROL_RISKNEUTRALCONSTRAINT_HPP
45 #define ROL_RISKNEUTRALCONSTRAINT_HPP
46 
47 #include "ROL_Ptr.hpp"
48 #include "ROL_Vector.hpp"
49 #include "ROL_Constraint.hpp"
50 #include "ROL_SampleGenerator.hpp"
51 
52 namespace ROL {
53 
54 template<class Real>
55 class RiskNeutralConstraint : public Constraint<Real> {
56 private:
57  const Ptr<Constraint<Real>> con_;
58  const Ptr<SampleGenerator<Real>> xsampler_;
59  const Ptr<BatchManager<Real>> cbman_;
60 
61  Ptr<Vector<Real>> conVec_;
62  Ptr<Vector<Real>> optVec_;
63 
65 
66  void init(const Vector<Real> &c, const Vector<Real> &x) {
67  if (!initialized_) {
68  conVec_ = c.clone();
69  optVec_ = x.dual().clone();
70  initialized_ = true;
71  }
72  }
73 
74 public:
76  const Ptr<SampleGenerator<Real>> &xsampler,
77  const Ptr<BatchManager<Real>> &cbman)
78  : con_(con), xsampler_(xsampler), cbman_(cbman), initialized_(false) {}
79 
80  void update( const Vector<Real> &x, bool flag = true, int iter = -1 ) {
81  con_->update(x,flag,iter);
82  }
83 
84  void update( const Vector<Real> &x, UpdateType type, int iter = -1 ) {
85  con_->update(x,type,iter);
86  }
87 
88  void value(Vector<Real> &c, const Vector<Real> &x, Real &tol ) {
89  init(c,x);
90  conVec_->zero();
91  for ( int i = 0; i < xsampler_->numMySamples(); ++i ) {
92  con_->setParameter(xsampler_->getMyPoint(i));
93  con_->value(c,x,tol);
94  conVec_->axpy(xsampler_->getMyWeight(i),c);
95  }
96  c.zero();
97  cbman_->sumAll(*conVec_,c);
98  }
99 
100  void applyJacobian(Vector<Real> &jv, const Vector<Real> &v, const Vector<Real> &x, Real &tol) {
101  init(jv,x);
102  conVec_->zero();
103  for ( int i = 0; i < xsampler_->numMySamples(); ++i ) {
104  con_->setParameter(xsampler_->getMyPoint(i));
105  con_->applyJacobian(jv,v,x,tol);
106  conVec_->axpy(xsampler_->getMyWeight(i),jv);
107  }
108  jv.zero();
109  cbman_->sumAll(*conVec_,jv);
110  }
111 
112  void applyAdjointJacobian(Vector<Real> &ajv, const Vector<Real> &v, const Vector<Real> &x, Real &tol) {
113  init(v.dual(),x);
114  optVec_->zero();
115  for ( int i = 0; i < xsampler_->numMySamples(); ++i ) {
116  con_->setParameter(xsampler_->getMyPoint(i));
117  con_->applyAdjointJacobian(ajv,v,x,tol);
118  optVec_->axpy(xsampler_->getMyWeight(i),ajv);
119  }
120  ajv.zero();
121  xsampler_->sumAll(*optVec_,ajv);
122  }
123 
124  void applyAdjointHessian(Vector<Real> &ahuv, const Vector<Real> &u, const Vector<Real> &v, const Vector<Real> &x, Real &tol) {
125  init(u.dual(),x);
126  optVec_->zero();
127  for ( int i = 0; i < xsampler_->numMySamples(); ++i ) {
128  con_->setParameter(xsampler_->getMyPoint(i));
129  con_->applyAdjointHessian(ahuv,u,v,x,tol);
130  optVec_->axpy(xsampler_->getMyWeight(i),ahuv);
131  }
132  ahuv.zero();
133  xsampler_->sumAll(*optVec_,ahuv);
134  }
135 
136 };
137 
138 }
139 
140 #endif
RiskNeutralConstraint(const Ptr< Constraint< Real >> &con, const Ptr< SampleGenerator< Real >> &xsampler, const Ptr< BatchManager< Real >> &cbman)
virtual const Vector & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
Definition: ROL_Vector.hpp:226
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
void init(const Vector< Real > &c, const Vector< Real > &x)
virtual void zero()
Set to zero vector.
Definition: ROL_Vector.hpp:167
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:80
void applyJacobian(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply the constraint Jacobian at , , to vector .
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_
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 ...
void value(Vector< Real > &c, const Vector< Real > &x, Real &tol)
Evaluate the constraint operator at .
const Ptr< Constraint< Real > > con_
const Ptr< BatchManager< Real > > cbman_
Defines the general constraint operator interface.
void update(const Vector< Real > &x, UpdateType type, int iter=-1)
Update constraint function.