ROL
test_19.cpp
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 
49 #include "ROL_StdVector.hpp"
52 #include "ROL_Stream.hpp"
53 #include "Teuchos_GlobalMPISession.hpp"
54 #include <cassert>
55 
56 template<typename Real>
57 class constraint1 : public ROL::Constraint_SimOpt<Real> {
58 public:
60  void value(ROL::Vector<Real> &c, const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
61  assert(c.dimension()==1);
62  assert(u.dimension()==1);
63  assert(z.dimension()==2);
64  ROL::StdVector<Real> cs = dynamic_cast<ROL::StdVector<Real>&>(c);
65  const ROL::StdVector<Real> us = dynamic_cast<const ROL::StdVector<Real>&>(u);
66  const ROL::StdVector<Real> zs = dynamic_cast<const ROL::StdVector<Real>&>(z);
67  Real u1 = (*(us.getVector()))[0];
68  Real z1 = (*(zs.getVector()))[0];
69  Real z2 = (*(zs.getVector()))[1];
70  (*(cs.getVector()))[0] = std::exp(z1*u1)-z2*z2;
71  }
72  void solve(ROL::Vector<Real> &c, ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
73  assert(c.dimension()==1);
74  assert(u.dimension()==1);
75  assert(z.dimension()==2);
76  ROL::StdVector<Real> us = dynamic_cast<ROL::StdVector<Real>&>(u);
77  const ROL::StdVector<Real> zs = dynamic_cast<const ROL::StdVector<Real>&>(z);
78  Real z1 = (*(zs.getVector()))[0];
79  Real z2 = (*(zs.getVector()))[1];
80  (*(us.getVector()))[0] = static_cast<Real>(2)*std::log(std::abs(z2)) / z1;
81  constraint1<Real>::value(c,u,z,tol);
82  }
83  void applyJacobian_1(ROL::Vector<Real> &jv, const ROL::Vector<Real> &v, const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
84  assert(jv.dimension()==1);
85  assert(v.dimension()==1);
86  assert(u.dimension()==1);
87  assert(z.dimension()==2);
88  ROL::StdVector<Real> jvs = dynamic_cast<ROL::StdVector<Real>&>(jv);
89  const ROL::StdVector<Real> vs = dynamic_cast<const ROL::StdVector<Real>&>(v);
90  const ROL::StdVector<Real> us = dynamic_cast<const ROL::StdVector<Real>&>(u);
91  const ROL::StdVector<Real> zs = dynamic_cast<const ROL::StdVector<Real>&>(z);
92  Real v1 = (*(vs.getVector()))[0];
93  Real u1 = (*(us.getVector()))[0];
94  Real z1 = (*(zs.getVector()))[0];
95  (*(jvs.getVector()))[0] = z1*std::exp(z1*u1)*v1;
96  }
97  void applyJacobian_2(ROL::Vector<Real> &jv, const ROL::Vector<Real> &v, const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
98  assert(jv.dimension()==1);
99  assert(v.dimension()==2);
100  assert(u.dimension()==1);
101  assert(z.dimension()==2);
102  ROL::StdVector<Real> jvs = dynamic_cast<ROL::StdVector<Real>&>(jv);
103  const ROL::StdVector<Real> vs = dynamic_cast<const ROL::StdVector<Real>&>(v);
104  const ROL::StdVector<Real> us = dynamic_cast<const ROL::StdVector<Real>&>(u);
105  const ROL::StdVector<Real> zs = dynamic_cast<const ROL::StdVector<Real>&>(z);
106  Real v1 = (*(vs.getVector()))[0];
107  Real v2 = (*(vs.getVector()))[1];
108  Real u1 = (*(us.getVector()))[0];
109  Real z1 = (*(zs.getVector()))[0];
110  Real z2 = (*(zs.getVector()))[1];
111  (*(jvs.getVector()))[0] = u1*std::exp(z1*u1)*v1 - static_cast<Real>(2)*z2*v2;
112  }
114  assert(ijv.dimension()==1);
115  assert(v.dimension()==1);
116  assert(u.dimension()==1);
117  assert(z.dimension()==2);
118  ROL::StdVector<Real> ijvs = dynamic_cast<ROL::StdVector<Real>&>(ijv);
119  const ROL::StdVector<Real> vs = dynamic_cast<const ROL::StdVector<Real>&>(v);
120  const ROL::StdVector<Real> us = dynamic_cast<const ROL::StdVector<Real>&>(u);
121  const ROL::StdVector<Real> zs = dynamic_cast<const ROL::StdVector<Real>&>(z);
122  Real v1 = (*(vs.getVector()))[0];
123  Real u1 = (*(us.getVector()))[0];
124  Real z1 = (*(zs.getVector()))[0];
125  (*(ijvs.getVector()))[0] = v1 / (z1*std::exp(z1*u1));
126  }
128  constraint1<Real>::applyJacobian_1(ajv,v,u,z,tol);
129  }
131  assert(ajv.dimension()==2);
132  assert(v.dimension()==1);
133  assert(u.dimension()==1);
134  assert(z.dimension()==2);
135  ROL::StdVector<Real> ajvs = dynamic_cast<ROL::StdVector<Real>&>(ajv);
136  const ROL::StdVector<Real> vs = dynamic_cast<const ROL::StdVector<Real>&>(v);
137  const ROL::StdVector<Real> us = dynamic_cast<const ROL::StdVector<Real>&>(u);
138  const ROL::StdVector<Real> zs = dynamic_cast<const ROL::StdVector<Real>&>(z);
139  Real v1 = (*(vs.getVector()))[0];
140  Real u1 = (*(us.getVector()))[0];
141  Real z1 = (*(zs.getVector()))[0];
142  Real z2 = (*(zs.getVector()))[1];
143  (*(ajvs.getVector()))[0] = u1*std::exp(z1*u1)*v1;
144  (*(ajvs.getVector()))[1] = -static_cast<Real>(2)*z2*v1;
145  }
148  }
150  assert(ahwv.dimension()==1);
151  assert(w.dimension()==1);
152  assert(v.dimension()==1);
153  assert(u.dimension()==1);
154  assert(z.dimension()==2);
155  ROL::StdVector<Real> ahwvs = dynamic_cast<ROL::StdVector<Real>&>(ahwv);
156  const ROL::StdVector<Real> ws = dynamic_cast<const ROL::StdVector<Real>&>(w);
157  const ROL::StdVector<Real> vs = dynamic_cast<const ROL::StdVector<Real>&>(v);
158  const ROL::StdVector<Real> us = dynamic_cast<const ROL::StdVector<Real>&>(u);
159  const ROL::StdVector<Real> zs = dynamic_cast<const ROL::StdVector<Real>&>(z);
160  Real w1 = (*(ws.getVector()))[0];
161  Real v1 = (*(vs.getVector()))[0];
162  Real u1 = (*(us.getVector()))[0];
163  Real z1 = (*(zs.getVector()))[0];
164  (*(ahwvs.getVector()))[0] = z1*z1*std::exp(z1*u1)*v1*w1;
165  }
167  assert(ahwv.dimension()==2);
168  assert(w.dimension()==1);
169  assert(v.dimension()==1);
170  assert(u.dimension()==1);
171  assert(z.dimension()==2);
172  ROL::StdVector<Real> ahwvs = dynamic_cast<ROL::StdVector<Real>&>(ahwv);
173  const ROL::StdVector<Real> ws = dynamic_cast<const ROL::StdVector<Real>&>(w);
174  const ROL::StdVector<Real> vs = dynamic_cast<const ROL::StdVector<Real>&>(v);
175  const ROL::StdVector<Real> us = dynamic_cast<const ROL::StdVector<Real>&>(u);
176  const ROL::StdVector<Real> zs = dynamic_cast<const ROL::StdVector<Real>&>(z);
177  Real w1 = (*(ws.getVector()))[0];
178  Real v1 = (*(vs.getVector()))[0];
179  Real u1 = (*(us.getVector()))[0];
180  Real z1 = (*(zs.getVector()))[0];
181  (*(ahwvs.getVector()))[0] = std::exp(z1*u1)*(static_cast<Real>(1)+u1*z1)*v1*w1;
182  (*(ahwvs.getVector()))[1] = static_cast<Real>(0);
183  }
185  assert(ahwv.dimension()==1);
186  assert(w.dimension()==1);
187  assert(v.dimension()==2);
188  assert(u.dimension()==1);
189  assert(z.dimension()==2);
190  ROL::StdVector<Real> ahwvs = dynamic_cast<ROL::StdVector<Real>&>(ahwv);
191  const ROL::StdVector<Real> ws = dynamic_cast<const ROL::StdVector<Real>&>(w);
192  const ROL::StdVector<Real> vs = dynamic_cast<const ROL::StdVector<Real>&>(v);
193  const ROL::StdVector<Real> us = dynamic_cast<const ROL::StdVector<Real>&>(u);
194  const ROL::StdVector<Real> zs = dynamic_cast<const ROL::StdVector<Real>&>(z);
195  Real w1 = (*(ws.getVector()))[0];
196  Real v1 = (*(vs.getVector()))[0];
197  Real u1 = (*(us.getVector()))[0];
198  Real z1 = (*(zs.getVector()))[0];
199  (*(ahwvs.getVector()))[0] = std::exp(z1*u1)*(static_cast<Real>(1)+u1*z1)*v1*w1;
200  }
202  assert(ahwv.dimension()==2);
203  assert(w.dimension()==1);
204  assert(v.dimension()==2);
205  assert(u.dimension()==1);
206  assert(z.dimension()==2);
207  ROL::StdVector<Real> ahwvs = dynamic_cast<ROL::StdVector<Real>&>(ahwv);
208  const ROL::StdVector<Real> ws = dynamic_cast<const ROL::StdVector<Real>&>(w);
209  const ROL::StdVector<Real> vs = dynamic_cast<const ROL::StdVector<Real>&>(v);
210  const ROL::StdVector<Real> us = dynamic_cast<const ROL::StdVector<Real>&>(u);
211  const ROL::StdVector<Real> zs = dynamic_cast<const ROL::StdVector<Real>&>(z);
212  Real w1 = (*(ws.getVector()))[0];
213  Real v1 = (*(vs.getVector()))[0];
214  Real v2 = (*(vs.getVector()))[1];
215  Real u1 = (*(us.getVector()))[0];
216  Real z1 = (*(zs.getVector()))[0];
217  (*(ahwvs.getVector()))[0] = u1*u1*std::exp(z1*u1)*v1*w1;
218  (*(ahwvs.getVector()))[1] = -static_cast<Real>(2)*v2*w1;
219  }
220 };
221 
222 
223 template<typename Real>
224 class constraint2 : public ROL::Constraint_SimOpt<Real> {
225 public:
227  void value(ROL::Vector<Real> &c, const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
228  assert(c.dimension()==3);
229  assert(u.dimension()==1);
230  assert(z.dimension()==2);
231  ROL::StdVector<Real> cs = dynamic_cast<ROL::StdVector<Real>&>(c);
232  const ROL::StdVector<Real> us = dynamic_cast<const ROL::StdVector<Real>&>(u);
233  const ROL::StdVector<Real> zs = dynamic_cast<const ROL::StdVector<Real>&>(z);
234  Real u1 = (*(us.getVector()))[0];
235  Real z1 = (*(zs.getVector()))[0];
236  Real z2 = (*(zs.getVector()))[1];
237  (*(cs.getVector()))[0] = z1*z2*u1;
238  (*(cs.getVector()))[1] = (z1-z2)*u1;
239  (*(cs.getVector()))[2] = u1*u1;
240  }
241  void applyJacobian_1(ROL::Vector<Real> &jv, const ROL::Vector<Real> &v, const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
242  assert(jv.dimension()==3);
243  assert(v.dimension()==1);
244  assert(u.dimension()==1);
245  assert(z.dimension()==2);
246  ROL::StdVector<Real> jvs = dynamic_cast<ROL::StdVector<Real>&>(jv);
247  const ROL::StdVector<Real> vs = dynamic_cast<const ROL::StdVector<Real>&>(v);
248  const ROL::StdVector<Real> us = dynamic_cast<const ROL::StdVector<Real>&>(u);
249  const ROL::StdVector<Real> zs = dynamic_cast<const ROL::StdVector<Real>&>(z);
250  const Real two(2);
251  Real v1 = (*(vs.getVector()))[0];
252  Real u1 = (*(us.getVector()))[0];
253  Real z1 = (*(zs.getVector()))[0];
254  Real z2 = (*(zs.getVector()))[1];
255  (*(jvs.getVector()))[0] = z1*z2*v1;
256  (*(jvs.getVector()))[1] = (z1-z2)*v1;
257  (*(jvs.getVector()))[2] = two*u1*v1;
258  }
259  void applyJacobian_2(ROL::Vector<Real> &jv, const ROL::Vector<Real> &v, const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
260  assert(jv.dimension()==3);
261  assert(v.dimension()==2);
262  assert(u.dimension()==1);
263  assert(z.dimension()==2);
264  ROL::StdVector<Real> jvs = dynamic_cast<ROL::StdVector<Real>&>(jv);
265  const ROL::StdVector<Real> vs = dynamic_cast<const ROL::StdVector<Real>&>(v);
266  const ROL::StdVector<Real> us = dynamic_cast<const ROL::StdVector<Real>&>(u);
267  const ROL::StdVector<Real> zs = dynamic_cast<const ROL::StdVector<Real>&>(z);
268  Real v1 = (*(vs.getVector()))[0];
269  Real v2 = (*(vs.getVector()))[1];
270  Real u1 = (*(us.getVector()))[0];
271  Real z1 = (*(zs.getVector()))[0];
272  Real z2 = (*(zs.getVector()))[1];
273  (*(jvs.getVector()))[0] = z2*u1*v1 + z1*u1*v2;
274  (*(jvs.getVector()))[1] = (v1-v2)*u1;
275  (*(jvs.getVector()))[2] = static_cast<Real>(0);
276  }
278  assert(ajv.dimension()==1);
279  assert(v.dimension()==3);
280  assert(u.dimension()==1);
281  assert(z.dimension()==2);
282  ROL::StdVector<Real> ajvs = dynamic_cast<ROL::StdVector<Real>&>(ajv);
283  const ROL::StdVector<Real> vs = dynamic_cast<const ROL::StdVector<Real>&>(v);
284  const ROL::StdVector<Real> us = dynamic_cast<const ROL::StdVector<Real>&>(u);
285  const ROL::StdVector<Real> zs = dynamic_cast<const ROL::StdVector<Real>&>(z);
286  const Real two(2);
287  Real v1 = (*(vs.getVector()))[0];
288  Real v2 = (*(vs.getVector()))[1];
289  Real v3 = (*(vs.getVector()))[2];
290  Real u1 = (*(us.getVector()))[0];
291  Real z1 = (*(zs.getVector()))[0];
292  Real z2 = (*(zs.getVector()))[1];
293  (*(ajvs.getVector()))[0] = z1*z2*v1 + (z1-z2)*v2 + two*u1*v3;
294  }
296  assert(ajv.dimension()==2);
297  assert(v.dimension()==3);
298  assert(u.dimension()==1);
299  assert(z.dimension()==2);
300  ROL::StdVector<Real> ajvs = dynamic_cast<ROL::StdVector<Real>&>(ajv);
301  const ROL::StdVector<Real> vs = dynamic_cast<const ROL::StdVector<Real>&>(v);
302  const ROL::StdVector<Real> us = dynamic_cast<const ROL::StdVector<Real>&>(u);
303  const ROL::StdVector<Real> zs = dynamic_cast<const ROL::StdVector<Real>&>(z);
304  Real v1 = (*(vs.getVector()))[0];
305  Real v2 = (*(vs.getVector()))[1];
306  Real u1 = (*(us.getVector()))[0];
307  Real z1 = (*(zs.getVector()))[0];
308  Real z2 = (*(zs.getVector()))[1];
309  (*(ajvs.getVector()))[0] = (z2*u1*v1 + u1*v2);
310  (*(ajvs.getVector()))[1] = (z1*u1*v1 - u1*v2);
311  }
313  assert(ahwv.dimension()==1);
314  assert(w.dimension()==3);
315  assert(v.dimension()==1);
316  assert(u.dimension()==1);
317  assert(z.dimension()==2);
318  ROL::StdVector<Real> ahwvs = dynamic_cast<ROL::StdVector<Real>&>(ahwv);
319  const ROL::StdVector<Real> ws = dynamic_cast<const ROL::StdVector<Real>&>(w);
320  const ROL::StdVector<Real> vs = dynamic_cast<const ROL::StdVector<Real>&>(v);
321  const ROL::StdVector<Real> us = dynamic_cast<const ROL::StdVector<Real>&>(u);
322  const ROL::StdVector<Real> zs = dynamic_cast<const ROL::StdVector<Real>&>(z);
323  const Real two(2);
324  Real w3 = (*(ws.getVector()))[2];
325  Real v1 = (*(vs.getVector()))[0];
326  (*(ahwvs.getVector()))[0] = two*v1*w3;
327  }
329  assert(ahwv.dimension()==2);
330  assert(w.dimension()==3);
331  assert(v.dimension()==1);
332  assert(u.dimension()==1);
333  assert(z.dimension()==2);
334  ROL::StdVector<Real> ahwvs = dynamic_cast<ROL::StdVector<Real>&>(ahwv);
335  const ROL::StdVector<Real> ws = dynamic_cast<const ROL::StdVector<Real>&>(w);
336  const ROL::StdVector<Real> vs = dynamic_cast<const ROL::StdVector<Real>&>(v);
337  const ROL::StdVector<Real> us = dynamic_cast<const ROL::StdVector<Real>&>(u);
338  const ROL::StdVector<Real> zs = dynamic_cast<const ROL::StdVector<Real>&>(z);
339  Real w1 = (*(ws.getVector()))[0];
340  Real w2 = (*(ws.getVector()))[1];
341  Real v1 = (*(vs.getVector()))[0];
342  Real z1 = (*(zs.getVector()))[0];
343  Real z2 = (*(zs.getVector()))[1];
344  (*(ahwvs.getVector()))[0] = (z2*v1*w1 + v1*w2);
345  (*(ahwvs.getVector()))[1] = (z1*v1*w1 - v1*w2);
346  }
348  assert(ahwv.dimension()==1);
349  assert(w.dimension()==3);
350  assert(v.dimension()==2);
351  assert(u.dimension()==1);
352  assert(z.dimension()==2);
353  ROL::StdVector<Real> ahwvs = dynamic_cast<ROL::StdVector<Real>&>(ahwv);
354  const ROL::StdVector<Real> ws = dynamic_cast<const ROL::StdVector<Real>&>(w);
355  const ROL::StdVector<Real> vs = dynamic_cast<const ROL::StdVector<Real>&>(v);
356  const ROL::StdVector<Real> us = dynamic_cast<const ROL::StdVector<Real>&>(u);
357  const ROL::StdVector<Real> zs = dynamic_cast<const ROL::StdVector<Real>&>(z);
358  Real w1 = (*(ws.getVector()))[0];
359  Real w2 = (*(ws.getVector()))[1];
360  Real v1 = (*(vs.getVector()))[0];
361  Real v2 = (*(vs.getVector()))[1];
362  Real z1 = (*(zs.getVector()))[0];
363  Real z2 = (*(zs.getVector()))[1];
364  (*(ahwvs.getVector()))[0] = (v1*z2+z1*v2)*w1 + (v1-v2)*w2;
365  }
367  assert(ahwv.dimension()==2);
368  assert(w.dimension()==3);
369  assert(v.dimension()==2);
370  assert(u.dimension()==1);
371  assert(z.dimension()==2);
372  ROL::StdVector<Real> ahwvs = dynamic_cast<ROL::StdVector<Real>&>(ahwv);
373  const ROL::StdVector<Real> ws = dynamic_cast<const ROL::StdVector<Real>&>(w);
374  const ROL::StdVector<Real> vs = dynamic_cast<const ROL::StdVector<Real>&>(v);
375  const ROL::StdVector<Real> us = dynamic_cast<const ROL::StdVector<Real>&>(u);
376  const ROL::StdVector<Real> zs = dynamic_cast<const ROL::StdVector<Real>&>(z);
377  Real w1 = (*(ws.getVector()))[0];
378  Real v1 = (*(vs.getVector()))[0];
379  Real v2 = (*(vs.getVector()))[1];
380  Real u1 = (*(us.getVector()))[0];
381  (*(ahwvs.getVector()))[0] = v2*u1*w1;
382  (*(ahwvs.getVector()))[1] = v1*u1*w1;
383  }
384 };
385 
386 int main(int argc, char *argv[]) {
387  using RealT = double;
388 
389  Teuchos::GlobalMPISession mpiSession(&argc, &argv);
390 
391  // This little trick lets us print to std::cout only if a (dummy) command-line argument is provided.
392  int iprint = argc - 1;
393  ROL::Ptr<std::ostream> outStream;
394  ROL::nullstream bhs; // outputs nothing
395  if (iprint > 0)
396  outStream = ROL::makePtrFromRef(std::cout);
397  else
398  outStream = ROL::makePtrFromRef(bhs);
399 
400  // Save the format state of the original std::cout.
401  ROL::nullstream oldFormatState;
402  oldFormatState.copyfmt(std::cout);
403 
404 // RealT errtol = std::sqrt(ROL::ROL_THRESHOLD<RealT>());
405 
406  int errorFlag = 0;
407 
408  // *** Test body.
409 
410  try {
411 
412  unsigned c1_dim = 1; // Constraint1 dimension
413  unsigned c2_dim = 3; // Constraint1 dimension
414  unsigned u_dim = 1; // State dimension
415  unsigned z_dim = 2; // Control dimension
416 
417  auto c1 = ROL::makePtr<ROL::StdVector<RealT>>(c1_dim);
418  auto c2 = ROL::makePtr<ROL::StdVector<RealT>>(c2_dim);
419  auto u = ROL::makePtr<ROL::StdVector<RealT>>(u_dim);
420  auto z = ROL::makePtr<ROL::StdVector<RealT>>(z_dim);
421  auto vc1 = ROL::makePtr<ROL::StdVector<RealT>>(c1_dim);
422  auto vc2 = ROL::makePtr<ROL::StdVector<RealT>>(c2_dim);
423  auto vu = ROL::makePtr<ROL::StdVector<RealT>>(u_dim);
424  auto vz = ROL::makePtr<ROL::StdVector<RealT>>(z_dim);
425  auto du = ROL::makePtr<ROL::StdVector<RealT>>(u_dim);
426  auto dz = ROL::makePtr<ROL::StdVector<RealT>>(z_dim);
427  c1->randomize(static_cast<RealT>(-1),static_cast<RealT>(1));
428  c2->randomize(static_cast<RealT>(-1),static_cast<RealT>(1));
429  u->randomize(static_cast<RealT>(-1),static_cast<RealT>(1));
430  z->randomize(static_cast<RealT>(-1),static_cast<RealT>(1));
431  vc1->randomize(static_cast<RealT>(-1),static_cast<RealT>(1));
432  vc2->randomize(static_cast<RealT>(-1),static_cast<RealT>(1));
433  vu->randomize(static_cast<RealT>(-1),static_cast<RealT>(1));
434  vz->randomize(static_cast<RealT>(-1),static_cast<RealT>(1));
435  du->randomize(static_cast<RealT>(-1),static_cast<RealT>(1));
436  dz->randomize(static_cast<RealT>(-1),static_cast<RealT>(1));
437 
438  auto con1 = ROL::makePtr<constraint1<RealT>>();
439  auto con2 = ROL::makePtr<constraint2<RealT>>();
440  auto stateStore = ROL::makePtr<ROL::VectorController<RealT>>();
441  auto rcon = ROL::makePtr<ROL::Reduced_Constraint_SimOpt<RealT>>(con2,con1,stateStore,u,z,vc1,c2,true,false);
442 
443  con1->checkSolve(*u,*z,*c1,true,*outStream);
444  con1->checkAdjointConsistencyJacobian_1(*vc1,*vu,*u,*z,true,*outStream);
445  con1->checkAdjointConsistencyJacobian_2(*vc1,*vz,*u,*z,true,*outStream);
446  con1->checkInverseJacobian_1(*c1,*vu,*u,*z,true,*outStream);
447  con1->checkInverseAdjointJacobian_1(*c1,*vu,*u,*z,true,*outStream);
448  con1->checkApplyJacobian_1(*u,*z,*vu,*vc1,true,*outStream);
449  con1->checkApplyJacobian_2(*u,*z,*vz,*vc1,true,*outStream);
450  con1->checkApplyAdjointHessian_11(*u,*z,*vc1,*vu,*du,true,*outStream);
451  con1->checkApplyAdjointHessian_12(*u,*z,*vc1,*vu,*dz,true,*outStream);
452  con1->checkApplyAdjointHessian_21(*u,*z,*vc1,*vz,*du,true,*outStream);
453  con1->checkApplyAdjointHessian_22(*u,*z,*vc1,*vz,*dz,true,*outStream);
454 
455  con2->checkAdjointConsistencyJacobian_1(*vc2,*vu,*u,*z,true,*outStream);
456  con2->checkAdjointConsistencyJacobian_2(*vc2,*vz,*u,*z,true,*outStream);
457  con2->checkApplyJacobian_1(*u,*z,*vu,*vc2,true,*outStream);
458  con2->checkApplyJacobian_2(*u,*z,*vz,*vc2,true,*outStream);
459  con2->checkApplyAdjointHessian_11(*u,*z,*vc2,*vu,*du,true,*outStream);
460  con2->checkApplyAdjointHessian_12(*u,*z,*vc2,*vu,*dz,true,*outStream);
461  con2->checkApplyAdjointHessian_21(*u,*z,*vc2,*vz,*du,true,*outStream);
462  con2->checkApplyAdjointHessian_22(*u,*z,*vc2,*vz,*dz,true,*outStream);
463 
464  rcon->checkAdjointConsistencyJacobian(*vc2,*vz,*z,true,*outStream);
465  rcon->checkApplyJacobian(*z,*vz,*vc2,true,*outStream);
466  rcon->checkApplyAdjointHessian(*z,*vc2,*vz,*dz,true,*outStream);
467  }
468  catch (std::logic_error& err) {
469  *outStream << err.what() << "\n";
470  errorFlag = -1000;
471  }; // end try
472 
473  if (errorFlag != 0)
474  std::cout << "End Result: TEST FAILED\n";
475  else
476  std::cout << "End Result: TEST PASSED\n";
477 
478  return 0;
479 
480 
481 }
482 
void applyAdjointHessian_21(ROL::Vector< Real > &ahwv, const ROL::Vector< Real > &w, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply the simulation-space derivative of the adjoint of the constraint optimization-space Jacobian at...
Definition: test_19.cpp:184
virtual int dimension() const
Return dimension of the vector space.
Definition: ROL_Vector.hpp:162
void applyAdjointJacobian_2(ROL::Vector< Real > &ajv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply the adjoint of the partial constraint Jacobian at , , to vector . This is the primary interface...
Definition: test_19.cpp:130
void solve(ROL::Vector< Real > &c, ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Given , solve for .
Definition: test_19.cpp:72
Ptr< const std::vector< Element > > getVector() const
void applyInverseAdjointJacobian_1(ROL::Vector< Real > &iajv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply the inverse of the adjoint of the partial constraint Jacobian at , , to the vector ...
Definition: test_19.cpp:146
void applyAdjointHessian_22(ROL::Vector< Real > &ahwv, const ROL::Vector< Real > &w, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply the optimization-space derivative of the adjoint of the constraint optimization-space Jacobian ...
Definition: test_19.cpp:366
void applyAdjointHessian_12(ROL::Vector< Real > &ahwv, const ROL::Vector< Real > &w, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply the optimization-space derivative of the adjoint of the constraint simulation-space Jacobian at...
Definition: test_19.cpp:166
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:46
void applyJacobian_2(ROL::Vector< Real > &jv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply the partial constraint Jacobian at , , to the vector .
Definition: test_19.cpp:97
Defines a no-output stream class ROL::NullStream and a function makeStreamPtr which either wraps a re...
basic_nullstream< char, std::char_traits< char >> nullstream
Definition: ROL_Stream.hpp:36
void applyAdjointJacobian_1(ROL::Vector< Real > &ajv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply the adjoint of the partial constraint Jacobian at , , to the vector . This is the primary inter...
Definition: test_19.cpp:277
void applyAdjointHessian_22(ROL::Vector< Real > &ahwv, const ROL::Vector< Real > &w, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply the optimization-space derivative of the adjoint of the constraint optimization-space Jacobian ...
Definition: test_19.cpp:201
void value(ROL::Vector< Real > &c, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Evaluate the constraint operator at .
Definition: test_19.cpp:227
void applyJacobian_1(ROL::Vector< Real > &jv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply the partial constraint Jacobian at , , to the vector .
Definition: test_19.cpp:241
void applyJacobian_1(ROL::Vector< Real > &jv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply the partial constraint Jacobian at , , to the vector .
Definition: test_19.cpp:83
void applyJacobian_2(ROL::Vector< Real > &jv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply the partial constraint Jacobian at , , to the vector .
Definition: test_19.cpp:259
void applyAdjointHessian_12(ROL::Vector< Real > &ahwv, const ROL::Vector< Real > &w, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply the optimization-space derivative of the adjoint of the constraint simulation-space Jacobian at...
Definition: test_19.cpp:328
void applyAdjointHessian_21(ROL::Vector< Real > &ahwv, const ROL::Vector< Real > &w, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply the simulation-space derivative of the adjoint of the constraint optimization-space Jacobian at...
Definition: test_19.cpp:347
void applyInverseJacobian_1(ROL::Vector< Real > &ijv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply the inverse partial constraint Jacobian at , , to the vector .
Definition: test_19.cpp:113
int main(int argc, char *argv[])
Defines the constraint operator interface for simulation-based optimization.
void applyAdjointJacobian_2(ROL::Vector< Real > &ajv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply the adjoint of the partial constraint Jacobian at , , to vector . This is the primary interface...
Definition: test_19.cpp:295
void applyAdjointHessian_11(ROL::Vector< Real > &ahwv, const ROL::Vector< Real > &w, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply the simulation-space derivative of the adjoint of the constraint simulation-space Jacobian at ...
Definition: test_19.cpp:149
void value(ROL::Vector< Real > &c, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Evaluate the constraint operator at .
Definition: test_19.cpp:60
void applyAdjointHessian_11(ROL::Vector< Real > &ahwv, const ROL::Vector< Real > &w, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply the simulation-space derivative of the adjoint of the constraint simulation-space Jacobian at ...
Definition: test_19.cpp:312
void applyAdjointJacobian_1(ROL::Vector< Real > &ajv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply the adjoint of the partial constraint Jacobian at , , to the vector . This is the primary inter...
Definition: test_19.cpp:127