10 #ifndef ROL_CONSTRAINT_SIMOPT_H
11 #define ROL_CONSTRAINT_SIMOPT_H
78 Ptr<Vector<Real>>
jv_;
199 Real cnorm = c.
norm();
200 const Real ctol = std::min(
atol_,
rtol_*cnorm);
208 Real alpha(1), tmp(0);
210 *stream << std::endl;
211 *stream <<
" Default Constraint_SimOpt::solve" << std::endl;
213 *stream << std::setw(6) << std::left <<
"iter";
214 *stream << std::setw(15) << std::left <<
"rnorm";
215 *stream << std::setw(15) << std::left <<
"alpha";
216 *stream << std::endl;
217 for (cnt = 0; cnt <
maxit_; ++cnt) {
226 while ( tmp > (one-
decr_*alpha)*cnorm &&
236 *stream << std::setw(6) << std::left << cnt;
237 *stream << std::scientific << std::setprecision(6);
238 *stream << std::setw(15) << std::left << tmp;
239 *stream << std::scientific << std::setprecision(6);
240 *stream << std::setw(15) << std::left << alpha;
241 *stream << std::endl;
246 if (cnorm <= ctol)
break;
251 Ptr<Constraint<Real>> con = makePtr<SimConstraint<Real>>(makePtrFromRef(*
this),makePtrFromRef(z),
true);
252 Ptr<Objective<Real>> obj = makePtr<NonlinearLeastSquaresObjective<Real>>(con,u,c,
true);
253 ParameterList parlist;
254 parlist.sublist(
"General").set(
"Output Level",1);
255 parlist.sublist(
"General").sublist(
"Krylov").set(
"Iteration Limit",100);
256 parlist.sublist(
"Step").sublist(
"Trust Region").set(
"Subproblem Solver",
"Truncated CG");
257 parlist.sublist(
"Status Test").set(
"Gradient Tolerance",ctol);
258 parlist.sublist(
"Status Test").set(
"Step Tolerance",
stol_);
259 parlist.sublist(
"Status Test").set(
"Iteration Limit",
maxit_);
260 Ptr<TypeU::Algorithm<Real>> algo = makePtr<TypeU::TrustRegionAlgorithm<Real>>(parlist);
261 algo->run(u,*obj,*stream);
265 Ptr<Constraint<Real>> con = makePtr<SimConstraint<Real>>(makePtrFromRef(*
this),makePtrFromRef(z),
true);
266 Ptr<Objective<Real>> obj = makePtr<Objective_FSsolver<Real>>();
267 Ptr<Vector<Real>> l = c.
dual().clone();
268 ParameterList parlist;
269 parlist.sublist(
"General").set(
"Output Level",1);
270 parlist.sublist(
"Status Test").set(
"Constraint Tolerance",ctol);
271 parlist.sublist(
"Status Test").set(
"Step Tolerance",
stol_);
272 parlist.sublist(
"Status Test").set(
"Iteration Limit",
maxit_);
273 Ptr<TypeE::Algorithm<Real>> algo = makePtr<TypeE::CompositeStepAlgorithm<Real>>(parlist);
274 algo->run(u,*obj,*con,*l,*stream);
278 ROL_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
279 ">>> ERROR (ROL:Constraint_SimOpt:solve): Invalid solver type!");
304 ParameterList & list = parlist.sublist(
"SimOpt").sublist(
"Solve");
336 Real ctol = std::sqrt(ROL_EPSILON<Real>());
339 if (v.
norm() > std::sqrt(ROL_EPSILON<Real>())) {
340 h = std::max(1.0,u.
norm()/v.
norm())*tol;
343 Ptr<Vector<Real>> unew = u.
clone();
348 value(jv,*unew,z,ctol);
350 Ptr<Vector<Real>> cold = jv.
clone();
352 value(*cold,u,z,ctol);
379 Real ctol = std::sqrt(ROL_EPSILON<Real>());
382 if (v.
norm() > std::sqrt(ROL_EPSILON<Real>())) {
383 h = std::max(1.0,u.
norm()/v.
norm())*tol;
386 Ptr<Vector<Real>> znew = z.
clone();
391 value(jv,u,*znew,ctol);
393 Ptr<Vector<Real>> cold = jv.
clone();
395 value(*cold,u,z,ctol);
421 ROL_TEST_FOR_EXCEPTION(
true, std::logic_error,
422 "The method applyInverseJacobian_1 is used but not implemented!\n");
472 Real ctol = std::sqrt(ROL_EPSILON<Real>());
474 if (v.
norm() > std::sqrt(ROL_EPSILON<Real>())) {
475 h = std::max(1.0,u.
norm()/v.
norm())*tol;
477 Ptr<Vector<Real>> cold = dualv.
clone();
478 Ptr<Vector<Real>> cnew = dualv.
clone();
480 value(*cold,u,z,ctol);
481 Ptr<Vector<Real>> unew = u.
clone();
483 for (
int i = 0; i < u.
dimension(); i++) {
485 unew->axpy(h,*(u.
basis(i)));
487 value(*cnew,*unew,z,ctol);
488 cnew->axpy(-1.0,*cold);
490 ajv.
axpy(cnew->dot(v),*((u.
dual()).basis(i)));
543 Real ctol = std::sqrt(ROL_EPSILON<Real>());
545 if (v.
norm() > std::sqrt(ROL_EPSILON<Real>())) {
546 h = std::max(1.0,u.
norm()/v.
norm())*tol;
548 Ptr<Vector<Real>> cold = dualv.
clone();
549 Ptr<Vector<Real>> cnew = dualv.
clone();
551 value(*cold,u,z,ctol);
552 Ptr<Vector<Real>> znew = z.
clone();
554 for (
int i = 0; i < z.
dimension(); i++) {
556 znew->axpy(h,*(z.
basis(i)));
558 value(*cnew,u,*znew,ctol);
559 cnew->axpy(-1.0,*cold);
561 ajv.
axpy(cnew->dot(v),*((z.
dual()).basis(i)));
586 ROL_TEST_FOR_EXCEPTION(
true, std::logic_error,
587 "The method applyInverseAdjointJacobian_1 is used but not implemented!\n");
613 Real jtol = std::sqrt(ROL_EPSILON<Real>());
616 if (v.
norm() > std::sqrt(ROL_EPSILON<Real>())) {
617 h = std::max(1.0,u.
norm()/v.
norm())*tol;
620 Ptr<Vector<Real>> unew = u.
clone();
626 Ptr<Vector<Real>> jv = ahwv.
clone();
658 Real jtol = std::sqrt(ROL_EPSILON<Real>());
661 if (v.
norm() > std::sqrt(ROL_EPSILON<Real>())) {
662 h = std::max(1.0,u.
norm()/v.
norm())*tol;
665 Ptr<Vector<Real>> unew = u.
clone();
671 Ptr<Vector<Real>> jv = ahwv.
clone();
703 Real jtol = std::sqrt(ROL_EPSILON<Real>());
706 if (v.
norm() > std::sqrt(ROL_EPSILON<Real>())) {
707 h = std::max(1.0,u.
norm()/v.
norm())*tol;
710 Ptr<Vector<Real>> znew = z.
clone();
716 Ptr<Vector<Real>> jv = ahwv.
clone();
747 Real jtol = std::sqrt(ROL_EPSILON<Real>());
750 if (v.
norm() > std::sqrt(ROL_EPSILON<Real>())) {
751 h = std::max(1.0,u.
norm()/v.
norm())*tol;
754 Ptr<Vector<Real>> znew = z.
clone();
760 Ptr<Vector<Real>> jv = ahwv.
clone();
840 Ptr<Vector<Real>> ijv = (xs.
get_1())->clone();
845 catch (
const std::logic_error &e) {
851 Ptr<Vector<Real>> ijv_dual = (gs.
get_1())->clone();
857 catch (
const std::logic_error &e) {
898 Ptr<Vector<Real>> jv2 = jv.
clone();
913 Ptr<Vector<Real>> ajv1 = (ajvs.
get_1())->clone();
916 Ptr<Vector<Real>> ajv2 = (ajvs.
get_2())->clone();
934 Ptr<Vector<Real>> C11 = (ahwvs.
get_1())->clone();
935 Ptr<Vector<Real>> C21 = (ahwvs.
get_1())->clone();
941 Ptr<Vector<Real>> C12 = (ahwvs.
get_2())->clone();
942 Ptr<Vector<Real>> C22 = (ahwvs.
get_2())->clone();
954 const bool printToStream =
true,
955 std::ostream & outStream = std::cout) {
957 Real tol = ROL_EPSILON<Real>();
958 Ptr<Vector<Real>> r = c.
clone();
959 Ptr<Vector<Real>> s = u.
clone();
962 Ptr<Vector<Real>> cs = c.
clone();
966 Real rnorm = r->norm();
967 Real cnorm = cs->norm();
968 if ( printToStream ) {
969 std::stringstream hist;
970 hist << std::scientific << std::setprecision(8);
971 hist <<
"\nTest SimOpt solve at feasible (u,z):\n";
972 hist <<
" Solver Residual = " << rnorm <<
"\n";
973 hist <<
" ||c(u,z)|| = " << cnorm <<
"\n";
974 outStream << hist.str();
996 const bool printToStream =
true,
997 std::ostream & outStream = std::cout) {
1023 const bool printToStream =
true,
1024 std::ostream & outStream = std::cout) {
1025 Real tol = ROL_EPSILON<Real>();
1026 Ptr<Vector<Real>> Jv = dualw.
clone();
1030 Real wJv = w.
apply(*Jv);
1031 Ptr<Vector<Real>> Jw = dualv.
clone();
1035 Real vJw = v.
apply(*Jw);
1036 Real diff = std::abs(wJv-vJw);
1037 if ( printToStream ) {
1038 std::stringstream hist;
1039 hist << std::scientific << std::setprecision(8);
1040 hist <<
"\nTest SimOpt consistency of Jacobian_1 and its adjoint: \n |<w,Jv> - <adj(J)w,v>| = "
1042 hist <<
" |<w,Jv>| = " << std::abs(wJv) <<
"\n";
1043 hist <<
" Relative Error = " << diff / (std::abs(wJv)+ROL_UNDERFLOW<Real>()) <<
"\n";
1044 outStream << hist.str();
1066 const bool printToStream =
true,
1067 std::ostream & outStream = std::cout) {
1092 const bool printToStream =
true,
1093 std::ostream & outStream = std::cout) {
1094 Real tol = ROL_EPSILON<Real>();
1095 Ptr<Vector<Real>> Jv = dualw.
clone();
1099 Real wJv = w.
apply(*Jv);
1100 Ptr<Vector<Real>> Jw = dualv.
clone();
1104 Real vJw = v.
apply(*Jw);
1105 Real diff = std::abs(wJv-vJw);
1106 if ( printToStream ) {
1107 std::stringstream hist;
1108 hist << std::scientific << std::setprecision(8);
1109 hist <<
"\nTest SimOpt consistency of Jacobian_2 and its adjoint: \n |<w,Jv> - <adj(J)w,v>| = "
1111 hist <<
" |<w,Jv>| = " << std::abs(wJv) <<
"\n";
1112 hist <<
" Relative Error = " << diff / (std::abs(wJv)+ROL_UNDERFLOW<Real>()) <<
"\n";
1113 outStream << hist.str();
1122 const bool printToStream =
true,
1123 std::ostream & outStream = std::cout) {
1124 Real tol = ROL_EPSILON<Real>();
1125 Ptr<Vector<Real>> Jv = jv.
clone();
1128 Ptr<Vector<Real>> iJJv = u.
clone();
1131 Ptr<Vector<Real>> diff = v.
clone();
1133 diff->axpy(-1.0,*iJJv);
1134 Real dnorm = diff->norm();
1135 Real vnorm = v.
norm();
1136 if ( printToStream ) {
1137 std::stringstream hist;
1138 hist << std::scientific << std::setprecision(8);
1139 hist <<
"\nTest SimOpt consistency of inverse Jacobian_1: \n ||v-inv(J)Jv|| = "
1141 hist <<
" ||v|| = " << vnorm <<
"\n";
1142 hist <<
" Relative Error = " << dnorm / (vnorm+ROL_UNDERFLOW<Real>()) <<
"\n";
1143 outStream << hist.str();
1152 const bool printToStream =
true,
1153 std::ostream & outStream = std::cout) {
1154 Real tol = ROL_EPSILON<Real>();
1155 Ptr<Vector<Real>> Jv = jv.
clone();
1158 Ptr<Vector<Real>> iJJv = v.
clone();
1161 Ptr<Vector<Real>> diff = v.
clone();
1163 diff->axpy(-1.0,*iJJv);
1164 Real dnorm = diff->norm();
1165 Real vnorm = v.
norm();
1166 if ( printToStream ) {
1167 std::stringstream hist;
1168 hist << std::scientific << std::setprecision(8);
1169 hist <<
"\nTest SimOpt consistency of inverse adjoint Jacobian_1: \n ||v-inv(adj(J))adj(J)v|| = "
1171 hist <<
" ||v|| = " << vnorm <<
"\n";
1172 hist <<
" Relative Error = " << dnorm / (vnorm+ROL_UNDERFLOW<Real>()) <<
"\n";
1173 outStream << hist.str();
1184 const bool printToStream =
true,
1185 std::ostream & outStream = std::cout,
1187 const int order = 1) {
1188 std::vector<Real> steps(numSteps);
1189 for(
int i=0;i<numSteps;++i) {
1190 steps[i] = pow(10,-i);
1203 const std::vector<Real> &steps,
1204 const bool printToStream =
true,
1205 std::ostream & outStream = std::cout,
1206 const int order = 1) {
1208 ROL_TEST_FOR_EXCEPTION( order<1 || order>4, std::invalid_argument,
1209 "Error: finite difference order must be 1,2,3, or 4" );
1216 Real tol = std::sqrt(ROL_EPSILON<Real>());
1218 int numSteps = steps.size();
1220 std::vector<Real> tmp(numVals);
1221 std::vector<std::vector<Real>> jvCheck(numSteps, tmp);
1225 oldFormatState.copyfmt(outStream);
1228 Ptr<Vector<Real>> c = jv.
clone();
1230 this->
value(*c, u, z, tol);
1233 Ptr<Vector<Real>> Jv = jv.
clone();
1235 Real normJv = Jv->norm();
1238 Ptr<Vector<Real>> cdif = jv.
clone();
1239 Ptr<Vector<Real>> cnew = jv.
clone();
1240 Ptr<Vector<Real>> unew = u.
clone();
1242 for (
int i=0; i<numSteps; i++) {
1244 Real eta = steps[i];
1249 cdif->scale(
weights[order-1][0]);
1251 for(
int j=0; j<order; ++j) {
1253 unew->axpy(eta*
shifts[order-1][j], v);
1255 if(
weights[order-1][j+1] != 0 ) {
1257 this->
value(*cnew,*unew,z,tol);
1258 cdif->axpy(
weights[order-1][j+1],*cnew);
1263 cdif->scale(one/eta);
1266 jvCheck[i][0] = eta;
1267 jvCheck[i][1] = normJv;
1268 jvCheck[i][2] = cdif->norm();
1269 cdif->axpy(-one, *Jv);
1270 jvCheck[i][3] = cdif->norm();
1272 if (printToStream) {
1273 std::stringstream hist;
1276 << std::setw(20) <<
"Step size"
1277 << std::setw(20) <<
"norm(Jac*vec)"
1278 << std::setw(20) <<
"norm(FD approx)"
1279 << std::setw(20) <<
"norm(abs error)"
1281 << std::setw(20) <<
"---------"
1282 << std::setw(20) <<
"-------------"
1283 << std::setw(20) <<
"---------------"
1284 << std::setw(20) <<
"---------------"
1287 hist << std::scientific << std::setprecision(11) << std::right
1288 << std::setw(20) << jvCheck[i][0]
1289 << std::setw(20) << jvCheck[i][1]
1290 << std::setw(20) << jvCheck[i][2]
1291 << std::setw(20) << jvCheck[i][3]
1293 outStream << hist.str();
1299 outStream.copyfmt(oldFormatState);
1309 const bool printToStream =
true,
1310 std::ostream & outStream = std::cout,
1312 const int order = 1) {
1313 std::vector<Real> steps(numSteps);
1314 for(
int i=0;i<numSteps;++i) {
1315 steps[i] = pow(10,-i);
1328 const std::vector<Real> &steps,
1329 const bool printToStream =
true,
1330 std::ostream & outStream = std::cout,
1331 const int order = 1) {
1333 ROL_TEST_FOR_EXCEPTION( order<1 || order>4, std::invalid_argument,
1334 "Error: finite difference order must be 1,2,3, or 4" );
1341 Real tol = std::sqrt(ROL_EPSILON<Real>());
1343 int numSteps = steps.size();
1345 std::vector<Real> tmp(numVals);
1346 std::vector<std::vector<Real>> jvCheck(numSteps, tmp);
1350 oldFormatState.copyfmt(outStream);
1353 Ptr<Vector<Real>> c = jv.
clone();
1355 this->
value(*c, u, z, tol);
1358 Ptr<Vector<Real>> Jv = jv.
clone();
1360 Real normJv = Jv->norm();
1363 Ptr<Vector<Real>> cdif = jv.
clone();
1364 Ptr<Vector<Real>> cnew = jv.
clone();
1365 Ptr<Vector<Real>> znew = z.
clone();
1367 for (
int i=0; i<numSteps; i++) {
1369 Real eta = steps[i];
1374 cdif->scale(
weights[order-1][0]);
1376 for(
int j=0; j<order; ++j) {
1378 znew->axpy(eta*
shifts[order-1][j], v);
1380 if(
weights[order-1][j+1] != 0 ) {
1382 this->
value(*cnew,u,*znew,tol);
1383 cdif->axpy(
weights[order-1][j+1],*cnew);
1388 cdif->scale(one/eta);
1391 jvCheck[i][0] = eta;
1392 jvCheck[i][1] = normJv;
1393 jvCheck[i][2] = cdif->norm();
1394 cdif->axpy(-one, *Jv);
1395 jvCheck[i][3] = cdif->norm();
1397 if (printToStream) {
1398 std::stringstream hist;
1401 << std::setw(20) <<
"Step size"
1402 << std::setw(20) <<
"norm(Jac*vec)"
1403 << std::setw(20) <<
"norm(FD approx)"
1404 << std::setw(20) <<
"norm(abs error)"
1406 << std::setw(20) <<
"---------"
1407 << std::setw(20) <<
"-------------"
1408 << std::setw(20) <<
"---------------"
1409 << std::setw(20) <<
"---------------"
1412 hist << std::scientific << std::setprecision(11) << std::right
1413 << std::setw(20) << jvCheck[i][0]
1414 << std::setw(20) << jvCheck[i][1]
1415 << std::setw(20) << jvCheck[i][2]
1416 << std::setw(20) << jvCheck[i][3]
1418 outStream << hist.str();
1424 outStream.copyfmt(oldFormatState);
1436 const bool printToStream =
true,
1437 std::ostream & outStream = std::cout,
1439 const int order = 1 ) {
1440 std::vector<Real> steps(numSteps);
1441 for(
int i=0;i<numSteps;++i) {
1442 steps[i] = pow(10,-i);
1454 const std::vector<Real> &steps,
1455 const bool printToStream =
true,
1456 std::ostream & outStream = std::cout,
1457 const int order = 1 ) {
1463 Real tol = std::sqrt(ROL_EPSILON<Real>());
1465 int numSteps = steps.size();
1467 std::vector<Real> tmp(numVals);
1468 std::vector<std::vector<Real>> ahpvCheck(numSteps, tmp);
1471 Ptr<Vector<Real>> AJdif = hv.
clone();
1472 Ptr<Vector<Real>> AJp = hv.
clone();
1473 Ptr<Vector<Real>> AHpv = hv.
clone();
1474 Ptr<Vector<Real>> AJnew = hv.
clone();
1475 Ptr<Vector<Real>> unew = u.
clone();
1479 oldFormatState.copyfmt(outStream);
1487 Real normAHpv = AHpv->norm();
1489 for (
int i=0; i<numSteps; i++) {
1491 Real eta = steps[i];
1497 AJdif->scale(
weights[order-1][0]);
1499 for(
int j=0; j<order; ++j) {
1501 unew->axpy(eta*
shifts[order-1][j],v);
1503 if(
weights[order-1][j+1] != 0 ) {
1506 AJdif->axpy(
weights[order-1][j+1],*AJnew);
1510 AJdif->scale(one/eta);
1513 ahpvCheck[i][0] = eta;
1514 ahpvCheck[i][1] = normAHpv;
1515 ahpvCheck[i][2] = AJdif->norm();
1516 AJdif->axpy(-one, *AHpv);
1517 ahpvCheck[i][3] = AJdif->norm();
1519 if (printToStream) {
1520 std::stringstream hist;
1523 << std::setw(20) <<
"Step size"
1524 << std::setw(20) <<
"norm(adj(H)(u,v))"
1525 << std::setw(20) <<
"norm(FD approx)"
1526 << std::setw(20) <<
"norm(abs error)"
1528 << std::setw(20) <<
"---------"
1529 << std::setw(20) <<
"-----------------"
1530 << std::setw(20) <<
"---------------"
1531 << std::setw(20) <<
"---------------"
1534 hist << std::scientific << std::setprecision(11) << std::right
1535 << std::setw(20) << ahpvCheck[i][0]
1536 << std::setw(20) << ahpvCheck[i][1]
1537 << std::setw(20) << ahpvCheck[i][2]
1538 << std::setw(20) << ahpvCheck[i][3]
1540 outStream << hist.str();
1546 outStream.copyfmt(oldFormatState);
1559 const bool printToStream =
true,
1560 std::ostream & outStream = std::cout,
1562 const int order = 1 ) {
1563 std::vector<Real> steps(numSteps);
1564 for(
int i=0;i<numSteps;++i) {
1565 steps[i] = pow(10,-i);
1580 const std::vector<Real> &steps,
1581 const bool printToStream =
true,
1582 std::ostream & outStream = std::cout,
1583 const int order = 1 ) {
1589 Real tol = std::sqrt(ROL_EPSILON<Real>());
1591 int numSteps = steps.size();
1593 std::vector<Real> tmp(numVals);
1594 std::vector<std::vector<Real>> ahpvCheck(numSteps, tmp);
1597 Ptr<Vector<Real>> AJdif = hv.
clone();
1598 Ptr<Vector<Real>> AJp = hv.
clone();
1599 Ptr<Vector<Real>> AHpv = hv.
clone();
1600 Ptr<Vector<Real>> AJnew = hv.
clone();
1601 Ptr<Vector<Real>> znew = z.
clone();
1605 oldFormatState.copyfmt(outStream);
1613 Real normAHpv = AHpv->norm();
1615 for (
int i=0; i<numSteps; i++) {
1617 Real eta = steps[i];
1623 AJdif->scale(
weights[order-1][0]);
1625 for(
int j=0; j<order; ++j) {
1627 znew->axpy(eta*
shifts[order-1][j],v);
1629 if(
weights[order-1][j+1] != 0 ) {
1632 AJdif->axpy(
weights[order-1][j+1],*AJnew);
1636 AJdif->scale(one/eta);
1639 ahpvCheck[i][0] = eta;
1640 ahpvCheck[i][1] = normAHpv;
1641 ahpvCheck[i][2] = AJdif->norm();
1642 AJdif->axpy(-one, *AHpv);
1643 ahpvCheck[i][3] = AJdif->norm();
1645 if (printToStream) {
1646 std::stringstream hist;
1649 << std::setw(20) <<
"Step size"
1650 << std::setw(20) <<
"norm(adj(H)(u,v))"
1651 << std::setw(20) <<
"norm(FD approx)"
1652 << std::setw(20) <<
"norm(abs error)"
1654 << std::setw(20) <<
"---------"
1655 << std::setw(20) <<
"-----------------"
1656 << std::setw(20) <<
"---------------"
1657 << std::setw(20) <<
"---------------"
1660 hist << std::scientific << std::setprecision(11) << std::right
1661 << std::setw(20) << ahpvCheck[i][0]
1662 << std::setw(20) << ahpvCheck[i][1]
1663 << std::setw(20) << ahpvCheck[i][2]
1664 << std::setw(20) << ahpvCheck[i][3]
1666 outStream << hist.str();
1672 outStream.copyfmt(oldFormatState);
1685 const bool printToStream =
true,
1686 std::ostream & outStream = std::cout,
1688 const int order = 1 ) {
1689 std::vector<Real> steps(numSteps);
1690 for(
int i=0;i<numSteps;++i) {
1691 steps[i] = pow(10,-i);
1704 const std::vector<Real> &steps,
1705 const bool printToStream =
true,
1706 std::ostream & outStream = std::cout,
1707 const int order = 1 ) {
1713 Real tol = std::sqrt(ROL_EPSILON<Real>());
1715 int numSteps = steps.size();
1717 std::vector<Real> tmp(numVals);
1718 std::vector<std::vector<Real>> ahpvCheck(numSteps, tmp);
1721 Ptr<Vector<Real>> AJdif = hv.
clone();
1722 Ptr<Vector<Real>> AJp = hv.
clone();
1723 Ptr<Vector<Real>> AHpv = hv.
clone();
1724 Ptr<Vector<Real>> AJnew = hv.
clone();
1725 Ptr<Vector<Real>> unew = u.
clone();
1729 oldFormatState.copyfmt(outStream);
1737 Real normAHpv = AHpv->norm();
1739 for (
int i=0; i<numSteps; i++) {
1741 Real eta = steps[i];
1747 AJdif->scale(
weights[order-1][0]);
1749 for(
int j=0; j<order; ++j) {
1751 unew->axpy(eta*
shifts[order-1][j],v);
1753 if(
weights[order-1][j+1] != 0 ) {
1756 AJdif->axpy(
weights[order-1][j+1],*AJnew);
1760 AJdif->scale(one/eta);
1763 ahpvCheck[i][0] = eta;
1764 ahpvCheck[i][1] = normAHpv;
1765 ahpvCheck[i][2] = AJdif->norm();
1766 AJdif->axpy(-one, *AHpv);
1767 ahpvCheck[i][3] = AJdif->norm();
1769 if (printToStream) {
1770 std::stringstream hist;
1773 << std::setw(20) <<
"Step size"
1774 << std::setw(20) <<
"norm(adj(H)(u,v))"
1775 << std::setw(20) <<
"norm(FD approx)"
1776 << std::setw(20) <<
"norm(abs error)"
1778 << std::setw(20) <<
"---------"
1779 << std::setw(20) <<
"-----------------"
1780 << std::setw(20) <<
"---------------"
1781 << std::setw(20) <<
"---------------"
1784 hist << std::scientific << std::setprecision(11) << std::right
1785 << std::setw(20) << ahpvCheck[i][0]
1786 << std::setw(20) << ahpvCheck[i][1]
1787 << std::setw(20) << ahpvCheck[i][2]
1788 << std::setw(20) << ahpvCheck[i][3]
1790 outStream << hist.str();
1796 outStream.copyfmt(oldFormatState);
1806 const bool printToStream =
true,
1807 std::ostream & outStream = std::cout,
1809 const int order = 1 ) {
1810 std::vector<Real> steps(numSteps);
1811 for(
int i=0;i<numSteps;++i) {
1812 steps[i] = pow(10,-i);
1824 const std::vector<Real> &steps,
1825 const bool printToStream =
true,
1826 std::ostream & outStream = std::cout,
1827 const int order = 1 ) {
1833 Real tol = std::sqrt(ROL_EPSILON<Real>());
1835 int numSteps = steps.size();
1837 std::vector<Real> tmp(numVals);
1838 std::vector<std::vector<Real>> ahpvCheck(numSteps, tmp);
1841 Ptr<Vector<Real>> AJdif = hv.
clone();
1842 Ptr<Vector<Real>> AJp = hv.
clone();
1843 Ptr<Vector<Real>> AHpv = hv.
clone();
1844 Ptr<Vector<Real>> AJnew = hv.
clone();
1845 Ptr<Vector<Real>> znew = z.
clone();
1849 oldFormatState.copyfmt(outStream);
1857 Real normAHpv = AHpv->norm();
1859 for (
int i=0; i<numSteps; i++) {
1861 Real eta = steps[i];
1867 AJdif->scale(
weights[order-1][0]);
1869 for(
int j=0; j<order; ++j) {
1871 znew->axpy(eta*
shifts[order-1][j],v);
1873 if(
weights[order-1][j+1] != 0 ) {
1876 AJdif->axpy(
weights[order-1][j+1],*AJnew);
1880 AJdif->scale(one/eta);
1883 ahpvCheck[i][0] = eta;
1884 ahpvCheck[i][1] = normAHpv;
1885 ahpvCheck[i][2] = AJdif->norm();
1886 AJdif->axpy(-one, *AHpv);
1887 ahpvCheck[i][3] = AJdif->norm();
1889 if (printToStream) {
1890 std::stringstream hist;
1893 << std::setw(20) <<
"Step size"
1894 << std::setw(20) <<
"norm(adj(H)(u,v))"
1895 << std::setw(20) <<
"norm(FD approx)"
1896 << std::setw(20) <<
"norm(abs error)"
1898 << std::setw(20) <<
"---------"
1899 << std::setw(20) <<
"-----------------"
1900 << std::setw(20) <<
"---------------"
1901 << std::setw(20) <<
"---------------"
1904 hist << std::scientific << std::setprecision(11) << std::right
1905 << std::setw(20) << ahpvCheck[i][0]
1906 << std::setw(20) << ahpvCheck[i][1]
1907 << std::setw(20) << ahpvCheck[i][2]
1908 << std::setw(20) << ahpvCheck[i][3]
1910 outStream << hist.str();
1916 outStream.copyfmt(oldFormatState);
virtual void value(Vector< Real > &c, const Vector< Real > &x, Real &tol)
Evaluate the constraint operator at .
virtual const Vector & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
virtual void applyPreconditioner(Vector< Real > &pv, const Vector< Real > &v, const Vector< Real > &x, const Vector< Real > &g, Real &tol)
Apply a constraint preconditioner at , , to vector . Ideally, this preconditioner satisfies the follo...
virtual void scale(const Real alpha)=0
Compute where .
virtual void applyJacobian_2(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Apply the partial constraint Jacobian at , , to the vector .
virtual void applyAdjointHessian(Vector< Real > &ahwv, const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply the derivative of the adjoint of the constraint Jacobian at to vector in direction ...
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
ROL::Ptr< const Vector< Real > > get_2() const
virtual int dimension() const
Return dimension of the vector space.
virtual Real apply(const Vector< Real > &x) const
Apply to a dual vector. This is equivalent to the call .
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 ROL::Ptr< Vector > basis(const int i) const
Return i-th basis vector.
virtual void plus(const Vector &x)=0
Compute , where .
const double weights[4][5]
virtual void solve_update(const Vector< Real > &u, const Vector< Real > &z, UpdateType type, int iter=-1)
Update SimOpt constraint during solve (disconnected from optimization updates).
virtual void axpy(const Real alpha, const Vector &x)
Compute where .
virtual void applyJacobian_1(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &u, const 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 std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
virtual Real checkSolve(const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &c, const bool printToStream=true, std::ostream &outStream=std::cout)
virtual void update(const Vector< Real > &x, UpdateType type, int iter=-1)
Update constraint function.
Defines the linear algebra or vector space interface for simulation-based optimization.
virtual void update_1(const Vector< Real > &u, UpdateType type, int iter=-1)
virtual void applyAdjointHessian_11(Vector< Real > &ahwv, const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Apply the simulation-space derivative of the adjoint of the constraint simulation-space Jacobian at ...
virtual Real checkAdjointConsistencyJacobian_2(const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &dualw, const Vector< Real > &dualv, const bool printToStream=true, std::ostream &outStream=std::cout)
Check the consistency of the Jacobian and its adjoint. This is the secondary interface, for use with dual spaces where the user does not define the dual() operation.
virtual void update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update constraint functions. x is the optimization variable, flag = true if optimization variable is ...
virtual void update_1(const Vector< Real > &u, bool flag=true, int iter=-1)
Update constraint functions with respect to Sim variable. x is the optimization variable, flag = true if optimization variable is changed, iter is the outer algorithm iterations count.
virtual void applyAdjointHessian_22(Vector< Real > &ahwv, const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Apply the optimization-space derivative of the adjoint of the constraint optimization-space Jacobian ...
Contains definitions of custom data types in ROL.
void set_1(const Vector< Real > &vec)
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)
virtual void applyJacobian(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply the constraint Jacobian at , , to vector .
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.
virtual void update(const Vector< Real > &u, const Vector< Real > &z, UpdateType type, int iter=-1)
virtual void applyPreconditioner(Vector< Real > &pv, const Vector< Real > &v, const Vector< Real > &x, const Vector< Real > &g, Real &tol)
Apply a constraint preconditioner at , , to vector . In general, this preconditioner satisfies the fo...
Ptr< std::ostream > makeStreamPtr(std::ostream &os, bool noSuppressOutput=true)
virtual void zero()
Set to zero vector.
Defines the linear algebra or vector space interface.
virtual void update(const Vector< Real > &u, const Vector< Real > &z, bool flag=true, int iter=-1)
Update constraint functions. x is the optimization variable, flag = true if optimization variable is ...
basic_nullstream< char, std::char_traits< char >> nullstream
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)
Ptr< Vector< Real > > jv_
virtual Real checkInverseJacobian_1(const Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const bool printToStream=true, std::ostream &outStream=std::cout)
virtual void applyInverseJacobian_1(Vector< Real > &ijv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Apply the inverse partial constraint Jacobian at , , to the vector .
const int DEFAULT_solverType_
virtual void applyAdjointHessian_12(Vector< Real > &ahwv, const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Apply the optimization-space derivative of the adjoint of the constraint simulation-space Jacobian at...
virtual std::vector< Real > solveAugmentedSystem(Vector< Real > &v1, Vector< Real > &v2, const Vector< Real > &b1, const Vector< Real > &b2, const Vector< Real > &x, Real &tol)
Approximately solves the augmented system where , , , , is an identity operator, and is a zero operator.
virtual void applyAdjointJacobian_1(Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Apply the adjoint of the partial constraint Jacobian at , , to the vector . This is the primary inter...
virtual void value(Vector< Real > &c, const Vector< Real > &u, const Vector< Real > &z, Real &tol)=0
Evaluate the constraint operator at .
virtual void applyAdjointJacobian_2(Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &dualv, Real &tol)
Apply the adjoint of the partial constraint Jacobian at , , to vector . This is the secondary interfa...
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 std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, 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 std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
virtual void applyAdjointJacobian_1(Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &dualv, Real &tol)
Apply the adjoint of the partial constraint Jacobian at , , to the vector . This is the secondary int...
virtual void applyAdjointJacobian_2(Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Apply the adjoint of the partial constraint Jacobian at , , to vector . This is the primary interface...
const Real DEFAULT_factor_
#define ROL_NUM_CHECKDERIV_STEPS
Number of steps for derivative checks.
virtual void applyAdjointHessian_21(Vector< Real > &ahwv, const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Apply the simulation-space derivative of the adjoint of the constraint optimization-space Jacobian at...
virtual void solve(Vector< Real > &c, Vector< Real > &u, const Vector< Real > &z, Real &tol)
Given , solve for .
virtual std::vector< Real > solveAugmentedSystem(Vector< Real > &v1, Vector< Real > &v2, const Vector< Real > &b1, const Vector< Real > &b2, const Vector< Real > &x, Real &tol)
Approximately solves the augmented system where , , , , is an identity or Riesz operator...
virtual void update_2(const Vector< Real > &z, UpdateType type, int iter=-1)
std::vector< std::vector< Real > > checkApplyJacobian_1(const Vector< Real > &u, const Vector< Real > &z, 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)
virtual Real checkAdjointConsistencyJacobian_1(const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &dualw, const Vector< Real > &dualv, const bool printToStream=true, std::ostream &outStream=std::cout)
Check the consistency of the Jacobian and its adjoint. This is the secondary interface, for use with dual spaces where the user does not define the dual() operation.
virtual void setSolveParameters(ParameterList &parlist)
Set solve parameters.
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.
virtual void set(const Vector &x)
Set where .
virtual Real norm() const =0
Returns where .
virtual void update_2(const Vector< Real > &z, bool flag=true, int iter=-1)
Update constraint functions with respect to Opt variable. x is the optimization variable, flag = true if optimization variable is changed, iter is the outer algorithm iterations count.
Ptr< Vector< Real > > unew_
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)
Real ROL_EPSILON(void)
Platform-dependent machine epsilon.
Defines the constraint operator interface for simulation-based optimization.
virtual void applyAdjointJacobian(Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply the adjoint of the the constraint Jacobian at , , to vector .
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_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)
, , , ,
virtual void applyInverseAdjointJacobian_1(Vector< Real > &iajv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Apply the inverse of the adjoint of the partial constraint Jacobian at , , to the vector ...
void set_2(const Vector< Real > &vec)
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 std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
, , , ,
Defines the general constraint operator interface.
const bool DEFAULT_print_
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 std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
virtual Real checkInverseAdjointJacobian_1(const Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const bool printToStream=true, std::ostream &outStream=std::cout)
const Vector< Real > & dual(void) const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
ROL::Ptr< const Vector< Real > > get_1() const