51 #include "Teuchos_GlobalMPISession.hpp"
61 ROL::Ptr<std::vector<Real> > cp
63 ROL::Ptr<const std::vector<Real> > up
65 ROL::Ptr<const std::vector<Real> > zp
68 Real half(0.5), two(2);
70 (*cp)[0] = (*up)[0]-(*zp)[0];
72 (*cp)[1] = half*std::pow((*up)[0]+(*up)[1]-(*zp)[0],two);
77 ROL::Ptr<std::vector<Real> > jvp
79 ROL::Ptr<const std::vector<Real> > vp
81 ROL::Ptr<const std::vector<Real> > up
83 ROL::Ptr<const std::vector<Real> > zp
86 (*jvp)[1] = ((*up)[0] + (*up)[1] - (*zp)[0]) * ((*vp)[0] + (*vp)[1]);
91 ROL::Ptr<std::vector<Real> > jvp
93 ROL::Ptr<const std::vector<Real> > vp
95 ROL::Ptr<const std::vector<Real> > up
97 ROL::Ptr<const std::vector<Real> > zp
99 (*jvp)[0] = -(*vp)[0];
100 (*jvp)[1] = ((*zp)[0] - (*up)[0] - (*up)[1]) * (*vp)[0];
105 ROL::Ptr<std::vector<Real> > ajvp
107 ROL::Ptr<const std::vector<Real> > vp
109 ROL::Ptr<const std::vector<Real> > up
111 ROL::Ptr<const std::vector<Real> > zp
113 (*ajvp)[0] = (*vp)[0] + ((*up)[0] + (*up)[1] - (*zp)[0]) * (*vp)[1];
114 (*ajvp)[1] = ((*up)[0] + (*up)[1] - (*zp)[0]) * (*vp)[1];
119 ROL::Ptr<std::vector<Real> > ajvp
121 ROL::Ptr<const std::vector<Real> > vp
123 ROL::Ptr<const std::vector<Real> > up
125 ROL::Ptr<const std::vector<Real> > zp
127 (*ajvp)[0] = ((*zp)[0] - (*up)[0] - (*up)[1]) * (*vp)[1] - (*vp)[0];
132 ROL::Ptr<std::vector<Real> > ahwvp
134 ROL::Ptr<const std::vector<Real> > wp
136 ROL::Ptr<const std::vector<Real> > vp
138 ROL::Ptr<const std::vector<Real> > up
140 ROL::Ptr<const std::vector<Real> > zp
142 (*ahwvp)[0] = (*wp)[1] * ((*vp)[0] + (*vp)[1]);
143 (*ahwvp)[1] = (*wp)[1] * ((*vp)[0] + (*vp)[1]);
148 ROL::Ptr<std::vector<Real> > ahwvp
150 ROL::Ptr<const std::vector<Real> > wp
152 ROL::Ptr<const std::vector<Real> > vp
154 ROL::Ptr<const std::vector<Real> > up
156 ROL::Ptr<const std::vector<Real> > zp
158 (*ahwvp)[0] = -(*wp)[1] * ((*vp)[0] + (*vp)[1]);
163 ROL::Ptr<std::vector<Real> > ahwvp
165 ROL::Ptr<const std::vector<Real> > wp
167 ROL::Ptr<const std::vector<Real> > vp
169 ROL::Ptr<const std::vector<Real> > up
171 ROL::Ptr<const std::vector<Real> > zp
173 (*ahwvp)[0] = -(*wp)[1] * (*vp)[0];
174 (*ahwvp)[1] = -(*wp)[1] * (*vp)[0];
179 ROL::Ptr<std::vector<Real> > ahwvp
181 ROL::Ptr<const std::vector<Real> > wp
183 ROL::Ptr<const std::vector<Real> > vp
185 ROL::Ptr<const std::vector<Real> > up
187 ROL::Ptr<const std::vector<Real> > zp
189 (*ahwvp)[0] = (*wp)[1] * (*vp)[0];
199 ROL::Ptr<std::vector<Real> > cp
201 ROL::Ptr<const std::vector<Real> > up
203 ROL::Ptr<const std::vector<Real> > zp
206 const Real one(1), two(2);
208 (*cp)[0] = std::exp((*up)[0])-(std::pow((*zp)[0],two) + one);
213 ROL::Ptr<std::vector<Real> > jvp
215 ROL::Ptr<const std::vector<Real> > vp
217 ROL::Ptr<const std::vector<Real> > up
219 ROL::Ptr<const std::vector<Real> > zp
221 (*jvp)[0] = std::exp((*up)[0]) * (*vp)[0];
226 ROL::Ptr<std::vector<Real> > jvp
228 ROL::Ptr<const std::vector<Real> > vp
230 ROL::Ptr<const std::vector<Real> > up
232 ROL::Ptr<const std::vector<Real> > zp
236 (*jvp)[0] = -two * (*zp)[0] * (*vp)[0];
241 ROL::Ptr<std::vector<Real> > ajvp
243 ROL::Ptr<const std::vector<Real> > vp
245 ROL::Ptr<const std::vector<Real> > up
247 ROL::Ptr<const std::vector<Real> > zp
249 (*ajvp)[0] = std::exp((*up)[0]) * (*vp)[0];
254 ROL::Ptr<std::vector<Real> > ajvp
256 ROL::Ptr<const std::vector<Real> > vp
258 ROL::Ptr<const std::vector<Real> > up
260 ROL::Ptr<const std::vector<Real> > zp
264 (*ajvp)[0] = -two * (*zp)[0] * (*vp)[0];
269 ROL::Ptr<std::vector<Real> > ijvp
271 ROL::Ptr<const std::vector<Real> > vp
273 ROL::Ptr<const std::vector<Real> > up
275 ROL::Ptr<const std::vector<Real> > zp
277 (*ijvp)[0] = (*vp)[0] / std::exp((*up)[0]);
282 ROL::Ptr<std::vector<Real> > ijvp
284 ROL::Ptr<const std::vector<Real> > vp
286 ROL::Ptr<const std::vector<Real> > up
288 ROL::Ptr<const std::vector<Real> > zp
290 (*ijvp)[0] = (*vp)[0] / std::exp((*up)[0]);
295 ROL::Ptr<std::vector<Real> > ahwvp
297 ROL::Ptr<const std::vector<Real> > wp
299 ROL::Ptr<const std::vector<Real> > vp
301 ROL::Ptr<const std::vector<Real> > up
303 ROL::Ptr<const std::vector<Real> > zp
305 (*ahwvp)[0] = std::exp((*up)[0]) * (*wp)[0] * (*vp)[0];
310 ROL::Ptr<std::vector<Real> > ahwvp
312 ROL::Ptr<const std::vector<Real> > wp
314 ROL::Ptr<const std::vector<Real> > vp
316 ROL::Ptr<const std::vector<Real> > up
318 ROL::Ptr<const std::vector<Real> > zp
320 (*ahwvp)[0] =
static_cast<Real
>(0);
325 ROL::Ptr<std::vector<Real> > ahwvp
327 ROL::Ptr<const std::vector<Real> > wp
329 ROL::Ptr<const std::vector<Real> > vp
331 ROL::Ptr<const std::vector<Real> > up
333 ROL::Ptr<const std::vector<Real> > zp
335 (*ahwvp)[0] =
static_cast<Real
>(0);
340 ROL::Ptr<std::vector<Real> > ahwvp
342 ROL::Ptr<const std::vector<Real> > wp
344 ROL::Ptr<const std::vector<Real> > vp
346 ROL::Ptr<const std::vector<Real> > up
348 ROL::Ptr<const std::vector<Real> > zp
350 (*ahwvp)[0] =
static_cast<Real
>(-2) * (*wp)[0] * (*vp)[0];
356 int main(
int argc,
char *argv[]) {
358 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
361 int iprint = argc - 1;
362 ROL::Ptr<std::ostream> outStream;
365 outStream = ROL::makePtrFromRef(std::cout);
367 outStream = ROL::makePtrFromRef(bhs);
377 ROL::Ptr<std::vector<RealT> > ustd = ROL::makePtr<std::vector<RealT>>(dim);
378 ROL::Ptr<std::vector<RealT> > dustd = ROL::makePtr<std::vector<RealT>>(dim);
379 ROL::Ptr<std::vector<RealT> > zstd = ROL::makePtr<std::vector<RealT>>(dimz);
380 ROL::Ptr<std::vector<RealT> > dzstd = ROL::makePtr<std::vector<RealT>>(dimz);
381 ROL::Ptr<std::vector<RealT> > cstd = ROL::makePtr<std::vector<RealT>>(dim);
382 ROL::Ptr<std::vector<RealT> > czstd = ROL::makePtr<std::vector<RealT>>(dimz);
383 ROL::Ptr<std::vector<RealT> > sstd = ROL::makePtr<std::vector<RealT>>(dimz);
384 ROL::Ptr<std::vector<RealT> > dsstd = ROL::makePtr<std::vector<RealT>>(dimz);
386 (*ustd)[0] =
static_cast<RealT>(rand())/static_cast<RealT>(RAND_MAX);
387 (*ustd)[1] =
static_cast<RealT>(rand())/static_cast<RealT>(RAND_MAX);
388 (*dustd)[0] =
static_cast<RealT>(rand())/static_cast<RealT>(RAND_MAX);
389 (*dustd)[1] =
static_cast<RealT>(rand())/static_cast<RealT>(RAND_MAX);
390 (*zstd)[0] =
static_cast<RealT>(rand())/static_cast<RealT>(RAND_MAX);
391 (*dzstd)[0] =
static_cast<RealT>(rand())/static_cast<RealT>(RAND_MAX);
392 (*cstd)[0] =
static_cast<RealT>(rand())/static_cast<RealT>(RAND_MAX);
393 (*cstd)[1] =
static_cast<RealT>(rand())/static_cast<RealT>(RAND_MAX);
394 (*czstd)[0] =
static_cast<RealT>(rand())/static_cast<RealT>(RAND_MAX);
395 (*sstd)[0] =
static_cast<RealT>(rand())/static_cast<RealT>(RAND_MAX);
396 (*dsstd)[0] =
static_cast<RealT>(rand())/static_cast<RealT>(RAND_MAX);
398 ROL::Ptr<ROL::Vector<RealT> > u = ROL::makePtr<ROL::StdVector<RealT>>(ustd);
399 ROL::Ptr<ROL::Vector<RealT> > du = ROL::makePtr<ROL::StdVector<RealT>>(dustd);
400 ROL::Ptr<ROL::Vector<RealT> > z = ROL::makePtr<ROL::StdVector<RealT>>(zstd);
401 ROL::Ptr<ROL::Vector<RealT> > dz = ROL::makePtr<ROL::StdVector<RealT>>(dzstd);
402 ROL::Ptr<ROL::Vector<RealT> > c = ROL::makePtr<ROL::StdVector<RealT>>(cstd);
403 ROL::Ptr<ROL::Vector<RealT> > cz = ROL::makePtr<ROL::StdVector<RealT>>(czstd);
404 ROL::Ptr<ROL::Vector<RealT> > s = ROL::makePtr<ROL::StdVector<RealT>>(sstd);
405 ROL::Ptr<ROL::Vector<RealT> > ds = ROL::makePtr<ROL::StdVector<RealT>>(dsstd);
414 ROL::Ptr<ROL::Constraint_SimOpt<RealT> > valCon = ROL::makePtr<valConstraint<RealT>>();
415 valCon->checkAdjointConsistencyJacobian_1(*c,*du,*u,*s,
true,*outStream);
416 valCon->checkAdjointConsistencyJacobian_2(*c,*dz,*u,*s,
true,*outStream);
417 valCon->checkApplyJacobian_1(*u,*s,*du,*c,
true,*outStream);
418 valCon->checkApplyJacobian_2(*u,*s,*ds,*c,
true,*outStream);
419 valCon->checkApplyJacobian(x,dx,*c,
true,*outStream);
420 valCon->checkApplyAdjointHessian_11(*u,*s,*c,*du,*u,
true,*outStream);
421 valCon->checkApplyAdjointHessian_12(*u,*s,*c,*du,*s,
true,*outStream);
422 valCon->checkApplyAdjointHessian_21(*u,*s,*c,*ds,*u,
true,*outStream);
423 valCon->checkApplyAdjointHessian_22(*u,*s,*c,*ds,*s,
true,*outStream);
424 valCon->checkApplyAdjointHessian(x,*c,dx,x,
true,*outStream);
426 ROL::Ptr<ROL::Constraint_SimOpt<RealT> > redCon = ROL::makePtr<redConstraint<RealT>>();
427 redCon->checkAdjointConsistencyJacobian_1(*cz,*ds,*s,*z,
true,*outStream);
428 redCon->checkAdjointConsistencyJacobian_2(*cz,*dz,*s,*z,
true,*outStream);
429 redCon->checkInverseJacobian_1(*cz,*ds,*s,*z,
true,*outStream);
430 redCon->checkInverseAdjointJacobian_1(*ds,*cz,*s,*z,
true,*outStream);
431 redCon->checkApplyJacobian_1(*s,*z,*ds,*cz,
true,*outStream);
432 redCon->checkApplyJacobian_2(*s,*z,*dz,*cz,
true,*outStream);
433 redCon->checkApplyJacobian(y,dy,*cz,
true,*outStream);
434 redCon->checkApplyAdjointHessian_11(*s,*z,*cz,*ds,*s,
true,*outStream);
435 redCon->checkApplyAdjointHessian_12(*s,*z,*cz,*ds,*z,
true,*outStream);
436 redCon->checkApplyAdjointHessian_21(*s,*z,*cz,*dz,*s,
true,*outStream);
437 redCon->checkApplyAdjointHessian_22(*s,*z,*cz,*dz,*z,
true,*outStream);
438 redCon->checkApplyAdjointHessian(y,*cz,dy,y,
true,*outStream);
452 catch (std::logic_error err) {
453 *outStream << err.what() <<
"\n";
458 std::cout <<
"End Result: TEST FAILED\n";
460 std::cout <<
"End Result: TEST PASSED\n";
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 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 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...
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 > > 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)
void value(ROL::Vector< Real > &c, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Evaluate the constraint operator at .
Defines the linear algebra or vector space interface for simulation-based optimization.
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 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 .
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)
, , , ,
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 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 ...
Defines the linear algebra or vector space interface.
Defines a no-output stream class ROL::NullStream and a function makeStreamPtr which either wraps a re...
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)
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 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 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 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...
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)
, , , ,
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 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_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...
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. ...
Defines a composite equality constraint operator interface for simulation-based optimization.
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.
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...
basic_nullstream< char, char_traits< char >> nullstream
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 .
int main(int argc, char *argv[])
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)
Defines the constraint operator interface for simulation-based optimization.
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 .
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_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)
void value(ROL::Vector< Real > &c, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Evaluate the constraint operator at .