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