ROL
function/test_10.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"
16#include "ROL_Stream.hpp"
17#include "Teuchos_GlobalMPISession.hpp"
18
19#include <iostream>
20
21template<class Real>
23public:
25
26 void value(ROL::Vector<Real> &c, const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
27 ROL::Ptr<std::vector<Real> > cp
28 = dynamic_cast<ROL::StdVector<Real>&>(c).getVector();
29 ROL::Ptr<const std::vector<Real> > up
30 = dynamic_cast<const ROL::StdVector<Real>&>(u).getVector();
31 ROL::Ptr<const std::vector<Real> > zp
32 = dynamic_cast<const ROL::StdVector<Real>&>(z).getVector();
33
34 Real half(0.5), two(2);
35 // C(0) = U(0) - Z(0)
36 (*cp)[0] = (*up)[0]-(*zp)[0];
37 // C(1) = 0.5 * (U(0) + U(1) - Z(0))^2
38 (*cp)[1] = half*std::pow((*up)[0]+(*up)[1]-(*zp)[0],two);
39 }
40
42 const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
43 ROL::Ptr<std::vector<Real> > jvp
44 = dynamic_cast<ROL::StdVector<Real>&>(jv).getVector();
45 ROL::Ptr<const std::vector<Real> > vp
46 = dynamic_cast<const ROL::StdVector<Real>&>(v).getVector();
47 ROL::Ptr<const std::vector<Real> > up
48 = dynamic_cast<const ROL::StdVector<Real>&>(u).getVector();
49 ROL::Ptr<const std::vector<Real> > zp
50 = dynamic_cast<const ROL::StdVector<Real>&>(z).getVector();
51 (*jvp)[0] = (*vp)[0];
52 (*jvp)[1] = ((*up)[0] + (*up)[1] - (*zp)[0]) * ((*vp)[0] + (*vp)[1]);
53 }
54
56 const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
57 ROL::Ptr<std::vector<Real> > jvp
58 = dynamic_cast<ROL::StdVector<Real>&>(jv).getVector();
59 ROL::Ptr<const std::vector<Real> > vp
60 = dynamic_cast<const ROL::StdVector<Real>&>(v).getVector();
61 ROL::Ptr<const std::vector<Real> > up
62 = dynamic_cast<const ROL::StdVector<Real>&>(u).getVector();
63 ROL::Ptr<const std::vector<Real> > zp
64 = dynamic_cast<const ROL::StdVector<Real>&>(z).getVector();
65 (*jvp)[0] = -(*vp)[0];
66 (*jvp)[1] = ((*zp)[0] - (*up)[0] - (*up)[1]) * (*vp)[0];
67 }
68
70 const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
71 ROL::Ptr<std::vector<Real> > ajvp
72 = dynamic_cast<ROL::StdVector<Real>&>(ajv).getVector();
73 ROL::Ptr<const std::vector<Real> > vp
74 = dynamic_cast<const ROL::StdVector<Real>&>(v).getVector();
75 ROL::Ptr<const std::vector<Real> > up
76 = dynamic_cast<const ROL::StdVector<Real>&>(u).getVector();
77 ROL::Ptr<const std::vector<Real> > zp
78 = dynamic_cast<const ROL::StdVector<Real>&>(z).getVector();
79 (*ajvp)[0] = (*vp)[0] + ((*up)[0] + (*up)[1] - (*zp)[0]) * (*vp)[1];
80 (*ajvp)[1] = ((*up)[0] + (*up)[1] - (*zp)[0]) * (*vp)[1];
81 }
82
84 const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
85 ROL::Ptr<std::vector<Real> > ajvp
86 = dynamic_cast<ROL::StdVector<Real>&>(ajv).getVector();
87 ROL::Ptr<const std::vector<Real> > vp
88 = dynamic_cast<const ROL::StdVector<Real>&>(v).getVector();
89 ROL::Ptr<const std::vector<Real> > up
90 = dynamic_cast<const ROL::StdVector<Real>&>(u).getVector();
91 ROL::Ptr<const std::vector<Real> > zp
92 = dynamic_cast<const ROL::StdVector<Real>&>(z).getVector();
93 (*ajvp)[0] = ((*zp)[0] - (*up)[0] - (*up)[1]) * (*vp)[1] - (*vp)[0];
94 }
95
97 const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
98 ROL::Ptr<std::vector<Real> > ahwvp
99 = dynamic_cast<ROL::StdVector<Real>&>(ahwv).getVector();
100 ROL::Ptr<const std::vector<Real> > wp
101 = dynamic_cast<const ROL::StdVector<Real>&>(w).getVector();
102 ROL::Ptr<const std::vector<Real> > vp
103 = dynamic_cast<const ROL::StdVector<Real>&>(v).getVector();
104 ROL::Ptr<const std::vector<Real> > up
105 = dynamic_cast<const ROL::StdVector<Real>&>(u).getVector();
106 ROL::Ptr<const std::vector<Real> > zp
107 = dynamic_cast<const ROL::StdVector<Real>&>(z).getVector();
108 (*ahwvp)[0] = (*wp)[1] * ((*vp)[0] + (*vp)[1]);
109 (*ahwvp)[1] = (*wp)[1] * ((*vp)[0] + (*vp)[1]);
110 }
111
113 const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
114 ROL::Ptr<std::vector<Real> > ahwvp
115 = dynamic_cast<ROL::StdVector<Real>&>(ahwv).getVector();
116 ROL::Ptr<const std::vector<Real> > wp
117 = dynamic_cast<const ROL::StdVector<Real>&>(w).getVector();
118 ROL::Ptr<const std::vector<Real> > vp
119 = dynamic_cast<const ROL::StdVector<Real>&>(v).getVector();
120 ROL::Ptr<const std::vector<Real> > up
121 = dynamic_cast<const ROL::StdVector<Real>&>(u).getVector();
122 ROL::Ptr<const std::vector<Real> > zp
123 = dynamic_cast<const ROL::StdVector<Real>&>(z).getVector();
124 (*ahwvp)[0] = -(*wp)[1] * ((*vp)[0] + (*vp)[1]);
125 }
126
128 const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
129 ROL::Ptr<std::vector<Real> > ahwvp
130 = dynamic_cast<ROL::StdVector<Real>&>(ahwv).getVector();
131 ROL::Ptr<const std::vector<Real> > wp
132 = dynamic_cast<const ROL::StdVector<Real>&>(w).getVector();
133 ROL::Ptr<const std::vector<Real> > vp
134 = dynamic_cast<const ROL::StdVector<Real>&>(v).getVector();
135 ROL::Ptr<const std::vector<Real> > up
136 = dynamic_cast<const ROL::StdVector<Real>&>(u).getVector();
137 ROL::Ptr<const std::vector<Real> > zp
138 = dynamic_cast<const ROL::StdVector<Real>&>(z).getVector();
139 (*ahwvp)[0] = -(*wp)[1] * (*vp)[0];
140 (*ahwvp)[1] = -(*wp)[1] * (*vp)[0];
141 }
142
144 const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
145 ROL::Ptr<std::vector<Real> > ahwvp
146 = dynamic_cast<ROL::StdVector<Real>&>(ahwv).getVector();
147 ROL::Ptr<const std::vector<Real> > wp
148 = dynamic_cast<const ROL::StdVector<Real>&>(w).getVector();
149 ROL::Ptr<const std::vector<Real> > vp
150 = dynamic_cast<const ROL::StdVector<Real>&>(v).getVector();
151 ROL::Ptr<const std::vector<Real> > up
152 = dynamic_cast<const ROL::StdVector<Real>&>(u).getVector();
153 ROL::Ptr<const std::vector<Real> > zp
154 = dynamic_cast<const ROL::StdVector<Real>&>(z).getVector();
155 (*ahwvp)[0] = (*wp)[1] * (*vp)[0];
156 }
157};
158
159template<class Real>
161public:
163
164 void value(ROL::Vector<Real> &c, const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
165 ROL::Ptr<std::vector<Real> > cp
166 = dynamic_cast<ROL::StdVector<Real>&>(c).getVector();
167 ROL::Ptr<const std::vector<Real> > up
168 = dynamic_cast<const ROL::StdVector<Real>&>(u).getVector();
169 ROL::Ptr<const std::vector<Real> > zp
170 = dynamic_cast<const ROL::StdVector<Real>&>(z).getVector();
171
172 const Real one(1), two(2);
173 // C = exp(U) - (Z^2 + 1)
174 (*cp)[0] = std::exp((*up)[0])-(std::pow((*zp)[0],two) + one);
175 }
176
178 const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
179 ROL::Ptr<std::vector<Real> > jvp
180 = dynamic_cast<ROL::StdVector<Real>&>(jv).getVector();
181 ROL::Ptr<const std::vector<Real> > vp
182 = dynamic_cast<const ROL::StdVector<Real>&>(v).getVector();
183 ROL::Ptr<const std::vector<Real> > up
184 = dynamic_cast<const ROL::StdVector<Real>&>(u).getVector();
185 ROL::Ptr<const std::vector<Real> > zp
186 = dynamic_cast<const ROL::StdVector<Real>&>(z).getVector();
187 (*jvp)[0] = std::exp((*up)[0]) * (*vp)[0];
188 }
189
191 const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
192 ROL::Ptr<std::vector<Real> > jvp
193 = dynamic_cast<ROL::StdVector<Real>&>(jv).getVector();
194 ROL::Ptr<const std::vector<Real> > vp
195 = dynamic_cast<const ROL::StdVector<Real>&>(v).getVector();
196 ROL::Ptr<const std::vector<Real> > up
197 = dynamic_cast<const ROL::StdVector<Real>&>(u).getVector();
198 ROL::Ptr<const std::vector<Real> > zp
199 = dynamic_cast<const ROL::StdVector<Real>&>(z).getVector();
200
201 const Real two(2);
202 (*jvp)[0] = -two * (*zp)[0] * (*vp)[0];
203 }
204
206 const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
207 ROL::Ptr<std::vector<Real> > ajvp
208 = dynamic_cast<ROL::StdVector<Real>&>(ajv).getVector();
209 ROL::Ptr<const std::vector<Real> > vp
210 = dynamic_cast<const ROL::StdVector<Real>&>(v).getVector();
211 ROL::Ptr<const std::vector<Real> > up
212 = dynamic_cast<const ROL::StdVector<Real>&>(u).getVector();
213 ROL::Ptr<const std::vector<Real> > zp
214 = dynamic_cast<const ROL::StdVector<Real>&>(z).getVector();
215 (*ajvp)[0] = std::exp((*up)[0]) * (*vp)[0];
216 }
217
219 const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
220 ROL::Ptr<std::vector<Real> > ajvp
221 = dynamic_cast<ROL::StdVector<Real>&>(ajv).getVector();
222 ROL::Ptr<const std::vector<Real> > vp
223 = dynamic_cast<const ROL::StdVector<Real>&>(v).getVector();
224 ROL::Ptr<const std::vector<Real> > up
225 = dynamic_cast<const ROL::StdVector<Real>&>(u).getVector();
226 ROL::Ptr<const std::vector<Real> > zp
227 = dynamic_cast<const ROL::StdVector<Real>&>(z).getVector();
228
229 const Real two(2);
230 (*ajvp)[0] = -two * (*zp)[0] * (*vp)[0];
231 }
232
234 const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
235 ROL::Ptr<std::vector<Real> > ijvp
236 = dynamic_cast<ROL::StdVector<Real>&>(ijv).getVector();
237 ROL::Ptr<const std::vector<Real> > vp
238 = dynamic_cast<const ROL::StdVector<Real>&>(v).getVector();
239 ROL::Ptr<const std::vector<Real> > up
240 = dynamic_cast<const ROL::StdVector<Real>&>(u).getVector();
241 ROL::Ptr<const std::vector<Real> > zp
242 = dynamic_cast<const ROL::StdVector<Real>&>(z).getVector();
243 (*ijvp)[0] = (*vp)[0] / std::exp((*up)[0]);
244 }
245
247 const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
248 ROL::Ptr<std::vector<Real> > ijvp
249 = dynamic_cast<ROL::StdVector<Real>&>(ijv).getVector();
250 ROL::Ptr<const std::vector<Real> > vp
251 = dynamic_cast<const ROL::StdVector<Real>&>(v).getVector();
252 ROL::Ptr<const std::vector<Real> > up
253 = dynamic_cast<const ROL::StdVector<Real>&>(u).getVector();
254 ROL::Ptr<const std::vector<Real> > zp
255 = dynamic_cast<const ROL::StdVector<Real>&>(z).getVector();
256 (*ijvp)[0] = (*vp)[0] / std::exp((*up)[0]);
257 }
258
260 const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
261 ROL::Ptr<std::vector<Real> > ahwvp
262 = dynamic_cast<ROL::StdVector<Real>&>(ahwv).getVector();
263 ROL::Ptr<const std::vector<Real> > wp
264 = dynamic_cast<const ROL::StdVector<Real>&>(w).getVector();
265 ROL::Ptr<const std::vector<Real> > vp
266 = dynamic_cast<const ROL::StdVector<Real>&>(v).getVector();
267 ROL::Ptr<const std::vector<Real> > up
268 = dynamic_cast<const ROL::StdVector<Real>&>(u).getVector();
269 ROL::Ptr<const std::vector<Real> > zp
270 = dynamic_cast<const ROL::StdVector<Real>&>(z).getVector();
271 (*ahwvp)[0] = std::exp((*up)[0]) * (*wp)[0] * (*vp)[0];
272 }
273
275 const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
276 ROL::Ptr<std::vector<Real> > ahwvp
277 = dynamic_cast<ROL::StdVector<Real>&>(ahwv).getVector();
278 ROL::Ptr<const std::vector<Real> > wp
279 = dynamic_cast<const ROL::StdVector<Real>&>(w).getVector();
280 ROL::Ptr<const std::vector<Real> > vp
281 = dynamic_cast<const ROL::StdVector<Real>&>(v).getVector();
282 ROL::Ptr<const std::vector<Real> > up
283 = dynamic_cast<const ROL::StdVector<Real>&>(u).getVector();
284 ROL::Ptr<const std::vector<Real> > zp
285 = dynamic_cast<const ROL::StdVector<Real>&>(z).getVector();
286 (*ahwvp)[0] = static_cast<Real>(0);
287 }
288
290 const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
291 ROL::Ptr<std::vector<Real> > ahwvp
292 = dynamic_cast<ROL::StdVector<Real>&>(ahwv).getVector();
293 ROL::Ptr<const std::vector<Real> > wp
294 = dynamic_cast<const ROL::StdVector<Real>&>(w).getVector();
295 ROL::Ptr<const std::vector<Real> > vp
296 = dynamic_cast<const ROL::StdVector<Real>&>(v).getVector();
297 ROL::Ptr<const std::vector<Real> > up
298 = dynamic_cast<const ROL::StdVector<Real>&>(u).getVector();
299 ROL::Ptr<const std::vector<Real> > zp
300 = dynamic_cast<const ROL::StdVector<Real>&>(z).getVector();
301 (*ahwvp)[0] = static_cast<Real>(0);
302 }
303
305 const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
306 ROL::Ptr<std::vector<Real> > ahwvp
307 = dynamic_cast<ROL::StdVector<Real>&>(ahwv).getVector();
308 ROL::Ptr<const std::vector<Real> > wp
309 = dynamic_cast<const ROL::StdVector<Real>&>(w).getVector();
310 ROL::Ptr<const std::vector<Real> > vp
311 = dynamic_cast<const ROL::StdVector<Real>&>(v).getVector();
312 ROL::Ptr<const std::vector<Real> > up
313 = dynamic_cast<const ROL::StdVector<Real>&>(u).getVector();
314 ROL::Ptr<const std::vector<Real> > zp
315 = dynamic_cast<const ROL::StdVector<Real>&>(z).getVector();
316 (*ahwvp)[0] = static_cast<Real>(-2) * (*wp)[0] * (*vp)[0];
317 }
318};
319
320typedef double RealT;
321
322int main(int argc, char *argv[]) {
323
324 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
325
326 // This little trick lets us print to std::cout only if a (dummy) command-line argument is provided.
327 int iprint = argc - 1;
328 ROL::Ptr<std::ostream> outStream;
329 ROL::nullstream bhs; // outputs nothing
330 if (iprint > 0)
331 outStream = ROL::makePtrFromRef(std::cout);
332 else
333 outStream = ROL::makePtrFromRef(bhs);
334
335 int errorFlag = 0;
336
337 // *** Example body.
338
339 try {
340
341 int dim = 2;
342 int dimz = 1;
343 ROL::Ptr<std::vector<RealT> > ustd = ROL::makePtr<std::vector<RealT>>(dim);
344 ROL::Ptr<std::vector<RealT> > dustd = ROL::makePtr<std::vector<RealT>>(dim);
345 ROL::Ptr<std::vector<RealT> > zstd = ROL::makePtr<std::vector<RealT>>(dimz);
346 ROL::Ptr<std::vector<RealT> > dzstd = ROL::makePtr<std::vector<RealT>>(dimz);
347 ROL::Ptr<std::vector<RealT> > cstd = ROL::makePtr<std::vector<RealT>>(dim);
348 ROL::Ptr<std::vector<RealT> > czstd = ROL::makePtr<std::vector<RealT>>(dimz);
349 ROL::Ptr<std::vector<RealT> > sstd = ROL::makePtr<std::vector<RealT>>(dimz);
350 ROL::Ptr<std::vector<RealT> > dsstd = ROL::makePtr<std::vector<RealT>>(dimz);
351
352 (*ustd)[0] = static_cast<RealT>(rand())/static_cast<RealT>(RAND_MAX);
353 (*ustd)[1] = static_cast<RealT>(rand())/static_cast<RealT>(RAND_MAX);
354 (*dustd)[0] = static_cast<RealT>(rand())/static_cast<RealT>(RAND_MAX);
355 (*dustd)[1] = static_cast<RealT>(rand())/static_cast<RealT>(RAND_MAX);
356 (*zstd)[0] = static_cast<RealT>(rand())/static_cast<RealT>(RAND_MAX);
357 (*dzstd)[0] = static_cast<RealT>(rand())/static_cast<RealT>(RAND_MAX);
358 (*cstd)[0] = static_cast<RealT>(rand())/static_cast<RealT>(RAND_MAX);
359 (*cstd)[1] = static_cast<RealT>(rand())/static_cast<RealT>(RAND_MAX);
360 (*czstd)[0] = static_cast<RealT>(rand())/static_cast<RealT>(RAND_MAX);
361 (*sstd)[0] = static_cast<RealT>(rand())/static_cast<RealT>(RAND_MAX);
362 (*dsstd)[0] = static_cast<RealT>(rand())/static_cast<RealT>(RAND_MAX);
363
364 ROL::Ptr<ROL::Vector<RealT> > u = ROL::makePtr<ROL::StdVector<RealT>>(ustd);
365 ROL::Ptr<ROL::Vector<RealT> > du = ROL::makePtr<ROL::StdVector<RealT>>(dustd);
366 ROL::Ptr<ROL::Vector<RealT> > z = ROL::makePtr<ROL::StdVector<RealT>>(zstd);
367 ROL::Ptr<ROL::Vector<RealT> > dz = ROL::makePtr<ROL::StdVector<RealT>>(dzstd);
368 ROL::Ptr<ROL::Vector<RealT> > c = ROL::makePtr<ROL::StdVector<RealT>>(cstd);
369 ROL::Ptr<ROL::Vector<RealT> > cz = ROL::makePtr<ROL::StdVector<RealT>>(czstd);
370 ROL::Ptr<ROL::Vector<RealT> > s = ROL::makePtr<ROL::StdVector<RealT>>(sstd);
371 ROL::Ptr<ROL::Vector<RealT> > ds = ROL::makePtr<ROL::StdVector<RealT>>(dsstd);
372
379
380 ROL::Ptr<ROL::Constraint_SimOpt<RealT> > valCon = ROL::makePtr<valConstraint<RealT>>();
381 valCon->checkAdjointConsistencyJacobian_1(*c,*du,*u,*s,true,*outStream);
382 valCon->checkAdjointConsistencyJacobian_2(*c,*dz,*u,*s,true,*outStream);
383 valCon->checkApplyJacobian_1(*u,*s,*du,*c,true,*outStream);
384 valCon->checkApplyJacobian_2(*u,*s,*ds,*c,true,*outStream);
385 valCon->checkApplyJacobian(x,dx,*c,true,*outStream);
386 valCon->checkApplyAdjointHessian_11(*u,*s,*c,*du,*u,true,*outStream);
387 valCon->checkApplyAdjointHessian_12(*u,*s,*c,*du,*s,true,*outStream);
388 valCon->checkApplyAdjointHessian_21(*u,*s,*c,*ds,*u,true,*outStream);
389 valCon->checkApplyAdjointHessian_22(*u,*s,*c,*ds,*s,true,*outStream);
390 valCon->checkApplyAdjointHessian(x,*c,dx,x,true,*outStream);
391
392 ROL::Ptr<ROL::Constraint_SimOpt<RealT> > redCon = ROL::makePtr<redConstraint<RealT>>();
393 redCon->checkAdjointConsistencyJacobian_1(*cz,*ds,*s,*z,true,*outStream);
394 redCon->checkAdjointConsistencyJacobian_2(*cz,*dz,*s,*z,true,*outStream);
395 redCon->checkInverseJacobian_1(*cz,*ds,*s,*z,true,*outStream);
396 redCon->checkInverseAdjointJacobian_1(*ds,*cz,*s,*z,true,*outStream);
397 redCon->checkApplyJacobian_1(*s,*z,*ds,*cz,true,*outStream);
398 redCon->checkApplyJacobian_2(*s,*z,*dz,*cz,true,*outStream);
399 redCon->checkApplyJacobian(y,dy,*cz,true,*outStream);
400 redCon->checkApplyAdjointHessian_11(*s,*z,*cz,*ds,*s,true,*outStream);
401 redCon->checkApplyAdjointHessian_12(*s,*z,*cz,*ds,*z,true,*outStream);
402 redCon->checkApplyAdjointHessian_21(*s,*z,*cz,*dz,*s,true,*outStream);
403 redCon->checkApplyAdjointHessian_22(*s,*z,*cz,*dz,*z,true,*outStream);
404 redCon->checkApplyAdjointHessian(y,*cz,dy,y,true,*outStream);
405
406 ROL::CompositeConstraint_SimOpt<RealT> con(valCon,redCon,*c,*cz,*u,*s,*z);
407 con.checkAdjointConsistencyJacobian_1(*c,*du,*u,*z,true,*outStream);
408 con.checkAdjointConsistencyJacobian_2(*c,*dz,*u,*z,true,*outStream);
409 con.checkApplyJacobian_1(*u,*z,*du,*c,true,*outStream);
410 con.checkApplyJacobian_2(*u,*z,*dz,*c,true,*outStream);
411 con.checkApplyJacobian(w,dw,*c,true,*outStream);
412 con.checkApplyAdjointHessian_11(*u,*z,*c,*du,*u,true,*outStream);
413 con.checkApplyAdjointHessian_12(*u,*z,*c,*du,*z,true,*outStream);
414 con.checkApplyAdjointHessian_21(*u,*z,*c,*dz,*u,true,*outStream);
415 con.checkApplyAdjointHessian_22(*u,*z,*c,*dz,*z,true,*outStream);
416 con.checkApplyAdjointHessian(w,*c,dw,w,true,*outStream);
417 }
418 catch (std::logic_error& err) {
419 *outStream << err.what() << "\n";
420 errorFlag = -1000;
421 }; // end try
422
423 if (errorFlag != 0)
424 std::cout << "End Result: TEST FAILED\n";
425 else
426 std::cout << "End Result: TEST PASSED\n";
427
428 return 0;
429
430}
431
Defines a no-output stream class ROL::NullStream and a function makeStreamPtr which either wraps a re...
Defines a composite equality constraint operator interface for simulation-based optimization.
Defines the constraint operator interface for simulation-based optimization.
virtual Real checkAdjointConsistencyJacobian_1(const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const bool printToStream=true, std::ostream &outStream=std::cout)
Check the consistency of the Jacobian and its adjoint. This is the primary interface.
std::vector< std::vector< Real > > checkApplyJacobian_1(const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &v, const Vector< Real > &jv, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
std::vector< std::vector< Real > > checkApplyJacobian_2(const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &v, const Vector< Real > &jv, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
virtual Real checkAdjointConsistencyJacobian_2(const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const bool printToStream=true, std::ostream &outStream=std::cout)
Check the consistency of the Jacobian and its adjoint. This is the primary interface.
std::vector< std::vector< Real > > checkApplyAdjointHessian_21(const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &p, const Vector< Real > &v, const Vector< Real > &hv, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
, , , ,
std::vector< std::vector< Real > > checkApplyAdjointHessian_11(const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &p, const Vector< Real > &v, const Vector< Real > &hv, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
std::vector< std::vector< Real > > checkApplyAdjointHessian_22(const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &p, const Vector< Real > &v, const Vector< Real > &hv, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
std::vector< std::vector< Real > > checkApplyAdjointHessian_12(const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &p, const Vector< Real > &v, const Vector< Real > &hv, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
, , , ,
virtual std::vector< std::vector< Real > > checkApplyJacobian(const Vector< Real > &x, const Vector< Real > &v, const Vector< Real > &jv, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
Finite-difference check for the constraint Jacobian application.
virtual std::vector< std::vector< Real > > checkApplyAdjointHessian(const Vector< Real > &x, const Vector< Real > &u, const Vector< Real > &v, const Vector< Real > &hv, const std::vector< Real > &step, const bool printToScreen=true, std::ostream &outStream=std::cout, const int order=1)
Finite-difference check for the application of the adjoint of constraint Hessian.
Provides the ROL::Vector interface for scalar values, to be used, for example, with scalar constraint...
Defines the linear algebra or vector space interface for simulation-based optimization.
Defines the linear algebra or vector space interface.
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...
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...
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...
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 .
void applyInverseAdjointJacobian_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 of the adjoint of the partial constraint Jacobian at , , to the vector .
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 .
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 ...
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...
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 .
void value(ROL::Vector< Real > &c, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Evaluate the constraint operator at .
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 ...
void value(ROL::Vector< Real > &c, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Evaluate the constraint operator at .
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...
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 ...
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...
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 .
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...
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 ...
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 .
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...
int main(int argc, char *argv[])
double RealT
constexpr auto dim