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 "Teuchos_oblackholestream.hpp"
52 #include "Teuchos_GlobalMPISession.hpp"
53 #include "Teuchos_XMLParameterListHelpers.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 
97  void update(std::vector<Real> &u, const std::vector<Real> &s, const Real alpha=1.0) {
98  for (unsigned i=0; i<u.size(); i++) {
99  u[i] += alpha*s[i];
100  }
101  }
102 
103  void scale(std::vector<Real> &u, const Real alpha=0.0) {
104  for (unsigned i=0; i<u.size(); i++) {
105  u[i] *= alpha;
106  }
107  }
108 
109  void compute_residual(std::vector<Real> &r, const std::vector<Real> &u,
110  const std::vector<Real> &z) {
111  r.clear();
112  r.resize(this->nx_,0.0);
113  for (int i=0; i<this->nx_; i++) {
114  // Contribution from stiffness term
115  if ( i==0 ) {
116  r[i] = this->nu_/this->dx_*(2.0*u[i]-u[i+1]);
117  }
118  else if (i==this->nx_-1) {
119  r[i] = this->nu_/this->dx_*(2.0*u[i]-u[i-1]);
120  }
121  else {
122  r[i] = this->nu_/this->dx_*(2.0*u[i]-u[i-1]-u[i+1]);
123  }
124  // Contribution from nonlinear term
125  if (i<this->nx_-1){
126  r[i] += u[i+1]*(u[i]+u[i+1])/6.0;
127  }
128  if (i>0) {
129  r[i] -= u[i-1]*(u[i-1]+u[i])/6.0;
130  }
131  // Contribution from control
132  r[i] -= this->dx_/6.0*(z[i]+4.0*z[i+1]+z[i+2]);
133  // Contribution from right-hand side
134  r[i] -= this->dx_*this->f_;
135  }
136  // Contribution from Dirichlet boundary terms
137  r[0] -= this->u0_*u[ 0]/6.0 + this->u0_*this->u0_/6.0 + this->nu_*this->u0_/this->dx_;
138  r[this->nx_-1] += this->u1_*u[this->nx_-1]/6.0 + this->u1_*this->u1_/6.0 - this->nu_*this->u1_/this->dx_;
139  }
140 
141  void compute_pde_jacobian(std::vector<Real> &dl, std::vector<Real> &d, std::vector<Real> &du,
142  const std::vector<Real> &u) {
143  // Get Diagonal and Off-Diagonal Entries of linear PDE Jacobian
144  d.clear();
145  d.resize(this->nx_,this->nu_*2.0/this->dx_);
146  dl.clear();
147  dl.resize(this->nx_-1,-this->nu_/this->dx_);
148  du.clear();
149  du.resize(this->nx_-1,-this->nu_/this->dx_);
150  // Contribution from nonlinearity
151  for (int i=0; i<this->nx_; i++) {
152  if (i<this->nx_-1) {
153  dl[i] += (-2.0*u[i]-u[i+1])/6.0;
154  d[i] += u[i+1]/6.0;
155  }
156  if (i>0) {
157  d[i] += -u[i-1]/6.0;
158  du[i-1] += (u[i-1]+2.0*u[i])/6.0;
159  }
160  }
161  // Contribution from Dirichlet boundary conditions
162  d[0] -= this->u0_/6.0;
163  d[this->nx_-1] += this->u1_/6.0;
164  }
165 
166  void linear_solve(std::vector<Real> &u, std::vector<Real> &dl, std::vector<Real> &d, std::vector<Real> &du,
167  const std::vector<Real> &r, const bool transpose = false) {
168  u.assign(r.begin(),r.end());
169  // Perform LDL factorization
170  Teuchos::LAPACK<int,Real> lp;
171  std::vector<Real> du2(this->nx_-2,0.0);
172  std::vector<int> ipiv(this->nx_,0);
173  int info;
174  int ldb = this->nx_;
175  int nhrs = 1;
176  lp.GTTRF(this->nx_,&dl[0],&d[0],&du[0],&du2[0],&ipiv[0],&info);
177  char trans = 'N';
178  if ( transpose ) {
179  trans = 'T';
180  }
181  lp.GTTRS(trans,this->nx_,nhrs,&dl[0],&d[0],&du[0],&du2[0],&ipiv[0],&u[0],ldb,&info);
182  }
183 
184 public:
185 
186  EqualityConstraint_BurgersControl(int nx = 128, Real nu = 1.e-2, Real u0 = 1.0, Real u1 = 0.0, Real f = 0.0)
187  : nx_(nx), nu_(nu), u0_(u0), u1_(u1), f_(f) {
188  dx_ = 1.0/((Real)nx+1.0);
189  }
190 
192  const ROL::Vector<Real> &z, Real &tol) {
193  Teuchos::RCP<std::vector<Real> > cp =
194  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(c)).getVector());
195  Teuchos::RCP<const std::vector<Real> > up =
196  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(u))).getVector();
197  Teuchos::RCP<const std::vector<Real> > zp =
198  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(z))).getVector();
199  this->compute_residual(*cp,*up,*zp);
200  }
201 
202  void solve(ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
203  Teuchos::RCP<std::vector<Real> > up =
204  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(u)).getVector());
205  up->assign(up->size(),z.norm()/up->size());
206  Teuchos::RCP<const std::vector<Real> > zp =
207  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(z))).getVector();
208  // Compute residual and residual norm
209  std::vector<Real> r(up->size(),0.0);
210  this->compute_residual(r,*up,*zp);
211  Real rnorm = this->compute_norm(r);
212  // Define tolerances
213  Real rtol = 1.e2*ROL::ROL_EPSILON;
214  Real maxit = 500;
215  // Initialize Jacobian storage
216  std::vector<Real> d(this->nx_,0.0);
217  std::vector<Real> dl(this->nx_-1,0.0);
218  std::vector<Real> du(this->nx_-1,0.0);
219  // Iterate Newton's method
220  Real alpha = 1.0, tmp = 0.0;
221  std::vector<Real> s(this->nx_,0.0);
222  std::vector<Real> utmp(this->nx_,0.0);
223  for (int i=0; i<maxit; i++) {
224  //std::cout << i << " " << rnorm << "\n";
225  // Get Jacobian
226  this->compute_pde_jacobian(dl,d,du,*up);
227  // Solve Newton system
228  this->linear_solve(s,dl,d,du,r);
229  // Perform line search
230  tmp = rnorm;
231  alpha = 1.0;
232  utmp.assign(up->begin(),up->end());
233  this->update(utmp,s,-alpha);
234  this->compute_residual(r,utmp,*zp);
235  rnorm = this->compute_norm(r);
236  while ( rnorm > (1.0-1.e-4*alpha)*tmp && alpha > std::sqrt(ROL::ROL_EPSILON) ) {
237  alpha /= 2.0;
238  utmp.assign(up->begin(),up->end());
239  this->update(utmp,s,-alpha);
240  this->compute_residual(r,utmp,*zp);
241  rnorm = this->compute_norm(r);
242  }
243  // Update iterate
244  up->assign(utmp.begin(),utmp.end());
245  if ( rnorm < rtol ) {
246  break;
247  }
248  }
249  }
250 
251 
252 
254  const ROL::Vector<Real> &z, Real &tol) {
255  Teuchos::RCP<std::vector<Real> > jvp =
256  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(jv)).getVector());
257  Teuchos::RCP<const std::vector<Real> > vp =
258  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(v))).getVector();
259  Teuchos::RCP<const std::vector<Real> > up =
260  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(u))).getVector();
261  Teuchos::RCP<const std::vector<Real> > zp =
262  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(z))).getVector();
263  // Fill jvp
264  for (int i = 0; i < this->nx_; i++) {
265  (*jvp)[i] = this->nu_/this->dx_*2.0*(*vp)[i];
266  if ( i > 0 ) {
267  (*jvp)[i] += -this->nu_/this->dx_*(*vp)[i-1]
268  -(*up)[i-1]/6.0*(*vp)[i]
269  -((*up)[i]+2.0*(*up)[i-1])/6.0*(*vp)[i-1];
270  }
271  if ( i < this->nx_-1 ) {
272  (*jvp)[i] += -this->nu_/this->dx_*(*vp)[i+1]
273  +(*up)[i+1]/6.0*(*vp)[i]
274  +((*up)[i]+2.0*(*up)[i+1])/6.0*(*vp)[i+1];
275  }
276  }
277  (*jvp)[0] -= this->u0_/6.0*(*vp)[0];
278  (*jvp)[this->nx_-1] += this->u1_/6.0*(*vp)[this->nx_-1];
279  }
280 
282  const ROL::Vector<Real> &z, Real &tol) {
283  Teuchos::RCP<std::vector<Real> > jvp =
284  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(jv)).getVector());
285  Teuchos::RCP<const std::vector<Real> > vp =
286  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(v))).getVector();
287  Teuchos::RCP<const std::vector<Real> > up =
288  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(u))).getVector();
289  Teuchos::RCP<const std::vector<Real> > zp =
290  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(z))).getVector();
291  for (int i=0; i<this->nx_; i++) {
292  // Contribution from control
293  (*jvp)[i] = -this->dx_/6.0*((*vp)[i]+4.0*(*vp)[i+1]+(*vp)[i+2]);
294  }
295  }
296 
298  const ROL::Vector<Real> &z, Real &tol) {
299  Teuchos::RCP<std::vector<Real> > ijvp =
300  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(ijv)).getVector());
301  Teuchos::RCP<const std::vector<Real> > vp =
302  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(v))).getVector();
303  Teuchos::RCP<const std::vector<Real> > up =
304  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(u))).getVector();
305  Teuchos::RCP<const std::vector<Real> > zp =
306  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(z))).getVector();
307  // Get PDE Jacobian
308  std::vector<Real> d(this->nx_,0.0);
309  std::vector<Real> dl(this->nx_-1,0.0);
310  std::vector<Real> du(this->nx_-1,0.0);
311  this->compute_pde_jacobian(dl,d,du,*up);
312  // Solve solve state sensitivity system at current time step
313  this->linear_solve(*ijvp,dl,d,du,*vp);
314  }
315 
317  const ROL::Vector<Real> &z, Real &tol) {
318  Teuchos::RCP<std::vector<Real> > jvp =
319  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(ajv)).getVector());
320  Teuchos::RCP<const std::vector<Real> > vp =
321  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(v))).getVector();
322  Teuchos::RCP<const std::vector<Real> > up =
323  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(u))).getVector();
324  Teuchos::RCP<const std::vector<Real> > zp =
325  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(z))).getVector();
326  // Fill jvp
327  for (int i = 0; i < this->nx_; i++) {
328  (*jvp)[i] = this->nu_/this->dx_*2.0*(*vp)[i];
329  if ( i > 0 ) {
330  (*jvp)[i] += -this->nu_/this->dx_*(*vp)[i-1]
331  -(*up)[i-1]/6.0*(*vp)[i]
332  +((*up)[i-1]+2.0*(*up)[i])/6.0*(*vp)[i-1];
333  }
334  if ( i < this->nx_-1 ) {
335  (*jvp)[i] += -this->nu_/this->dx_*(*vp)[i+1]
336  +(*up)[i+1]/6.0*(*vp)[i]
337  -((*up)[i+1]+2.0*(*up)[i])/6.0*(*vp)[i+1];
338  }
339  }
340  (*jvp)[0] -= this->u0_/6.0*(*vp)[0];
341  (*jvp)[this->nx_-1] += this->u1_/6.0*(*vp)[this->nx_-1];
342  }
343 
345  const ROL::Vector<Real> &z, Real &tol) {
346  Teuchos::RCP<std::vector<Real> > jvp =
347  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(jv)).getVector());
348  Teuchos::RCP<const std::vector<Real> > vp =
349  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(v))).getVector();
350  Teuchos::RCP<const std::vector<Real> > up =
351  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(u))).getVector();
352  Teuchos::RCP<const std::vector<Real> > zp =
353  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(z))).getVector();
354  for (int i=0; i<this->nx_+2; i++) {
355  if ( i == 0 ) {
356  (*jvp)[i] = -this->dx_/6.0*(*vp)[i];
357  }
358  else if ( i == 1 ) {
359  (*jvp)[i] = -this->dx_/6.0*(4.0*(*vp)[i-1]+(*vp)[i]);
360  }
361  else if ( i == this->nx_ ) {
362  (*jvp)[i] = -this->dx_/6.0*(4.0*(*vp)[i-1]+(*vp)[i-2]);
363  }
364  else if ( i == this->nx_+1 ) {
365  (*jvp)[i] = -this->dx_/6.0*(*vp)[i-2];
366  }
367  else {
368  (*jvp)[i] = -this->dx_/6.0*((*vp)[i-2]+4.0*(*vp)[i-1]+(*vp)[i]);
369  }
370  }
371  }
372 
374  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
375  Teuchos::RCP<std::vector<Real> > iajvp =
376  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(iajv)).getVector());
377  Teuchos::RCP<const std::vector<Real> > vp =
378  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(v))).getVector();
379  Teuchos::RCP<const std::vector<Real> > up =
380  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(u))).getVector();
381  // Get PDE Jacobian
382  std::vector<Real> d(this->nx_,0.0);
383  std::vector<Real> du(this->nx_-1,0.0);
384  std::vector<Real> dl(this->nx_-1,0.0);
385  this->compute_pde_jacobian(dl,d,du,*up);
386  // Solve solve adjoint system at current time step
387  this->linear_solve(*iajvp,dl,d,du,*vp,true);
388  }
389 
391  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
392  Teuchos::RCP<std::vector<Real> > ahwvp =
393  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(ahwv)).getVector());
394  Teuchos::RCP<const std::vector<Real> > wp =
395  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(w))).getVector();
396  Teuchos::RCP<const std::vector<Real> > vp =
397  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(v))).getVector();
398  Teuchos::RCP<const std::vector<Real> > up =
399  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(u))).getVector();
400  Teuchos::RCP<const std::vector<Real> > zp =
401  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(z))).getVector();
402  for (int i=0; i<this->nx_; i++) {
403  // Contribution from nonlinear term
404  (*ahwvp)[i] = 0.0;
405  if (i<this->nx_-1){
406  (*ahwvp)[i] += ((*wp)[i]*(*vp)[i+1] - (*wp)[i+1]*(2.0*(*vp)[i]+(*vp)[i+1]))/6.0;
407  }
408  if (i>0) {
409  (*ahwvp)[i] += ((*wp)[i-1]*((*vp)[i-1]+2.0*(*vp)[i]) - (*wp)[i]*(*vp)[i-1])/6.0;
410  }
411  }
412  }
413 
415  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
416  ahwv.zero();
417  }
419  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
420  ahwv.zero();
421  }
423  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
424  ahwv.zero();
425  }
426 
427 // void solveAugmentedSystem(ROL::Vector<Real> &v1, ROL::Vector<Real> &v2, const ROL::Vector<Real> &b1,
428 // const ROL::Vector<Real> &b2, const ROL::Vector<Real> &x, Real &tol) {}
429 };
430 
431 template<class Real>
432 class Objective_BurgersControl : public ROL::Objective_SimOpt<Real> {
433 private:
434  Real alpha_; // Penalty Parameter
435 
436  int nx_; // Number of interior nodes
437  Real dx_; // Mesh spacing (i.e. 1/(nx+1))
438 
439 /***************************************************************/
440 /********** BEGIN PRIVATE MEMBER FUNCTION DECLARATION **********/
441 /***************************************************************/
442  Real evaluate_target(Real x) {
443  Real val = 0.0;
444  int example = 2;
445  switch (example) {
446  case 1: val = ((x<0.5) ? 1.0 : 0.0); break;
447  case 2: val = 1.0; break;
448  case 3: val = std::abs(std::sin(8.0*M_PI*x)); break;
449  case 4: val = std::exp(-0.5*(x-0.5)*(x-0.5)); break;
450  }
451  return val;
452  }
453 
454  Real dot(const std::vector<Real> &x, const std::vector<Real> &y) {
455  Real ip = 0.0;
456  Real c = (((int)x.size()==this->nx_) ? 4.0 : 2.0);
457  for (unsigned i=0; i<x.size(); i++) {
458  if ( i == 0 ) {
459  ip += this->dx_/6.0*(c*x[i] + x[i+1])*y[i];
460  }
461  else if ( i == x.size()-1 ) {
462  ip += this->dx_/6.0*(x[i-1] + c*x[i])*y[i];
463  }
464  else {
465  ip += this->dx_/6.0*(x[i-1] + 4.0*x[i] + x[i+1])*y[i];
466  }
467  }
468  return ip;
469  }
470 
471  void apply_mass(std::vector<Real> &Mu, const std::vector<Real> &u ) {
472  Mu.resize(u.size(),0.0);
473  Real c = (((int)u.size()==this->nx_) ? 4.0 : 2.0);
474  for (unsigned i=0; i<u.size(); i++) {
475  if ( i == 0 ) {
476  Mu[i] = this->dx_/6.0*(c*u[i] + u[i+1]);
477  }
478  else if ( i == u.size()-1 ) {
479  Mu[i] = this->dx_/6.0*(u[i-1] + c*u[i]);
480  }
481  else {
482  Mu[i] = this->dx_/6.0*(u[i-1] + 4.0*u[i] + u[i+1]);
483  }
484  }
485  }
486 /*************************************************************/
487 /********** END PRIVATE MEMBER FUNCTION DECLARATION **********/
488 /*************************************************************/
489 
490 public:
491 
492  Objective_BurgersControl(Real alpha = 1.e-4, int nx = 128) : alpha_(alpha), nx_(nx) {
493  dx_ = 1.0/((Real)nx+1.0);
494  }
495 
496  Real value( const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
497  Teuchos::RCP<const std::vector<Real> > up =
498  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(u))).getVector();
499  Teuchos::RCP<const std::vector<Real> > zp =
500  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(z))).getVector();
501  // COMPUTE RESIDUAL
502  Real res1 = 0.0, res2 = 0.0, res3 = 0.0;
503  Real valu = 0.0, valz = this->dot(*zp,*zp);
504  for (int i=0; i<this->nx_; i++) {
505  if ( i == 0 ) {
506  res1 = (*up)[i]-evaluate_target((Real)(i+1)*this->dx_);
507  res2 = (*up)[i+1]-evaluate_target((Real)(i+2)*this->dx_);
508  valu += this->dx_/6.0*(4.0*res1 + res2)*res1;
509  }
510  else if ( i == this->nx_-1 ) {
511  res1 = (*up)[i-1]-evaluate_target((Real)i*this->dx_);
512  res2 = (*up)[i]-evaluate_target((Real)(i+1)*this->dx_);
513  valu += this->dx_/6.0*(res1 + 4.0*res2)*res2;
514  }
515  else {
516  res1 = (*up)[i-1]-evaluate_target((Real)i*this->dx_);
517  res2 = (*up)[i]-evaluate_target((Real)(i+1)*this->dx_);
518  res3 = (*up)[i+1]-evaluate_target((Real)(i+2)*this->dx_);
519  valu += this->dx_/6.0*(res1 + 4.0*res2 + res3)*res2;
520  }
521  }
522  return 0.5*(valu + this->alpha_*valz);
523  }
524 
525  void gradient_1( ROL::Vector<Real> &g, const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
526  // Unwrap g
527  Teuchos::RCP<std::vector<Real> > gup = Teuchos::rcp_const_cast<std::vector<Real> >(
528  (Teuchos::dyn_cast<const ROL::StdVector<Real> >(g)).getVector());
529  // Unwrap x
530  Teuchos::RCP<const std::vector<Real> > up =
531  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(u))).getVector();
532  Teuchos::RCP<const std::vector<Real> > zp =
533  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(z))).getVector();
534  // COMPUTE GRADIENT WRT U
535  std::vector<Real> diff(this->nx_,0.0);
536  for (int i=0; i<this->nx_; i++) {
537  diff[i] = ((*up)[i]-this->evaluate_target((Real)(i+1)*this->dx_));
538  }
539  this->apply_mass(*gup,diff);
540  }
541 
542  void gradient_2( ROL::Vector<Real> &g, const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
543  // Unwrap g
544  Teuchos::RCP<std::vector<Real> > gzp = Teuchos::rcp_const_cast<std::vector<Real> >(
545  (Teuchos::dyn_cast<const ROL::StdVector<Real> >(g)).getVector());
546  // Unwrap x
547  Teuchos::RCP<const std::vector<Real> > up =
548  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(u))).getVector();
549  Teuchos::RCP<const std::vector<Real> > zp =
550  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(z))).getVector();
551  // COMPUTE GRADIENT WRT Z
552  for (int i=0; i<this->nx_+2; i++) {
553  if (i==0) {
554  (*gzp)[i] = this->alpha_*this->dx_/6.0*(2.0*(*zp)[i]+(*zp)[i+1]);
555  }
556  else if (i==this->nx_+1) {
557  (*gzp)[i] = this->alpha_*this->dx_/6.0*(2.0*(*zp)[i]+(*zp)[i-1]);
558  }
559  else {
560  (*gzp)[i] = this->alpha_*this->dx_/6.0*((*zp)[i-1]+4.0*(*zp)[i]+(*zp)[i+1]);
561  }
562  }
563  }
564 
566  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
567  Teuchos::RCP<std::vector<Real> > hvup = Teuchos::rcp_const_cast<std::vector<Real> >(
568  (Teuchos::dyn_cast<const ROL::StdVector<Real> >(hv)).getVector());
569  // Unwrap v
570  Teuchos::RCP<const std::vector<Real> > vup =
571  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(v))).getVector();
572  // COMPUTE GRADIENT WRT U
573  this->apply_mass(*hvup,*vup);
574  }
575 
577  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
578  hv.zero();
579  }
580 
582  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
583  hv.zero();
584  }
585 
587  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
588  Teuchos::RCP<std::vector<Real> > hvzp = Teuchos::rcp_const_cast<std::vector<Real> >(
589  (Teuchos::dyn_cast<const ROL::StdVector<Real> >(hv)).getVector());
590  // Unwrap v
591  Teuchos::RCP<const std::vector<Real> > vzp =
592  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(v))).getVector();
593  // COMPUTE GRADIENT WRT Z
594  for (int i=0; i<this->nx_+2; i++) {
595  if (i==0) {
596  (*hvzp)[i] = this->alpha_*this->dx_/6.0*(2.0*(*vzp)[i]+(*vzp)[i+1]);
597  }
598  else if (i==this->nx_+1) {
599  (*hvzp)[i] = this->alpha_*this->dx_/6.0*(2.0*(*vzp)[i]+(*vzp)[i-1]);
600  }
601  else {
602  (*hvzp)[i] = this->alpha_*this->dx_/6.0*((*vzp)[i-1]+4.0*(*vzp)[i]+(*vzp)[i+1]);
603  }
604  }
605  }
606 };
Provides the interface to evaluate simulation-based objective functions.
void compute_pde_jacobian(std::vector< Real > &dl, std::vector< Real > &d, std::vector< Real > &du, const std::vector< Real > &u)
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 adjoint of the partial constraint Hessian at , , to vector in direction ...
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)
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 adjoint of the partial constraint Hessian at , , to vector in direction ...
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 .
Real dot(const std::vector< Real > &x, const std::vector< Real > &y)
Contains definitions of custom data types in ROL.
EqualityConstraint_BurgersControl(int nx=128, Real nu=1.e-2, Real u0=1.0, Real u1=0.0, Real f=0.0)
Real compute_norm(const std::vector< Real > &r)
virtual void zero()
Set to zero vector.
Definition: ROL_Vector.hpp:157
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:74
void compute_residual(std::vector< Real > &r, const std::vector< Real > &u, const std::vector< Real > &z)
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 adjoint of the partial constraint Hessian at , , to vector in direction ...
Defines the equality constraint operator interface for simulation-based optimization.
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 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 ...
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 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...
Provides the std::vector implementation of the ROL::Vector interface.
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 .
Real dot(const std::vector< Real > &x, const std::vector< Real > &y)
Objective_BurgersControl(Real alpha=1.e-4, int nx=128)
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)
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.
void solve(ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Given , solve for .
void scale(std::vector< Real > &u, const Real alpha=0.0)
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 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 value(ROL::Vector< Real > &c, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Evaluate the constraint operator at .
void apply_mass(std::vector< Real > &Mu, const std::vector< Real > &u)
virtual Real norm() const =0
Returns where .
void update(std::vector< Real > &u, const std::vector< Real > &s, const Real alpha=1.0)
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 adjoint of the partial constraint Hessian at , , to vector in direction ...
static const double ROL_EPSILON
Platform-dependent machine epsilon.
Definition: ROL_Types.hpp:118