ROL
ROL_SimulatedConstraint.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_SIMULATED_CONSTRAINT_H
45 #define ROL_SIMULATED_CONSTRAINT_H
46 
47 #include "ROL_SimulatedVector.hpp"
48 #include "ROL_RiskVector.hpp"
50 
51 namespace ROL {
52 
53 template <class Real>
54 class SimulatedConstraint : public Constraint<Real> {
55 private:
56  const ROL::Ptr<SampleGenerator<Real> > sampler_;
57  const ROL::Ptr<Constraint_SimOpt<Real> > pcon_;
58  const bool useWeights_;
59 
60 public:
61 
62  virtual ~SimulatedConstraint() {}
63 
64  SimulatedConstraint(const ROL::Ptr<SampleGenerator<Real> > & sampler,
65  const ROL::Ptr<Constraint_SimOpt<Real> > & pcon,
66  const bool useWeights = true)
67  : sampler_(sampler), pcon_(pcon), useWeights_(useWeights) {}
68 
69  void update( const Vector<Real> &x, bool flag = true, int iter = -1 ) {
70  pcon_->update(x,flag,iter);
71  }
72  void update( const Vector<Real> &x, UpdateType type, int iter = -1 ) {
73  pcon_->update(x,type,iter);
74  }
75 
77  const Vector<Real> &x,
78  Real &tol) {
79  c.zero();
80  SimulatedVector<Real> &pc = dynamic_cast<SimulatedVector<Real>&>(c);
81  const Vector_SimOpt<Real> &uz = dynamic_cast<const Vector_SimOpt<Real>&>(x);
82  ROL::Ptr<const Vector<Real> > uptr = uz.get_1();
83  ROL::Ptr<const Vector<Real> > zptr = uz.get_2();
84  try {
85  const RiskVector<Real> &rz = dynamic_cast<const RiskVector<Real>&>(*zptr);
86  zptr = rz.getVector();
87  }
88  catch (const std::bad_cast &e) {}
89  const SimulatedVector<Real> &pu = dynamic_cast<const SimulatedVector<Real>&>(*uptr);
90 
91  std::vector<Real> param;
92  Real weight(0), one(1);
93  for (typename std::vector<SimulatedVector<Real> >::size_type i=0; i<pu.numVectors(); ++i) {
94  param = sampler_->getMyPoint(static_cast<int>(i));
95  weight = sampler_->getMyWeight(static_cast<int>(i));
96  pcon_->setParameter(param);
97  pcon_->update(*(pu.get(i)), *zptr);
98  pcon_->value(*(pc.get(i)), *(pu.get(i)), *zptr, tol);
99  weight = (useWeights_) ? weight : one;
100  pc.get(i)->scale(weight);
101  }
102 
103  }
104 
105 
106  virtual void applyJacobian(Vector<Real> &jv,
107  const Vector<Real> &v,
108  const Vector<Real> &x,
109  Real &tol) {
110  jv.zero();
111  // cast jv
112  SimulatedVector<Real> &pjv = dynamic_cast<SimulatedVector<Real>&>(jv);
113  // split x
114  const Vector_SimOpt<Real> &xuz = dynamic_cast<const Vector_SimOpt<Real>&>(x);
115  ROL::Ptr<const Vector<Real> > xuptr = xuz.get_1();
116  ROL::Ptr<const Vector<Real> > xzptr = xuz.get_2();
117  try {
118  const RiskVector<Real> &rxz = dynamic_cast<const RiskVector<Real>&>(*xzptr);
119  xzptr = rxz.getVector();
120  }
121  catch (const std::bad_cast &e) {}
122  const SimulatedVector<Real> &pxu = dynamic_cast<const SimulatedVector<Real>&>(*xuptr);
123  // split v
124  const Vector_SimOpt<Real> &vuz = dynamic_cast<const Vector_SimOpt<Real>&>(v);
125  ROL::Ptr<const Vector<Real> > vuptr = vuz.get_1();
126  ROL::Ptr<const Vector<Real> > vzptr = vuz.get_2();
127  try {
128  const RiskVector<Real> &rvz = dynamic_cast<const RiskVector<Real>&>(*vzptr);
129  vzptr = rvz.getVector();
130  }
131  catch (const std::bad_cast &e) {}
132  const SimulatedVector<Real> &pvu = dynamic_cast<const SimulatedVector<Real>&>(*vuptr);
133 
134  std::vector<Real> param;
135  Real weight(0), one(1);
136  for (typename std::vector<SimulatedVector<Real> >::size_type i=0; i<pvu.numVectors(); ++i) {
137  param = sampler_->getMyPoint(static_cast<int>(i));
138  weight = sampler_->getMyWeight(static_cast<int>(i));
139  pcon_->setParameter(param);
140  Vector_SimOpt<Real> xi(ROL::constPtrCast<Vector<Real> >(pxu.get(i)), ROL::constPtrCast<Vector<Real> >(xzptr));
141  Vector_SimOpt<Real> vi(ROL::constPtrCast<Vector<Real> >(pvu.get(i)), ROL::constPtrCast<Vector<Real> >(vzptr));
142  pcon_->update(xi);
143  pcon_->applyJacobian(*(pjv.get(i)), vi, xi, tol);
144  weight = (useWeights_) ? weight : one;
145  pjv.get(i)->scale(weight);
146  }
147  }
148 
149 
151  const Vector<Real> &v,
152  const Vector<Real> &x,
153  Real &tol) {
154  ajv.zero();
155  // split ajv
156  Vector_SimOpt<Real> &ajvuz = dynamic_cast<Vector_SimOpt<Real>&>(ajv);
157  ROL::Ptr<Vector<Real> > ajvuptr = ajvuz.get_1();
158  ROL::Ptr<Vector<Real> > ajvzptr = ajvuz.get_2();
159  try {
160  RiskVector<Real> &rajvz = dynamic_cast<RiskVector<Real>&>(*ajvzptr);
161  ajvzptr = rajvz.getVector();
162  }
163  catch (const std::bad_cast &e) {}
164  SimulatedVector<Real> &pajvu = dynamic_cast<SimulatedVector<Real>&>(*ajvuptr);
165  // split x
166  const Vector_SimOpt<Real> &xuz = dynamic_cast<const Vector_SimOpt<Real>&>(x);
167  ROL::Ptr<const Vector<Real> > xuptr = xuz.get_1();
168  ROL::Ptr<const Vector<Real> > xzptr = xuz.get_2();
169  try {
170  const RiskVector<Real> &rxz = dynamic_cast<const RiskVector<Real>&>(*xzptr);
171  xzptr = rxz.getVector();
172  }
173  catch (const std::bad_cast &e) {}
174  const SimulatedVector<Real> &pxu = dynamic_cast<const SimulatedVector<Real>&>(*xuptr);
175  // cast v
176  const SimulatedVector<Real> &pv = dynamic_cast<const SimulatedVector<Real>&>(v);
177 
178  std::vector<Real> param;
179  Real weight(0), one(1);
180  ROL::Ptr<Vector<Real> > tmp1 = ajvzptr->clone();
181  ROL::Ptr<Vector<Real> > tmp2 = ajvzptr->clone();
182  for (typename std::vector<SimulatedVector<Real> >::size_type i=0; i<pv.numVectors(); ++i) {
183  param = sampler_->getMyPoint(static_cast<int>(i));
184  weight = sampler_->getMyWeight(static_cast<int>(i));
185  pcon_->setParameter(param);
186  Vector_SimOpt<Real> xi(ROL::constPtrCast<Vector<Real> >(pxu.get(i)), ROL::constPtrCast<Vector<Real> >(xzptr));
187  Vector_SimOpt<Real> ajvi(pajvu.get(i), tmp1);
188  pcon_->update(xi);
189  pcon_->applyAdjointJacobian(ajvi, *(pv.get(i)), xi, tol);
190  weight = (useWeights_) ? weight : one;
191  ajvi.scale(weight);
192  tmp2->plus(*tmp1);
193  }
194  sampler_->sumAll(*tmp2, *ajvzptr);
195 
196  }
197 
198 
199  virtual void applyAdjointHessian(Vector<Real> &ahuv,
200  const Vector<Real> &u,
201  const Vector<Real> &v,
202  const Vector<Real> &x,
203  Real &tol) {
204  ahuv.zero();
205  // split ahuv
206  Vector_SimOpt<Real> &ahuvuz = dynamic_cast<Vector_SimOpt<Real>&>(ahuv);
207  ROL::Ptr<Vector<Real> > ahuvuptr = ahuvuz.get_1();
208  ROL::Ptr<Vector<Real> > ahuvzptr = ahuvuz.get_2();
209  try {
210  RiskVector<Real> &rahuvz = dynamic_cast<RiskVector<Real>&>(*ahuvzptr);
211  ahuvzptr = rahuvz.getVector();
212  }
213  catch (const std::bad_cast &e) {}
214  SimulatedVector<Real> &pahuvu = dynamic_cast<SimulatedVector<Real>&>(*ahuvuptr);
215  // cast u
216  const SimulatedVector<Real> &pu = dynamic_cast<const SimulatedVector<Real>&>(u);
217  // split v
218  const Vector_SimOpt<Real> &vuz = dynamic_cast<const Vector_SimOpt<Real>&>(v);
219  ROL::Ptr<const Vector<Real> > vuptr = vuz.get_1();
220  ROL::Ptr<const Vector<Real> > vzptr = vuz.get_2();
221  try {
222  const RiskVector<Real> &rvz = dynamic_cast<const RiskVector<Real>&>(*vzptr);
223  vzptr = rvz.getVector();
224  }
225  catch (const std::bad_cast &e) {}
226  const SimulatedVector<Real> &pvu = dynamic_cast<const SimulatedVector<Real>&>(*vuptr);
227  // split x
228  const Vector_SimOpt<Real> &xuz = dynamic_cast<const Vector_SimOpt<Real>&>(x);
229  ROL::Ptr<const Vector<Real> > xuptr = xuz.get_1();
230  ROL::Ptr<const Vector<Real> > xzptr = xuz.get_2();
231  try {
232  const RiskVector<Real> &rxz = dynamic_cast<const RiskVector<Real>&>(*xzptr);
233  xzptr = rxz.getVector();
234  }
235  catch (const std::bad_cast &e) {}
236  const SimulatedVector<Real> &pxu = dynamic_cast<const SimulatedVector<Real>&>(*xuptr);
237 
238  std::vector<Real> param;
239  Real weight(0), one(1);
240  ROL::Ptr<Vector<Real> > tmp1 = ahuvzptr->clone();
241  ROL::Ptr<Vector<Real> > tmp2 = ahuvzptr->clone();
242  for (typename std::vector<SimulatedVector<Real> >::size_type i=0; i<pxu.numVectors(); ++i) {
243  param = sampler_->getMyPoint(static_cast<int>(i));
244  weight = sampler_->getMyWeight(static_cast<int>(i));
245  pcon_->setParameter(param);
246  Vector_SimOpt<Real> xi(ROL::constPtrCast<Vector<Real> >(pxu.get(i)), ROL::constPtrCast<Vector<Real> >(xzptr));
247  Vector_SimOpt<Real> vi(ROL::constPtrCast<Vector<Real> >(pvu.get(i)), ROL::constPtrCast<Vector<Real> >(vzptr));
248  Vector_SimOpt<Real> ahuvi(pahuvu.get(i), tmp1);
249  pcon_->update(xi);
250  pcon_->applyAdjointHessian(ahuvi, *(pu.get(i)), vi, xi, tol);
251  weight = (useWeights_) ? weight : one;
252  ahuvi.scale(weight);
253  tmp2->plus(*tmp1);
254  }
255  sampler_->sumAll(*tmp2, *ahuvzptr);
256 
257  }
258 
260  const Vector<Real> &v,
261  const Vector<Real> &x,
262  const Vector<Real> &g,
263  Real &tol) {
264  Pv.zero();
265  // cast Pv
266  SimulatedVector<Real> &ppv = dynamic_cast<SimulatedVector<Real>&>(Pv);
267  // split x
268  const Vector_SimOpt<Real> &xuz = dynamic_cast<const Vector_SimOpt<Real>&>(x);
269  ROL::Ptr<const Vector<Real> > xuptr = xuz.get_1();
270  ROL::Ptr<const Vector<Real> > xzptr = xuz.get_2();
271  try {
272  const RiskVector<Real> &rxz = dynamic_cast<const RiskVector<Real>&>(*xzptr);
273  xzptr = rxz.getVector();
274  }
275  catch (const std::bad_cast &e) {}
276  const SimulatedVector<Real> &pxu = dynamic_cast<const SimulatedVector<Real>&>(*xuptr);
277  // split g
278  const Vector_SimOpt<Real> &guz = dynamic_cast<const Vector_SimOpt<Real>&>(g);
279  ROL::Ptr<const Vector<Real> > guptr = guz.get_1();
280  ROL::Ptr<const Vector<Real> > gzptr = guz.get_2();
281  try {
282  const RiskVector<Real> &rgz = dynamic_cast<const RiskVector<Real>&>(*gzptr);
283  gzptr = rgz.getVector();
284  }
285  catch (const std::bad_cast &e) {}
286  const SimulatedVector<Real> &pgu = dynamic_cast<const SimulatedVector<Real>&>(*guptr);
287  // cast v
288  const SimulatedVector<Real> &pv = dynamic_cast<const SimulatedVector<Real>&>(v);
289 
290  std::vector<Real> param;
291  Real weight(0), one(1);
292  for (typename std::vector<SimulatedVector<Real> >::size_type i=0; i<pv.numVectors(); ++i) {
293  param = sampler_->getMyPoint(static_cast<int>(i));
294  weight = sampler_->getMyWeight(static_cast<int>(i));
295  pcon_->setParameter(param);
296  Vector_SimOpt<Real> xi(ROL::constPtrCast<Vector<Real> >(pxu.get(i)), ROL::constPtrCast<Vector<Real> >(xzptr));
297  Vector_SimOpt<Real> gi(ROL::constPtrCast<Vector<Real> >(pgu.get(i)), ROL::constPtrCast<Vector<Real> >(gzptr));
298  pcon_->update(xi);
299  pcon_->applyPreconditioner(*(ppv.get(i)), *(pv.get(i)), xi, gi, tol);
300  weight = (useWeights_) ? weight : one;
301  ppv.get(i)->scale(one/(weight*weight));
302  }
303 
304  }
305 
306 
307 }; // class SimulatedConstraint
308 
309 } // namespace ROL
310 
311 #endif
typename PV< Real >::size_type size_type
virtual void applyJacobian(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply the constraint Jacobian at , , to vector .
ROL::Ptr< const Vector< Real > > get_2() const
Defines the linear algebra or vector space interface for simulation-based optimization.
virtual 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 ROL::Ptr< Constraint_SimOpt< Real > > pcon_
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
const ROL::Ptr< SampleGenerator< Real > > sampler_
Defines the linear algebra of a vector space on a generic partitioned vector where the individual vec...
ROL::Ptr< const Vector< Real > > get(size_type i) const
virtual void applyPreconditioner(Vector< Real > &Pv, const Vector< Real > &v, const Vector< Real > &x, const Vector< Real > &g, Real &tol)
Apply a constraint preconditioner at , , to vector . Ideally, this preconditioner satisfies the follo...
size_type numVectors() const
virtual 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 .
SimulatedConstraint(const ROL::Ptr< SampleGenerator< Real > > &sampler, const ROL::Ptr< Constraint_SimOpt< Real > > &pcon, const bool useWeights=true)
void update(const Vector< Real > &x, UpdateType type, int iter=-1)
Update constraint function.
void value(Vector< Real > &c, const Vector< Real > &x, Real &tol)
Evaluate the constraint operator at .
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 ...
Ptr< const Vector< Real > > getVector(void) const
Defines the constraint operator interface for simulation-based optimization.
Defines the general constraint operator interface.
ROL::Ptr< const Vector< Real > > get_1() const