ROL
ROL_EqualityConstraint_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_EQUALITY_CONSTRAINT_SIMOPT_H
45 #define ROL_EQUALITY_CONSTRAINT_SIMOPT_H
46 
48 #include "ROL_Vector_SimOpt.hpp"
49 #include "ROL_Types.hpp"
50 #include <iostream>
51 
92 namespace ROL {
93 
94 template <class Real>
96 public:
102  virtual void update( const Vector<Real> &u, const Vector<Real> &z, bool flag = true, int iter = -1 ) {}
103 
104 
118  virtual void value(Vector<Real> &c,
119  const Vector<Real> &u,
120  const Vector<Real> &z,
121  Real &tol) = 0;
122 
131  virtual void solve(Vector<Real> &u,
132  const Vector<Real> &z,
133  Real &tol) {
134  TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error,
135  "The method solve is used but not implemented!\n");
136  }
137 
153  virtual void applyJacobian_1(Vector<Real> &jv,
154  const Vector<Real> &v,
155  const Vector<Real> &u,
156  const Vector<Real> &z,
157  Real &tol) {
158  Real ctol = std::sqrt(ROL_EPSILON);
159  // Compute step length
160  Real h = tol;
161  if (v.norm() > std::sqrt(ROL_EPSILON)) {
162  h = std::max(1.0,u.norm()/v.norm())*tol;
163  }
164  // Update state vector to u + hv
165  Teuchos::RCP<Vector<Real> > unew = u.clone();
166  unew->set(u);
167  unew->axpy(h,v);
168  // Compute new constraint value
169  update(*unew,z);
170  value(jv,*unew,z,ctol);
171  // Compute current constraint value
172  Teuchos::RCP<Vector<Real> > cold = jv.clone();
173  update(u,z);
174  value(*cold,u,z,ctol);
175  // Compute Newton quotient
176  jv.axpy(-1.0,*cold);
177  jv.scale(1.0/h);
178  }
179 
180 
196  virtual void applyJacobian_2(Vector<Real> &jv,
197  const Vector<Real> &v,
198  const Vector<Real> &u,
199  const Vector<Real> &z,
200  Real &tol) {
201  Real ctol = std::sqrt(ROL_EPSILON);
202  // Compute step length
203  Real h = tol;
204  if (v.norm() > std::sqrt(ROL_EPSILON)) {
205  h = std::max(1.0,u.norm()/v.norm())*tol;
206  }
207  // Update state vector to u + hv
208  Teuchos::RCP<Vector<Real> > znew = z.clone();
209  znew->set(z);
210  znew->axpy(h,v);
211  // Compute new constraint value
212  update(u,*znew);
213  value(jv,u,*znew,ctol);
214  // Compute current constraint value
215  Teuchos::RCP<Vector<Real> > cold = jv.clone();
216  update(u,z);
217  value(*cold,u,z,ctol);
218  // Compute Newton quotient
219  jv.axpy(-1.0,*cold);
220  jv.scale(1.0/h);
221  }
222 
239  const Vector<Real> &v,
240  const Vector<Real> &u,
241  const Vector<Real> &z,
242  Real &tol) {
243  TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error,
244  "The method applyInverseJacobian_1 is used but not implemented!\n");
245  }
246 
263  const Vector<Real> &v,
264  const Vector<Real> &u,
265  const Vector<Real> &z,
266  Real &tol) {
267  applyAdjointJacobian_1(ajv, v, u, z, v.dual(), tol);
268  }
269 
270 
289  const Vector<Real> &v,
290  const Vector<Real> &u,
291  const Vector<Real> &z,
292  const Vector<Real> &dualv,
293  Real &tol) {
294  Real ctol = std::sqrt(ROL_EPSILON);
295  Real h = tol;
296  if (v.norm() > std::sqrt(ROL_EPSILON)) {
297  h = std::max(1.0,u.norm()/v.norm())*tol;
298  }
299  Teuchos::RCP<Vector<Real> > cold = dualv.clone();
300  Teuchos::RCP<Vector<Real> > cnew = dualv.clone();
301  update(u,z);
302  value(*cold,u,z,ctol);
303  Teuchos::RCP<Vector<Real> > unew = u.clone();
304  ajv.zero();
305  for (int i = 0; i < u.dimension(); i++) {
306  unew->set(u);
307  unew->axpy(h,*(u.basis(i)));
308  update(*unew,z);
309  value(*cnew,*unew,z,ctol);
310  cnew->axpy(-1.0,*cold);
311  cnew->scale(1.0/h);
312  ajv.axpy(cnew->dot(v),*((u.dual()).basis(i)));
313  }
314  update(u,z);
315  }
316 
317 
334  const Vector<Real> &v,
335  const Vector<Real> &u,
336  const Vector<Real> &z,
337  Real &tol) {
338  applyAdjointJacobian_2(ajv, v, u, z, v.dual(), tol);
339  }
340 
341 
360  const Vector<Real> &v,
361  const Vector<Real> &u,
362  const Vector<Real> &z,
363  const Vector<Real> &dualv,
364  Real &tol) {
365  Real ctol = std::sqrt(ROL_EPSILON);
366  Real h = tol;
367  if (v.norm() > std::sqrt(ROL_EPSILON)) {
368  h = std::max(1.0,u.norm()/v.norm())*tol;
369  }
370  Teuchos::RCP<Vector<Real> > cold = dualv.clone();
371  Teuchos::RCP<Vector<Real> > cnew = dualv.clone();
372  update(u,z);
373  value(*cold,u,z,ctol);
374  Teuchos::RCP<Vector<Real> > znew = z.clone();
375  ajv.zero();
376  for (int i = 0; i < z.dimension(); i++) {
377  znew->set(z);
378  znew->axpy(h,*(z.basis(i)));
379  update(u,*znew);
380  value(*cnew,u,*znew,ctol);
381  cnew->axpy(-1.0,*cold);
382  cnew->scale(1.0/h);
383  ajv.axpy(cnew->dot(v),*((z.dual()).basis(i)));
384  }
385  update(u,z);
386  }
387 
404  const Vector<Real> &v,
405  const Vector<Real> &u,
406  const Vector<Real> &z,
407  Real &tol) {
408  TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error,
409  "The method applyInverseAdjointJacobian_1 is used but not implemented!\n");
410  };
411 
429  const Vector<Real> &w,
430  const Vector<Real> &v,
431  const Vector<Real> &u,
432  const Vector<Real> &z,
433  Real &tol) {
434  Real jtol = std::sqrt(ROL_EPSILON);
435  // Compute step size
436  Real h = tol;
437  if (v.norm() > std::sqrt(ROL_EPSILON)) {
438  h = std::max(1.0,u.norm()/v.norm())*tol;
439  }
440  // Evaluate Jacobian at new state
441  Teuchos::RCP<Vector<Real> > unew = u.clone();
442  unew->set(u);
443  unew->axpy(h,v);
444  update(*unew,z);
445  applyAdjointJacobian_1(ahwv,w,*unew,z,jtol);
446  // Evaluate Jacobian at old state
447  Teuchos::RCP<Vector<Real> > jv = ahwv.clone();
448  update(u,z);
449  applyAdjointJacobian_1(*jv,w,u,z,jtol);
450  // Compute Newton quotient
451  ahwv.axpy(-1.0,*jv);
452  ahwv.scale(1.0/h);
453  }
454 
455 
473  const Vector<Real> &w,
474  const Vector<Real> &v,
475  const Vector<Real> &u,
476  const Vector<Real> &z,
477  Real &tol) {
478  Real jtol = std::sqrt(ROL_EPSILON);
479  // Compute step size
480  Real h = tol;
481  if (v.norm() > std::sqrt(ROL_EPSILON)) {
482  h = std::max(1.0,u.norm()/v.norm())*tol;
483  }
484  // Evaluate Jacobian at new state
485  Teuchos::RCP<Vector<Real> > unew = u.clone();
486  unew->set(u);
487  unew->axpy(h,v);
488  update(*unew,z);
489  applyAdjointJacobian_2(ahwv,w,*unew,z,jtol);
490  // Evaluate Jacobian at old state
491  Teuchos::RCP<Vector<Real> > jv = ahwv.clone();
492  update(u,z);
493  applyAdjointJacobian_2(*jv,w,u,z,jtol);
494  // Compute Newton quotient
495  ahwv.axpy(-1.0,*jv);
496  ahwv.scale(1.0/h);
497  }
498 
499 
517  const Vector<Real> &w,
518  const Vector<Real> &v,
519  const Vector<Real> &u,
520  const Vector<Real> &z,
521  Real &tol) {
522  Real jtol = std::sqrt(ROL_EPSILON);
523  // Compute step size
524  Real h = tol;
525  if (v.norm() > std::sqrt(ROL_EPSILON)) {
526  h = std::max(1.0,u.norm()/v.norm())*tol;
527  }
528  // Evaluate Jacobian at new control
529  Teuchos::RCP<Vector<Real> > znew = z.clone();
530  znew->set(z);
531  znew->axpy(h,v);
532  update(u,*znew);
533  applyAdjointJacobian_1(ahwv,w,u,*znew,jtol);
534  // Evaluate Jacobian at old control
535  Teuchos::RCP<Vector<Real> > jv = ahwv.clone();
536  update(u,z);
537  applyAdjointJacobian_1(*jv,w,u,z,jtol);
538  // Compute Newton quotient
539  ahwv.axpy(-1.0,*jv);
540  ahwv.scale(1.0/h);
541  }
542 
560  const Vector<Real> &w,
561  const Vector<Real> &v,
562  const Vector<Real> &u,
563  const Vector<Real> &z,
564  Real &tol) {
565  Real jtol = std::sqrt(ROL_EPSILON);
566  // Compute step size
567  Real h = tol;
568  if (v.norm() > std::sqrt(ROL_EPSILON)) {
569  h = std::max(1.0,u.norm()/v.norm())*tol;
570  }
571  // Evaluate Jacobian at new control
572  Teuchos::RCP<Vector<Real> > znew = z.clone();
573  znew->set(z);
574  znew->axpy(h,v);
575  update(u,*znew);
576  applyAdjointJacobian_2(ahwv,w,u,*znew,jtol);
577  // Evaluate Jacobian at old control
578  Teuchos::RCP<Vector<Real> > jv = ahwv.clone();
579  update(u,z);
580  applyAdjointJacobian_2(*jv,w,u,z,jtol);
581  // Compute Newton quotient
582  ahwv.axpy(-1.0,*jv);
583  ahwv.scale(1.0/h);
584 }
585 
624  virtual std::vector<Real> solveAugmentedSystem(Vector<Real> &v1,
625  Vector<Real> &v2,
626  const Vector<Real> &b1,
627  const Vector<Real> &b2,
628  const Vector<Real> &x,
629  Real &tol) {
630  return EqualityConstraint<Real>::solveAugmentedSystem(v1,v2,b1,b2,x,tol);
631  }
632 
633 
653  const Vector<Real> &v,
654  const Vector<Real> &x,
655  const Vector<Real> &g,
656  Real &tol) {
657  const Vector_SimOpt<Real> &xs = Teuchos::dyn_cast<const Vector_SimOpt<Real> >(x);
658  Teuchos::RCP<ROL::Vector<Real> > ijv = (xs.get_1())->clone();
659 
660  try {
661  applyInverseJacobian_1(*ijv, v, *(xs.get_1()), *(xs.get_2()), tol);
662  }
663  catch (const std::logic_error &e) {
665  return;
666  }
667 
668  const Vector_SimOpt<Real> &gs = Teuchos::dyn_cast<const Vector_SimOpt<Real> >(g);
669  Teuchos::RCP<ROL::Vector<Real> > ijv_dual = (gs.get_1())->clone();
670  ijv_dual->set(ijv->dual());
671 
672  try {
673  applyInverseAdjointJacobian_1(pv, *ijv_dual, *(xs.get_1()), *(xs.get_2()), tol);
674  }
675  catch (const std::logic_error &e) {
677  return;
678  }
679 
680  }
681 
682 
684 
690  virtual void update( const Vector<Real> &x, bool flag = true, int iter = -1 ) {
691  const Vector_SimOpt<Real> &xs = Teuchos::dyn_cast<const Vector_SimOpt<Real> >(
692  Teuchos::dyn_cast<const Vector<Real> >(x));
693  update(*(xs.get_1()),*(xs.get_2()),flag,iter);
694  }
695 
698  virtual bool isFeasible( const Vector<Real> &v ) { return true; }
699 
700  virtual void value(Vector<Real> &c,
701  const Vector<Real> &x,
702  Real &tol) {
703  const Vector_SimOpt<Real> &xs = Teuchos::dyn_cast<const Vector_SimOpt<Real> >(
704  Teuchos::dyn_cast<const Vector<Real> >(x));
705  value(c,*(xs.get_1()),*(xs.get_2()),tol);
706  }
707 
708 
709  virtual void applyJacobian(Vector<Real> &jv,
710  const Vector<Real> &v,
711  const Vector<Real> &x,
712  Real &tol) {
713  const Vector_SimOpt<Real> &xs = Teuchos::dyn_cast<const Vector_SimOpt<Real> >(
714  Teuchos::dyn_cast<const Vector<Real> >(x));
715  const Vector_SimOpt<Real> &vs = Teuchos::dyn_cast<const Vector_SimOpt<Real> >(
716  Teuchos::dyn_cast<const Vector<Real> >(v));
717  applyJacobian_1(jv,*(vs.get_1()),*(xs.get_1()),*(xs.get_2()),tol);
718  Teuchos::RCP<Vector<Real> > jv2 = jv.clone();
719  applyJacobian_2(*jv2,*(vs.get_2()),*(xs.get_1()),*(xs.get_2()),tol);
720  jv.plus(*jv2);
721  }
722 
723 
725  const Vector<Real> &v,
726  const Vector<Real> &x,
727  Real &tol) {
728  Vector_SimOpt<Real> &ajvs = Teuchos::dyn_cast<Vector_SimOpt<Real> >(
729  Teuchos::dyn_cast<Vector<Real> >(ajv));
730  const Vector_SimOpt<Real> &xs = Teuchos::dyn_cast<const Vector_SimOpt<Real> >(
731  Teuchos::dyn_cast<const Vector<Real> >(x));
732  Teuchos::RCP<Vector<Real> > ajv1 = (ajvs.get_1())->clone();
733  applyAdjointJacobian_1(*ajv1,v,*(xs.get_1()),*(xs.get_2()),tol);
734  ajvs.set_1(*ajv1);
735  Teuchos::RCP<Vector<Real> > ajv2 = (ajvs.get_2())->clone();
736  applyAdjointJacobian_2(*ajv2,v,*(xs.get_1()),*(xs.get_2()),tol);
737  ajvs.set_2(*ajv2);
738  }
739 
740 
741  virtual void applyAdjointHessian(Vector<Real> &ahwv,
742  const Vector<Real> &w,
743  const Vector<Real> &v,
744  const Vector<Real> &x,
745  Real &tol) {
746  Vector_SimOpt<Real> &ahwvs = Teuchos::dyn_cast<Vector_SimOpt<Real> >(
747  Teuchos::dyn_cast<Vector<Real> >(ahwv));
748  const Vector_SimOpt<Real> &xs = Teuchos::dyn_cast<const Vector_SimOpt<Real> >(
749  Teuchos::dyn_cast<const Vector<Real> >(x));
750  const Vector_SimOpt<Real> &vs = Teuchos::dyn_cast<const Vector_SimOpt<Real> >(
751  Teuchos::dyn_cast<const Vector<Real> >(v));
752  // Block-row 1
753  Teuchos::RCP<Vector<Real> > C11 = (ahwvs.get_1())->clone();
754  Teuchos::RCP<Vector<Real> > C21 = (ahwvs.get_1())->clone();
755  applyAdjointHessian_11(*C11,w,*(vs.get_1()),*(xs.get_1()),*(xs.get_2()),tol);
756  applyAdjointHessian_21(*C21,w,*(vs.get_2()),*(xs.get_1()),*(xs.get_2()),tol);
757  C11->plus(*C21);
758  ahwvs.set_1(*C11);
759  // Block-row 2
760  Teuchos::RCP<Vector<Real> > C12 = (ahwvs.get_2())->clone();
761  Teuchos::RCP<Vector<Real> > C22 = (ahwvs.get_2())->clone();
762  applyAdjointHessian_12(*C12,w,*(vs.get_1()),*(xs.get_1()),*(xs.get_2()),tol);
763  applyAdjointHessian_22(*C22,w,*(vs.get_2()),*(xs.get_1()),*(xs.get_2()),tol);
764  C22->plus(*C12);
765  ahwvs.set_2(*C22);
766  }
767 
768 
769 
770  virtual Real checkSolve(const ROL::Vector<Real> &u,
771  const ROL::Vector<Real> &z,
772  const ROL::Vector<Real> &c,
773  const bool printToStream = true,
774  std::ostream & outStream = std::cout) {
775  // Solve equality constraint for u.
776  Real tol = ROL_EPSILON;
777  Teuchos::RCP<ROL::Vector<Real> > s = u.clone();
778  solve(*s,z,tol);
779  // Evaluate equality constraint residual at (u,z).
780  Teuchos::RCP<ROL::Vector<Real> > cs = c.clone();
781  value(*cs,*s,z,tol);
782  // Output norm of residual.
783  Real cnorm = cs->norm();
784  if ( printToStream ) {
785  std::stringstream hist;
786  hist << std::scientific << std::setprecision(8);
787  hist << "\nTest SimOpt solve at feasible (u,z): \n ||c(u,z)|| = " << cnorm << "\n";
788  outStream << hist.str();
789  }
790  return cnorm;
791  }
792 
793 
807  const Vector<Real> &v,
808  const Vector<Real> &u,
809  const Vector<Real> &z,
810  const bool printToStream = true,
811  std::ostream & outStream = std::cout) {
812  return checkAdjointConsistencyJacobian_1(w, v, u, z, w.dual(), v.dual(), printToStream, outStream);
813  }
814 
815 
832  const Vector<Real> &v,
833  const Vector<Real> &u,
834  const Vector<Real> &z,
835  const Vector<Real> &dualw,
836  const Vector<Real> &dualv,
837  const bool printToStream = true,
838  std::ostream & outStream = std::cout) {
839  Real tol = ROL_EPSILON;
840  Teuchos::RCP<Vector<Real> > Jv = dualw.clone();
841  applyJacobian_1(*Jv,v,u,z,tol);
842  Real wJv = w.dot(Jv->dual());
843  Teuchos::RCP<Vector<Real> > Jw = dualv.clone();
844  applyAdjointJacobian_1(*Jw,w,u,z,tol);
845  Real vJw = v.dot(Jw->dual());
846  Real diff = std::abs(wJv-vJw);
847  if ( printToStream ) {
848  std::stringstream hist;
849  hist << std::scientific << std::setprecision(8);
850  hist << "\nTest SimOpt consistency of Jacobian_1 and its adjoint: \n |<w,Jv> - <adj(J)w,v>| = "
851  << diff << "\n";
852  hist << " |<w,Jv>| = " << std::abs(wJv) << "\n";
853  hist << " Relative Error = " << diff / (std::abs(wJv)+ROL_UNDERFLOW) << "\n";
854  outStream << hist.str();
855  }
856  return diff;
857  }
858 
859 
873  const Vector<Real> &v,
874  const Vector<Real> &u,
875  const Vector<Real> &z,
876  const bool printToStream = true,
877  std::ostream & outStream = std::cout) {
878  return checkAdjointConsistencyJacobian_2(w, v, u, z, w.dual(), v.dual(), printToStream, outStream);
879  }
880 
897  const Vector<Real> &v,
898  const Vector<Real> &u,
899  const Vector<Real> &z,
900  const Vector<Real> &dualw,
901  const Vector<Real> &dualv,
902  const bool printToStream = true,
903  std::ostream & outStream = std::cout) {
904  Real tol = ROL_EPSILON;
905  Teuchos::RCP<Vector<Real> > Jv = dualw.clone();
906  applyJacobian_2(*Jv,v,u,z,tol);
907  Real wJv = w.dot(Jv->dual());
908  Teuchos::RCP<Vector<Real> > Jw = dualv.clone();
909  applyAdjointJacobian_2(*Jw,w,u,z,tol);
910  Real vJw = v.dot(Jw->dual());
911  Real diff = std::abs(wJv-vJw);
912  if ( printToStream ) {
913  std::stringstream hist;
914  hist << std::scientific << std::setprecision(8);
915  hist << "\nTest SimOpt consistency of Jacobian_2 and its adjoint: \n |<w,Jv> - <adj(J)w,v>| = "
916  << diff << "\n";
917  hist << " |<w,Jv>| = " << std::abs(wJv) << "\n";
918  hist << " Relative Error = " << diff / (std::abs(wJv)+ROL_UNDERFLOW) << "\n";
919  outStream << hist.str();
920  }
921  return diff;
922  }
923 
924  virtual Real checkInverseJacobian_1(const Vector<Real> &jv,
925  const Vector<Real> &v,
926  const Vector<Real> &u,
927  const Vector<Real> &z,
928  const bool printToStream = true,
929  std::ostream & outStream = std::cout) {
930  Real tol = ROL_EPSILON;
931  Teuchos::RCP<Vector<Real> > Jv = jv.clone();
932  applyJacobian_1(*Jv,v,u,z,tol);
933  Teuchos::RCP<Vector<Real> > iJJv = u.clone();
934  applyInverseJacobian_1(*iJJv,*Jv,u,z,tol);
935  Teuchos::RCP<Vector<Real> > diff = v.clone();
936  diff->set(v);
937  diff->axpy(-1.0,*iJJv);
938  Real dnorm = diff->norm();
939  Real vnorm = v.norm();
940  if ( printToStream ) {
941  std::stringstream hist;
942  hist << std::scientific << std::setprecision(8);
943  hist << "\nTest SimOpt consistency of inverse Jacobian_1: \n ||v-inv(J)Jv|| = "
944  << dnorm << "\n";
945  hist << " ||v|| = " << vnorm << "\n";
946  hist << " Relative Error = " << dnorm / (vnorm+ROL_UNDERFLOW) << "\n";
947  outStream << hist.str();
948  }
949  return dnorm;
950  }
951 
953  const Vector<Real> &v,
954  const Vector<Real> &u,
955  const Vector<Real> &z,
956  const bool printToStream = true,
957  std::ostream & outStream = std::cout) {
958  Real tol = ROL_EPSILON;
959  Teuchos::RCP<Vector<Real> > Jv = jv.clone();
960  applyAdjointJacobian_1(*Jv,v,u,z,tol);
961  Teuchos::RCP<Vector<Real> > iJJv = v.clone();
962  applyInverseAdjointJacobian_1(*iJJv,*Jv,u,z,tol);
963  Teuchos::RCP<Vector<Real> > diff = v.clone();
964  diff->set(v);
965  diff->axpy(-1.0,*iJJv);
966  Real dnorm = diff->norm();
967  Real vnorm = v.norm();
968  if ( printToStream ) {
969  std::stringstream hist;
970  hist << std::scientific << std::setprecision(8);
971  hist << "\nTest SimOpt consistency of inverse adjoint Jacobian_1: \n ||v-inv(adj(J))adj(J)v|| = "
972  << dnorm << "\n";
973  hist << " ||v|| = " << vnorm << "\n";
974  hist << " Relative Error = " << dnorm / (vnorm+ROL_UNDERFLOW) << "\n";
975  outStream << hist.str();
976  }
977  return dnorm;
978  }
979 
980 }; // class EqualityConstraint_SimOpt
981 
982 } // namespace ROL
983 
984 #endif
virtual Real checkInverseJacobian_1(const Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const bool printToStream=true, std::ostream &outStream=std::cout)
virtual void applyAdjointJacobian_2(Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &dualv, Real &tol)
Apply the adjoint of the partial constraint Jacobian at , , to vector . This is the secondary interfa...
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 applyAdjointJacobian_1(Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Apply the adjoint of the partial constraint Jacobian at , , to the vector . This is the primary inter...
virtual void applyAdjointJacobian_2(Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Apply the adjoint of the partial constraint Jacobian at , , to vector . This is the primary interface...
virtual void scale(const Real alpha)=0
Compute where .
Teuchos::RCP< const Vector< Real > > get_1() const
virtual int dimension() const
Return dimension of the vector space.
Definition: ROL_Vector.hpp:181
virtual void applyJacobian_2(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Apply the partial constraint Jacobian at , , to the vector .
virtual void plus(const Vector &x)=0
Compute , where .
virtual void axpy(const Real alpha, const Vector &x)
Compute where .
Definition: ROL_Vector.hpp:141
Defines the linear algebra or vector space interface for simulation-based optimization.
virtual void value(Vector< Real > &c, const Vector< Real > &u, const Vector< Real > &z, Real &tol)=0
Evaluate the constraint operator at .
virtual void applyInverseAdjointJacobian_1(Vector< Real > &iajv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Apply the inverse of the adjoint of the partial constraint Jacobian at , , to the vector ...
virtual void applyAdjointHessian_11(Vector< Real > &ahwv, const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Apply the adjoint of the partial constraint Hessian at , , to vector in direction ...
virtual void applyAdjointHessian_22(Vector< Real > &ahwv, const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Apply the adjoint of the partial constraint Hessian at , , to vector in direction ...
Contains definitions of custom data types in ROL.
void set_1(const Vector< Real > &vec)
virtual void applyAdjointHessian_21(Vector< Real > &ahwv, const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Apply the adjoint of the partial constraint Hessian at , , to vector in direction ...
virtual Teuchos::RCP< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
virtual void applyAdjointHessian(Vector< Real > &ahwv, const Vector< Real > &w, 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 void zero()
Set to zero vector.
Definition: ROL_Vector.hpp:155
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:72
Defines the equality constraint operator interface for simulation-based optimization.
virtual Real dot(const Vector &x) const =0
Compute where .
virtual void value(Vector< Real > &c, const Vector< Real > &x, Real &tol)
Evaluate the constraint operator at .
virtual void applyPreconditioner(Vector< Real > &pv, const Vector< Real > &v, const Vector< Real > &x, const Vector< Real > &g, Real &tol)
Apply a constraint preconditioner at , , to vector . In general, this preconditioner satisfies the fo...
virtual void update(const Vector< Real > &u, const Vector< Real > &z, bool flag=true, int iter=-1)
Update constraint functions. x is the optimization variable, flag = true if optimization variable is ...
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.
Defines the equality constraint operator interface.
virtual Real checkSolve(const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, const ROL::Vector< Real > &c, const bool printToStream=true, std::ostream &outStream=std::cout)
virtual void applyInverseJacobian_1(Vector< Real > &ijv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Apply the inverse partial constraint Jacobian at , , to the vector .
virtual void applyAdjointHessian_12(Vector< Real > &ahwv, const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Apply the adjoint of the partial constraint Hessian at , , to vector in direction ...
virtual bool isFeasible(const Vector< Real > &v)
Check if the vector, v, is feasible.
virtual void applyJacobian_1(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Apply the partial constraint Jacobian at , , to the vector .
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 Real checkAdjointConsistencyJacobian_2(const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &dualw, const Vector< Real > &dualv, const bool printToStream=true, std::ostream &outStream=std::cout)
Check the consistency of the Jacobian and its adjoint. This is the secondary interface, for use with dual spaces where the user does not define the dual() operation.
virtual void solve(Vector< Real > &u, const Vector< Real > &z, Real &tol)
Given , solve for .
virtual void applyPreconditioner(Vector< Real > &pv, const Vector< Real > &v, const Vector< Real > &x, const Vector< Real > &g, Real &tol)
Apply a constraint preconditioner at , , to vector . Ideally, this preconditioner satisfies the follo...
virtual void applyJacobian(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply the constraint Jacobian at , , to vector .
virtual Real checkAdjointConsistencyJacobian_1(const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &dualw, const Vector< Real > &dualv, const bool printToStream=true, std::ostream &outStream=std::cout)
Check the consistency of the Jacobian and its adjoint. This is the secondary interface, for use with dual spaces where the user does not define the dual() operation.
virtual void update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update constraint functions. x is the optimization variable, flag = true if optimization variable is ...
virtual Real checkAdjointConsistencyJacobian_1(const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const bool printToStream=true, std::ostream &outStream=std::cout)
Check the consistency of the Jacobian and its adjoint. This is the primary interface.
virtual void set(const Vector &x)
Set where .
Definition: ROL_Vector.hpp:194
virtual void applyAdjointJacobian_1(Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const Vector< Real > &dualv, Real &tol)
Apply the adjoint of the partial constraint Jacobian at , , to the vector . This is the secondary int...
virtual Real checkInverseAdjointJacobian_1(const Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const bool printToStream=true, std::ostream &outStream=std::cout)
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< 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.
Teuchos::RCP< const Vector< Real > > get_2() const
void set_2(const Vector< Real > &vec)
virtual Real checkAdjointConsistencyJacobian_2(const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, const bool printToStream=true, std::ostream &outStream=std::cout)
Check the consistency of the Jacobian and its adjoint. This is the primary interface.
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