ROL
burgers-control/example_02.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 
49 #include "ROL_Algorithm.hpp"
50 #include "ROL_Types.hpp"
51 #include "ROL_ParameterList.hpp"
52 #include "ROL_Stream.hpp"
53 #include "Teuchos_GlobalMPISession.hpp"
54 #include "Teuchos_LAPACK.hpp"
55 
56 #include <iostream>
57 #include <algorithm>
58 
59 #include "ROL_StdVector.hpp"
60 #include "ROL_Vector_SimOpt.hpp"
62 #include "ROL_Objective_SimOpt.hpp"
64 
65 template<class Real>
67 private:
68  int nx_;
69  Real dx_;
70  Real nu_;
71  Real u0_;
72  Real u1_;
73  Real f_;
74 
75 private:
76  Real compute_norm(const std::vector<Real> &r) {
77  return std::sqrt(this->dot(r,r));
78  }
79 
80  Real dot(const std::vector<Real> &x, const std::vector<Real> &y) {
81  Real ip = 0.0;
82  Real c = (((int)x.size()==this->nx_) ? 4.0 : 2.0);
83  for (unsigned i=0; i<x.size(); i++) {
84  if ( i == 0 ) {
85  ip += this->dx_/6.0*(c*x[i] + x[i+1])*y[i];
86  }
87  else if ( i == x.size()-1 ) {
88  ip += this->dx_/6.0*(x[i-1] + c*x[i])*y[i];
89  }
90  else {
91  ip += this->dx_/6.0*(x[i-1] + 4.0*x[i] + x[i+1])*y[i];
92  }
93  }
94  return ip;
95  }
96 
98 
99  void update(std::vector<Real> &u, const std::vector<Real> &s, const Real alpha=1.0) {
100  for (unsigned i=0; i<u.size(); i++) {
101  u[i] += alpha*s[i];
102  }
103  }
104 
105  void scale(std::vector<Real> &u, const Real alpha=0.0) {
106  for (unsigned i=0; i<u.size(); i++) {
107  u[i] *= alpha;
108  }
109  }
110 
111  void compute_residual(std::vector<Real> &r, const std::vector<Real> &u,
112  const std::vector<Real> &z) {
113  r.clear();
114  r.resize(this->nx_,0.0);
115  for (int i=0; i<this->nx_; i++) {
116  // Contribution from stiffness term
117  if ( i==0 ) {
118  r[i] = this->nu_/this->dx_*(2.0*u[i]-u[i+1]);
119  }
120  else if (i==this->nx_-1) {
121  r[i] = this->nu_/this->dx_*(2.0*u[i]-u[i-1]);
122  }
123  else {
124  r[i] = this->nu_/this->dx_*(2.0*u[i]-u[i-1]-u[i+1]);
125  }
126  // Contribution from nonlinear term
127  if (i<this->nx_-1){
128  r[i] += u[i+1]*(u[i]+u[i+1])/6.0;
129  }
130  if (i>0) {
131  r[i] -= u[i-1]*(u[i-1]+u[i])/6.0;
132  }
133  // Contribution from control
134  r[i] -= this->dx_/6.0*(z[i]+4.0*z[i+1]+z[i+2]);
135  // Contribution from right-hand side
136  r[i] -= this->dx_*this->f_;
137  }
138  // Contribution from Dirichlet boundary terms
139  r[0] -= this->u0_*u[ 0]/6.0 + this->u0_*this->u0_/6.0 + this->nu_*this->u0_/this->dx_;
140  r[this->nx_-1] += this->u1_*u[this->nx_-1]/6.0 + this->u1_*this->u1_/6.0 - this->nu_*this->u1_/this->dx_;
141  }
142 
143  void compute_pde_jacobian(std::vector<Real> &dl, std::vector<Real> &d, std::vector<Real> &du,
144  const std::vector<Real> &u) {
145  // Get Diagonal and Off-Diagonal Entries of linear PDE Jacobian
146  d.clear();
147  d.resize(this->nx_,this->nu_*2.0/this->dx_);
148  dl.clear();
149  dl.resize(this->nx_-1,-this->nu_/this->dx_);
150  du.clear();
151  du.resize(this->nx_-1,-this->nu_/this->dx_);
152  // Contribution from nonlinearity
153  for (int i=0; i<this->nx_; i++) {
154  if (i<this->nx_-1) {
155  dl[i] += (-2.0*u[i]-u[i+1])/6.0;
156  d[i] += u[i+1]/6.0;
157  }
158  if (i>0) {
159  d[i] += -u[i-1]/6.0;
160  du[i-1] += (u[i-1]+2.0*u[i])/6.0;
161  }
162  }
163  // Contribution from Dirichlet boundary conditions
164  d[0] -= this->u0_/6.0;
165  d[this->nx_-1] += this->u1_/6.0;
166  }
167 
168  void linear_solve(std::vector<Real> &u, std::vector<Real> &dl, std::vector<Real> &d, std::vector<Real> &du,
169  const std::vector<Real> &r, const bool transpose = false) {
170  u.assign(r.begin(),r.end());
171  // Perform LDL factorization
172  Teuchos::LAPACK<int,Real> lp;
173  std::vector<Real> du2(this->nx_-2,0.0);
174  std::vector<int> ipiv(this->nx_,0);
175  int info;
176  int ldb = this->nx_;
177  int nhrs = 1;
178  lp.GTTRF(this->nx_,&dl[0],&d[0],&du[0],&du2[0],&ipiv[0],&info);
179  char trans = 'N';
180  if ( transpose ) {
181  trans = 'T';
182  }
183  lp.GTTRS(trans,this->nx_,nhrs,&dl[0],&d[0],&du[0],&du2[0],&ipiv[0],&u[0],ldb,&info);
184  }
185 
186 public:
187 
188  Constraint_BurgersControl(int nx = 128, Real nu = 1.e-2, Real u0 = 1.0, Real u1 = 0.0, Real f = 0.0)
189  : nx_(nx), nu_(nu), u0_(u0), u1_(u1), f_(f) {
190  dx_ = 1.0/((Real)nx+1.0);
191  }
192 
194  const ROL::Vector<Real> &z, Real &tol) {
195  ROL::Ptr<std::vector<Real> > cp =
196  dynamic_cast<ROL::StdVector<Real>&>(c).getVector();
197  ROL::Ptr<const std::vector<Real> > up =
198  dynamic_cast<const ROL::StdVector<Real>&>(u).getVector();
199  ROL::Ptr<const std::vector<Real> > zp =
200  dynamic_cast<const ROL::StdVector<Real>&>(z).getVector();
201  this->compute_residual(*cp,*up,*zp);
202  }
203 
205  const ROL::Vector<Real> &z, Real &tol) {
206  ROL::Ptr<std::vector<Real> > jvp =
207  dynamic_cast<ROL::StdVector<Real>&>(jv).getVector();
208  ROL::Ptr<const std::vector<Real> > vp =
209  dynamic_cast<const ROL::StdVector<Real>&>(v).getVector();
210  ROL::Ptr<const std::vector<Real> > up =
211  dynamic_cast<const ROL::StdVector<Real>&>(u).getVector();
212  ROL::Ptr<const std::vector<Real> > zp =
213  dynamic_cast<const ROL::StdVector<Real>&>(z).getVector();
214  // Fill jvp
215  for (int i = 0; i < this->nx_; i++) {
216  (*jvp)[i] = this->nu_/this->dx_*2.0*(*vp)[i];
217  if ( i > 0 ) {
218  (*jvp)[i] += -this->nu_/this->dx_*(*vp)[i-1]
219  -(*up)[i-1]/6.0*(*vp)[i]
220  -((*up)[i]+2.0*(*up)[i-1])/6.0*(*vp)[i-1];
221  }
222  if ( i < this->nx_-1 ) {
223  (*jvp)[i] += -this->nu_/this->dx_*(*vp)[i+1]
224  +(*up)[i+1]/6.0*(*vp)[i]
225  +((*up)[i]+2.0*(*up)[i+1])/6.0*(*vp)[i+1];
226  }
227  }
228  (*jvp)[0] -= this->u0_/6.0*(*vp)[0];
229  (*jvp)[this->nx_-1] += this->u1_/6.0*(*vp)[this->nx_-1];
230  }
231 
233  const ROL::Vector<Real> &z, Real &tol) {
234  ROL::Ptr<std::vector<Real> > jvp =
235  dynamic_cast<ROL::StdVector<Real>&>(jv).getVector();
236  ROL::Ptr<const std::vector<Real> > vp =
237  dynamic_cast<const ROL::StdVector<Real>&>(v).getVector();
238  ROL::Ptr<const std::vector<Real> > up =
239  dynamic_cast<const ROL::StdVector<Real>&>(u).getVector();
240  ROL::Ptr<const std::vector<Real> > zp =
241  dynamic_cast<const ROL::StdVector<Real>&>(z).getVector();
242  for (int i=0; i<this->nx_; i++) {
243  // Contribution from control
244  (*jvp)[i] = -this->dx_/6.0*((*vp)[i]+4.0*(*vp)[i+1]+(*vp)[i+2]);
245  }
246  }
247 
249  const ROL::Vector<Real> &z, Real &tol) {
250  ROL::Ptr<std::vector<Real> > ijvp =
251  dynamic_cast<ROL::StdVector<Real>&>(ijv).getVector();
252  ROL::Ptr<const std::vector<Real> > vp =
253  dynamic_cast<const ROL::StdVector<Real>&>(v).getVector();
254  ROL::Ptr<const std::vector<Real> > up =
255  dynamic_cast<const ROL::StdVector<Real>&>(u).getVector();
256  ROL::Ptr<const std::vector<Real> > zp =
257  dynamic_cast<const ROL::StdVector<Real>&>(z).getVector();
258  // Get PDE Jacobian
259  std::vector<Real> d(this->nx_,0.0);
260  std::vector<Real> dl(this->nx_-1,0.0);
261  std::vector<Real> du(this->nx_-1,0.0);
262  this->compute_pde_jacobian(dl,d,du,*up);
263  // Solve solve state sensitivity system at current time step
264  this->linear_solve(*ijvp,dl,d,du,*vp);
265  }
266 
268  const ROL::Vector<Real> &z, Real &tol) {
269  ROL::Ptr<std::vector<Real> > jvp =
270  dynamic_cast<ROL::StdVector<Real>&>(ajv).getVector();
271  ROL::Ptr<const std::vector<Real> > vp =
272  dynamic_cast<const ROL::StdVector<Real>&>(v).getVector();
273  ROL::Ptr<const std::vector<Real> > up =
274  dynamic_cast<const ROL::StdVector<Real>&>(u).getVector();
275  ROL::Ptr<const std::vector<Real> > zp =
276  dynamic_cast<const ROL::StdVector<Real>&>(z).getVector();
277  // Fill jvp
278  for (int i = 0; i < this->nx_; i++) {
279  (*jvp)[i] = this->nu_/this->dx_*2.0*(*vp)[i];
280  if ( i > 0 ) {
281  (*jvp)[i] += -this->nu_/this->dx_*(*vp)[i-1]
282  -(*up)[i-1]/6.0*(*vp)[i]
283  +((*up)[i-1]+2.0*(*up)[i])/6.0*(*vp)[i-1];
284  }
285  if ( i < this->nx_-1 ) {
286  (*jvp)[i] += -this->nu_/this->dx_*(*vp)[i+1]
287  +(*up)[i+1]/6.0*(*vp)[i]
288  -((*up)[i+1]+2.0*(*up)[i])/6.0*(*vp)[i+1];
289  }
290  }
291  (*jvp)[0] -= this->u0_/6.0*(*vp)[0];
292  (*jvp)[this->nx_-1] += this->u1_/6.0*(*vp)[this->nx_-1];
293  }
294 
296  const ROL::Vector<Real> &z, Real &tol) {
297  ROL::Ptr<std::vector<Real> > jvp =
298  dynamic_cast<ROL::StdVector<Real>&>(jv).getVector();
299  ROL::Ptr<const std::vector<Real> > vp =
300  dynamic_cast<const ROL::StdVector<Real>&>(v).getVector();
301  ROL::Ptr<const std::vector<Real> > up =
302  dynamic_cast<const ROL::StdVector<Real>&>(u).getVector();
303  ROL::Ptr<const std::vector<Real> > zp =
304  dynamic_cast<const ROL::StdVector<Real>&>(z).getVector();
305  for (int i=0; i<this->nx_+2; i++) {
306  if ( i == 0 ) {
307  (*jvp)[i] = -this->dx_/6.0*(*vp)[i];
308  }
309  else if ( i == 1 ) {
310  (*jvp)[i] = -this->dx_/6.0*(4.0*(*vp)[i-1]+(*vp)[i]);
311  }
312  else if ( i == this->nx_ ) {
313  (*jvp)[i] = -this->dx_/6.0*(4.0*(*vp)[i-1]+(*vp)[i-2]);
314  }
315  else if ( i == this->nx_+1 ) {
316  (*jvp)[i] = -this->dx_/6.0*(*vp)[i-2];
317  }
318  else {
319  (*jvp)[i] = -this->dx_/6.0*((*vp)[i-2]+4.0*(*vp)[i-1]+(*vp)[i]);
320  }
321  }
322  }
323 
325  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
326  ROL::Ptr<std::vector<Real> > iajvp =
327  dynamic_cast<ROL::StdVector<Real>&>(iajv).getVector();
328  ROL::Ptr<const std::vector<Real> > vp =
329  dynamic_cast<const ROL::StdVector<Real>&>(v).getVector();
330  ROL::Ptr<const std::vector<Real> > up =
331  dynamic_cast<const ROL::StdVector<Real>&>(u).getVector();
332  // Get PDE Jacobian
333  std::vector<Real> d(this->nx_,0.0);
334  std::vector<Real> du(this->nx_-1,0.0);
335  std::vector<Real> dl(this->nx_-1,0.0);
336  this->compute_pde_jacobian(dl,d,du,*up);
337  // Solve solve adjoint system at current time step
338  this->linear_solve(*iajvp,dl,d,du,*vp,true);
339  }
340 
342  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
343  ROL::Ptr<std::vector<Real> > ahwvp =
344  dynamic_cast<ROL::StdVector<Real>&>(ahwv).getVector();
345  ROL::Ptr<const std::vector<Real> > wp =
346  dynamic_cast<const ROL::StdVector<Real>&>(w).getVector();
347  ROL::Ptr<const std::vector<Real> > vp =
348  dynamic_cast<const ROL::StdVector<Real>&>(v).getVector();
349  ROL::Ptr<const std::vector<Real> > up =
350  dynamic_cast<const ROL::StdVector<Real>&>(u).getVector();
351  ROL::Ptr<const std::vector<Real> > zp =
352  dynamic_cast<const ROL::StdVector<Real>&>(z).getVector();
353  for (int i=0; i<this->nx_; i++) {
354  // Contribution from nonlinear term
355  (*ahwvp)[i] = 0.0;
356  if (i<this->nx_-1){
357  (*ahwvp)[i] += ((*wp)[i]*(*vp)[i+1] - (*wp)[i+1]*(2.0*(*vp)[i]+(*vp)[i+1]))/6.0;
358  }
359  if (i>0) {
360  (*ahwvp)[i] += ((*wp)[i-1]*((*vp)[i-1]+2.0*(*vp)[i]) - (*wp)[i]*(*vp)[i-1])/6.0;
361  }
362  }
363  }
364 
366  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
367  ahwv.zero();
368  }
370  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
371  ahwv.zero();
372  }
374  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
375  ahwv.zero();
376  }
377 
378 // void solveAugmentedSystem(ROL::Vector<Real> &v1, ROL::Vector<Real> &v2, const ROL::Vector<Real> &b1,
379 // const ROL::Vector<Real> &b2, const ROL::Vector<Real> &x, Real &tol) {}
380 };
381 
382 template<class Real>
383 class Objective_BurgersControl : public ROL::Objective_SimOpt<Real> {
384 private:
385  Real alpha_; // Penalty Parameter
386 
387  int nx_; // Number of interior nodes
388  Real dx_; // Mesh spacing (i.e. 1/(nx+1))
389 
390 /***************************************************************/
391 /********** BEGIN PRIVATE MEMBER FUNCTION DECLARATION **********/
392 /***************************************************************/
393  Real evaluate_target(Real x) {
394  Real val = 0.0;
395  int example = 2;
396  switch (example) {
397  case 1: val = ((x<0.5) ? 1.0 : 0.0); break;
398  case 2: val = 1.0; break;
399  case 3: val = std::abs(std::sin(8.0*M_PI*x)); break;
400  case 4: val = std::exp(-0.5*(x-0.5)*(x-0.5)); break;
401  }
402  return val;
403  }
404 
405  Real dot(const std::vector<Real> &x, const std::vector<Real> &y) {
406  Real ip = 0.0;
407  Real c = (((int)x.size()==this->nx_) ? 4.0 : 2.0);
408  for (unsigned i=0; i<x.size(); i++) {
409  if ( i == 0 ) {
410  ip += this->dx_/6.0*(c*x[i] + x[i+1])*y[i];
411  }
412  else if ( i == x.size()-1 ) {
413  ip += this->dx_/6.0*(x[i-1] + c*x[i])*y[i];
414  }
415  else {
416  ip += this->dx_/6.0*(x[i-1] + 4.0*x[i] + x[i+1])*y[i];
417  }
418  }
419  return ip;
420  }
421 
422  void apply_mass(std::vector<Real> &Mu, const std::vector<Real> &u ) {
423  Mu.resize(u.size(),0.0);
424  Real c = (((int)u.size()==this->nx_) ? 4.0 : 2.0);
425  for (unsigned i=0; i<u.size(); i++) {
426  if ( i == 0 ) {
427  Mu[i] = this->dx_/6.0*(c*u[i] + u[i+1]);
428  }
429  else if ( i == u.size()-1 ) {
430  Mu[i] = this->dx_/6.0*(u[i-1] + c*u[i]);
431  }
432  else {
433  Mu[i] = this->dx_/6.0*(u[i-1] + 4.0*u[i] + u[i+1]);
434  }
435  }
436  }
437 /*************************************************************/
438 /********** END PRIVATE MEMBER FUNCTION DECLARATION **********/
439 /*************************************************************/
440 
441 public:
442 
443  Objective_BurgersControl(Real alpha = 1.e-4, int nx = 128) : alpha_(alpha), nx_(nx) {
444  dx_ = 1.0/((Real)nx+1.0);
445  }
446 
447  Real value( const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
448  ROL::Ptr<const std::vector<Real> > up =
449  dynamic_cast<const ROL::StdVector<Real>&>(u).getVector();
450  ROL::Ptr<const std::vector<Real> > zp =
451  dynamic_cast<const ROL::StdVector<Real>&>(z).getVector();
452  // COMPUTE RESIDUAL
453  Real res1 = 0.0, res2 = 0.0, res3 = 0.0;
454  Real valu = 0.0, valz = this->dot(*zp,*zp);
455  for (int i=0; i<this->nx_; i++) {
456  if ( i == 0 ) {
457  res1 = (*up)[i]-evaluate_target((Real)(i+1)*this->dx_);
458  res2 = (*up)[i+1]-evaluate_target((Real)(i+2)*this->dx_);
459  valu += this->dx_/6.0*(4.0*res1 + res2)*res1;
460  }
461  else if ( i == this->nx_-1 ) {
462  res1 = (*up)[i-1]-evaluate_target((Real)i*this->dx_);
463  res2 = (*up)[i]-evaluate_target((Real)(i+1)*this->dx_);
464  valu += this->dx_/6.0*(res1 + 4.0*res2)*res2;
465  }
466  else {
467  res1 = (*up)[i-1]-evaluate_target((Real)i*this->dx_);
468  res2 = (*up)[i]-evaluate_target((Real)(i+1)*this->dx_);
469  res3 = (*up)[i+1]-evaluate_target((Real)(i+2)*this->dx_);
470  valu += this->dx_/6.0*(res1 + 4.0*res2 + res3)*res2;
471  }
472  }
473  return 0.5*(valu + this->alpha_*valz);
474  }
475 
476  void gradient_1( ROL::Vector<Real> &g, const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
477  // Unwrap g
478  ROL::Ptr<std::vector<Real> > gup = ROL::constPtrCast<std::vector<Real> >(
479  (dynamic_cast<const ROL::StdVector<Real>&>(g)).getVector());
480  // Unwrap x
481  ROL::Ptr<const std::vector<Real> > up =
482  dynamic_cast<const ROL::StdVector<Real>&>(u).getVector();
483  ROL::Ptr<const std::vector<Real> > zp =
484  dynamic_cast<const ROL::StdVector<Real>&>(z).getVector();
485  // COMPUTE GRADIENT WRT U
486  std::vector<Real> diff(this->nx_,0.0);
487  for (int i=0; i<this->nx_; i++) {
488  diff[i] = ((*up)[i]-this->evaluate_target((Real)(i+1)*this->dx_));
489  }
490  this->apply_mass(*gup,diff);
491  }
492 
493  void gradient_2( ROL::Vector<Real> &g, const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
494  // Unwrap g
495  ROL::Ptr<std::vector<Real> > gzp = ROL::constPtrCast<std::vector<Real> >(
496  (dynamic_cast<const ROL::StdVector<Real>&>(g)).getVector());
497  // Unwrap x
498  ROL::Ptr<const std::vector<Real> > up =
499  dynamic_cast<const ROL::StdVector<Real>&>(u).getVector();
500  ROL::Ptr<const std::vector<Real> > zp =
501  dynamic_cast<const ROL::StdVector<Real>&>(z).getVector();
502  // COMPUTE GRADIENT WRT Z
503  for (int i=0; i<this->nx_+2; i++) {
504  if (i==0) {
505  (*gzp)[i] = this->alpha_*this->dx_/6.0*(2.0*(*zp)[i]+(*zp)[i+1]);
506  }
507  else if (i==this->nx_+1) {
508  (*gzp)[i] = this->alpha_*this->dx_/6.0*(2.0*(*zp)[i]+(*zp)[i-1]);
509  }
510  else {
511  (*gzp)[i] = this->alpha_*this->dx_/6.0*((*zp)[i-1]+4.0*(*zp)[i]+(*zp)[i+1]);
512  }
513  }
514  }
515 
517  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
518  ROL::Ptr<std::vector<Real> > hvup =
519  dynamic_cast<ROL::StdVector<Real>&>(hv).getVector();
520  // Unwrap v
521  ROL::Ptr<const std::vector<Real> > vup =
522  dynamic_cast<const ROL::StdVector<Real>&>(v).getVector();
523  // COMPUTE GRADIENT WRT U
524  this->apply_mass(*hvup,*vup);
525  }
526 
528  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
529  hv.zero();
530  }
531 
533  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
534  hv.zero();
535  }
536 
538  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
539  ROL::Ptr<std::vector<Real> > hvzp =
540  dynamic_cast<ROL::StdVector<Real>&>(hv).getVector();
541  // Unwrap v
542  ROL::Ptr<const std::vector<Real> > vzp =
543  dynamic_cast<const ROL::StdVector<Real>&>(v).getVector();
544  // COMPUTE GRADIENT WRT Z
545  for (int i=0; i<this->nx_+2; i++) {
546  if (i==0) {
547  (*hvzp)[i] = this->alpha_*this->dx_/6.0*(2.0*(*vzp)[i]+(*vzp)[i+1]);
548  }
549  else if (i==this->nx_+1) {
550  (*hvzp)[i] = this->alpha_*this->dx_/6.0*(2.0*(*vzp)[i]+(*vzp)[i-1]);
551  }
552  else {
553  (*hvzp)[i] = this->alpha_*this->dx_/6.0*((*vzp)[i-1]+4.0*(*vzp)[i]+(*vzp)[i+1]);
554  }
555  }
556  }
557 };
Provides the interface to evaluate simulation-based objective functions.
void applyAdjointJacobian_1(ROL::Vector< Real > &ajv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply the adjoint of the partial constraint Jacobian at , , to the vector . This is the primary inter...
void applyAdjointHessian_12(ROL::Vector< Real > &ahwv, const ROL::Vector< Real > &w, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply the optimization-space derivative of the adjoint of the constraint simulation-space Jacobian at...
Contains definitions of custom data types in ROL.
void applyAdjointHessian_11(ROL::Vector< Real > &ahwv, const ROL::Vector< Real > &w, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply the simulation-space derivative of the adjoint of the constraint simulation-space Jacobian at ...
void applyInverseAdjointJacobian_1(ROL::Vector< Real > &iajv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply the inverse of the adjoint of the partial constraint Jacobian at , , to the vector ...
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
Defines a no-output stream class ROL::NullStream and a function makeStreamPtr which either wraps a re...
void hessVec_22(ROL::Vector< Real > &hv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Real value(const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Compute value.
void gradient_1(ROL::Vector< Real > &g, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Compute gradient with respect to first component.
void applyAdjointHessian_21(ROL::Vector< Real > &ahwv, const ROL::Vector< Real > &w, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply the simulation-space derivative of the adjoint of the constraint optimization-space Jacobian at...
void compute_pde_jacobian(std::vector< Real > &dl, std::vector< Real > &d, std::vector< Real > &du, const std::vector< Real > &u)
Real dot(const std::vector< Real > &x, const std::vector< Real > &y)
Objective_BurgersControl(Real alpha=1.e-4, int nx=128)
void applyJacobian_1(ROL::Vector< Real > &jv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply the partial constraint Jacobian at , , to the vector .
void hessVec_21(ROL::Vector< Real > &hv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
void hessVec_12(ROL::Vector< Real > &hv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Constraint_BurgersControl(int nx=128, Real nu=1.e-2, Real u0=1.0, Real u1=0.0, Real f=0.0)
void gradient_2(ROL::Vector< Real > &g, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Compute gradient with respect to second component.
Real dot(const std::vector< Real > &x, const std::vector< Real > &y)
void applyInverseJacobian_1(ROL::Vector< Real > &ijv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply the inverse partial constraint Jacobian at , , to the vector .
void compute_residual(std::vector< Real > &r, const std::vector< Real > &u, const std::vector< Real > &z)
void scale(std::vector< Real > &u, const Real alpha=0.0)
void hessVec_11(ROL::Vector< Real > &hv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply Hessian approximation to vector.
void applyJacobian_2(ROL::Vector< Real > &jv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply the partial constraint Jacobian at , , to the vector .
void apply_mass(std::vector< Real > &Mu, const std::vector< Real > &u)
void applyAdjointHessian_22(ROL::Vector< Real > &ahwv, const ROL::Vector< Real > &w, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply the optimization-space derivative of the adjoint of the constraint optimization-space Jacobian ...
Real compute_norm(const std::vector< Real > &r)
void update(std::vector< Real > &u, const std::vector< Real > &s, const Real alpha=1.0)
void applyAdjointJacobian_2(ROL::Vector< Real > &jv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Apply the adjoint of the partial constraint Jacobian at , , to vector . This is the primary interface...
Defines the constraint operator interface for simulation-based optimization.
void value(ROL::Vector< Real > &c, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Evaluate the constraint operator at .
void linear_solve(std::vector< Real > &u, std::vector< Real > &dl, std::vector< Real > &d, std::vector< Real > &du, const std::vector< Real > &r, const bool transpose=false)