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