ROL
ROL_Constraint_TimeSimOpt.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Rapid Optimization Library (ROL) Package
4 //
5 // Copyright 2014 NTESS and the ROL contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef ROL_CONSTRAINT_TIMESIMOPT_H
11 #define ROL_CONSTRAINT_TIMESIMOPT_H
12 
14 #include "ROL_VectorWorkspace.hpp"
15 
63 namespace ROL {
64 
65 template <class Real>
67 private:
68 
69  // Get the end point of the time intervals vector
71  {
72  PartitionedVector<Real> & xpv = dynamic_cast<PartitionedVector<Real>&>(x);
73  return *xpv.get(1);
74  }
75 
76  const Vector<Real> & getNewVector(const Vector<Real> & x) const
77  {
78  const PartitionedVector<Real> & xpv = dynamic_cast<const PartitionedVector<Real>&>(x);
79  return *xpv.get(1);
80  }
81 
82  // Get the start point of the time intervals vector
84  {
85  PartitionedVector<Real> & xpv = dynamic_cast<PartitionedVector<Real>&>(x);
86  return *xpv.get(0);
87  }
88 
89  const Vector<Real> & getOldVector(const Vector<Real> & x) const
90  {
91  const PartitionedVector<Real> & xpv = dynamic_cast<const PartitionedVector<Real>&>(x);
92  return *xpv.get(0);
93  }
94 
96 
97 protected:
98 
100 
101 public:
103  : Constraint_SimOpt<Real>()
104  { }
105 
106  // Interface functions (to be overloaded)
107 
115  virtual void update( const Vector<Real> & u_old,
116  const Vector<Real> & u_new,
117  const Vector<Real> &z,
118  bool flag = true, int iter = -1 ) {
119  update_1_old(u_old,flag,iter);
120  update_1_new(u_new,flag,iter);
121  update_2(z,flag,iter);
122  }
123 
129  virtual void update_1_old( const Vector<Real> &u_old, bool flag = true, int iter = -1 ) {}
130 
136  virtual void update_1_new( const Vector<Real> &u_new, bool flag = true, int iter = -1 ) {}
137 
143  virtual void update_2( const Vector<Real> &z, bool flag = true, int iter = -1 ) override {}
144 
161  virtual void value(Vector<Real> &c,
162  const Vector<Real> &u_old,
163  const Vector<Real> &u_new,
164  const Vector<Real> &z,
165  Real &tol) = 0;
166 
167  virtual void solve(Vector<Real> &c,
168  const Vector<Real> &u_old,
169  Vector<Real> &u_new,
170  const Vector<Real> &z,
171  Real &tol) = 0;
172 
173  virtual void applyJacobian_1_old(Vector<Real> &jv,
174  const Vector<Real> &v_old,
175  const Vector<Real> &u_old, const Vector<Real> &u_new,
176  const Vector<Real> &z,
177  Real &tol) = 0;
178 
179  virtual void applyJacobian_1_new(Vector<Real> &jv,
180  const Vector<Real> &v_new,
181  const Vector<Real> &u_old, const Vector<Real> &u_new,
182  const Vector<Real> &z,
183  Real &tol) = 0;
184 
185  virtual void applyInverseJacobian_1_new(Vector<Real> &ijv,
186  const Vector<Real> &v_new,
187  const Vector<Real> &u_old, const Vector<Real> &u_new,
188  const Vector<Real> &z,
189  Real &tol) = 0;
190 
191 
192  virtual void applyJacobian_2(Vector<Real> &jv,
193  const Vector<Real> &v,
194  const Vector<Real> &u_old, const Vector<Real> &u_new,
195  const Vector<Real> &z,
196  Real &tol) = 0;
197 
198  virtual void applyAdjointJacobian_1_old(Vector<Real> &ajv_old,
199  const Vector<Real> &dualv,
200  const Vector<Real> &u_old, const Vector<Real> &u_new,
201  const Vector<Real> &z,
202  Real &tol) = 0;
203 
204  virtual void applyAdjointJacobian_1_new(Vector<Real> &ajv_new,
205  const Vector<Real> &dualv,
206  const Vector<Real> &u_old, const Vector<Real> &u_new,
207  const Vector<Real> &z,
208  Real &tol) = 0;
209 
211  const Vector<Real> &v_new,
212  const Vector<Real> &u_old, const Vector<Real> &u_new,
213  const Vector<Real> &z,
214  Real &tol) = 0;
215 
216  virtual void applyAdjointJacobian_2_time(Vector<Real> &ajv,
217  const Vector<Real> &dualv,
218  const Vector<Real> &u_old, const Vector<Real> &u_new,
219  const Vector<Real> &z,
220  Real &tol) = 0;
221 
222  virtual void applyAdjointHessian_11_old(Vector<Real> &ahwv_old,
223  const Vector<Real> &w,
224  const Vector<Real> &v_new,
225  const Vector<Real> &u_old, const Vector<Real> &u_new,
226  const Vector<Real> &z,
227  Real &tol) = 0;
228 
229  virtual void applyAdjointHessian_11_new(Vector<Real> &ahwv_new,
230  const Vector<Real> &w,
231  const Vector<Real> &v_new,
232  const Vector<Real> &u_old, const Vector<Real> &u_new,
233  const Vector<Real> &z,
234  Real &tol) = 0;
235 
236  // Functions from SimOpt that are overriden
238 
239  virtual void update( const Vector<Real> & u,
240  const Vector<Real> & z,
241  bool flag = true, int iter = -1 ) override {
242  update(getOldVector(u),
243  getNewVector(u),
244  z,
245  flag,iter);
246  }
247 
248  virtual void value(Vector<Real> &c,
249  const Vector<Real> &u,
250  const Vector<Real> &z,
251  Real &tol) override {
252 
253  value(c,
254  getOldVector(u),
255  getNewVector(u),
256  z,
257  tol);
258  }
259 
260  virtual void solve(Vector<Real> &c,
261  Vector<Real> &u,
262  const Vector<Real> &z,
263  Real &tol) override {
264  solve(c,
265  getOldVector(u),
266  getNewVector(u),
267  z,
268  tol);
269  }
270 
271  virtual void applyJacobian_1(Vector<Real> &jv,
272  const Vector<Real> &v,
273  const Vector<Real> &u,
274  const Vector<Real> &z,
275  Real &tol) override {
276  const Vector<Real> & v_old = getOldVector(v);
277  const Vector<Real> & v_new = getNewVector(v);
278  const Vector<Real> & u_old = getOldVector(u);
279  const Vector<Real> & u_new = getNewVector(u);
280 
281  // evaluate derivative against "old" time variable
282  applyJacobian_1_old(jv,v_old,
283  u_old,u_new,
284  z,
285  tol);
286 
287  ROL::Ptr<Vector<Real> > jv_new = workspace_.clone(jv);
288 
289  // evaluate derivative against "new" time variable
290  applyJacobian_1_new(*jv_new,v_new,
291  u_old,u_new,
292  z,
293  tol);
294 
295  jv.axpy(1.0,*jv_new);
296  }
297 
298  virtual void applyJacobian_2(Vector<Real> &jv,
299  const Vector<Real> &v,
300  const Vector<Real> &u,
301  const Vector<Real> &z,
302  Real &tol) override {
303  const Vector<Real> & u_old = getOldVector(u);
304  const Vector<Real> & u_new = getNewVector(u);
305 
306  // evaluate derivative against "old" time variable
307  applyJacobian_2(jv,v,
308  u_old,u_new,
309  z,
310  tol);
311  }
312 
314  const Vector<Real> &v,
315  const Vector<Real> &u,
316  const Vector<Real> &z,
317  Real &tol) override final {
318  ROL_TEST_FOR_EXCEPTION(true, std::logic_error,
319  "The method applyInverseJacobian_1 is used but not implemented!\n");
320  }
321 
323  const Vector<Real> &v,
324  const Vector<Real> &u,
325  const Vector<Real> &z,
326  Real &tol) override {
327  Vector<Real> & ajv_old = getOldVector(ajv);
328  Vector<Real> & ajv_new = getNewVector(ajv);
329  const Vector<Real> & u_old = getOldVector(u);
330  const Vector<Real> & u_new = getNewVector(u);
331 
332  applyAdjointJacobian_1_old(ajv_old,v,u_old,u_new,z,tol);
333  applyAdjointJacobian_1_new(ajv_new,v,u_old,u_new,z,tol);
334  }
335 
337  const Vector<Real> &v,
338  const Vector<Real> &u,
339  const Vector<Real> &z,
340  Real &tol) override {
341  const Vector<Real> & u_old = getOldVector(u);
342  const Vector<Real> & u_new = getNewVector(u);
343 
344  applyAdjointJacobian_2_time(ajv,v,u_old,u_new,z,tol);
345  }
346 
348  const Vector<Real> &v,
349  const Vector<Real> &u,
350  const Vector<Real> &z,
351  Real &tol) override final {
352  ROL_TEST_FOR_EXCEPTION(true, std::logic_error,
353  "The method applyInverseAdjointJacobian_1 is used but not implemented!\n");
354  };
355 
374  const Vector<Real> &w,
375  const Vector<Real> &v,
376  const Vector<Real> &u,
377  const Vector<Real> &z,
378  Real &tol) override
379  {
380  Vector<Real> & ahwv_old = getOldVector(ahwv);
381  Vector<Real> & ahwv_new = getNewVector(ahwv);
382  const Vector<Real> & v_old = getOldVector(v);
383  const Vector<Real> & v_new = getNewVector(v);
384  const Vector<Real> & u_old = getOldVector(u);
385  const Vector<Real> & u_new = getNewVector(u);
386 
387  // this implicitly assumes that there is no cross coupling
388  // between the old state and the new state. Is that true? For
389  // simple (Euler, Theta method) integrators yes.
390  applyAdjointHessian_11_old(ahwv_old,w,v_old,u_old,u_new,z,tol);
391  applyAdjointHessian_11_new(ahwv_new,w,v_new,u_old,u_new,z,tol);
392  }
393 
412  const Vector<Real> &w,
413  const Vector<Real> &v,
414  const Vector<Real> &u,
415  const Vector<Real> &z,
416  Real &tol) override
417  {
418  ahwv.zero();
419  }
420 
439  const Vector<Real> &w,
440  const Vector<Real> &v,
441  const Vector<Real> &u,
442  const Vector<Real> &z,
443  Real &tol) override {
444  ahwv.zero();
445  }
446 
465  const Vector<Real> &w,
466  const Vector<Real> &v,
467  const Vector<Real> &u,
468  const Vector<Real> &z,
469  Real &tol) override {
470  ahwv.zero();
471  }
472 
473  // We override the check solve routine because we are abusing SimOpt
474  virtual Real checkSolve(const ROL::Vector<Real> &u,
475  const ROL::Vector<Real> &z,
476  const ROL::Vector<Real> &c,
477  const bool printToStream = true,
478  std::ostream & outStream = std::cout) override {
479  // Solve constraint for u.
480  Real tol = ROL_EPSILON<Real>();
481  ROL::Ptr<ROL::Vector<Real> > r = workspace_.clone(c);
482  ROL::Ptr<ROL::Vector<Real> > s = workspace_.clone(u);
483  s->set(u);
484  solve(*r,*s,z,tol);
485  // Evaluate constraint residual at (u,z).
486  ROL::Ptr<ROL::Vector<Real> > cs = workspace_.clone(c);
487  update(*s,z);
488  value(*cs,*s,z,tol);
489  // Output norm of residual.
490  Real rnorm = r->norm();
491  Real cnorm = cs->norm();
492  if ( printToStream ) {
493  std::stringstream hist;
494  hist << std::scientific << std::setprecision(8);
495  hist << "\nTest SimOpt solve at feasible (u,z):\n";
496  hist << " Solver Residual = " << rnorm << "\n";
497  hist << " ||c(u,z)|| = " << cnorm << "\n";
498  outStream << hist.str();
499  }
500  return cnorm;
501  }
502 
503  // Verify that ||v-Jinv*J*v|| < tol
505  const ROL::Vector<Real> &u_new,
506  const ROL::Vector<Real> &u_old,
507  const ROL::Vector<Real> &z,
508  const ROL::Vector<Real> &v_new,
509  const bool printToStream = true,
510  std::ostream & outStream = std::cout) {
511  Real tol = ROL_EPSILON<Real>();
512  auto Jv = workspace_.clone(c);
513  update( u_new, u_old, z );
514  applyJacobian_1_new( *Jv, v_new, u_old, u_new, z, tol );
515  auto iJJv = workspace_.clone(u_new);
516  update( u_new, u_old, z );
517  applyInverseJacobian_1_new( *iJJv, *Jv, u_old, u_new, z, tol );
518  auto diff = workspace_.clone(v_new);
519  diff->set(v_new);
520  diff->axpy(-1.0,*iJJv);
521  Real dnorm = diff->norm();
522  Real vnorm = v_new.norm();
523  if ( printToStream ) {
524  std::stringstream hist;
525  hist << std::scientific << std::setprecision(8);
526  hist << "\nTest TimeSimOpt consistency of inverse Jacobian_1_new: \n ||v-inv(J)Jv|| = "
527  << dnorm << "\n";
528  hist << " ||v|| = " << vnorm << "\n";
529  hist << " Relative Error = " << dnorm / (vnorm+ROL_UNDERFLOW<Real>()) << "\n";
530  outStream << hist.str();
531  }
532  return dnorm;
533  }
534 
536  const ROL::Vector<Real> &u_new,
537  const ROL::Vector<Real> &u_old,
538  const ROL::Vector<Real> &z,
539  const ROL::Vector<Real> &v_new,
540  const bool printToStream = true,
541  std::ostream & outStream = std::cout) {
542  Real tol = ROL_EPSILON<Real>();
543  auto Jv = workspace_.clone(c);
544  update( u_new, u_old, z );
545  applyAdjointJacobian_1_new( *Jv, v_new, u_old, u_new, z, tol );
546  auto iJJv = workspace_.clone(u_new);
547  update( u_new, u_old, z );
548  applyInverseAdjointJacobian_1_new( *iJJv, *Jv, u_old, u_new, z, tol );
549  auto diff = workspace_.clone(v_new);
550  diff->set(v_new);
551  diff->axpy(-1.0,*iJJv);
552  Real dnorm = diff->norm();
553  Real vnorm = v_new.norm();
554  if ( printToStream ) {
555  std::stringstream hist;
556  hist << std::scientific << std::setprecision(8);
557  hist << "\nTest TimeSimOpt consistency of inverse adjoint Jacobian_1_new: \n ||v-inv(adj(J))adj(J)v|| = "
558  << dnorm << "\n";
559  hist << " ||v|| = " << vnorm << "\n";
560  hist << " Relative Error = " << dnorm / (vnorm+ROL_UNDERFLOW<Real>()) << "\n";
561  outStream << hist.str();
562  }
563  return dnorm;
564  }
565 
566  std::vector<std::vector<Real> > checkApplyJacobian_1_new(const Vector<Real> &u_new,
567  const Vector<Real> &u_old,
568  const Vector<Real> &z,
569  const Vector<Real> &v,
570  const Vector<Real> &jv,
571  const bool printToStream = true,
572  std::ostream & outStream = std::cout,
573  const int numSteps = ROL_NUM_CHECKDERIV_STEPS,
574  const int order = 1) {
575  std::vector<Real> steps(numSteps);
576  for(int i=0;i<numSteps;++i) {
577  steps[i] = pow(10,-i);
578  }
579 
580  return checkApplyJacobian_1_new(u_new,u_old,z,v,jv,steps,printToStream,outStream,order);
581  }
582 
583  std::vector<std::vector<Real> > checkApplyJacobian_1_new(const Vector<Real> &u_new,
584  const Vector<Real> &u_old,
585  const Vector<Real> &z,
586  const Vector<Real> &v,
587  const Vector<Real> &jv,
588  const std::vector<Real> &steps,
589  const bool printToStream = true,
590  std::ostream & outStream = std::cout,
591  const int order = 1) {
592 
593  ROL_TEST_FOR_EXCEPTION( order<1 || order>4, std::invalid_argument,
594  "Error: finite difference order must be 1,2,3, or 4" );
595 
596  Real one(1.0);
597 
600 
601  Real tol = std::sqrt(ROL_EPSILON<Real>());
602 
603  int numSteps = steps.size();
604  int numVals = 4;
605  std::vector<Real> tmp(numVals);
606  std::vector<std::vector<Real> > jvCheck(numSteps, tmp);
607 
608  // Save the format state of the original outStream.
609  ROL::nullstream oldFormatState;
610  oldFormatState.copyfmt(outStream);
611 
612  // Compute constraint value at x.
613  ROL::Ptr<Vector<Real> > c = workspace_.clone(jv);
614  this->update(u_new, u_old, z);
615  this->value(*c, u_new, u_old, z, tol);
616 
617  // Compute (Jacobian at x) times (vector v).
618  ROL::Ptr<Vector<Real> > Jv = workspace_.clone(jv);
619  this->applyJacobian_1_new(*Jv, v, u_new, u_old, z, tol);
620  Real normJv = Jv->norm();
621 
622  // Temporary vectors.
623  ROL::Ptr<Vector<Real> > cdif = workspace_.clone(jv);
624  ROL::Ptr<Vector<Real> > cnew = workspace_.clone(jv);
625  ROL::Ptr<Vector<Real> > u_2 = workspace_.clone(u_new);
626 
627  for (int i=0; i<numSteps; i++) {
628 
629  Real eta = steps[i];
630 
631  u_2->set(u_new);
632 
633  cdif->set(*c);
634  cdif->scale(weights[order-1][0]);
635 
636  for(int j=0; j<order; ++j) {
637 
638  u_2->axpy(eta*shifts[order-1][j], v);
639 
640  if( weights[order-1][j+1] != 0 ) {
641  this->update(*u_2,u_old,z);
642  this->value(*cnew,*u_2,u_old,z,tol);
643  cdif->axpy(weights[order-1][j+1],*cnew);
644  }
645 
646  }
647 
648  cdif->scale(one/eta);
649 
650  // Compute norms of Jacobian-vector products, finite-difference approximations, and error.
651  jvCheck[i][0] = eta;
652  jvCheck[i][1] = normJv;
653  jvCheck[i][2] = cdif->norm();
654  cdif->axpy(-one, *Jv);
655  jvCheck[i][3] = cdif->norm();
656 
657  if (printToStream) {
658  std::stringstream hist;
659  if (i==0) {
660  hist << std::right
661  << std::setw(20) << "Step size"
662  << std::setw(20) << "norm(Jac*vec)"
663  << std::setw(20) << "norm(FD approx)"
664  << std::setw(20) << "norm(abs error)"
665  << "\n"
666  << std::setw(20) << "---------"
667  << std::setw(20) << "-------------"
668  << std::setw(20) << "---------------"
669  << std::setw(20) << "---------------"
670  << "\n";
671  }
672  hist << std::scientific << std::setprecision(11) << std::right
673  << std::setw(20) << jvCheck[i][0]
674  << std::setw(20) << jvCheck[i][1]
675  << std::setw(20) << jvCheck[i][2]
676  << std::setw(20) << jvCheck[i][3]
677  << "\n";
678  outStream << hist.str();
679  }
680 
681  }
682 
683  // Reset format state of outStream.
684  outStream.copyfmt(oldFormatState);
685 
686  return jvCheck;
687  } // checkApplyJacobian_1_new
688 
689 
690 }; // class Constraint_SimOpt
691 
692 } // namespace ROL
693 
694 #endif
virtual void applyAdjointHessian_11_old(Vector< Real > &ahwv_old, const Vector< Real > &w, const Vector< Real > &v_new, const Vector< Real > &u_old, const Vector< Real > &u_new, const Vector< Real > &z, Real &tol)=0
virtual void applyJacobian_2(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &u_old, const Vector< Real > &u_new, const Vector< Real > &z, Real &tol)=0
virtual void update(const Vector< Real > &u_old, const Vector< Real > &u_new, const Vector< Real > &z, bool flag=true, int iter=-1)
Update constraint functions. u_old Is the state from the end of the previous time step...
virtual void applyInverseAdjointJacobian_1(Vector< Real > &iajv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) overridefinal
Apply the inverse of the adjoint of the partial constraint Jacobian at , , to the vector ...
ROL::Ptr< const Vector< Real > > get(size_type i) const
const double weights[4][5]
Definition: ROL_Types.hpp:838
virtual void axpy(const Real alpha, const Vector &x)
Compute where .
Definition: ROL_Vector.hpp:119
Defines the linear algebra of vector space on a generic partitioned vector.
std::vector< std::vector< Real > > checkApplyJacobian_1_new(const Vector< Real > &u_new, const Vector< Real > &u_old, const Vector< Real > &z, 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)
virtual void applyAdjointJacobian_1_old(Vector< Real > &ajv_old, const Vector< Real > &dualv, const Vector< Real > &u_old, const Vector< Real > &u_new, const Vector< Real > &z, Real &tol)=0
virtual void solve(Vector< Real > &c, const Vector< Real > &u_old, Vector< Real > &u_new, const Vector< Real > &z, Real &tol)=0
virtual void zero()
Set to zero vector.
Definition: ROL_Vector.hpp:133
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:46
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) override
Apply the optimization-space derivative of the adjoint of the constraint simulation-space Jacobian at...
Defines the time dependent constraint operator interface for simulation-based optimization.
virtual Real checkInverseJacobian_1_new(const ROL::Vector< Real > &c, const ROL::Vector< Real > &u_new, const ROL::Vector< Real > &u_old, const ROL::Vector< Real > &z, const ROL::Vector< Real > &v_new, const bool printToStream=true, std::ostream &outStream=std::cout)
basic_nullstream< char, std::char_traits< char >> nullstream
Definition: ROL_Stream.hpp:36
virtual void applyAdjointJacobian_2(Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
Apply the adjoint of the partial constraint Jacobian at , , to vector . This is the primary interface...
Vector< Real > & getNewVector(Vector< Real > &x) const
virtual void applyJacobian_1_old(Vector< Real > &jv, const Vector< Real > &v_old, const Vector< Real > &u_old, const Vector< Real > &u_new, const Vector< Real > &z, Real &tol)=0
virtual void applyInverseAdjointJacobian_1_new(Vector< Real > &iajv, const Vector< Real > &v_new, const Vector< Real > &u_old, const Vector< Real > &u_new, const Vector< Real > &z, Real &tol)=0
virtual void update_1_old(const Vector< Real > &u_old, bool flag=true, int iter=-1)
Update constraint functions with respect to Sim variable. u_old is the state variable flag = true if ...
virtual void applyJacobian_1(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
Apply the partial constraint Jacobian at , , to the vector .
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) override
Apply the optimization-space derivative of the adjoint of the constraint optimization-space Jacobian ...
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) override
Apply the simulation-space derivative of the adjoint of the constraint simulation-space Jacobian at ...
virtual void solve(Vector< Real > &c, Vector< Real > &u, const Vector< Real > &z, Real &tol) override
Given , solve for .
virtual void value(Vector< Real > &c, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
Evaluate the constraint operator at .
const Vector< Real > & getNewVector(const Vector< Real > &x) const
virtual void applyAdjointJacobian_1_new(Vector< Real > &ajv_new, const Vector< Real > &dualv, const Vector< Real > &u_old, const Vector< Real > &u_new, const Vector< Real > &z, Real &tol)=0
virtual void applyAdjointJacobian_1(Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
Apply the adjoint of the partial constraint Jacobian at , , to the vector . This is the primary inter...
virtual void update(const Vector< Real > &u, const Vector< Real > &z, bool flag=true, int iter=-1) override
Update constraint functions. x is the optimization variable, flag = true if optimization variable is ...
#define ROL_NUM_CHECKDERIV_STEPS
Number of steps for derivative checks.
Definition: ROL_Types.hpp:40
virtual void applyJacobian_2(Vector< Real > &jv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) override
Apply the partial constraint Jacobian at , , to the vector .
virtual void update_1_new(const Vector< Real > &u_new, bool flag=true, int iter=-1)
Update constraint functions with respect to Sim variable. u_new is the state variable flag = true if ...
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) override
virtual void applyInverseJacobian_1(Vector< Real > &ijv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol) overridefinal
Apply the inverse partial constraint Jacobian at , , to the vector .
std::vector< std::vector< Real > > checkApplyJacobian_1_new(const Vector< Real > &u_new, const Vector< Real > &u_old, const Vector< Real > &z, const Vector< Real > &v, const Vector< Real > &jv, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
virtual Real norm() const =0
Returns where .
virtual void update_2(const Vector< Real > &z, bool flag=true, int iter=-1) override
Update constraint functions with respect to Opt variable. z is the control variable, flag = true if optimization variable is changed, iter is the outer algorithm iterations count.
Defines the constraint operator interface for simulation-based optimization.
VectorWorkspace< Real > & getVectorWorkspace() const
virtual void value(Vector< Real > &c, const Vector< Real > &u_old, const Vector< Real > &u_new, const Vector< Real > &z, Real &tol)=0
Evaluate the constraint operator at .
virtual void applyInverseJacobian_1_new(Vector< Real > &ijv, const Vector< Real > &v_new, const Vector< Real > &u_old, const Vector< Real > &u_new, const Vector< Real > &z, Real &tol)=0
virtual void applyJacobian_1_new(Vector< Real > &jv, const Vector< Real > &v_new, const Vector< Real > &u_old, const Vector< Real > &u_new, const Vector< Real > &z, Real &tol)=0
virtual void applyAdjointJacobian_2_time(Vector< Real > &ajv, const Vector< Real > &dualv, const Vector< Real > &u_old, const Vector< Real > &u_new, const Vector< Real > &z, Real &tol)=0
virtual void applyAdjointHessian_11_new(Vector< Real > &ahwv_new, const Vector< Real > &w, const Vector< Real > &v_new, const Vector< Real > &u_old, const Vector< Real > &u_new, const Vector< Real > &z, Real &tol)=0
const Vector< Real > & getOldVector(const Vector< Real > &x) const
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) override
Apply the simulation-space derivative of the adjoint of the constraint optimization-space Jacobian at...
virtual Real checkInverseAdjointJacobian_1_new(const ROL::Vector< Real > &c, const ROL::Vector< Real > &u_new, const ROL::Vector< Real > &u_old, const ROL::Vector< Real > &z, const ROL::Vector< Real > &v_new, const bool printToStream=true, std::ostream &outStream=std::cout)
Vector< Real > & getOldVector(Vector< Real > &x) const