10 #ifndef ROL_OBJECTIVE_SIMOPT_H
11 #define ROL_OBJECTIVE_SIMOPT_H
65 Real ftol = std::sqrt(ROL_EPSILON<Real>());
68 Real v = this->
value(u,z,ftol);
70 ROL::Ptr<Vector<Real> > unew = u.
clone();
75 unew->axpy(h,*(u.
basis(i)));
77 deriv = (this->
value(*unew,z,ftol) - v)/h;
85 Real ftol = std::sqrt(ROL_EPSILON<Real>());
88 Real v = this->
value(u,z,ftol);
90 ROL::Ptr<Vector<Real> > znew = z.
clone();
95 znew->axpy(h,*(z.
basis(i)));
97 deriv = (this->
value(u,*znew,ftol) - v)/h;
108 ROL::Ptr<Vector<Real> > g1 = gs.
get_1()->clone();
109 ROL::Ptr<Vector<Real> > g2 = gs.
get_2()->clone();
121 Real gtol = std::sqrt(ROL_EPSILON<Real>());
124 if (v.
norm() > std::sqrt(ROL_EPSILON<Real>())) {
125 h = std::max(1.0,u.
norm()/v.
norm())*tol;
128 ROL::Ptr<Vector<Real> > unew = u.
clone();
135 ROL::Ptr<Vector<Real> > g = hv.
clone();
145 Real gtol = std::sqrt(ROL_EPSILON<Real>());
148 if (v.
norm() > std::sqrt(ROL_EPSILON<Real>())) {
149 h = std::max(1.0,u.
norm()/v.
norm())*tol;
152 ROL::Ptr<Vector<Real> > znew = z.
clone();
159 ROL::Ptr<Vector<Real> > g = hv.
clone();
169 Real gtol = std::sqrt(ROL_EPSILON<Real>());
172 if (v.
norm() > std::sqrt(ROL_EPSILON<Real>())) {
173 h = std::max(1.0,u.
norm()/v.
norm())*tol;
176 ROL::Ptr<Vector<Real> > unew = u.
clone();
183 ROL::Ptr<Vector<Real> > g = hv.
clone();
193 Real gtol = std::sqrt(ROL_EPSILON<Real>());
196 if (v.
norm() > std::sqrt(ROL_EPSILON<Real>())) {
197 h = std::max(1.0,u.
norm()/v.
norm())*tol;
200 ROL::Ptr<Vector<Real> > znew = z.
clone();
207 ROL::Ptr<Vector<Real> > g = hv.
clone();
222 ROL::Ptr<Vector<Real> > h11 = (hvs.
get_1())->clone();
224 ROL::Ptr<Vector<Real> > h12 = (hvs.
get_1())->clone();
226 ROL::Ptr<Vector<Real> > h21 = (hvs.
get_2())->clone();
228 ROL::Ptr<Vector<Real> > h22 = (hvs.
get_2())->clone();
239 const bool printToStream =
true,
240 std::ostream & outStream = std::cout,
242 const int order = 1 ) {
250 const bool printToStream,
251 std::ostream & outStream,
254 std::vector<Real> steps(numSteps);
255 for(
int i=0;i<numSteps;++i) {
256 steps[i] = pow(10,-i);
266 const std::vector<Real> &steps,
267 const bool printToStream,
268 std::ostream & outStream,
270 ROL_TEST_FOR_EXCEPTION( order<1 || order>4, std::invalid_argument,
271 "Error: finite difference order must be 1,2,3, or 4" );
276 Real tol = std::sqrt(ROL_EPSILON<Real>());
278 int numSteps = steps.size();
280 std::vector<Real> tmp(numVals);
281 std::vector<std::vector<Real> > gCheck(numSteps, tmp);
285 oldFormatState.copyfmt(outStream);
289 Real val = this->
value(u,z,tol);
292 ROL::Ptr<Vector<Real> > gtmp = g.
clone();
295 Real dtg = d.
apply(*gtmp);
298 ROL::Ptr<Vector<Real> > unew = u.
clone();
300 for (
int i=0; i<numSteps; i++) {
310 gCheck[i][2] =
weights[order-1][0] * val;
312 for(
int j=0; j<order; ++j) {
314 unew->axpy(eta*
shifts[order-1][j], d);
317 if(
weights[order-1][j+1] != 0 ) {
319 gCheck[i][2] +=
weights[order-1][j+1] * this->
value(*unew,z,tol);
324 gCheck[i][3] = std::abs(gCheck[i][2] - gCheck[i][1]);
328 outStream << std::right
329 << std::setw(20) <<
"Step size"
330 << std::setw(20) <<
"grad'*dir"
331 << std::setw(20) <<
"FD approx"
332 << std::setw(20) <<
"abs error"
334 << std::setw(20) <<
"---------"
335 << std::setw(20) <<
"---------"
336 << std::setw(20) <<
"---------"
337 << std::setw(20) <<
"---------"
340 outStream << std::scientific << std::setprecision(11) << std::right
341 << std::setw(20) << gCheck[i][0]
342 << std::setw(20) << gCheck[i][1]
343 << std::setw(20) << gCheck[i][2]
344 << std::setw(20) << gCheck[i][3]
351 outStream.copyfmt(oldFormatState);
360 const bool printToStream =
true,
361 std::ostream & outStream = std::cout,
363 const int order = 1 ) {
371 const bool printToStream,
372 std::ostream & outStream,
375 std::vector<Real> steps(numSteps);
376 for(
int i=0;i<numSteps;++i) {
377 steps[i] = pow(10,-i);
387 const std::vector<Real> &steps,
388 const bool printToStream,
389 std::ostream & outStream,
391 ROL_TEST_FOR_EXCEPTION( order<1 || order>4, std::invalid_argument,
392 "Error: finite difference order must be 1,2,3, or 4" );
397 Real tol = std::sqrt(ROL_EPSILON<Real>());
399 int numSteps = steps.size();
401 std::vector<Real> tmp(numVals);
402 std::vector<std::vector<Real> > gCheck(numSteps, tmp);
406 oldFormatState.copyfmt(outStream);
410 Real val = this->
value(u,z,tol);
413 ROL::Ptr<Vector<Real> > gtmp = g.
clone();
416 Real dtg = d.
apply(*gtmp);
419 ROL::Ptr<Vector<Real> > znew = z.
clone();
421 for (
int i=0; i<numSteps; i++) {
431 gCheck[i][2] =
weights[order-1][0] * val;
433 for(
int j=0; j<order; ++j) {
435 znew->axpy(eta*
shifts[order-1][j], d);
438 if(
weights[order-1][j+1] != 0 ) {
440 gCheck[i][2] +=
weights[order-1][j+1] * this->
value(u,*znew,tol);
445 gCheck[i][3] = std::abs(gCheck[i][2] - gCheck[i][1]);
449 outStream << std::right
450 << std::setw(20) <<
"Step size"
451 << std::setw(20) <<
"grad'*dir"
452 << std::setw(20) <<
"FD approx"
453 << std::setw(20) <<
"abs error"
455 << std::setw(20) <<
"---------"
456 << std::setw(20) <<
"---------"
457 << std::setw(20) <<
"---------"
458 << std::setw(20) <<
"---------"
461 outStream << std::scientific << std::setprecision(11) << std::right
462 << std::setw(20) << gCheck[i][0]
463 << std::setw(20) << gCheck[i][1]
464 << std::setw(20) << gCheck[i][2]
465 << std::setw(20) << gCheck[i][3]
472 outStream.copyfmt(oldFormatState);
481 const bool printToStream =
true,
482 std::ostream & outStream = std::cout,
484 const int order = 1 ) {
493 const std::vector<Real> &steps,
494 const bool printToStream =
true,
495 std::ostream & outStream = std::cout,
496 const int order = 1 ) {
506 const bool printToStream,
507 std::ostream & outStream,
510 std::vector<Real> steps(numSteps);
511 for(
int i=0;i<numSteps;++i) {
512 steps[i] = pow(10,-i);
523 const std::vector<Real> &steps,
524 const bool printToStream,
525 std::ostream & outStream,
528 ROL_TEST_FOR_EXCEPTION( order<1 || order>4, std::invalid_argument,
529 "Error: finite difference order must be 1,2,3, or 4" );
535 Real tol = std::sqrt(ROL_EPSILON<Real>());
537 int numSteps = steps.size();
539 std::vector<Real> tmp(numVals);
540 std::vector<std::vector<Real> > hvCheck(numSteps, tmp);
544 oldFormatState.copyfmt(outStream);
547 ROL::Ptr<Vector<Real> > g = hv.
clone();
552 ROL::Ptr<Vector<Real> > Hv = hv.
clone();
554 Real normHv = Hv->norm();
557 ROL::Ptr<Vector<Real> > gdif = hv.
clone();
558 ROL::Ptr<Vector<Real> > gnew = hv.
clone();
559 ROL::Ptr<Vector<Real> > unew = u.
clone();
561 for (
int i=0; i<numSteps; i++) {
569 gdif->scale(
weights[order-1][0]);
571 for(
int j=0; j<order; ++j) {
574 unew->axpy(eta*
shifts[order-1][j], v);
577 if(
weights[order-1][j+1] != 0 ) {
580 gdif->axpy(
weights[order-1][j+1],*gnew);
585 gdif->scale(1.0/eta);
589 hvCheck[i][1] = normHv;
590 hvCheck[i][2] = gdif->norm();
591 gdif->axpy(-1.0, *Hv);
592 hvCheck[i][3] = gdif->norm();
596 outStream << std::right
597 << std::setw(20) <<
"Step size"
598 << std::setw(20) <<
"norm(Hess*vec)"
599 << std::setw(20) <<
"norm(FD approx)"
600 << std::setw(20) <<
"norm(abs error)"
602 << std::setw(20) <<
"---------"
603 << std::setw(20) <<
"--------------"
604 << std::setw(20) <<
"---------------"
605 << std::setw(20) <<
"---------------"
608 outStream << std::scientific << std::setprecision(11) << std::right
609 << std::setw(20) << hvCheck[i][0]
610 << std::setw(20) << hvCheck[i][1]
611 << std::setw(20) << hvCheck[i][2]
612 << std::setw(20) << hvCheck[i][3]
619 outStream.copyfmt(oldFormatState);
628 const bool printToStream =
true,
629 std::ostream & outStream = std::cout,
631 const int order = 1 ) {
638 const std::vector<Real> &steps,
639 const bool printToStream =
true,
640 std::ostream & outStream = std::cout,
641 const int order = 1 ) {
650 const bool printToStream,
651 std::ostream & outStream,
654 std::vector<Real> steps(numSteps);
655 for(
int i=0;i<numSteps;++i) {
656 steps[i] = pow(10,-i);
667 const std::vector<Real> &steps,
668 const bool printToStream,
669 std::ostream & outStream,
672 ROL_TEST_FOR_EXCEPTION( order<1 || order>4, std::invalid_argument,
673 "Error: finite difference order must be 1,2,3, or 4" );
679 Real tol = std::sqrt(ROL_EPSILON<Real>());
681 int numSteps = steps.size();
683 std::vector<Real> tmp(numVals);
684 std::vector<std::vector<Real> > hvCheck(numSteps, tmp);
688 oldFormatState.copyfmt(outStream);
691 ROL::Ptr<Vector<Real> > g = hv.
clone();
696 ROL::Ptr<Vector<Real> > Hv = hv.
clone();
698 Real normHv = Hv->norm();
701 ROL::Ptr<Vector<Real> > gdif = hv.
clone();
702 ROL::Ptr<Vector<Real> > gnew = hv.
clone();
703 ROL::Ptr<Vector<Real> > znew = z.
clone();
705 for (
int i=0; i<numSteps; i++) {
713 gdif->scale(
weights[order-1][0]);
715 for(
int j=0; j<order; ++j) {
718 znew->axpy(eta*
shifts[order-1][j], v);
721 if(
weights[order-1][j+1] != 0 ) {
724 gdif->axpy(
weights[order-1][j+1],*gnew);
729 gdif->scale(1.0/eta);
733 hvCheck[i][1] = normHv;
734 hvCheck[i][2] = gdif->norm();
735 gdif->axpy(-1.0, *Hv);
736 hvCheck[i][3] = gdif->norm();
740 outStream << std::right
741 << std::setw(20) <<
"Step size"
742 << std::setw(20) <<
"norm(Hess*vec)"
743 << std::setw(20) <<
"norm(FD approx)"
744 << std::setw(20) <<
"norm(abs error)"
746 << std::setw(20) <<
"---------"
747 << std::setw(20) <<
"--------------"
748 << std::setw(20) <<
"---------------"
749 << std::setw(20) <<
"---------------"
752 outStream << std::scientific << std::setprecision(11) << std::right
753 << std::setw(20) << hvCheck[i][0]
754 << std::setw(20) << hvCheck[i][1]
755 << std::setw(20) << hvCheck[i][2]
756 << std::setw(20) << hvCheck[i][3]
763 outStream.copyfmt(oldFormatState);
772 const bool printToStream =
true,
773 std::ostream & outStream = std::cout,
775 const int order = 1 ) {
784 const std::vector<Real> &steps,
785 const bool printToStream =
true,
786 std::ostream & outStream = std::cout,
787 const int order = 1 ) {
797 const bool printToStream,
798 std::ostream & outStream,
801 std::vector<Real> steps(numSteps);
802 for(
int i=0;i<numSteps;++i) {
803 steps[i] = pow(10,-i);
814 const std::vector<Real> &steps,
815 const bool printToStream,
816 std::ostream & outStream,
819 ROL_TEST_FOR_EXCEPTION( order<1 || order>4, std::invalid_argument,
820 "Error: finite difference order must be 1,2,3, or 4" );
826 Real tol = std::sqrt(ROL_EPSILON<Real>());
828 int numSteps = steps.size();
830 std::vector<Real> tmp(numVals);
831 std::vector<std::vector<Real> > hvCheck(numSteps, tmp);
835 oldFormatState.copyfmt(outStream);
838 ROL::Ptr<Vector<Real> > g = hv.
clone();
843 ROL::Ptr<Vector<Real> > Hv = hv.
clone();
845 Real normHv = Hv->norm();
848 ROL::Ptr<Vector<Real> > gdif = hv.
clone();
849 ROL::Ptr<Vector<Real> > gnew = hv.
clone();
850 ROL::Ptr<Vector<Real> > unew = u.
clone();
852 for (
int i=0; i<numSteps; i++) {
860 gdif->scale(
weights[order-1][0]);
862 for(
int j=0; j<order; ++j) {
865 unew->axpy(eta*
shifts[order-1][j], v);
868 if(
weights[order-1][j+1] != 0 ) {
871 gdif->axpy(
weights[order-1][j+1],*gnew);
876 gdif->scale(1.0/eta);
880 hvCheck[i][1] = normHv;
881 hvCheck[i][2] = gdif->norm();
882 gdif->axpy(-1.0, *Hv);
883 hvCheck[i][3] = gdif->norm();
887 outStream << std::right
888 << std::setw(20) <<
"Step size"
889 << std::setw(20) <<
"norm(Hess*vec)"
890 << std::setw(20) <<
"norm(FD approx)"
891 << std::setw(20) <<
"norm(abs error)"
893 << std::setw(20) <<
"---------"
894 << std::setw(20) <<
"--------------"
895 << std::setw(20) <<
"---------------"
896 << std::setw(20) <<
"---------------"
899 outStream << std::scientific << std::setprecision(11) << std::right
900 << std::setw(20) << hvCheck[i][0]
901 << std::setw(20) << hvCheck[i][1]
902 << std::setw(20) << hvCheck[i][2]
903 << std::setw(20) << hvCheck[i][3]
910 outStream.copyfmt(oldFormatState);
919 const bool printToStream =
true,
920 std::ostream & outStream = std::cout,
922 const int order = 1 ) {
931 const std::vector<Real> &steps,
932 const bool printToStream =
true,
933 std::ostream & outStream = std::cout,
934 const int order = 1 ) {
944 const bool printToStream,
945 std::ostream & outStream,
948 std::vector<Real> steps(numSteps);
949 for(
int i=0;i<numSteps;++i) {
950 steps[i] = pow(10,-i);
961 const std::vector<Real> &steps,
962 const bool printToStream,
963 std::ostream & outStream,
966 ROL_TEST_FOR_EXCEPTION( order<1 || order>4, std::invalid_argument,
967 "Error: finite difference order must be 1,2,3, or 4" );
973 Real tol = std::sqrt(ROL_EPSILON<Real>());
975 int numSteps = steps.size();
977 std::vector<Real> tmp(numVals);
978 std::vector<std::vector<Real> > hvCheck(numSteps, tmp);
982 oldFormatState.copyfmt(outStream);
985 ROL::Ptr<Vector<Real> > g = hv.
clone();
990 ROL::Ptr<Vector<Real> > Hv = hv.
clone();
992 Real normHv = Hv->norm();
995 ROL::Ptr<Vector<Real> > gdif = hv.
clone();
996 ROL::Ptr<Vector<Real> > gnew = hv.
clone();
997 ROL::Ptr<Vector<Real> > znew = z.
clone();
999 for (
int i=0; i<numSteps; i++) {
1001 Real eta = steps[i];
1007 gdif->scale(
weights[order-1][0]);
1009 for(
int j=0; j<order; ++j) {
1012 znew->axpy(eta*
shifts[order-1][j], v);
1015 if(
weights[order-1][j+1] != 0 ) {
1018 gdif->axpy(
weights[order-1][j+1],*gnew);
1023 gdif->scale(1.0/eta);
1026 hvCheck[i][0] = eta;
1027 hvCheck[i][1] = normHv;
1028 hvCheck[i][2] = gdif->norm();
1029 gdif->axpy(-1.0, *Hv);
1030 hvCheck[i][3] = gdif->norm();
1032 if (printToStream) {
1034 outStream << std::right
1035 << std::setw(20) <<
"Step size"
1036 << std::setw(20) <<
"norm(Hess*vec)"
1037 << std::setw(20) <<
"norm(FD approx)"
1038 << std::setw(20) <<
"norm(abs error)"
1040 << std::setw(20) <<
"---------"
1041 << std::setw(20) <<
"--------------"
1042 << std::setw(20) <<
"---------------"
1043 << std::setw(20) <<
"---------------"
1046 outStream << std::scientific << std::setprecision(11) << std::right
1047 << std::setw(20) << hvCheck[i][0]
1048 << std::setw(20) << hvCheck[i][1]
1049 << std::setw(20) << hvCheck[i][2]
1050 << std::setw(20) << hvCheck[i][3]
1057 outStream.copyfmt(oldFormatState);
Provides the interface to evaluate objective functions.
Provides the interface to evaluate simulation-based objective functions.
std::vector< std::vector< Real > > checkGradient_2(const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &g, const Vector< Real > &d, const std::vector< Real > &steps, const bool printToStream, std::ostream &outStream, const int order)
virtual const Vector & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
virtual void scale(const Real alpha)=0
Compute where .
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
ROL::Ptr< const Vector< Real > > get_2() const
std::vector< std::vector< Real > > checkGradient_1(const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &d, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
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 .
virtual ROL::Ptr< Vector > basis(const int i) const
Return i-th basis vector.
const double weights[4][5]
virtual void axpy(const Real alpha, const Vector &x)
Compute where .
Defines the linear algebra or vector space interface for simulation-based optimization.
void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
std::vector< std::vector< Real > > checkHessVec_21(const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &v, 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 > > checkHessVec_12(const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &v, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
void update(const Vector< Real > &x, UpdateType type, int iter=-1)
Update objective function.
virtual void hessVec_12(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
virtual void update(const Vector< Real > &u, const Vector< Real > &z, UpdateType type, int iter=-1)
void set_1(const Vector< Real > &vec)
std::vector< std::vector< Real > > checkGradient_2(const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &d, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
virtual void gradient_2(Vector< Real > &g, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Compute gradient with respect to second component.
std::vector< std::vector< Real > > checkHessVec_21(const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &v, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
virtual void zero()
Set to zero vector.
Defines the linear algebra or vector space interface.
std::vector< std::vector< Real > > checkHessVec_22(const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &v, 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 > > checkHessVec_12(const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &hv, const Vector< Real > &v, const bool printToStream, std::ostream &outStream, const int numSteps, const int order)
virtual Real dot(const Vector &x) const =0
Compute where .
std::vector< std::vector< Real > > checkHessVec_12(const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &hv, const Vector< Real > &v, const std::vector< Real > &steps, const bool printToStream, std::ostream &outStream, const int order)
std::vector< std::vector< Real > > checkHessVec_11(const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &v, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
basic_nullstream< char, std::char_traits< char >> nullstream
virtual void hessVec_21(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
virtual Real value(const Vector< Real > &u, const Vector< Real > &z, Real &tol)=0
Compute value.
std::vector< std::vector< Real > > checkHessVec_11(const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &v, 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 > > checkHessVec_12(const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &v, 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 > > checkHessVec_22(const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &hv, const Vector< Real > &v, const bool printToStream, std::ostream &outStream, const int numSteps, const int order)
std::vector< std::vector< Real > > checkHessVec_22(const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &v, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
std::vector< std::vector< Real > > checkHessVec_11(const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &hv, const Vector< Real > &v, const std::vector< Real > &steps, const bool printToStream, std::ostream &outStream, const int order)
void update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update objective function.
void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply Hessian approximation to vector.
virtual void gradient_1(Vector< Real > &g, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Compute gradient with respect to first component.
virtual void hessVec_22(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
#define ROL_NUM_CHECKDERIV_STEPS
Number of steps for derivative checks.
std::vector< std::vector< Real > > checkGradient_1(const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &g, const Vector< Real > &d, const std::vector< Real > &steps, const bool printToStream, std::ostream &outStream, const int order)
virtual void update(const Vector< Real > &u, const Vector< Real > &z, bool flag=true, int iter=-1)
Update objective function. u is an iterate, z is an iterate, flag = true if the iterate has changed...
std::vector< std::vector< Real > > checkHessVec_21(const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &hv, const Vector< Real > &v, const std::vector< Real > &steps, const bool printToStream, std::ostream &outStream, const int order)
virtual Real norm() const =0
Returns where .
Real value(const Vector< Real > &x, Real &tol)
Compute value.
std::vector< std::vector< Real > > checkHessVec_22(const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &hv, const Vector< Real > &v, const std::vector< Real > &steps, const bool printToStream, std::ostream &outStream, const int order)
std::vector< std::vector< Real > > checkGradient_2(const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &g, const Vector< Real > &d, const bool printToStream, std::ostream &outStream, const int numSteps, const int order)
virtual void hessVec_11(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Apply Hessian approximation to vector.
std::vector< std::vector< Real > > checkHessVec_21(const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &hv, const Vector< Real > &v, const bool printToStream, std::ostream &outStream, const int numSteps, const int order)
std::vector< std::vector< Real > > checkGradient_1(const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &g, const Vector< Real > &d, const bool printToStream, std::ostream &outStream, const int numSteps, const int order)
std::vector< std::vector< Real > > checkHessVec_11(const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &hv, const Vector< Real > &v, const bool printToStream, std::ostream &outStream, const int numSteps, const int order)
void set_2(const Vector< Real > &vec)
ROL::Ptr< const Vector< Real > > get_1() const