ROL
ROL_Constraint_TimeSimOpt.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_CONSTRAINT_TIMESIMOPT_H
45 #define ROL_CONSTRAINT_TIMESIMOPT_H
46 
48 #include "ROL_VectorWorkspace.hpp"
49 
97 namespace ROL {
98 
99 template <class Real>
101 private:
102 
103  // Get the end point of the time intervals vector
105  {
106  PartitionedVector<Real> & xpv = dynamic_cast<PartitionedVector<Real>&>(x);
107  return *xpv.get(1);
108  }
109 
110  const Vector<Real> & getNewVector(const Vector<Real> & x) const
111  {
112  const PartitionedVector<Real> & xpv = dynamic_cast<const PartitionedVector<Real>&>(x);
113  return *xpv.get(1);
114  }
115 
116  // Get the start point of the time intervals vector
118  {
119  PartitionedVector<Real> & xpv = dynamic_cast<PartitionedVector<Real>&>(x);
120  return *xpv.get(0);
121  }
122 
123  const Vector<Real> & getOldVector(const Vector<Real> & x) const
124  {
125  const PartitionedVector<Real> & xpv = dynamic_cast<const PartitionedVector<Real>&>(x);
126  return *xpv.get(0);
127  }
128 
130 
131 protected:
132 
134 
135 public:
137  : Constraint_SimOpt<Real>()
138  { }
139 
140  // Interface functions (to be overloaded)
141 
149  virtual void update( const Vector<Real> & u_old,
150  const Vector<Real> & u_new,
151  const Vector<Real> &z,
152  bool flag = true, int iter = -1 ) {
153  update_1_old(u_old,flag,iter);
154  update_1_new(u_new,flag,iter);
155  update_2(z,flag,iter);
156  }
157 
163  virtual void update_1_old( const Vector<Real> &u_old, bool flag = true, int iter = -1 ) {}
164 
170  virtual void update_1_new( const Vector<Real> &u_new, bool flag = true, int iter = -1 ) {}
171 
177  virtual void update_2( const Vector<Real> &z, bool flag = true, int iter = -1 ) override {}
178 
195  virtual void value(Vector<Real> &c,
196  const Vector<Real> &u_old,
197  const Vector<Real> &u_new,
198  const Vector<Real> &z,
199  Real &tol) = 0;
200 
201  virtual void solve(Vector<Real> &c,
202  const Vector<Real> &u_old,
203  Vector<Real> &u_new,
204  const Vector<Real> &z,
205  Real &tol) = 0;
206 
207  virtual void applyJacobian_1_old(Vector<Real> &jv,
208  const Vector<Real> &v_old,
209  const Vector<Real> &u_old, const Vector<Real> &u_new,
210  const Vector<Real> &z,
211  Real &tol) = 0;
212 
213  virtual void applyJacobian_1_new(Vector<Real> &jv,
214  const Vector<Real> &v_new,
215  const Vector<Real> &u_old, const Vector<Real> &u_new,
216  const Vector<Real> &z,
217  Real &tol) = 0;
218 
219  virtual void applyInverseJacobian_1_new(Vector<Real> &ijv,
220  const Vector<Real> &v_new,
221  const Vector<Real> &u_old, const Vector<Real> &u_new,
222  const Vector<Real> &z,
223  Real &tol) = 0;
224 
225 
226  virtual void applyJacobian_2(Vector<Real> &jv,
227  const Vector<Real> &v,
228  const Vector<Real> &u_old, const Vector<Real> &u_new,
229  const Vector<Real> &z,
230  Real &tol) = 0;
231 
232  virtual void applyAdjointJacobian_1_old(Vector<Real> &ajv_old,
233  const Vector<Real> &dualv,
234  const Vector<Real> &u_old, const Vector<Real> &u_new,
235  const Vector<Real> &z,
236  Real &tol) = 0;
237 
238  virtual void applyAdjointJacobian_1_new(Vector<Real> &ajv_new,
239  const Vector<Real> &dualv,
240  const Vector<Real> &u_old, const Vector<Real> &u_new,
241  const Vector<Real> &z,
242  Real &tol) = 0;
243 
245  const Vector<Real> &v_new,
246  const Vector<Real> &u_old, const Vector<Real> &u_new,
247  const Vector<Real> &z,
248  Real &tol) = 0;
249 
250  virtual void applyAdjointJacobian_2_time(Vector<Real> &ajv,
251  const Vector<Real> &dualv,
252  const Vector<Real> &u_old, const Vector<Real> &u_new,
253  const Vector<Real> &z,
254  Real &tol) = 0;
255 
256  virtual void applyAdjointHessian_11_old(Vector<Real> &ahwv_old,
257  const Vector<Real> &w,
258  const Vector<Real> &v_new,
259  const Vector<Real> &u_old, const Vector<Real> &u_new,
260  const Vector<Real> &z,
261  Real &tol) = 0;
262 
263  virtual void applyAdjointHessian_11_new(Vector<Real> &ahwv_new,
264  const Vector<Real> &w,
265  const Vector<Real> &v_new,
266  const Vector<Real> &u_old, const Vector<Real> &u_new,
267  const Vector<Real> &z,
268  Real &tol) = 0;
269 
270  // Functions from SimOpt that are overriden
272 
273  virtual void update( const Vector<Real> & u,
274  const Vector<Real> & z,
275  bool flag = true, int iter = -1 ) override {
276  update(getOldVector(u),
277  getNewVector(u),
278  z,
279  flag,iter);
280  }
281 
282  virtual void value(Vector<Real> &c,
283  const Vector<Real> &u,
284  const Vector<Real> &z,
285  Real &tol) override {
286 
287  value(c,
288  getOldVector(u),
289  getNewVector(u),
290  z,
291  tol);
292  }
293 
294  virtual void solve(Vector<Real> &c,
295  Vector<Real> &u,
296  const Vector<Real> &z,
297  Real &tol) override {
298  solve(c,
299  getOldVector(u),
300  getNewVector(u),
301  z,
302  tol);
303  }
304 
305  virtual void applyJacobian_1(Vector<Real> &jv,
306  const Vector<Real> &v,
307  const Vector<Real> &u,
308  const Vector<Real> &z,
309  Real &tol) override {
310  const Vector<Real> & v_old = getOldVector(v);
311  const Vector<Real> & v_new = getNewVector(v);
312  const Vector<Real> & u_old = getOldVector(u);
313  const Vector<Real> & u_new = getNewVector(u);
314 
315  // evaluate derivative against "old" time variable
316  applyJacobian_1_old(jv,v_old,
317  u_old,u_new,
318  z,
319  tol);
320 
321  ROL::Ptr<Vector<Real> > jv_new = workspace_.clone(jv);
322 
323  // evaluate derivative against "new" time variable
324  applyJacobian_1_new(*jv_new,v_new,
325  u_old,u_new,
326  z,
327  tol);
328 
329  jv.axpy(1.0,*jv_new);
330  }
331 
332  virtual void applyJacobian_2(Vector<Real> &jv,
333  const Vector<Real> &v,
334  const Vector<Real> &u,
335  const Vector<Real> &z,
336  Real &tol) override {
337  const Vector<Real> & u_old = getOldVector(u);
338  const Vector<Real> & u_new = getNewVector(u);
339 
340  // evaluate derivative against "old" time variable
341  applyJacobian_2(jv,v,
342  u_old,u_new,
343  z,
344  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 applyInverseJacobian_1 is used but not implemented!\n");
354  }
355 
357  const Vector<Real> &v,
358  const Vector<Real> &u,
359  const Vector<Real> &z,
360  Real &tol) override {
361  Vector<Real> & ajv_old = getOldVector(ajv);
362  Vector<Real> & ajv_new = getNewVector(ajv);
363  const Vector<Real> & u_old = getOldVector(u);
364  const Vector<Real> & u_new = getNewVector(u);
365 
366  applyAdjointJacobian_1_old(ajv_old,v,u_old,u_new,z,tol);
367  applyAdjointJacobian_1_new(ajv_new,v,u_old,u_new,z,tol);
368  }
369 
371  const Vector<Real> &v,
372  const Vector<Real> &u,
373  const Vector<Real> &z,
374  Real &tol) override {
375  const Vector<Real> & u_old = getOldVector(u);
376  const Vector<Real> & u_new = getNewVector(u);
377 
378  applyAdjointJacobian_2_time(ajv,v,u_old,u_new,z,tol);
379  }
380 
382  const Vector<Real> &v,
383  const Vector<Real> &u,
384  const Vector<Real> &z,
385  Real &tol) override final {
386  ROL_TEST_FOR_EXCEPTION(true, std::logic_error,
387  "The method applyInverseAdjointJacobian_1 is used but not implemented!\n");
388  };
389 
408  const Vector<Real> &w,
409  const Vector<Real> &v,
410  const Vector<Real> &u,
411  const Vector<Real> &z,
412  Real &tol) override
413  {
414  Vector<Real> & ahwv_old = getOldVector(ahwv);
415  Vector<Real> & ahwv_new = getNewVector(ahwv);
416  const Vector<Real> & v_old = getOldVector(v);
417  const Vector<Real> & v_new = getNewVector(v);
418  const Vector<Real> & u_old = getOldVector(u);
419  const Vector<Real> & u_new = getNewVector(u);
420 
421  // this implicitly assumes that there is no cross coupling
422  // between the old state and the new state. Is that true? For
423  // simple (Euler, Theta method) integrators yes.
424  applyAdjointHessian_11_old(ahwv_old,w,v_old,u_old,u_new,z,tol);
425  applyAdjointHessian_11_new(ahwv_new,w,v_new,u_old,u_new,z,tol);
426  }
427 
446  const Vector<Real> &w,
447  const Vector<Real> &v,
448  const Vector<Real> &u,
449  const Vector<Real> &z,
450  Real &tol) override
451  {
452  ahwv.zero();
453  }
454 
473  const Vector<Real> &w,
474  const Vector<Real> &v,
475  const Vector<Real> &u,
476  const Vector<Real> &z,
477  Real &tol) override {
478  ahwv.zero();
479  }
480 
499  const Vector<Real> &w,
500  const Vector<Real> &v,
501  const Vector<Real> &u,
502  const Vector<Real> &z,
503  Real &tol) override {
504  ahwv.zero();
505  }
506 
507  // We override the check solve routine because we are abusing SimOpt
508  virtual Real checkSolve(const ROL::Vector<Real> &u,
509  const ROL::Vector<Real> &z,
510  const ROL::Vector<Real> &c,
511  const bool printToStream = true,
512  std::ostream & outStream = std::cout) override {
513  // Solve constraint for u.
514  Real tol = ROL_EPSILON<Real>();
515  ROL::Ptr<ROL::Vector<Real> > r = workspace_.clone(c);
516  ROL::Ptr<ROL::Vector<Real> > s = workspace_.clone(u);
517  s->set(u);
518  solve(*r,*s,z,tol);
519  // Evaluate constraint residual at (u,z).
520  ROL::Ptr<ROL::Vector<Real> > cs = workspace_.clone(c);
521  update(*s,z);
522  value(*cs,*s,z,tol);
523  // Output norm of residual.
524  Real rnorm = r->norm();
525  Real cnorm = cs->norm();
526  if ( printToStream ) {
527  std::stringstream hist;
528  hist << std::scientific << std::setprecision(8);
529  hist << "\nTest SimOpt solve at feasible (u,z):\n";
530  hist << " Solver Residual = " << rnorm << "\n";
531  hist << " ||c(u,z)|| = " << cnorm << "\n";
532  outStream << hist.str();
533  }
534  return cnorm;
535  }
536 
537  // Verify that ||v-Jinv*J*v|| < tol
539  const ROL::Vector<Real> &u_new,
540  const ROL::Vector<Real> &u_old,
541  const ROL::Vector<Real> &z,
542  const ROL::Vector<Real> &v_new,
543  const bool printToStream = true,
544  std::ostream & outStream = std::cout) {
545  Real tol = ROL_EPSILON<Real>();
546  auto Jv = workspace_.clone(c);
547  update( u_new, u_old, z );
548  applyJacobian_1_new( *Jv, v_new, u_old, u_new, z, tol );
549  auto iJJv = workspace_.clone(u_new);
550  update( u_new, u_old, z );
551  applyInverseJacobian_1_new( *iJJv, *Jv, u_old, u_new, z, tol );
552  auto diff = workspace_.clone(v_new);
553  diff->set(v_new);
554  diff->axpy(-1.0,*iJJv);
555  Real dnorm = diff->norm();
556  Real vnorm = v_new.norm();
557  if ( printToStream ) {
558  std::stringstream hist;
559  hist << std::scientific << std::setprecision(8);
560  hist << "\nTest TimeSimOpt consistency of inverse Jacobian_1_new: \n ||v-inv(J)Jv|| = "
561  << dnorm << "\n";
562  hist << " ||v|| = " << vnorm << "\n";
563  hist << " Relative Error = " << dnorm / (vnorm+ROL_UNDERFLOW<Real>()) << "\n";
564  outStream << hist.str();
565  }
566  return dnorm;
567  }
568 
570  const ROL::Vector<Real> &u_new,
571  const ROL::Vector<Real> &u_old,
572  const ROL::Vector<Real> &z,
573  const ROL::Vector<Real> &v_new,
574  const bool printToStream = true,
575  std::ostream & outStream = std::cout) {
576  Real tol = ROL_EPSILON<Real>();
577  auto Jv = workspace_.clone(c);
578  update( u_new, u_old, z );
579  applyAdjointJacobian_1_new( *Jv, v_new, u_old, u_new, z, tol );
580  auto iJJv = workspace_.clone(u_new);
581  update( u_new, u_old, z );
582  applyInverseAdjointJacobian_1_new( *iJJv, *Jv, u_old, u_new, z, tol );
583  auto diff = workspace_.clone(v_new);
584  diff->set(v_new);
585  diff->axpy(-1.0,*iJJv);
586  Real dnorm = diff->norm();
587  Real vnorm = v_new.norm();
588  if ( printToStream ) {
589  std::stringstream hist;
590  hist << std::scientific << std::setprecision(8);
591  hist << "\nTest TimeSimOpt consistency of inverse adjoint Jacobian_1_new: \n ||v-inv(adj(J))adj(J)v|| = "
592  << dnorm << "\n";
593  hist << " ||v|| = " << vnorm << "\n";
594  hist << " Relative Error = " << dnorm / (vnorm+ROL_UNDERFLOW<Real>()) << "\n";
595  outStream << hist.str();
596  }
597  return dnorm;
598  }
599 
600  std::vector<std::vector<Real> > checkApplyJacobian_1_new(const Vector<Real> &u_new,
601  const Vector<Real> &u_old,
602  const Vector<Real> &z,
603  const Vector<Real> &v,
604  const Vector<Real> &jv,
605  const bool printToStream = true,
606  std::ostream & outStream = std::cout,
607  const int numSteps = ROL_NUM_CHECKDERIV_STEPS,
608  const int order = 1) {
609  std::vector<Real> steps(numSteps);
610  for(int i=0;i<numSteps;++i) {
611  steps[i] = pow(10,-i);
612  }
613 
614  return checkApplyJacobian_1_new(u_new,u_old,z,v,jv,steps,printToStream,outStream,order);
615  }
616 
617  std::vector<std::vector<Real> > checkApplyJacobian_1_new(const Vector<Real> &u_new,
618  const Vector<Real> &u_old,
619  const Vector<Real> &z,
620  const Vector<Real> &v,
621  const Vector<Real> &jv,
622  const std::vector<Real> &steps,
623  const bool printToStream = true,
624  std::ostream & outStream = std::cout,
625  const int order = 1) {
626 
627  ROL_TEST_FOR_EXCEPTION( order<1 || order>4, std::invalid_argument,
628  "Error: finite difference order must be 1,2,3, or 4" );
629 
630  Real one(1.0);
631 
634 
635  Real tol = std::sqrt(ROL_EPSILON<Real>());
636 
637  int numSteps = steps.size();
638  int numVals = 4;
639  std::vector<Real> tmp(numVals);
640  std::vector<std::vector<Real> > jvCheck(numSteps, tmp);
641 
642  // Save the format state of the original outStream.
643  ROL::nullstream oldFormatState;
644  oldFormatState.copyfmt(outStream);
645 
646  // Compute constraint value at x.
647  ROL::Ptr<Vector<Real> > c = workspace_.clone(jv);
648  this->update(u_new, u_old, z);
649  this->value(*c, u_new, u_old, z, tol);
650 
651  // Compute (Jacobian at x) times (vector v).
652  ROL::Ptr<Vector<Real> > Jv = workspace_.clone(jv);
653  this->applyJacobian_1_new(*Jv, v, u_new, u_old, z, tol);
654  Real normJv = Jv->norm();
655 
656  // Temporary vectors.
657  ROL::Ptr<Vector<Real> > cdif = workspace_.clone(jv);
658  ROL::Ptr<Vector<Real> > cnew = workspace_.clone(jv);
659  ROL::Ptr<Vector<Real> > u_2 = workspace_.clone(u_new);
660 
661  for (int i=0; i<numSteps; i++) {
662 
663  Real eta = steps[i];
664 
665  u_2->set(u_new);
666 
667  cdif->set(*c);
668  cdif->scale(weights[order-1][0]);
669 
670  for(int j=0; j<order; ++j) {
671 
672  u_2->axpy(eta*shifts[order-1][j], v);
673 
674  if( weights[order-1][j+1] != 0 ) {
675  this->update(*u_2,u_old,z);
676  this->value(*cnew,*u_2,u_old,z,tol);
677  cdif->axpy(weights[order-1][j+1],*cnew);
678  }
679 
680  }
681 
682  cdif->scale(one/eta);
683 
684  // Compute norms of Jacobian-vector products, finite-difference approximations, and error.
685  jvCheck[i][0] = eta;
686  jvCheck[i][1] = normJv;
687  jvCheck[i][2] = cdif->norm();
688  cdif->axpy(-one, *Jv);
689  jvCheck[i][3] = cdif->norm();
690 
691  if (printToStream) {
692  std::stringstream hist;
693  if (i==0) {
694  hist << std::right
695  << std::setw(20) << "Step size"
696  << std::setw(20) << "norm(Jac*vec)"
697  << std::setw(20) << "norm(FD approx)"
698  << std::setw(20) << "norm(abs error)"
699  << "\n"
700  << std::setw(20) << "---------"
701  << std::setw(20) << "-------------"
702  << std::setw(20) << "---------------"
703  << std::setw(20) << "---------------"
704  << "\n";
705  }
706  hist << std::scientific << std::setprecision(11) << std::right
707  << std::setw(20) << jvCheck[i][0]
708  << std::setw(20) << jvCheck[i][1]
709  << std::setw(20) << jvCheck[i][2]
710  << std::setw(20) << jvCheck[i][3]
711  << "\n";
712  outStream << hist.str();
713  }
714 
715  }
716 
717  // Reset format state of outStream.
718  outStream.copyfmt(oldFormatState);
719 
720  return jvCheck;
721  } // checkApplyJacobian_1_new
722 
723 
724 }; // class Constraint_SimOpt
725 
726 } // namespace ROL
727 
728 #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:866
virtual void axpy(const Real alpha, const Vector &x)
Compute where .
Definition: ROL_Vector.hpp:153
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:167
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:80
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)
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:74
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
basic_nullstream< char, char_traits< char >> nullstream
Definition: ROL_Stream.hpp:72
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