ROL
ROL_ConstraintDef.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 
45 #ifndef ROL_CONSTRAINT_DEF_H
46 #define ROL_CONSTRAINT_DEF_H
47 
48 #include "ROL_LinearAlgebra.hpp"
49 #include "ROL_LAPACK.hpp"
50 
51 namespace ROL {
52 
53 template <class Real>
55  const Vector<Real> &v,
56  const Vector<Real> &x,
57  Real &tol) {
58  // By default we compute the finite-difference approximation.
59  const Real one(1);
60  Real ctol = std::sqrt(ROL_EPSILON<Real>());
61 
62  // Get step length.
63  Real h = std::max(one,x.norm()/v.norm())*tol;
64  //Real h = 2.0/(v.norm()*v.norm())*tol;
65 
66  // Compute constraint at x.
67  ROL::Ptr<Vector<Real> > c = jv.clone();
68  this->value(*c,x,ctol);
69 
70  // Compute perturbation x + h*v.
71  ROL::Ptr<Vector<Real> > xnew = x.clone();
72  xnew->set(x);
73  xnew->axpy(h,v);
74  this->update(*xnew);
75 
76  // Compute constraint at x + h*v.
77  jv.zero();
78  this->value(jv,*xnew,ctol);
79 
80  // Compute Newton quotient.
81  jv.axpy(-one,*c);
82  jv.scale(one/h);
83 }
84 
85 
86 template <class Real>
88  const Vector<Real> &v,
89  const Vector<Real> &x,
90  Real &tol) {
91  applyAdjointJacobian(ajv,v,x,v.dual(),tol);
92 }
93 
94 
95 
96 
97 
98 template <class Real>
100  const Vector<Real> &v,
101  const Vector<Real> &x,
102  const Vector<Real> &dualv,
103  Real &tol) {
104  // By default we compute the finite-difference approximation.
105  // This requires the implementation of a vector-space basis for the optimization variables.
106  // The default implementation requires that the constraint space is equal to its dual.
107  const Real one(1);
108  Real h(0), ctol = std::sqrt(ROL_EPSILON<Real>());
109 
110  ROL::Ptr<Vector<Real> > xnew = x.clone();
111  ROL::Ptr<Vector<Real> > ex = x.clone();
112  ROL::Ptr<Vector<Real> > eajv = ajv.clone();
113  ROL::Ptr<Vector<Real> > cnew = dualv.clone(); // in general, should be in the constraint space
114  ROL::Ptr<Vector<Real> > c0 = dualv.clone(); // in general, should be in the constraint space
115  this->value(*c0,x,ctol);
116 
117  ajv.zero();
118  for ( int i = 0; i < ajv.dimension(); i++ ) {
119  ex = x.basis(i);
120  eajv = ajv.basis(i);
121  h = std::max(one,x.norm()/ex->norm())*tol;
122  xnew->set(x);
123  xnew->axpy(h,*ex);
124  this->update(*xnew);
125  this->value(*cnew,*xnew,ctol);
126  cnew->axpy(-one,*c0);
127  cnew->scale(one/h);
128  ajv.axpy(cnew->dot(v.dual()),*eajv);
129  }
130 }
131 
132 
133 /*template <class Real>
134 void Constraint<Real>::applyHessian(Vector<Real> &huv,
135  const Vector<Real> &u,
136  const Vector<Real> &v,
137  const Vector<Real> &x,
138  Real &tol ) {
139  Real jtol = std::sqrt(ROL_EPSILON<Real>());
140 
141  // Get step length.
142  Real h = std::max(1.0,x.norm()/v.norm())*tol;
143  //Real h = 2.0/(v.norm()*v.norm())*tol;
144 
145  // Compute constraint Jacobian at x.
146  ROL::Ptr<Vector<Real> > ju = huv.clone();
147  this->applyJacobian(*ju,u,x,jtol);
148 
149  // Compute new step x + h*v.
150  ROL::Ptr<Vector<Real> > xnew = x.clone();
151  xnew->set(x);
152  xnew->axpy(h,v);
153  this->update(*xnew);
154 
155  // Compute constraint Jacobian at x + h*v.
156  huv.zero();
157  this->applyJacobian(huv,u,*xnew,jtol);
158 
159  // Compute Newton quotient.
160  huv.axpy(-1.0,*ju);
161  huv.scale(1.0/h);
162 }*/
163 
164 
165 template <class Real>
167  const Vector<Real> &u,
168  const Vector<Real> &v,
169  const Vector<Real> &x,
170  Real &tol ) {
171  // Get step length.
172  Real h = std::max(static_cast<Real>(1),x.norm()/v.norm())*tol;
173 
174  // Compute constraint Jacobian at x.
175  ROL::Ptr<Vector<Real> > aju = huv.clone();
176  applyAdjointJacobian(*aju,u,x,tol);
177 
178  // Compute new step x + h*v.
179  ROL::Ptr<Vector<Real> > xnew = x.clone();
180  xnew->set(x);
181  xnew->axpy(h,v);
182  update(*xnew);
183 
184  // Compute constraint Jacobian at x + h*v.
185  huv.zero();
186  applyAdjointJacobian(huv,u,*xnew,tol);
187 
188  // Compute Newton quotient.
189  huv.axpy(static_cast<Real>(-1),*aju);
190  huv.scale(static_cast<Real>(1)/h);
191 }
192 
193 
194 template <class Real>
196  Vector<Real> &v2,
197  const Vector<Real> &b1,
198  const Vector<Real> &b2,
199  const Vector<Real> &x,
200  Real &tol) {
201 
202  /*** Initialization. ***/
203  const Real zero(0), one(1);
204  int m = 200; // Krylov space size.
205  Real zerotol = zero;
206  int i = 0;
207  int k = 0;
208  Real temp = zero;
209  Real resnrm = zero;
210 
211  //tol = std::sqrt(b1.dot(b1)+b2.dot(b2))*1e-8;
212  tol = std::sqrt(b1.dot(b1)+b2.dot(b2))*tol;
213 
214  // Set initial guess to zero.
215  v1.zero(); v2.zero();
216 
217  // Allocate static memory.
218  ROL::Ptr<Vector<Real> > r1 = b1.clone();
219  ROL::Ptr<Vector<Real> > r2 = b2.clone();
220  ROL::Ptr<Vector<Real> > z1 = v1.clone();
221  ROL::Ptr<Vector<Real> > z2 = v2.clone();
222  ROL::Ptr<Vector<Real> > w1 = b1.clone();
223  ROL::Ptr<Vector<Real> > w2 = b2.clone();
224  std::vector<ROL::Ptr<Vector<Real> > > V1;
225  std::vector<ROL::Ptr<Vector<Real> > > V2;
226  ROL::Ptr<Vector<Real> > V2temp = b2.clone();
227  std::vector<ROL::Ptr<Vector<Real> > > Z1;
228  std::vector<ROL::Ptr<Vector<Real> > > Z2;
229  ROL::Ptr<Vector<Real> > w1temp = b1.clone();
230  ROL::Ptr<Vector<Real> > Z2temp = v2.clone();
231 
232  std::vector<Real> res(m+1, zero);
233  LA::Matrix<Real> H(m+1,m);
234  LA::Vector<Real> cs(m);
235  LA::Vector<Real> sn(m);
236  LA::Vector<Real> s(m+1);
237  LA::Vector<Real> y(m+1);
238  LA::Vector<Real> cnorm(m);
239  ROL::LAPACK<int, Real> lapack;
240 
241  // Compute initial residual.
242  applyAdjointJacobian(*r1, v2, x, zerotol);
243  r1->scale(-one); r1->axpy(-one, v1.dual()); r1->plus(b1);
244  applyJacobian(*r2, v1, x, zerotol);
245  r2->scale(-one); r2->plus(b2);
246  res[0] = std::sqrt(r1->dot(*r1) + r2->dot(*r2));
247 
248  // Check if residual is identically zero.
249  if (res[0] == zero) {
250  res.resize(0);
251  return res;
252  }
253 
254  V1.push_back(b1.clone()); (V1[0])->set(*r1); (V1[0])->scale(one/res[0]);
255  V2.push_back(b2.clone()); (V2[0])->set(*r2); (V2[0])->scale(one/res[0]);
256 
257  s(0) = res[0];
258 
259  for (i=0; i<m; i++) {
260 
261  // Apply right preconditioner.
262  V2temp->set(*(V2[i]));
263  applyPreconditioner(*Z2temp, *V2temp, x, b1, zerotol);
264  Z2.push_back(v2.clone()); (Z2[i])->set(*Z2temp);
265  Z1.push_back(v1.clone()); (Z1[i])->set((V1[i])->dual());
266 
267  // Apply operator.
268  applyJacobian(*w2, *(Z1[i]), x, zerotol);
269  applyAdjointJacobian(*w1temp, *Z2temp, x, zerotol);
270  w1->set(*(V1[i])); w1->plus(*w1temp);
271 
272  // Evaluate coefficients and orthogonalize using Gram-Schmidt.
273  for (k=0; k<=i; k++) {
274  H(k,i) = w1->dot(*(V1[k])) + w2->dot(*(V2[k]));
275  w1->axpy(-H(k,i), *(V1[k]));
276  w2->axpy(-H(k,i), *(V2[k]));
277  }
278  H(i+1,i) = std::sqrt(w1->dot(*w1) + w2->dot(*w2));
279 
280  V1.push_back(b1.clone()); (V1[i+1])->set(*w1); (V1[i+1])->scale(one/H(i+1,i));
281  V2.push_back(b2.clone()); (V2[i+1])->set(*w2); (V2[i+1])->scale(one/H(i+1,i));
282 
283  // Apply Givens rotations.
284  for (k=0; k<=i-1; k++) {
285  temp = cs(k)*H(k,i) + sn(k)*H(k+1,i);
286  H(k+1,i) = -sn(k)*H(k,i) + cs(k)*H(k+1,i);
287  H(k,i) = temp;
288  }
289 
290  // Form i-th rotation matrix.
291  if ( H(i+1,i) == zero ) {
292  cs(i) = one;
293  sn(i) = zero;
294  }
295  else if ( std::abs(H(i+1,i)) > std::abs(H(i,i)) ) {
296  temp = H(i,i) / H(i+1,i);
297  sn(i) = one / std::sqrt( one + temp*temp );
298  cs(i) = temp * sn(i);
299  }
300  else {
301  temp = H(i+1,i) / H(i,i);
302  cs(i) = one / std::sqrt( one + temp*temp );
303  sn(i) = temp * cs(i);
304  }
305 
306  // Approximate residual norm.
307  temp = cs(i)*s(i);
308  s(i+1) = -sn(i)*s(i);
309  s(i) = temp;
310  H(i,i) = cs(i)*H(i,i) + sn(i)*H(i+1,i);
311  H(i+1,i) = zero;
312  resnrm = std::abs(s(i+1));
313  res[i+1] = resnrm;
314 
315  // Update solution approximation.
316  const char uplo = 'U';
317  const char trans = 'N';
318  const char diag = 'N';
319  const char normin = 'N';
320  Real scaling = zero;
321  int info = 0;
322  y = s;
323  lapack.LATRS(uplo, trans, diag, normin, i+1, H.values(), m+1, y.values(), &scaling, cnorm.values(), &info);
324  z1->zero();
325  z2->zero();
326  for (k=0; k<=i; k++) {
327  z1->axpy(y(k), *(Z1[k]));
328  z2->axpy(y(k), *(Z2[k]));
329  }
330 
331  // Evaluate special stopping condition.
332  //tol = ???;
333 
334 // std::cout << " " << i+1 << ": " << res[i+1]/res[0] << std::endl;
335  if (res[i+1] <= tol) {
336 // std::cout << " solved in " << i+1 << " iterations to " << res[i+1] << " (" << res[i+1]/res[0] << ")" << std::endl;
337  // Update solution vector.
338  v1.plus(*z1);
339  v2.plus(*z2);
340  break;
341  }
342 
343  } // for (int i=0; i++; i<m)
344 
345  res.resize(i+2);
346 
347  /*
348  std::stringstream hist;
349  hist << std::scientific << std::setprecision(8);
350  hist << "\n Augmented System Solver:\n";
351  hist << " Iter Residual\n";
352  for (unsigned j=0; j<res.size(); j++) {
353  hist << " " << std::left << std::setw(14) << res[j] << "\n";
354  }
355  hist << "\n";
356  std::cout << hist.str();
357  */
358 
359  return res;
360 }
361 
362 
363 template <class Real>
364 std::vector<std::vector<Real> > Constraint<Real>::checkApplyJacobian(const Vector<Real> &x,
365  const Vector<Real> &v,
366  const Vector<Real> &jv,
367  const bool printToStream,
368  std::ostream & outStream,
369  const int numSteps,
370  const int order) {
371  std::vector<Real> steps(numSteps);
372  for(int i=0;i<numSteps;++i) {
373  steps[i] = pow(10,-i);
374  }
375 
376  return checkApplyJacobian(x,v,jv,steps,printToStream,outStream,order);
377 }
378 
379 
380 
381 
382 template <class Real>
383 std::vector<std::vector<Real> > Constraint<Real>::checkApplyJacobian(const Vector<Real> &x,
384  const Vector<Real> &v,
385  const Vector<Real> &jv,
386  const std::vector<Real> &steps,
387  const bool printToStream,
388  std::ostream & outStream,
389  const int order) {
390  ROL_TEST_FOR_EXCEPTION( order<1 || order>4, std::invalid_argument,
391  "Error: finite difference order must be 1,2,3, or 4" );
392 
393  const Real one(1.0);
394 
397 
398  Real tol = std::sqrt(ROL_EPSILON<Real>());
399 
400  int numSteps = steps.size();
401  int numVals = 4;
402  std::vector<Real> tmp(numVals);
403  std::vector<std::vector<Real> > jvCheck(numSteps, tmp);
404 
405  // Save the format state of the original outStream.
406  ROL::nullstream oldFormatState;
407  oldFormatState.copyfmt(outStream);
408 
409  // Compute constraint value at x.
410  ROL::Ptr<Vector<Real> > c = jv.clone();
411  this->update(x);
412  this->value(*c, x, tol);
413 
414  // Compute (Jacobian at x) times (vector v).
415  ROL::Ptr<Vector<Real> > Jv = jv.clone();
416  this->applyJacobian(*Jv, v, x, tol);
417  Real normJv = Jv->norm();
418 
419  // Temporary vectors.
420  ROL::Ptr<Vector<Real> > cdif = jv.clone();
421  ROL::Ptr<Vector<Real> > cnew = jv.clone();
422  ROL::Ptr<Vector<Real> > xnew = x.clone();
423 
424  for (int i=0; i<numSteps; i++) {
425 
426  Real eta = steps[i];
427 
428  xnew->set(x);
429 
430  cdif->set(*c);
431  cdif->scale(weights[order-1][0]);
432 
433  for(int j=0; j<order; ++j) {
434 
435  xnew->axpy(eta*shifts[order-1][j], v);
436 
437  if( weights[order-1][j+1] != 0 ) {
438  this->update(*xnew);
439  this->value(*cnew,*xnew,tol);
440  cdif->axpy(weights[order-1][j+1],*cnew);
441  }
442 
443  }
444 
445  cdif->scale(one/eta);
446 
447  // Compute norms of Jacobian-vector products, finite-difference approximations, and error.
448  jvCheck[i][0] = eta;
449  jvCheck[i][1] = normJv;
450  jvCheck[i][2] = cdif->norm();
451  cdif->axpy(-one, *Jv);
452  jvCheck[i][3] = cdif->norm();
453 
454  if (printToStream) {
455  std::stringstream hist;
456  if (i==0) {
457  hist << std::right
458  << std::setw(20) << "Step size"
459  << std::setw(20) << "norm(Jac*vec)"
460  << std::setw(20) << "norm(FD approx)"
461  << std::setw(20) << "norm(abs error)"
462  << "\n"
463  << std::setw(20) << "---------"
464  << std::setw(20) << "-------------"
465  << std::setw(20) << "---------------"
466  << std::setw(20) << "---------------"
467  << "\n";
468  }
469  hist << std::scientific << std::setprecision(11) << std::right
470  << std::setw(20) << jvCheck[i][0]
471  << std::setw(20) << jvCheck[i][1]
472  << std::setw(20) << jvCheck[i][2]
473  << std::setw(20) << jvCheck[i][3]
474  << "\n";
475  outStream << hist.str();
476  }
477 
478  }
479 
480  // Reset format state of outStream.
481  outStream.copyfmt(oldFormatState);
482 
483  return jvCheck;
484 } // checkApplyJacobian
485 
486 
487 template <class Real>
488 std::vector<std::vector<Real> > Constraint<Real>::checkApplyAdjointJacobian(const Vector<Real> &x,
489  const Vector<Real> &v,
490  const Vector<Real> &c,
491  const Vector<Real> &ajv,
492  const bool printToStream,
493  std::ostream & outStream,
494  const int numSteps) {
495  Real tol = std::sqrt(ROL_EPSILON<Real>());
496  const Real one(1);
497 
498  int numVals = 4;
499  std::vector<Real> tmp(numVals);
500  std::vector<std::vector<Real> > ajvCheck(numSteps, tmp);
501  Real eta_factor = 1e-1;
502  Real eta = one;
503 
504  // Temporary vectors.
505  ROL::Ptr<Vector<Real> > c0 = c.clone();
506  ROL::Ptr<Vector<Real> > cnew = c.clone();
507  ROL::Ptr<Vector<Real> > xnew = x.clone();
508  ROL::Ptr<Vector<Real> > ajv0 = ajv.clone();
509  ROL::Ptr<Vector<Real> > ajv1 = ajv.clone();
510  ROL::Ptr<Vector<Real> > ex = x.clone();
511  ROL::Ptr<Vector<Real> > eajv = ajv.clone();
512 
513  // Save the format state of the original outStream.
514  ROL::nullstream oldFormatState;
515  oldFormatState.copyfmt(outStream);
516 
517  // Compute constraint value at x.
518  this->update(x);
519  this->value(*c0, x, tol);
520 
521  // Compute (Jacobian at x) times (vector v).
522  this->applyAdjointJacobian(*ajv0, v, x, tol);
523  Real normAJv = ajv0->norm();
524 
525  for (int i=0; i<numSteps; i++) {
526 
527  ajv1->zero();
528 
529  for ( int j = 0; j < ajv.dimension(); j++ ) {
530  ex = x.basis(j);
531  eajv = ajv.basis(j);
532  xnew->set(x);
533  xnew->axpy(eta,*ex);
534  this->update(*xnew);
535  this->value(*cnew,*xnew,tol);
536  cnew->axpy(-one,*c0);
537  cnew->scale(one/eta);
538  ajv1->axpy(cnew->dot(v.dual()),*eajv);
539  }
540 
541  // Compute norms of Jacobian-vector products, finite-difference approximations, and error.
542  ajvCheck[i][0] = eta;
543  ajvCheck[i][1] = normAJv;
544  ajvCheck[i][2] = ajv1->norm();
545  ajv1->axpy(-one, *ajv0);
546  ajvCheck[i][3] = ajv1->norm();
547 
548  if (printToStream) {
549  std::stringstream hist;
550  if (i==0) {
551  hist << std::right
552  << std::setw(20) << "Step size"
553  << std::setw(20) << "norm(adj(Jac)*vec)"
554  << std::setw(20) << "norm(FD approx)"
555  << std::setw(20) << "norm(abs error)"
556  << "\n"
557  << std::setw(20) << "---------"
558  << std::setw(20) << "------------------"
559  << std::setw(20) << "---------------"
560  << std::setw(20) << "---------------"
561  << "\n";
562  }
563  hist << std::scientific << std::setprecision(11) << std::right
564  << std::setw(20) << ajvCheck[i][0]
565  << std::setw(20) << ajvCheck[i][1]
566  << std::setw(20) << ajvCheck[i][2]
567  << std::setw(20) << ajvCheck[i][3]
568  << "\n";
569  outStream << hist.str();
570  }
571 
572  // Update eta.
573  eta = eta*eta_factor;
574  }
575 
576  // Reset format state of outStream.
577  outStream.copyfmt(oldFormatState);
578 
579  return ajvCheck;
580 } // checkApplyAdjointJacobian
581 
582 template <class Real>
584  const Vector<Real> &v,
585  const Vector<Real> &x,
586  const Vector<Real> &dualw,
587  const Vector<Real> &dualv,
588  const bool printToStream,
589  std::ostream & outStream) {
590  Real tol = ROL_EPSILON<Real>();
591 
592  ROL::Ptr<Vector<Real> > Jv = dualw.clone();
593  ROL::Ptr<Vector<Real> > Jw = dualv.clone();
594 
595  this->update(x);
596  applyJacobian(*Jv,v,x,tol);
597  applyAdjointJacobian(*Jw,w,x,tol);
598 
599  Real vJw = v.dot(Jw->dual());
600  Real wJv = w.dot(Jv->dual());
601 
602  Real diff = std::abs(wJv-vJw);
603 
604  if ( printToStream ) {
605  std::stringstream hist;
606  hist << std::scientific << std::setprecision(8);
607  hist << "\nTest Consistency of Jacobian and its adjoint: \n |<w,Jv> - <adj(J)w,v>| = "
608  << diff << "\n";
609  hist << " |<w,Jv>| = " << std::abs(wJv) << "\n";
610  hist << " Relative Error = " << diff / (std::abs(wJv)+ROL_UNDERFLOW<Real>()) << "\n";
611  outStream << hist.str();
612  }
613  return diff;
614 } // checkAdjointConsistencyJacobian
615 
616 template <class Real>
617 std::vector<std::vector<Real> > Constraint<Real>::checkApplyAdjointHessian(const Vector<Real> &x,
618  const Vector<Real> &u,
619  const Vector<Real> &v,
620  const Vector<Real> &hv,
621  const bool printToStream,
622  std::ostream & outStream,
623  const int numSteps,
624  const int order) {
625  std::vector<Real> steps(numSteps);
626  for(int i=0;i<numSteps;++i) {
627  steps[i] = pow(10,-i);
628  }
629 
630  return checkApplyAdjointHessian(x,u,v,hv,steps,printToStream,outStream,order);
631 }
632 
633 
634 template <class Real>
635 std::vector<std::vector<Real> > Constraint<Real>::checkApplyAdjointHessian(const Vector<Real> &x,
636  const Vector<Real> &u,
637  const Vector<Real> &v,
638  const Vector<Real> &hv,
639  const std::vector<Real> &steps,
640  const bool printToStream,
641  std::ostream & outStream,
642  const int order) {
645 
646  const Real one(1);
647  Real tol = std::sqrt(ROL_EPSILON<Real>());
648 
649  int numSteps = steps.size();
650  int numVals = 4;
651  std::vector<Real> tmp(numVals);
652  std::vector<std::vector<Real> > ahuvCheck(numSteps, tmp);
653 
654  // Temporary vectors.
655  ROL::Ptr<Vector<Real> > AJdif = hv.clone();
656  ROL::Ptr<Vector<Real> > AJu = hv.clone();
657  ROL::Ptr<Vector<Real> > AHuv = hv.clone();
658  ROL::Ptr<Vector<Real> > AJnew = hv.clone();
659  ROL::Ptr<Vector<Real> > xnew = x.clone();
660 
661  // Save the format state of the original outStream.
662  ROL::nullstream oldFormatState;
663  oldFormatState.copyfmt(outStream);
664 
665  // Apply adjoint Jacobian to u.
666  this->update(x);
667  this->applyAdjointJacobian(*AJu, u, x, tol);
668 
669  // Apply adjoint Hessian at x, in direction v, to u.
670  this->applyAdjointHessian(*AHuv, u, v, x, tol);
671  Real normAHuv = AHuv->norm();
672 
673  for (int i=0; i<numSteps; i++) {
674 
675  Real eta = steps[i];
676 
677  // Apply adjoint Jacobian to u at x+eta*v.
678  xnew->set(x);
679 
680  AJdif->set(*AJu);
681  AJdif->scale(weights[order-1][0]);
682 
683  for(int j=0; j<order; ++j) {
684 
685  xnew->axpy(eta*shifts[order-1][j],v);
686 
687  if( weights[order-1][j+1] != 0 ) {
688  this->update(*xnew);
689  this->applyAdjointJacobian(*AJnew, u, *xnew, tol);
690  AJdif->axpy(weights[order-1][j+1],*AJnew);
691  }
692  }
693 
694  AJdif->scale(one/eta);
695 
696  // Compute norms of Jacobian-vector products, finite-difference approximations, and error.
697  ahuvCheck[i][0] = eta;
698  ahuvCheck[i][1] = normAHuv;
699  ahuvCheck[i][2] = AJdif->norm();
700  AJdif->axpy(-one, *AHuv);
701  ahuvCheck[i][3] = AJdif->norm();
702 
703  if (printToStream) {
704  std::stringstream hist;
705  if (i==0) {
706  hist << std::right
707  << std::setw(20) << "Step size"
708  << std::setw(20) << "norm(adj(H)(u,v))"
709  << std::setw(20) << "norm(FD approx)"
710  << std::setw(20) << "norm(abs error)"
711  << "\n"
712  << std::setw(20) << "---------"
713  << std::setw(20) << "-----------------"
714  << std::setw(20) << "---------------"
715  << std::setw(20) << "---------------"
716  << "\n";
717  }
718  hist << std::scientific << std::setprecision(11) << std::right
719  << std::setw(20) << ahuvCheck[i][0]
720  << std::setw(20) << ahuvCheck[i][1]
721  << std::setw(20) << ahuvCheck[i][2]
722  << std::setw(20) << ahuvCheck[i][3]
723  << "\n";
724  outStream << hist.str();
725  }
726 
727  }
728 
729  // Reset format state of outStream.
730  outStream.copyfmt(oldFormatState);
731 
732  return ahuvCheck;
733 } // checkApplyAdjointHessian
734 
735 } // namespace ROL
736 
737 #endif
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 applyJacobian(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply the constraint Jacobian at , , to vector .
virtual void scale(const Real alpha)=0
Compute where .
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
virtual Real checkAdjointConsistencyJacobian(const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &x, const bool printToStream=true, std::ostream &outStream=std::cout)
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
virtual void plus(const Vector &x)=0
Compute , where .
const double weights[4][5]
Definition: ROL_Types.hpp:866
virtual void axpy(const Real alpha, const Vector &x)
Compute where .
Definition: ROL_Vector.hpp:153
virtual void update(const Vector< Real > &u, const Vector< Real > &z, bool flag=true, int iter=-1) override
virtual std::vector< std::vector< Real > > checkApplyAdjointJacobian(const Vector< Real > &x, const Vector< Real > &v, const Vector< Real > &c, const Vector< Real > &ajv, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS)
Finite-difference check for the application of the adjoint of constraint Jacobian.
ROL::Objective_SimOpt value
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
virtual Real dot(const Vector &x) const =0
Compute where .
void applyJacobian(ROL::Vector< Real > &jv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &sol)
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0_ zero()
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0_ zero()
virtual void applyAdjointHessian(Vector< Real > &ahuv, const Vector< Real > &u, 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 std::vector< std::vector< Real > > checkApplyAdjointHessian(const Vector< Real > &x, const Vector< Real > &u, const Vector< Real > &v, const Vector< Real > &hv, const std::vector< Real > &step, const bool printToScreen=true, std::ostream &outStream=std::cout, const int order=1)
Finite-difference check for the application of the adjoint of constraint Hessian. ...
virtual std::vector< std::vector< Real > > checkApplyJacobian(const Vector< Real > &x, const Vector< Real > &v, const Vector< Real > &jv, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
Finite-difference check for the constraint Jacobian application.
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 .
basic_nullstream< char, char_traits< char >> nullstream
Definition: ROL_Stream.hpp:72
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 Real norm() const =0
Returns where .