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