17 #include "Teuchos_GlobalMPISession.hpp"
27 ROL::Ptr<std::vector<Real> > cp
29 ROL::Ptr<const std::vector<Real> > up
31 ROL::Ptr<const std::vector<Real> > zp
34 Real half(0.5), two(2);
36 (*cp)[0] = (*up)[0]-(*zp)[0];
38 (*cp)[1] = half*std::pow((*up)[0]+(*up)[1]-(*zp)[0],two);
43 ROL::Ptr<std::vector<Real> > jvp
45 ROL::Ptr<const std::vector<Real> > vp
47 ROL::Ptr<const std::vector<Real> > up
49 ROL::Ptr<const std::vector<Real> > zp
52 (*jvp)[1] = ((*up)[0] + (*up)[1] - (*zp)[0]) * ((*vp)[0] + (*vp)[1]);
57 ROL::Ptr<std::vector<Real> > jvp
59 ROL::Ptr<const std::vector<Real> > vp
61 ROL::Ptr<const std::vector<Real> > up
63 ROL::Ptr<const std::vector<Real> > zp
65 (*jvp)[0] = -(*vp)[0];
66 (*jvp)[1] = ((*zp)[0] - (*up)[0] - (*up)[1]) * (*vp)[0];
71 ROL::Ptr<std::vector<Real> > ajvp
73 ROL::Ptr<const std::vector<Real> > vp
75 ROL::Ptr<const std::vector<Real> > up
77 ROL::Ptr<const std::vector<Real> > zp
79 (*ajvp)[0] = (*vp)[0] + ((*up)[0] + (*up)[1] - (*zp)[0]) * (*vp)[1];
80 (*ajvp)[1] = ((*up)[0] + (*up)[1] - (*zp)[0]) * (*vp)[1];
85 ROL::Ptr<std::vector<Real> > ajvp
87 ROL::Ptr<const std::vector<Real> > vp
89 ROL::Ptr<const std::vector<Real> > up
91 ROL::Ptr<const std::vector<Real> > zp
93 (*ajvp)[0] = ((*zp)[0] - (*up)[0] - (*up)[1]) * (*vp)[1] - (*vp)[0];
98 ROL::Ptr<std::vector<Real> > ahwvp
100 ROL::Ptr<const std::vector<Real> > wp
102 ROL::Ptr<const std::vector<Real> > vp
104 ROL::Ptr<const std::vector<Real> > up
106 ROL::Ptr<const std::vector<Real> > zp
108 (*ahwvp)[0] = (*wp)[1] * ((*vp)[0] + (*vp)[1]);
109 (*ahwvp)[1] = (*wp)[1] * ((*vp)[0] + (*vp)[1]);
114 ROL::Ptr<std::vector<Real> > ahwvp
116 ROL::Ptr<const std::vector<Real> > wp
118 ROL::Ptr<const std::vector<Real> > vp
120 ROL::Ptr<const std::vector<Real> > up
122 ROL::Ptr<const std::vector<Real> > zp
124 (*ahwvp)[0] = -(*wp)[1] * ((*vp)[0] + (*vp)[1]);
129 ROL::Ptr<std::vector<Real> > ahwvp
131 ROL::Ptr<const std::vector<Real> > wp
133 ROL::Ptr<const std::vector<Real> > vp
135 ROL::Ptr<const std::vector<Real> > up
137 ROL::Ptr<const std::vector<Real> > zp
139 (*ahwvp)[0] = -(*wp)[1] * (*vp)[0];
140 (*ahwvp)[1] = -(*wp)[1] * (*vp)[0];
145 ROL::Ptr<std::vector<Real> > ahwvp
147 ROL::Ptr<const std::vector<Real> > wp
149 ROL::Ptr<const std::vector<Real> > vp
151 ROL::Ptr<const std::vector<Real> > up
153 ROL::Ptr<const std::vector<Real> > zp
155 (*ahwvp)[0] = (*wp)[1] * (*vp)[0];
165 ROL::Ptr<std::vector<Real> > cp
167 ROL::Ptr<const std::vector<Real> > up
169 ROL::Ptr<const std::vector<Real> > zp
172 const Real one(1), two(2);
174 (*cp)[0] = std::exp((*up)[0])-(std::pow((*zp)[0],two) + one);
179 ROL::Ptr<std::vector<Real> > jvp
181 ROL::Ptr<const std::vector<Real> > vp
183 ROL::Ptr<const std::vector<Real> > up
185 ROL::Ptr<const std::vector<Real> > zp
187 (*jvp)[0] = std::exp((*up)[0]) * (*vp)[0];
192 ROL::Ptr<std::vector<Real> > jvp
194 ROL::Ptr<const std::vector<Real> > vp
196 ROL::Ptr<const std::vector<Real> > up
198 ROL::Ptr<const std::vector<Real> > zp
202 (*jvp)[0] = -two * (*zp)[0] * (*vp)[0];
207 ROL::Ptr<std::vector<Real> > ajvp
209 ROL::Ptr<const std::vector<Real> > vp
211 ROL::Ptr<const std::vector<Real> > up
213 ROL::Ptr<const std::vector<Real> > zp
215 (*ajvp)[0] = std::exp((*up)[0]) * (*vp)[0];
220 ROL::Ptr<std::vector<Real> > ajvp
222 ROL::Ptr<const std::vector<Real> > vp
224 ROL::Ptr<const std::vector<Real> > up
226 ROL::Ptr<const std::vector<Real> > zp
230 (*ajvp)[0] = -two * (*zp)[0] * (*vp)[0];
235 ROL::Ptr<std::vector<Real> > ijvp
237 ROL::Ptr<const std::vector<Real> > vp
239 ROL::Ptr<const std::vector<Real> > up
241 ROL::Ptr<const std::vector<Real> > zp
243 (*ijvp)[0] = (*vp)[0] / std::exp((*up)[0]);
248 ROL::Ptr<std::vector<Real> > ijvp
250 ROL::Ptr<const std::vector<Real> > vp
252 ROL::Ptr<const std::vector<Real> > up
254 ROL::Ptr<const std::vector<Real> > zp
256 (*ijvp)[0] = (*vp)[0] / std::exp((*up)[0]);
261 ROL::Ptr<std::vector<Real> > ahwvp
263 ROL::Ptr<const std::vector<Real> > wp
265 ROL::Ptr<const std::vector<Real> > vp
267 ROL::Ptr<const std::vector<Real> > up
269 ROL::Ptr<const std::vector<Real> > zp
271 (*ahwvp)[0] = std::exp((*up)[0]) * (*wp)[0] * (*vp)[0];
276 ROL::Ptr<std::vector<Real> > ahwvp
278 ROL::Ptr<const std::vector<Real> > wp
280 ROL::Ptr<const std::vector<Real> > vp
282 ROL::Ptr<const std::vector<Real> > up
284 ROL::Ptr<const std::vector<Real> > zp
286 (*ahwvp)[0] =
static_cast<Real
>(0);
291 ROL::Ptr<std::vector<Real> > ahwvp
293 ROL::Ptr<const std::vector<Real> > wp
295 ROL::Ptr<const std::vector<Real> > vp
297 ROL::Ptr<const std::vector<Real> > up
299 ROL::Ptr<const std::vector<Real> > zp
301 (*ahwvp)[0] =
static_cast<Real
>(0);
306 ROL::Ptr<std::vector<Real> > ahwvp
308 ROL::Ptr<const std::vector<Real> > wp
310 ROL::Ptr<const std::vector<Real> > vp
312 ROL::Ptr<const std::vector<Real> > up
314 ROL::Ptr<const std::vector<Real> > zp
316 (*ahwvp)[0] =
static_cast<Real
>(-2) * (*wp)[0] * (*vp)[0];
322 int main(
int argc,
char *argv[]) {
324 Teuchos::GlobalMPISession mpiSession(&argc, &argv);
327 int iprint = argc - 1;
328 ROL::Ptr<std::ostream> outStream;
331 outStream = ROL::makePtrFromRef(std::cout);
333 outStream = ROL::makePtrFromRef(bhs);
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);
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);
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);
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);
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);
418 catch (std::logic_error& err) {
419 *outStream << err.what() <<
"\n";
424 std::cout <<
"End Result: TEST FAILED\n";
426 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)
basic_nullstream< char, std::char_traits< char >> nullstream
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...
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 .