ROL
ROL_Objective_SimOpt.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ************************************************************************
3 //
4 // Rapid Optimization Library (ROL) Package
5 // Copyright (2014) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact lead developers:
38 // Drew Kouri (dpkouri@sandia.gov) and
39 // Denis Ridzal (dridzal@sandia.gov)
40 //
41 // ************************************************************************
42 // @HEADER
43 
44 #ifndef ROL_OBJECTIVE_SIMOPT_H
45 #define ROL_OBJECTIVE_SIMOPT_H
46 
47 #include "ROL_Objective.hpp"
48 #include "ROL_Vector_SimOpt.hpp"
49 
56 namespace ROL {
57 
58 template <class Real>
59 class Objective_SimOpt : public Objective<Real> {
60 public:
61 
68  virtual void update( const Vector<Real> &u, const Vector<Real> &z, bool flag = true, int iter = -1 ) {}
69 
70  void update( const Vector<Real> &x, bool flag = true, int iter = -1 ) {
71  const ROL::Vector_SimOpt<Real> &xs = dynamic_cast<const ROL::Vector_SimOpt<Real>&>(
72  dynamic_cast<const ROL::Vector<Real>&>(x));
73  this->update(*(xs.get_1()),*(xs.get_2()),flag,iter);
74  }
75 
76 
79  virtual Real value( const Vector<Real> &u, const Vector<Real> &z, Real &tol ) = 0;
80 
81  Real value( const Vector<Real> &x, Real &tol ) {
82  const ROL::Vector_SimOpt<Real> &xs = dynamic_cast<const ROL::Vector_SimOpt<Real>&>(
83  dynamic_cast<const ROL::Vector<Real>&>(x));
84  return this->value(*(xs.get_1()),*(xs.get_2()),tol);
85  }
86 
87 
90  virtual void gradient_1( Vector<Real> &g, const Vector<Real> &u, const Vector<Real> &z, Real &tol ) {
91  Real ftol = std::sqrt(ROL_EPSILON<Real>());
92  Real h = 0.0;
93  this->update(u,z);
94  Real v = this->value(u,z,ftol);
95  Real deriv = 0.0;
96  ROL::Ptr<Vector<Real> > unew = u.clone();
97  g.zero();
98  for (int i = 0; i < g.dimension(); i++) {
99  h = u.dot(*u.basis(i))*tol;
100  unew->set(u);
101  unew->axpy(h,*(u.basis(i)));
102  this->update(*unew,z);
103  deriv = (this->value(*unew,z,ftol) - v)/h;
104  g.axpy(deriv,*(g.basis(i)));
105  }
106  this->update(u,z);
107  }
110  virtual void gradient_2( Vector<Real> &g, const Vector<Real> &u, const Vector<Real> &z, Real &tol ) {
111  Real ftol = std::sqrt(ROL_EPSILON<Real>());
112  Real h = 0.0;
113  this->update(u,z);
114  Real v = this->value(u,z,ftol);
115  Real deriv = 0.0;
116  ROL::Ptr<Vector<Real> > znew = z.clone();
117  g.zero();
118  for (int i = 0; i < g.dimension(); i++) {
119  h = z.dot(*z.basis(i))*tol;
120  znew->set(z);
121  znew->axpy(h,*(z.basis(i)));
122  this->update(u,*znew);
123  deriv = (this->value(u,*znew,ftol) - v)/h;
124  g.axpy(deriv,*(g.basis(i)));
125  }
126  this->update(u,z);
127  }
128 
129  void gradient( Vector<Real> &g, const Vector<Real> &x, Real &tol ) {
131  dynamic_cast<ROL::Vector<Real>&>(g));
132  const ROL::Vector_SimOpt<Real> &xs = dynamic_cast<const ROL::Vector_SimOpt<Real>&>(
133  dynamic_cast<const ROL::Vector<Real>&>(x));
134  ROL::Ptr<Vector<Real> > g1 = gs.get_1()->clone();
135  ROL::Ptr<Vector<Real> > g2 = gs.get_2()->clone();
136  this->gradient_1(*g1,*(xs.get_1()),*(xs.get_2()),tol);
137  this->gradient_2(*g2,*(xs.get_1()),*(xs.get_2()),tol);
138  gs.set_1(*g1);
139  gs.set_2(*g2);
140  }
141 
142 
145  virtual void hessVec_11( Vector<Real> &hv, const Vector<Real> &v,
146  const Vector<Real> &u, const Vector<Real> &z, Real &tol ) {
147  Real gtol = std::sqrt(ROL_EPSILON<Real>());
148  // Compute step length
149  Real h = tol;
150  if (v.norm() > std::sqrt(ROL_EPSILON<Real>())) {
151  h = std::max(1.0,u.norm()/v.norm())*tol;
152  }
153  // Evaluate gradient of first component at (u+hv,z)
154  ROL::Ptr<Vector<Real> > unew = u.clone();
155  unew->set(u);
156  unew->axpy(h,v);
157  this->update(*unew,z);
158  hv.zero();
159  this->gradient_1(hv,*unew,z,gtol);
160  // Evaluate gradient of first component at (u,z)
161  ROL::Ptr<Vector<Real> > g = hv.clone();
162  this->update(u,z);
163  this->gradient_1(*g,u,z,gtol);
164  // Compute Newton quotient
165  hv.axpy(-1.0,*g);
166  hv.scale(1.0/h);
167  }
168 
169  virtual void hessVec_12( Vector<Real> &hv, const Vector<Real> &v,
170  const Vector<Real> &u, const Vector<Real> &z, Real &tol ) {
171  Real gtol = std::sqrt(ROL_EPSILON<Real>());
172  // Compute step length
173  Real h = tol;
174  if (v.norm() > std::sqrt(ROL_EPSILON<Real>())) {
175  h = std::max(1.0,u.norm()/v.norm())*tol;
176  }
177  // Evaluate gradient of first component at (u,z+hv)
178  ROL::Ptr<Vector<Real> > znew = z.clone();
179  znew->set(z);
180  znew->axpy(h,v);
181  this->update(u,*znew);
182  hv.zero();
183  this->gradient_1(hv,u,*znew,gtol);
184  // Evaluate gradient of first component at (u,z)
185  ROL::Ptr<Vector<Real> > g = hv.clone();
186  this->update(u,z);
187  this->gradient_1(*g,u,z,gtol);
188  // Compute Newton quotient
189  hv.axpy(-1.0,*g);
190  hv.scale(1.0/h);
191  }
192 
193  virtual void hessVec_21( Vector<Real> &hv, const Vector<Real> &v,
194  const Vector<Real> &u, const Vector<Real> &z, Real &tol ) {
195  Real gtol = std::sqrt(ROL_EPSILON<Real>());
196  // Compute step length
197  Real h = tol;
198  if (v.norm() > std::sqrt(ROL_EPSILON<Real>())) {
199  h = std::max(1.0,u.norm()/v.norm())*tol;
200  }
201  // Evaluate gradient of first component at (u+hv,z)
202  ROL::Ptr<Vector<Real> > unew = u.clone();
203  unew->set(u);
204  unew->axpy(h,v);
205  this->update(*unew,z);
206  hv.zero();
207  this->gradient_2(hv,*unew,z,gtol);
208  // Evaluate gradient of first component at (u,z)
209  ROL::Ptr<Vector<Real> > g = hv.clone();
210  this->update(u,z);
211  this->gradient_2(*g,u,z,gtol);
212  // Compute Newton quotient
213  hv.axpy(-1.0,*g);
214  hv.scale(1.0/h);
215  }
216 
217  virtual void hessVec_22( Vector<Real> &hv, const Vector<Real> &v,
218  const Vector<Real> &u, const Vector<Real> &z, Real &tol ) {
219  Real gtol = std::sqrt(ROL_EPSILON<Real>());
220  // Compute step length
221  Real h = tol;
222  if (v.norm() > std::sqrt(ROL_EPSILON<Real>())) {
223  h = std::max(1.0,u.norm()/v.norm())*tol;
224  }
225  // Evaluate gradient of first component at (u,z+hv)
226  ROL::Ptr<Vector<Real> > znew = z.clone();
227  znew->set(z);
228  znew->axpy(h,v);
229  this->update(u,*znew);
230  hv.zero();
231  this->gradient_2(hv,u,*znew,gtol);
232  // Evaluate gradient of first component at (u,z)
233  ROL::Ptr<Vector<Real> > g = hv.clone();
234  this->update(u,z);
235  this->gradient_2(*g,u,z,gtol);
236  // Compute Newton quotient
237  hv.axpy(-1.0,*g);
238  hv.scale(1.0/h);
239  }
240 
241  void hessVec( Vector<Real> &hv, const Vector<Real> &v, const Vector<Real> &x, Real &tol ) {
243  dynamic_cast<ROL::Vector<Real>&>(hv));
244  const ROL::Vector_SimOpt<Real> &vs = dynamic_cast<const ROL::Vector_SimOpt<Real>&>(
245  dynamic_cast<const ROL::Vector<Real>&>(v));
246  const ROL::Vector_SimOpt<Real> &xs = dynamic_cast<const ROL::Vector_SimOpt<Real>&>(
247  dynamic_cast<const ROL::Vector<Real>&>(x));
248  ROL::Ptr<Vector<Real> > h11 = (hvs.get_1())->clone();
249  this->hessVec_11(*h11,*(vs.get_1()),*(xs.get_1()),*(xs.get_2()),tol);
250  ROL::Ptr<Vector<Real> > h12 = (hvs.get_1())->clone();
251  this->hessVec_12(*h12,*(vs.get_2()),*(xs.get_1()),*(xs.get_2()),tol);
252  ROL::Ptr<Vector<Real> > h21 = (hvs.get_2())->clone();
253  this->hessVec_21(*h21,*(vs.get_1()),*(xs.get_1()),*(xs.get_2()),tol);
254  ROL::Ptr<Vector<Real> > h22 = (hvs.get_2())->clone();
255  this->hessVec_22(*h22,*(vs.get_2()),*(xs.get_1()),*(xs.get_2()),tol);
256  h11->plus(*h12);
257  hvs.set_1(*h11);
258  h22->plus(*h21);
259  hvs.set_2(*h22);
260  }
261 
262  std::vector<std::vector<Real> > checkGradient_1( const Vector<Real> &u,
263  const Vector<Real> &z,
264  const Vector<Real> &d,
265  const bool printToStream = true,
266  std::ostream & outStream = std::cout,
267  const int numSteps = ROL_NUM_CHECKDERIV_STEPS,
268  const int order = 1 ) {
269  return checkGradient_1(u, z, u.dual(), d, printToStream, outStream, numSteps, order);
270  }
271 
272  std::vector<std::vector<Real> > checkGradient_1( const Vector<Real> &u,
273  const Vector<Real> &z,
274  const Vector<Real> &g,
275  const Vector<Real> &d,
276  const bool printToStream,
277  std::ostream & outStream,
278  const int numSteps,
279  const int order ) {
280  std::vector<Real> steps(numSteps);
281  for(int i=0;i<numSteps;++i) {
282  steps[i] = pow(10,-i);
283  }
284 
285  return checkGradient_1(u,z,g,d,steps,printToStream,outStream,order);
286  } // checkGradient_1
287 
288  std::vector<std::vector<Real> > checkGradient_1( const Vector<Real> &u,
289  const Vector<Real> &z,
290  const Vector<Real> &g,
291  const Vector<Real> &d,
292  const std::vector<Real> &steps,
293  const bool printToStream,
294  std::ostream & outStream,
295  const int order ) {
296  ROL_TEST_FOR_EXCEPTION( order<1 || order>4, std::invalid_argument,
297  "Error: finite difference order must be 1,2,3, or 4" );
298 
301 
302  Real tol = std::sqrt(ROL_EPSILON<Real>());
303 
304  int numSteps = steps.size();
305  int numVals = 4;
306  std::vector<Real> tmp(numVals);
307  std::vector<std::vector<Real> > gCheck(numSteps, tmp);
308 
309  // Save the format state of the original outStream.
310  ROL::nullstream oldFormatState;
311  oldFormatState.copyfmt(outStream);
312 
313  // Evaluate objective value at x.
314  this->update(u,z);
315  Real val = this->value(u,z,tol);
316 
317  // Compute gradient at x.
318  ROL::Ptr<Vector<Real> > gtmp = g.clone();
319  this->gradient_1(*gtmp, u, z, tol);
320  Real dtg = d.dot(gtmp->dual());
321 
322  // Temporary vectors.
323  ROL::Ptr<Vector<Real> > unew = u.clone();
324 
325  for (int i=0; i<numSteps; i++) {
326 
327  Real eta = steps[i];
328 
329  unew->set(u);
330 
331  // Compute gradient, finite-difference gradient, and absolute error.
332  gCheck[i][0] = eta;
333  gCheck[i][1] = dtg;
334 
335  gCheck[i][2] = weights[order-1][0] * val;
336 
337  for(int j=0; j<order; ++j) {
338  // Evaluate at x <- x+eta*c_i*d.
339  unew->axpy(eta*shifts[order-1][j], d);
340 
341  // Only evaluate at shifts where the weight is nonzero
342  if( weights[order-1][j+1] != 0 ) {
343  this->update(*unew,z);
344  gCheck[i][2] += weights[order-1][j+1] * this->value(*unew,z,tol);
345  }
346  }
347  gCheck[i][2] /= eta;
348 
349  gCheck[i][3] = std::abs(gCheck[i][2] - gCheck[i][1]);
350 
351  if (printToStream) {
352  if (i==0) {
353  outStream << std::right
354  << std::setw(20) << "Step size"
355  << std::setw(20) << "grad'*dir"
356  << std::setw(20) << "FD approx"
357  << std::setw(20) << "abs error"
358  << "\n"
359  << std::setw(20) << "---------"
360  << std::setw(20) << "---------"
361  << std::setw(20) << "---------"
362  << std::setw(20) << "---------"
363  << "\n";
364  }
365  outStream << std::scientific << std::setprecision(11) << std::right
366  << std::setw(20) << gCheck[i][0]
367  << std::setw(20) << gCheck[i][1]
368  << std::setw(20) << gCheck[i][2]
369  << std::setw(20) << gCheck[i][3]
370  << "\n";
371  }
372 
373  }
374 
375  // Reset format state of outStream.
376  outStream.copyfmt(oldFormatState);
377 
378  return gCheck;
379  } // checkGradient_1
380 
381 
382  std::vector<std::vector<Real> > checkGradient_2( const Vector<Real> &u,
383  const Vector<Real> &z,
384  const Vector<Real> &d,
385  const bool printToStream = true,
386  std::ostream & outStream = std::cout,
387  const int numSteps = ROL_NUM_CHECKDERIV_STEPS,
388  const int order = 1 ) {
389  return checkGradient_2(u, z, z.dual(), d, printToStream, outStream, numSteps, order);
390  }
391 
392  std::vector<std::vector<Real> > checkGradient_2( const Vector<Real> &u,
393  const Vector<Real> &z,
394  const Vector<Real> &g,
395  const Vector<Real> &d,
396  const bool printToStream,
397  std::ostream & outStream,
398  const int numSteps,
399  const int order ) {
400  std::vector<Real> steps(numSteps);
401  for(int i=0;i<numSteps;++i) {
402  steps[i] = pow(10,-i);
403  }
404 
405  return checkGradient_2(u,z,g,d,steps,printToStream,outStream,order);
406  } // checkGradient_2
407 
408  std::vector<std::vector<Real> > checkGradient_2( const Vector<Real> &u,
409  const Vector<Real> &z,
410  const Vector<Real> &g,
411  const Vector<Real> &d,
412  const std::vector<Real> &steps,
413  const bool printToStream,
414  std::ostream & outStream,
415  const int order ) {
416  ROL_TEST_FOR_EXCEPTION( order<1 || order>4, std::invalid_argument,
417  "Error: finite difference order must be 1,2,3, or 4" );
418 
421 
422  Real tol = std::sqrt(ROL_EPSILON<Real>());
423 
424  int numSteps = steps.size();
425  int numVals = 4;
426  std::vector<Real> tmp(numVals);
427  std::vector<std::vector<Real> > gCheck(numSteps, tmp);
428 
429  // Save the format state of the original outStream.
430  ROL::nullstream oldFormatState;
431  oldFormatState.copyfmt(outStream);
432 
433  // Evaluate objective value at x.
434  this->update(u,z);
435  Real val = this->value(u,z,tol);
436 
437  // Compute gradient at x.
438  ROL::Ptr<Vector<Real> > gtmp = g.clone();
439  this->gradient_2(*gtmp, u, z, tol);
440  Real dtg = d.dot(gtmp->dual());
441 
442  // Temporary vectors.
443  ROL::Ptr<Vector<Real> > znew = z.clone();
444 
445  for (int i=0; i<numSteps; i++) {
446 
447  Real eta = steps[i];
448 
449  znew->set(z);
450 
451  // Compute gradient, finite-difference gradient, and absolute error.
452  gCheck[i][0] = eta;
453  gCheck[i][1] = dtg;
454 
455  gCheck[i][2] = weights[order-1][0] * val;
456 
457  for(int j=0; j<order; ++j) {
458  // Evaluate at x <- x+eta*c_i*d.
459  znew->axpy(eta*shifts[order-1][j], d);
460 
461  // Only evaluate at shifts where the weight is nonzero
462  if( weights[order-1][j+1] != 0 ) {
463  this->update(u,*znew);
464  gCheck[i][2] += weights[order-1][j+1] * this->value(u,*znew,tol);
465  }
466  }
467  gCheck[i][2] /= eta;
468 
469  gCheck[i][3] = std::abs(gCheck[i][2] - gCheck[i][1]);
470 
471  if (printToStream) {
472  if (i==0) {
473  outStream << std::right
474  << std::setw(20) << "Step size"
475  << std::setw(20) << "grad'*dir"
476  << std::setw(20) << "FD approx"
477  << std::setw(20) << "abs error"
478  << "\n"
479  << std::setw(20) << "---------"
480  << std::setw(20) << "---------"
481  << std::setw(20) << "---------"
482  << std::setw(20) << "---------"
483  << "\n";
484  }
485  outStream << std::scientific << std::setprecision(11) << std::right
486  << std::setw(20) << gCheck[i][0]
487  << std::setw(20) << gCheck[i][1]
488  << std::setw(20) << gCheck[i][2]
489  << std::setw(20) << gCheck[i][3]
490  << "\n";
491  }
492 
493  }
494 
495  // Reset format state of outStream.
496  outStream.copyfmt(oldFormatState);
497 
498  return gCheck;
499  } // checkGradient_2
500 
501 
502  std::vector<std::vector<Real> > checkHessVec_11( const Vector<Real> &u,
503  const Vector<Real> &z,
504  const Vector<Real> &v,
505  const bool printToStream = true,
506  std::ostream & outStream = std::cout,
507  const int numSteps = ROL_NUM_CHECKDERIV_STEPS,
508  const int order = 1 ) {
509 
510  return checkHessVec_11(u, z, u.dual(), v, printToStream, outStream, numSteps, order);
511 
512  }
513 
514  std::vector<std::vector<Real> > checkHessVec_11( const Vector<Real> &u,
515  const Vector<Real> &z,
516  const Vector<Real> &v,
517  const std::vector<Real> &steps,
518  const bool printToStream = true,
519  std::ostream & outStream = std::cout,
520  const int order = 1 ) {
521 
522  return checkHessVec_11(u, z, u.dual(), v, steps, printToStream, outStream, order);
523  }
524 
525 
526  std::vector<std::vector<Real> > checkHessVec_11( const Vector<Real> &u,
527  const Vector<Real> &z,
528  const Vector<Real> &hv,
529  const Vector<Real> &v,
530  const bool printToStream,
531  std::ostream & outStream,
532  const int numSteps,
533  const int order ) {
534  std::vector<Real> steps(numSteps);
535  for(int i=0;i<numSteps;++i) {
536  steps[i] = pow(10,-i);
537  }
538 
539  return checkHessVec_11(u,z,hv,v,steps,printToStream,outStream,order);
540  } // checkHessVec_11
541 
542 
543  std::vector<std::vector<Real> > checkHessVec_11( const Vector<Real> &u,
544  const Vector<Real> &z,
545  const Vector<Real> &hv,
546  const Vector<Real> &v,
547  const std::vector<Real> &steps,
548  const bool printToStream,
549  std::ostream & outStream,
550  const int order ) {
551 
552  ROL_TEST_FOR_EXCEPTION( order<1 || order>4, std::invalid_argument,
553  "Error: finite difference order must be 1,2,3, or 4" );
554 
557 
558 
559  Real tol = std::sqrt(ROL_EPSILON<Real>());
560 
561  int numSteps = steps.size();
562  int numVals = 4;
563  std::vector<Real> tmp(numVals);
564  std::vector<std::vector<Real> > hvCheck(numSteps, tmp);
565 
566  // Save the format state of the original outStream.
567  ROL::nullstream oldFormatState;
568  oldFormatState.copyfmt(outStream);
569 
570  // Compute gradient at x.
571  ROL::Ptr<Vector<Real> > g = hv.clone();
572  this->update(u,z);
573  this->gradient_1(*g, u, z, tol);
574 
575  // Compute (Hessian at x) times (vector v).
576  ROL::Ptr<Vector<Real> > Hv = hv.clone();
577  this->hessVec_11(*Hv, v, u, z, tol);
578  Real normHv = Hv->norm();
579 
580  // Temporary vectors.
581  ROL::Ptr<Vector<Real> > gdif = hv.clone();
582  ROL::Ptr<Vector<Real> > gnew = hv.clone();
583  ROL::Ptr<Vector<Real> > unew = u.clone();
584 
585  for (int i=0; i<numSteps; i++) {
586 
587  Real eta = steps[i];
588 
589  // Evaluate objective value at x+eta*d.
590  unew->set(u);
591 
592  gdif->set(*g);
593  gdif->scale(weights[order-1][0]);
594 
595  for(int j=0; j<order; ++j) {
596 
597  // Evaluate at x <- x+eta*c_i*d.
598  unew->axpy(eta*shifts[order-1][j], v);
599 
600  // Only evaluate at shifts where the weight is nonzero
601  if( weights[order-1][j+1] != 0 ) {
602  this->update(*unew,z);
603  this->gradient_1(*gnew, *unew, z, tol);
604  gdif->axpy(weights[order-1][j+1],*gnew);
605  }
606 
607  }
608 
609  gdif->scale(1.0/eta);
610 
611  // Compute norms of hessvec, finite-difference hessvec, and error.
612  hvCheck[i][0] = eta;
613  hvCheck[i][1] = normHv;
614  hvCheck[i][2] = gdif->norm();
615  gdif->axpy(-1.0, *Hv);
616  hvCheck[i][3] = gdif->norm();
617 
618  if (printToStream) {
619  if (i==0) {
620  outStream << std::right
621  << std::setw(20) << "Step size"
622  << std::setw(20) << "norm(Hess*vec)"
623  << std::setw(20) << "norm(FD approx)"
624  << std::setw(20) << "norm(abs error)"
625  << "\n"
626  << std::setw(20) << "---------"
627  << std::setw(20) << "--------------"
628  << std::setw(20) << "---------------"
629  << std::setw(20) << "---------------"
630  << "\n";
631  }
632  outStream << std::scientific << std::setprecision(11) << std::right
633  << std::setw(20) << hvCheck[i][0]
634  << std::setw(20) << hvCheck[i][1]
635  << std::setw(20) << hvCheck[i][2]
636  << std::setw(20) << hvCheck[i][3]
637  << "\n";
638  }
639 
640  }
641 
642  // Reset format state of outStream.
643  outStream.copyfmt(oldFormatState);
644 
645  return hvCheck;
646  } // checkHessVec_11
647 
648 
649  std::vector<std::vector<Real> > checkHessVec_12( const Vector<Real> &u,
650  const Vector<Real> &z,
651  const Vector<Real> &v,
652  const bool printToStream = true,
653  std::ostream & outStream = std::cout,
654  const int numSteps = ROL_NUM_CHECKDERIV_STEPS,
655  const int order = 1 ) {
656  return checkHessVec_12(u, z, u.dual(), v, printToStream, outStream, numSteps, order);
657  }
658 
659  std::vector<std::vector<Real> > checkHessVec_12( const Vector<Real> &u,
660  const Vector<Real> &z,
661  const Vector<Real> &v,
662  const std::vector<Real> &steps,
663  const bool printToStream = true,
664  std::ostream & outStream = std::cout,
665  const int order = 1 ) {
666  return checkHessVec_12(u, z, u.dual(), v, steps, printToStream, outStream, order);
667  }
668 
669 
670  std::vector<std::vector<Real> > checkHessVec_12( const Vector<Real> &u,
671  const Vector<Real> &z,
672  const Vector<Real> &hv,
673  const Vector<Real> &v,
674  const bool printToStream,
675  std::ostream & outStream,
676  const int numSteps,
677  const int order ) {
678  std::vector<Real> steps(numSteps);
679  for(int i=0;i<numSteps;++i) {
680  steps[i] = pow(10,-i);
681  }
682 
683  return checkHessVec_12(u,z,hv,v,steps,printToStream,outStream,order);
684  } // checkHessVec_12
685 
686 
687  std::vector<std::vector<Real> > checkHessVec_12( const Vector<Real> &u,
688  const Vector<Real> &z,
689  const Vector<Real> &hv,
690  const Vector<Real> &v,
691  const std::vector<Real> &steps,
692  const bool printToStream,
693  std::ostream & outStream,
694  const int order ) {
695 
696  ROL_TEST_FOR_EXCEPTION( order<1 || order>4, std::invalid_argument,
697  "Error: finite difference order must be 1,2,3, or 4" );
698 
701 
702 
703  Real tol = std::sqrt(ROL_EPSILON<Real>());
704 
705  int numSteps = steps.size();
706  int numVals = 4;
707  std::vector<Real> tmp(numVals);
708  std::vector<std::vector<Real> > hvCheck(numSteps, tmp);
709 
710  // Save the format state of the original outStream.
711  ROL::nullstream oldFormatState;
712  oldFormatState.copyfmt(outStream);
713 
714  // Compute gradient at x.
715  ROL::Ptr<Vector<Real> > g = hv.clone();
716  this->update(u,z);
717  this->gradient_1(*g, u, z, tol);
718 
719  // Compute (Hessian at x) times (vector v).
720  ROL::Ptr<Vector<Real> > Hv = hv.clone();
721  this->hessVec_12(*Hv, v, u, z, tol);
722  Real normHv = Hv->norm();
723 
724  // Temporary vectors.
725  ROL::Ptr<Vector<Real> > gdif = hv.clone();
726  ROL::Ptr<Vector<Real> > gnew = hv.clone();
727  ROL::Ptr<Vector<Real> > znew = z.clone();
728 
729  for (int i=0; i<numSteps; i++) {
730 
731  Real eta = steps[i];
732 
733  // Evaluate objective value at x+eta*d.
734  znew->set(z);
735 
736  gdif->set(*g);
737  gdif->scale(weights[order-1][0]);
738 
739  for(int j=0; j<order; ++j) {
740 
741  // Evaluate at x <- x+eta*c_i*d.
742  znew->axpy(eta*shifts[order-1][j], v);
743 
744  // Only evaluate at shifts where the weight is nonzero
745  if( weights[order-1][j+1] != 0 ) {
746  this->update(u,*znew);
747  this->gradient_1(*gnew, u, *znew, tol);
748  gdif->axpy(weights[order-1][j+1],*gnew);
749  }
750 
751  }
752 
753  gdif->scale(1.0/eta);
754 
755  // Compute norms of hessvec, finite-difference hessvec, and error.
756  hvCheck[i][0] = eta;
757  hvCheck[i][1] = normHv;
758  hvCheck[i][2] = gdif->norm();
759  gdif->axpy(-1.0, *Hv);
760  hvCheck[i][3] = gdif->norm();
761 
762  if (printToStream) {
763  if (i==0) {
764  outStream << std::right
765  << std::setw(20) << "Step size"
766  << std::setw(20) << "norm(Hess*vec)"
767  << std::setw(20) << "norm(FD approx)"
768  << std::setw(20) << "norm(abs error)"
769  << "\n"
770  << std::setw(20) << "---------"
771  << std::setw(20) << "--------------"
772  << std::setw(20) << "---------------"
773  << std::setw(20) << "---------------"
774  << "\n";
775  }
776  outStream << std::scientific << std::setprecision(11) << std::right
777  << std::setw(20) << hvCheck[i][0]
778  << std::setw(20) << hvCheck[i][1]
779  << std::setw(20) << hvCheck[i][2]
780  << std::setw(20) << hvCheck[i][3]
781  << "\n";
782  }
783 
784  }
785 
786  // Reset format state of outStream.
787  outStream.copyfmt(oldFormatState);
788 
789  return hvCheck;
790  } // checkHessVec_12
791 
792 
793  std::vector<std::vector<Real> > checkHessVec_21( const Vector<Real> &u,
794  const Vector<Real> &z,
795  const Vector<Real> &v,
796  const bool printToStream = true,
797  std::ostream & outStream = std::cout,
798  const int numSteps = ROL_NUM_CHECKDERIV_STEPS,
799  const int order = 1 ) {
800 
801  return checkHessVec_21(u, z, z.dual(), v, printToStream, outStream, numSteps, order);
802 
803  }
804 
805  std::vector<std::vector<Real> > checkHessVec_21( const Vector<Real> &u,
806  const Vector<Real> &z,
807  const Vector<Real> &v,
808  const std::vector<Real> &steps,
809  const bool printToStream = true,
810  std::ostream & outStream = std::cout,
811  const int order = 1 ) {
812 
813  return checkHessVec_21(u, z, z.dual(), v, steps, printToStream, outStream, order);
814  }
815 
816 
817  std::vector<std::vector<Real> > checkHessVec_21( const Vector<Real> &u,
818  const Vector<Real> &z,
819  const Vector<Real> &hv,
820  const Vector<Real> &v,
821  const bool printToStream,
822  std::ostream & outStream,
823  const int numSteps,
824  const int order ) {
825  std::vector<Real> steps(numSteps);
826  for(int i=0;i<numSteps;++i) {
827  steps[i] = pow(10,-i);
828  }
829 
830  return checkHessVec_21(u,z,hv,v,steps,printToStream,outStream,order);
831  } // checkHessVec_21
832 
833 
834  std::vector<std::vector<Real> > checkHessVec_21( const Vector<Real> &u,
835  const Vector<Real> &z,
836  const Vector<Real> &hv,
837  const Vector<Real> &v,
838  const std::vector<Real> &steps,
839  const bool printToStream,
840  std::ostream & outStream,
841  const int order ) {
842 
843  ROL_TEST_FOR_EXCEPTION( order<1 || order>4, std::invalid_argument,
844  "Error: finite difference order must be 1,2,3, or 4" );
845 
848 
849 
850  Real tol = std::sqrt(ROL_EPSILON<Real>());
851 
852  int numSteps = steps.size();
853  int numVals = 4;
854  std::vector<Real> tmp(numVals);
855  std::vector<std::vector<Real> > hvCheck(numSteps, tmp);
856 
857  // Save the format state of the original outStream.
858  ROL::nullstream oldFormatState;
859  oldFormatState.copyfmt(outStream);
860 
861  // Compute gradient at x.
862  ROL::Ptr<Vector<Real> > g = hv.clone();
863  this->update(u,z);
864  this->gradient_2(*g, u, z, tol);
865 
866  // Compute (Hessian at x) times (vector v).
867  ROL::Ptr<Vector<Real> > Hv = hv.clone();
868  this->hessVec_21(*Hv, v, u, z, tol);
869  Real normHv = Hv->norm();
870 
871  // Temporary vectors.
872  ROL::Ptr<Vector<Real> > gdif = hv.clone();
873  ROL::Ptr<Vector<Real> > gnew = hv.clone();
874  ROL::Ptr<Vector<Real> > unew = u.clone();
875 
876  for (int i=0; i<numSteps; i++) {
877 
878  Real eta = steps[i];
879 
880  // Evaluate objective value at x+eta*d.
881  unew->set(u);
882 
883  gdif->set(*g);
884  gdif->scale(weights[order-1][0]);
885 
886  for(int j=0; j<order; ++j) {
887 
888  // Evaluate at x <- x+eta*c_i*d.
889  unew->axpy(eta*shifts[order-1][j], v);
890 
891  // Only evaluate at shifts where the weight is nonzero
892  if( weights[order-1][j+1] != 0 ) {
893  this->update(*unew,z);
894  this->gradient_2(*gnew, *unew, z, tol);
895  gdif->axpy(weights[order-1][j+1],*gnew);
896  }
897 
898  }
899 
900  gdif->scale(1.0/eta);
901 
902  // Compute norms of hessvec, finite-difference hessvec, and error.
903  hvCheck[i][0] = eta;
904  hvCheck[i][1] = normHv;
905  hvCheck[i][2] = gdif->norm();
906  gdif->axpy(-1.0, *Hv);
907  hvCheck[i][3] = gdif->norm();
908 
909  if (printToStream) {
910  if (i==0) {
911  outStream << std::right
912  << std::setw(20) << "Step size"
913  << std::setw(20) << "norm(Hess*vec)"
914  << std::setw(20) << "norm(FD approx)"
915  << std::setw(20) << "norm(abs error)"
916  << "\n"
917  << std::setw(20) << "---------"
918  << std::setw(20) << "--------------"
919  << std::setw(20) << "---------------"
920  << std::setw(20) << "---------------"
921  << "\n";
922  }
923  outStream << std::scientific << std::setprecision(11) << std::right
924  << std::setw(20) << hvCheck[i][0]
925  << std::setw(20) << hvCheck[i][1]
926  << std::setw(20) << hvCheck[i][2]
927  << std::setw(20) << hvCheck[i][3]
928  << "\n";
929  }
930 
931  }
932 
933  // Reset format state of outStream.
934  outStream.copyfmt(oldFormatState);
935 
936  return hvCheck;
937  } // checkHessVec_21
938 
939 
940  std::vector<std::vector<Real> > checkHessVec_22( const Vector<Real> &u,
941  const Vector<Real> &z,
942  const Vector<Real> &v,
943  const bool printToStream = true,
944  std::ostream & outStream = std::cout,
945  const int numSteps = ROL_NUM_CHECKDERIV_STEPS,
946  const int order = 1 ) {
947 
948  return checkHessVec_22(u, z, z.dual(), v, printToStream, outStream, numSteps, order);
949 
950  }
951 
952  std::vector<std::vector<Real> > checkHessVec_22( const Vector<Real> &u,
953  const Vector<Real> &z,
954  const Vector<Real> &v,
955  const std::vector<Real> &steps,
956  const bool printToStream = true,
957  std::ostream & outStream = std::cout,
958  const int order = 1 ) {
959 
960  return checkHessVec_22(u, z, z.dual(), v, steps, printToStream, outStream, order);
961  }
962 
963 
964  std::vector<std::vector<Real> > checkHessVec_22( const Vector<Real> &u,
965  const Vector<Real> &z,
966  const Vector<Real> &hv,
967  const Vector<Real> &v,
968  const bool printToStream,
969  std::ostream & outStream,
970  const int numSteps,
971  const int order ) {
972  std::vector<Real> steps(numSteps);
973  for(int i=0;i<numSteps;++i) {
974  steps[i] = pow(10,-i);
975  }
976 
977  return checkHessVec_22(u,z,hv,v,steps,printToStream,outStream,order);
978  } // checkHessVec_22
979 
980 
981  std::vector<std::vector<Real> > checkHessVec_22( const Vector<Real> &u,
982  const Vector<Real> &z,
983  const Vector<Real> &hv,
984  const Vector<Real> &v,
985  const std::vector<Real> &steps,
986  const bool printToStream,
987  std::ostream & outStream,
988  const int order ) {
989 
990  ROL_TEST_FOR_EXCEPTION( order<1 || order>4, std::invalid_argument,
991  "Error: finite difference order must be 1,2,3, or 4" );
992 
995 
996 
997  Real tol = std::sqrt(ROL_EPSILON<Real>());
998 
999  int numSteps = steps.size();
1000  int numVals = 4;
1001  std::vector<Real> tmp(numVals);
1002  std::vector<std::vector<Real> > hvCheck(numSteps, tmp);
1003 
1004  // Save the format state of the original outStream.
1005  ROL::nullstream oldFormatState;
1006  oldFormatState.copyfmt(outStream);
1007 
1008  // Compute gradient at x.
1009  ROL::Ptr<Vector<Real> > g = hv.clone();
1010  this->update(u,z);
1011  this->gradient_2(*g, u, z, tol);
1012 
1013  // Compute (Hessian at x) times (vector v).
1014  ROL::Ptr<Vector<Real> > Hv = hv.clone();
1015  this->hessVec_22(*Hv, v, u, z, tol);
1016  Real normHv = Hv->norm();
1017 
1018  // Temporary vectors.
1019  ROL::Ptr<Vector<Real> > gdif = hv.clone();
1020  ROL::Ptr<Vector<Real> > gnew = hv.clone();
1021  ROL::Ptr<Vector<Real> > znew = z.clone();
1022 
1023  for (int i=0; i<numSteps; i++) {
1024 
1025  Real eta = steps[i];
1026 
1027  // Evaluate objective value at x+eta*d.
1028  znew->set(z);
1029 
1030  gdif->set(*g);
1031  gdif->scale(weights[order-1][0]);
1032 
1033  for(int j=0; j<order; ++j) {
1034 
1035  // Evaluate at x <- x+eta*c_i*d.
1036  znew->axpy(eta*shifts[order-1][j], v);
1037 
1038  // Only evaluate at shifts where the weight is nonzero
1039  if( weights[order-1][j+1] != 0 ) {
1040  this->update(u,*znew);
1041  this->gradient_2(*gnew, u, *znew, tol);
1042  gdif->axpy(weights[order-1][j+1],*gnew);
1043  }
1044 
1045  }
1046 
1047  gdif->scale(1.0/eta);
1048 
1049  // Compute norms of hessvec, finite-difference hessvec, and error.
1050  hvCheck[i][0] = eta;
1051  hvCheck[i][1] = normHv;
1052  hvCheck[i][2] = gdif->norm();
1053  gdif->axpy(-1.0, *Hv);
1054  hvCheck[i][3] = gdif->norm();
1055 
1056  if (printToStream) {
1057  if (i==0) {
1058  outStream << std::right
1059  << std::setw(20) << "Step size"
1060  << std::setw(20) << "norm(Hess*vec)"
1061  << std::setw(20) << "norm(FD approx)"
1062  << std::setw(20) << "norm(abs error)"
1063  << "\n"
1064  << std::setw(20) << "---------"
1065  << std::setw(20) << "--------------"
1066  << std::setw(20) << "---------------"
1067  << std::setw(20) << "---------------"
1068  << "\n";
1069  }
1070  outStream << std::scientific << std::setprecision(11) << std::right
1071  << std::setw(20) << hvCheck[i][0]
1072  << std::setw(20) << hvCheck[i][1]
1073  << std::setw(20) << hvCheck[i][2]
1074  << std::setw(20) << hvCheck[i][3]
1075  << "\n";
1076  }
1077 
1078  }
1079 
1080  // Reset format state of outStream.
1081  outStream.copyfmt(oldFormatState);
1082 
1083  return hvCheck;
1084  } // checkHessVec_22
1085 
1086 }; // class Objective_SimOpt
1087 
1088 } // namespace ROL
1089 
1090 #endif
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...
Definition: ROL_Vector.hpp:226
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.
Definition: ROL_Vector.hpp:196
virtual ROL::Ptr< Vector > basis(const int i) const
Return i-th basis vector.
Definition: ROL_Vector.hpp:182
const double weights[4][5]
Definition: ROL_Types.hpp:937
virtual void axpy(const Real alpha, const Vector &x)
Compute where .
Definition: ROL_Vector.hpp:153
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)
virtual void hessVec_12(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
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.
Definition: ROL_Vector.hpp:167
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:80
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.
Definition: ROL_Types.hpp:74
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
Definition: ROL_Stream.hpp:72
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