44 #ifndef ROL_OBJECTIVE_SIMOPT_H
45 #define ROL_OBJECTIVE_SIMOPT_H
99 Real ftol = std::sqrt(ROL_EPSILON<Real>());
102 Real v = this->
value(u,z,ftol);
104 ROL::Ptr<Vector<Real> > unew = u.
clone();
106 for (
int i = 0; i < g.
dimension(); i++) {
109 unew->axpy(h,*(u.
basis(i)));
111 deriv = (this->
value(*unew,z,ftol) - v)/h;
119 Real ftol = std::sqrt(ROL_EPSILON<Real>());
122 Real v = this->
value(u,z,ftol);
124 ROL::Ptr<Vector<Real> > znew = z.
clone();
126 for (
int i = 0; i < g.
dimension(); i++) {
129 znew->axpy(h,*(z.
basis(i)));
131 deriv = (this->
value(u,*znew,ftol) - v)/h;
142 ROL::Ptr<Vector<Real> > g1 = gs.
get_1()->clone();
143 ROL::Ptr<Vector<Real> > g2 = gs.
get_2()->clone();
155 Real gtol = std::sqrt(ROL_EPSILON<Real>());
158 if (v.
norm() > std::sqrt(ROL_EPSILON<Real>())) {
159 h = std::max(1.0,u.
norm()/v.
norm())*tol;
162 ROL::Ptr<Vector<Real> > unew = u.
clone();
169 ROL::Ptr<Vector<Real> > g = hv.
clone();
179 Real gtol = std::sqrt(ROL_EPSILON<Real>());
182 if (v.
norm() > std::sqrt(ROL_EPSILON<Real>())) {
183 h = std::max(1.0,u.
norm()/v.
norm())*tol;
186 ROL::Ptr<Vector<Real> > znew = z.
clone();
193 ROL::Ptr<Vector<Real> > g = hv.
clone();
203 Real gtol = std::sqrt(ROL_EPSILON<Real>());
206 if (v.
norm() > std::sqrt(ROL_EPSILON<Real>())) {
207 h = std::max(1.0,u.
norm()/v.
norm())*tol;
210 ROL::Ptr<Vector<Real> > unew = u.
clone();
217 ROL::Ptr<Vector<Real> > g = hv.
clone();
227 Real gtol = std::sqrt(ROL_EPSILON<Real>());
230 if (v.
norm() > std::sqrt(ROL_EPSILON<Real>())) {
231 h = std::max(1.0,u.
norm()/v.
norm())*tol;
234 ROL::Ptr<Vector<Real> > znew = z.
clone();
241 ROL::Ptr<Vector<Real> > g = hv.
clone();
256 ROL::Ptr<Vector<Real> > h11 = (hvs.
get_1())->clone();
258 ROL::Ptr<Vector<Real> > h12 = (hvs.
get_1())->clone();
260 ROL::Ptr<Vector<Real> > h21 = (hvs.
get_2())->clone();
262 ROL::Ptr<Vector<Real> > h22 = (hvs.
get_2())->clone();
273 const bool printToStream =
true,
274 std::ostream & outStream = std::cout,
276 const int order = 1 ) {
284 const bool printToStream,
285 std::ostream & outStream,
288 std::vector<Real> steps(numSteps);
289 for(
int i=0;i<numSteps;++i) {
290 steps[i] = pow(10,-i);
300 const std::vector<Real> &steps,
301 const bool printToStream,
302 std::ostream & outStream,
304 ROL_TEST_FOR_EXCEPTION( order<1 || order>4, std::invalid_argument,
305 "Error: finite difference order must be 1,2,3, or 4" );
310 Real tol = std::sqrt(ROL_EPSILON<Real>());
312 int numSteps = steps.size();
314 std::vector<Real> tmp(numVals);
315 std::vector<std::vector<Real> > gCheck(numSteps, tmp);
319 oldFormatState.copyfmt(outStream);
323 Real val = this->
value(u,z,tol);
326 ROL::Ptr<Vector<Real> > gtmp = g.
clone();
329 Real dtg = d.
apply(*gtmp);
332 ROL::Ptr<Vector<Real> > unew = u.
clone();
334 for (
int i=0; i<numSteps; i++) {
344 gCheck[i][2] =
weights[order-1][0] * val;
346 for(
int j=0; j<order; ++j) {
348 unew->axpy(eta*
shifts[order-1][j], d);
351 if(
weights[order-1][j+1] != 0 ) {
353 gCheck[i][2] +=
weights[order-1][j+1] * this->
value(*unew,z,tol);
358 gCheck[i][3] = std::abs(gCheck[i][2] - gCheck[i][1]);
362 outStream << std::right
363 << std::setw(20) <<
"Step size"
364 << std::setw(20) <<
"grad'*dir"
365 << std::setw(20) <<
"FD approx"
366 << std::setw(20) <<
"abs error"
368 << std::setw(20) <<
"---------"
369 << std::setw(20) <<
"---------"
370 << std::setw(20) <<
"---------"
371 << std::setw(20) <<
"---------"
374 outStream << std::scientific << std::setprecision(11) << std::right
375 << std::setw(20) << gCheck[i][0]
376 << std::setw(20) << gCheck[i][1]
377 << std::setw(20) << gCheck[i][2]
378 << std::setw(20) << gCheck[i][3]
385 outStream.copyfmt(oldFormatState);
394 const bool printToStream =
true,
395 std::ostream & outStream = std::cout,
397 const int order = 1 ) {
405 const bool printToStream,
406 std::ostream & outStream,
409 std::vector<Real> steps(numSteps);
410 for(
int i=0;i<numSteps;++i) {
411 steps[i] = pow(10,-i);
421 const std::vector<Real> &steps,
422 const bool printToStream,
423 std::ostream & outStream,
425 ROL_TEST_FOR_EXCEPTION( order<1 || order>4, std::invalid_argument,
426 "Error: finite difference order must be 1,2,3, or 4" );
431 Real tol = std::sqrt(ROL_EPSILON<Real>());
433 int numSteps = steps.size();
435 std::vector<Real> tmp(numVals);
436 std::vector<std::vector<Real> > gCheck(numSteps, tmp);
440 oldFormatState.copyfmt(outStream);
444 Real val = this->
value(u,z,tol);
447 ROL::Ptr<Vector<Real> > gtmp = g.
clone();
450 Real dtg = d.
apply(*gtmp);
453 ROL::Ptr<Vector<Real> > znew = z.
clone();
455 for (
int i=0; i<numSteps; i++) {
465 gCheck[i][2] =
weights[order-1][0] * val;
467 for(
int j=0; j<order; ++j) {
469 znew->axpy(eta*
shifts[order-1][j], d);
472 if(
weights[order-1][j+1] != 0 ) {
474 gCheck[i][2] +=
weights[order-1][j+1] * this->
value(u,*znew,tol);
479 gCheck[i][3] = std::abs(gCheck[i][2] - gCheck[i][1]);
483 outStream << std::right
484 << std::setw(20) <<
"Step size"
485 << std::setw(20) <<
"grad'*dir"
486 << std::setw(20) <<
"FD approx"
487 << std::setw(20) <<
"abs error"
489 << std::setw(20) <<
"---------"
490 << std::setw(20) <<
"---------"
491 << std::setw(20) <<
"---------"
492 << std::setw(20) <<
"---------"
495 outStream << std::scientific << std::setprecision(11) << std::right
496 << std::setw(20) << gCheck[i][0]
497 << std::setw(20) << gCheck[i][1]
498 << std::setw(20) << gCheck[i][2]
499 << std::setw(20) << gCheck[i][3]
506 outStream.copyfmt(oldFormatState);
515 const bool printToStream =
true,
516 std::ostream & outStream = std::cout,
518 const int order = 1 ) {
527 const std::vector<Real> &steps,
528 const bool printToStream =
true,
529 std::ostream & outStream = std::cout,
530 const int order = 1 ) {
540 const bool printToStream,
541 std::ostream & outStream,
544 std::vector<Real> steps(numSteps);
545 for(
int i=0;i<numSteps;++i) {
546 steps[i] = pow(10,-i);
557 const std::vector<Real> &steps,
558 const bool printToStream,
559 std::ostream & outStream,
562 ROL_TEST_FOR_EXCEPTION( order<1 || order>4, std::invalid_argument,
563 "Error: finite difference order must be 1,2,3, or 4" );
569 Real tol = std::sqrt(ROL_EPSILON<Real>());
571 int numSteps = steps.size();
573 std::vector<Real> tmp(numVals);
574 std::vector<std::vector<Real> > hvCheck(numSteps, tmp);
578 oldFormatState.copyfmt(outStream);
581 ROL::Ptr<Vector<Real> > g = hv.
clone();
586 ROL::Ptr<Vector<Real> > Hv = hv.
clone();
588 Real normHv = Hv->norm();
591 ROL::Ptr<Vector<Real> > gdif = hv.
clone();
592 ROL::Ptr<Vector<Real> > gnew = hv.
clone();
593 ROL::Ptr<Vector<Real> > unew = u.
clone();
595 for (
int i=0; i<numSteps; i++) {
603 gdif->scale(
weights[order-1][0]);
605 for(
int j=0; j<order; ++j) {
608 unew->axpy(eta*
shifts[order-1][j], v);
611 if(
weights[order-1][j+1] != 0 ) {
614 gdif->axpy(
weights[order-1][j+1],*gnew);
619 gdif->scale(1.0/eta);
623 hvCheck[i][1] = normHv;
624 hvCheck[i][2] = gdif->norm();
625 gdif->axpy(-1.0, *Hv);
626 hvCheck[i][3] = gdif->norm();
630 outStream << std::right
631 << std::setw(20) <<
"Step size"
632 << std::setw(20) <<
"norm(Hess*vec)"
633 << std::setw(20) <<
"norm(FD approx)"
634 << std::setw(20) <<
"norm(abs error)"
636 << std::setw(20) <<
"---------"
637 << std::setw(20) <<
"--------------"
638 << std::setw(20) <<
"---------------"
639 << std::setw(20) <<
"---------------"
642 outStream << std::scientific << std::setprecision(11) << std::right
643 << std::setw(20) << hvCheck[i][0]
644 << std::setw(20) << hvCheck[i][1]
645 << std::setw(20) << hvCheck[i][2]
646 << std::setw(20) << hvCheck[i][3]
653 outStream.copyfmt(oldFormatState);
662 const bool printToStream =
true,
663 std::ostream & outStream = std::cout,
665 const int order = 1 ) {
672 const std::vector<Real> &steps,
673 const bool printToStream =
true,
674 std::ostream & outStream = std::cout,
675 const int order = 1 ) {
684 const bool printToStream,
685 std::ostream & outStream,
688 std::vector<Real> steps(numSteps);
689 for(
int i=0;i<numSteps;++i) {
690 steps[i] = pow(10,-i);
701 const std::vector<Real> &steps,
702 const bool printToStream,
703 std::ostream & outStream,
706 ROL_TEST_FOR_EXCEPTION( order<1 || order>4, std::invalid_argument,
707 "Error: finite difference order must be 1,2,3, or 4" );
713 Real tol = std::sqrt(ROL_EPSILON<Real>());
715 int numSteps = steps.size();
717 std::vector<Real> tmp(numVals);
718 std::vector<std::vector<Real> > hvCheck(numSteps, tmp);
722 oldFormatState.copyfmt(outStream);
725 ROL::Ptr<Vector<Real> > g = hv.
clone();
730 ROL::Ptr<Vector<Real> > Hv = hv.
clone();
732 Real normHv = Hv->norm();
735 ROL::Ptr<Vector<Real> > gdif = hv.
clone();
736 ROL::Ptr<Vector<Real> > gnew = hv.
clone();
737 ROL::Ptr<Vector<Real> > znew = z.
clone();
739 for (
int i=0; i<numSteps; i++) {
747 gdif->scale(
weights[order-1][0]);
749 for(
int j=0; j<order; ++j) {
752 znew->axpy(eta*
shifts[order-1][j], v);
755 if(
weights[order-1][j+1] != 0 ) {
758 gdif->axpy(
weights[order-1][j+1],*gnew);
763 gdif->scale(1.0/eta);
767 hvCheck[i][1] = normHv;
768 hvCheck[i][2] = gdif->norm();
769 gdif->axpy(-1.0, *Hv);
770 hvCheck[i][3] = gdif->norm();
774 outStream << std::right
775 << std::setw(20) <<
"Step size"
776 << std::setw(20) <<
"norm(Hess*vec)"
777 << std::setw(20) <<
"norm(FD approx)"
778 << std::setw(20) <<
"norm(abs error)"
780 << std::setw(20) <<
"---------"
781 << std::setw(20) <<
"--------------"
782 << std::setw(20) <<
"---------------"
783 << std::setw(20) <<
"---------------"
786 outStream << std::scientific << std::setprecision(11) << std::right
787 << std::setw(20) << hvCheck[i][0]
788 << std::setw(20) << hvCheck[i][1]
789 << std::setw(20) << hvCheck[i][2]
790 << std::setw(20) << hvCheck[i][3]
797 outStream.copyfmt(oldFormatState);
806 const bool printToStream =
true,
807 std::ostream & outStream = std::cout,
809 const int order = 1 ) {
818 const std::vector<Real> &steps,
819 const bool printToStream =
true,
820 std::ostream & outStream = std::cout,
821 const int order = 1 ) {
831 const bool printToStream,
832 std::ostream & outStream,
835 std::vector<Real> steps(numSteps);
836 for(
int i=0;i<numSteps;++i) {
837 steps[i] = pow(10,-i);
848 const std::vector<Real> &steps,
849 const bool printToStream,
850 std::ostream & outStream,
853 ROL_TEST_FOR_EXCEPTION( order<1 || order>4, std::invalid_argument,
854 "Error: finite difference order must be 1,2,3, or 4" );
860 Real tol = std::sqrt(ROL_EPSILON<Real>());
862 int numSteps = steps.size();
864 std::vector<Real> tmp(numVals);
865 std::vector<std::vector<Real> > hvCheck(numSteps, tmp);
869 oldFormatState.copyfmt(outStream);
872 ROL::Ptr<Vector<Real> > g = hv.
clone();
877 ROL::Ptr<Vector<Real> > Hv = hv.
clone();
879 Real normHv = Hv->norm();
882 ROL::Ptr<Vector<Real> > gdif = hv.
clone();
883 ROL::Ptr<Vector<Real> > gnew = hv.
clone();
884 ROL::Ptr<Vector<Real> > unew = u.
clone();
886 for (
int i=0; i<numSteps; i++) {
894 gdif->scale(
weights[order-1][0]);
896 for(
int j=0; j<order; ++j) {
899 unew->axpy(eta*
shifts[order-1][j], v);
902 if(
weights[order-1][j+1] != 0 ) {
905 gdif->axpy(
weights[order-1][j+1],*gnew);
910 gdif->scale(1.0/eta);
914 hvCheck[i][1] = normHv;
915 hvCheck[i][2] = gdif->norm();
916 gdif->axpy(-1.0, *Hv);
917 hvCheck[i][3] = gdif->norm();
921 outStream << std::right
922 << std::setw(20) <<
"Step size"
923 << std::setw(20) <<
"norm(Hess*vec)"
924 << std::setw(20) <<
"norm(FD approx)"
925 << std::setw(20) <<
"norm(abs error)"
927 << std::setw(20) <<
"---------"
928 << std::setw(20) <<
"--------------"
929 << std::setw(20) <<
"---------------"
930 << std::setw(20) <<
"---------------"
933 outStream << std::scientific << std::setprecision(11) << std::right
934 << std::setw(20) << hvCheck[i][0]
935 << std::setw(20) << hvCheck[i][1]
936 << std::setw(20) << hvCheck[i][2]
937 << std::setw(20) << hvCheck[i][3]
944 outStream.copyfmt(oldFormatState);
953 const bool printToStream =
true,
954 std::ostream & outStream = std::cout,
956 const int order = 1 ) {
965 const std::vector<Real> &steps,
966 const bool printToStream =
true,
967 std::ostream & outStream = std::cout,
968 const int order = 1 ) {
978 const bool printToStream,
979 std::ostream & outStream,
982 std::vector<Real> steps(numSteps);
983 for(
int i=0;i<numSteps;++i) {
984 steps[i] = pow(10,-i);
995 const std::vector<Real> &steps,
996 const bool printToStream,
997 std::ostream & outStream,
1000 ROL_TEST_FOR_EXCEPTION( order<1 || order>4, std::invalid_argument,
1001 "Error: finite difference order must be 1,2,3, or 4" );
1007 Real tol = std::sqrt(ROL_EPSILON<Real>());
1009 int numSteps = steps.size();
1011 std::vector<Real> tmp(numVals);
1012 std::vector<std::vector<Real> > hvCheck(numSteps, tmp);
1016 oldFormatState.copyfmt(outStream);
1019 ROL::Ptr<Vector<Real> > g = hv.
clone();
1024 ROL::Ptr<Vector<Real> > Hv = hv.
clone();
1026 Real normHv = Hv->norm();
1029 ROL::Ptr<Vector<Real> > gdif = hv.
clone();
1030 ROL::Ptr<Vector<Real> > gnew = hv.
clone();
1031 ROL::Ptr<Vector<Real> > znew = z.
clone();
1033 for (
int i=0; i<numSteps; i++) {
1035 Real eta = steps[i];
1041 gdif->scale(
weights[order-1][0]);
1043 for(
int j=0; j<order; ++j) {
1046 znew->axpy(eta*
shifts[order-1][j], v);
1049 if(
weights[order-1][j+1] != 0 ) {
1052 gdif->axpy(
weights[order-1][j+1],*gnew);
1057 gdif->scale(1.0/eta);
1060 hvCheck[i][0] = eta;
1061 hvCheck[i][1] = normHv;
1062 hvCheck[i][2] = gdif->norm();
1063 gdif->axpy(-1.0, *Hv);
1064 hvCheck[i][3] = gdif->norm();
1066 if (printToStream) {
1068 outStream << std::right
1069 << std::setw(20) <<
"Step size"
1070 << std::setw(20) <<
"norm(Hess*vec)"
1071 << std::setw(20) <<
"norm(FD approx)"
1072 << std::setw(20) <<
"norm(abs error)"
1074 << std::setw(20) <<
"---------"
1075 << std::setw(20) <<
"--------------"
1076 << std::setw(20) <<
"---------------"
1077 << std::setw(20) <<
"---------------"
1080 outStream << std::scientific << std::setprecision(11) << std::right
1081 << std::setw(20) << hvCheck[i][0]
1082 << std::setw(20) << hvCheck[i][1]
1083 << std::setw(20) << hvCheck[i][2]
1084 << std::setw(20) << hvCheck[i][3]
1091 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)
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...
basic_nullstream< char, char_traits< char >> nullstream
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