ROL
burgers-control/example_02.cpp
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_CompositeStepSQP.hpp"
51 #include "ROL_TrustRegionStep.hpp"
52 #include "ROL_StatusTest.hpp"
53 #include "ROL_Types.hpp"
54 #include "Teuchos_oblackholestream.hpp"
55 #include "Teuchos_GlobalMPISession.hpp"
56 #include "Teuchos_XMLParameterListHelpers.hpp"
57 #include "Teuchos_LAPACK.hpp"
58 
59 #include <iostream>
60 #include <algorithm>
61 
62 #include "ROL_StdVector.hpp"
63 #include "ROL_Vector_SimOpt.hpp"
65 #include "ROL_Objective_SimOpt.hpp"
67 
68 template<class Real>
70 private:
71  int nx_;
72  Real dx_;
73  Real nu_;
74  Real u0_;
75  Real u1_;
76  Real f_;
77 
78 private:
79  Real compute_norm(const std::vector<Real> &r) {
80  return std::sqrt(this->dot(r,r));
81  }
82 
83  Real dot(const std::vector<Real> &x, const std::vector<Real> &y) {
84  Real ip = 0.0;
85  Real c = (((int)x.size()==this->nx_) ? 4.0 : 2.0);
86  for (unsigned i=0; i<x.size(); i++) {
87  if ( i == 0 ) {
88  ip += this->dx_/6.0*(c*x[i] + x[i+1])*y[i];
89  }
90  else if ( i == x.size()-1 ) {
91  ip += this->dx_/6.0*(x[i-1] + c*x[i])*y[i];
92  }
93  else {
94  ip += this->dx_/6.0*(x[i-1] + 4.0*x[i] + x[i+1])*y[i];
95  }
96  }
97  return ip;
98  }
99 
100  void update(std::vector<Real> &u, const std::vector<Real> &s, const Real alpha=1.0) {
101  for (unsigned i=0; i<u.size(); i++) {
102  u[i] += alpha*s[i];
103  }
104  }
105 
106  void scale(std::vector<Real> &u, const Real alpha=0.0) {
107  for (unsigned i=0; i<u.size(); i++) {
108  u[i] *= alpha;
109  }
110  }
111 
112  void compute_residual(std::vector<Real> &r, const std::vector<Real> &u,
113  const std::vector<Real> &z) {
114  r.clear();
115  r.resize(this->nx_,0.0);
116  for (int i=0; i<this->nx_; i++) {
117  // Contribution from stiffness term
118  if ( i==0 ) {
119  r[i] = this->nu_/this->dx_*(2.0*u[i]-u[i+1]);
120  }
121  else if (i==this->nx_-1) {
122  r[i] = this->nu_/this->dx_*(2.0*u[i]-u[i-1]);
123  }
124  else {
125  r[i] = this->nu_/this->dx_*(2.0*u[i]-u[i-1]-u[i+1]);
126  }
127  // Contribution from nonlinear term
128  if (i<this->nx_-1){
129  r[i] += u[i+1]*(u[i]+u[i+1])/6.0;
130  }
131  if (i>0) {
132  r[i] -= u[i-1]*(u[i-1]+u[i])/6.0;
133  }
134  // Contribution from control
135  r[i] -= this->dx_/6.0*(z[i]+4.0*z[i+1]+z[i+2]);
136  // Contribution from right-hand side
137  r[i] -= this->dx_*this->f_;
138  }
139  // Contribution from Dirichlet boundary terms
140  r[0] -= this->u0_*u[ 0]/6.0 + this->u0_*this->u0_/6.0 + this->nu_*this->u0_/this->dx_;
141  r[this->nx_-1] += this->u1_*u[this->nx_-1]/6.0 + this->u1_*this->u1_/6.0 - this->nu_*this->u1_/this->dx_;
142  }
143 
144  void compute_pde_jacobian(std::vector<Real> &dl, std::vector<Real> &d, std::vector<Real> &du,
145  const std::vector<Real> &u) {
146  // Get Diagonal and Off-Diagonal Entries of linear PDE Jacobian
147  d.clear();
148  d.resize(this->nx_,this->nu_*2.0/this->dx_);
149  dl.clear();
150  dl.resize(this->nx_-1,-this->nu_/this->dx_);
151  du.clear();
152  du.resize(this->nx_-1,-this->nu_/this->dx_);
153  // Contribution from nonlinearity
154  for (int i=0; i<this->nx_; i++) {
155  if (i<this->nx_-1) {
156  dl[i] += (-2.0*u[i]-u[i+1])/6.0;
157  d[i] += u[i+1]/6.0;
158  }
159  if (i>0) {
160  d[i] += -u[i-1]/6.0;
161  du[i-1] += (u[i-1]+2.0*u[i])/6.0;
162  }
163  }
164  // Contribution from Dirichlet boundary conditions
165  d[0] -= this->u0_/6.0;
166  d[this->nx_-1] += this->u1_/6.0;
167  }
168 
169  void linear_solve(std::vector<Real> &u, std::vector<Real> &dl, std::vector<Real> &d, std::vector<Real> &du,
170  const std::vector<Real> &r, const bool transpose = false) {
171  u.assign(r.begin(),r.end());
172  // Perform LDL factorization
173  Teuchos::LAPACK<int,Real> lp;
174  std::vector<Real> du2(this->nx_-2,0.0);
175  std::vector<int> ipiv(this->nx_,0);
176  int info;
177  int ldb = this->nx_;
178  int nhrs = 1;
179  lp.GTTRF(this->nx_,&dl[0],&d[0],&du[0],&du2[0],&ipiv[0],&info);
180  char trans = 'N';
181  if ( transpose ) {
182  trans = 'T';
183  }
184  lp.GTTRS(trans,this->nx_,nhrs,&dl[0],&d[0],&du[0],&du2[0],&ipiv[0],&u[0],ldb,&info);
185  }
186 
187 public:
188 
189  EqualityConstraint_BurgersControl(int nx = 128, Real nu = 1.e-2, Real u0 = 1.0, Real u1 = 0.0, Real f = 0.0)
190  : nx_(nx), nu_(nu), u0_(u0), u1_(u1), f_(f) {
191  dx_ = 1.0/((Real)nx+1.0);
192  }
193 
195  const ROL::Vector<Real> &z, Real &tol) {
196  Teuchos::RCP<std::vector<Real> > cp =
197  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(c)).getVector());
198  Teuchos::RCP<const std::vector<Real> > up =
199  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(u))).getVector();
200  Teuchos::RCP<const std::vector<Real> > zp =
201  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(z))).getVector();
202  this->compute_residual(*cp,*up,*zp);
203  }
204 
205  void solve(ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
206  Teuchos::RCP<std::vector<Real> > up =
207  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(u)).getVector());
208  up->assign(up->size(),z.norm()/up->size());
209  Teuchos::RCP<const std::vector<Real> > zp =
210  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(z))).getVector();
211  // Compute residual and residual norm
212  std::vector<Real> r(up->size(),0.0);
213  this->compute_residual(r,*up,*zp);
214  Real rnorm = this->compute_norm(r);
215  // Define tolerances
216  Real rtol = 1.e2*ROL::ROL_EPSILON;
217  Real maxit = 500;
218  // Initialize Jacobian storage
219  std::vector<Real> d(this->nx_,0.0);
220  std::vector<Real> dl(this->nx_-1,0.0);
221  std::vector<Real> du(this->nx_-1,0.0);
222  // Iterate Newton's method
223  Real alpha = 1.0, tmp = 0.0;
224  std::vector<Real> s(this->nx_,0.0);
225  std::vector<Real> utmp(this->nx_,0.0);
226  for (int i=0; i<maxit; i++) {
227  //std::cout << i << " " << rnorm << "\n";
228  // Get Jacobian
229  this->compute_pde_jacobian(dl,d,du,*up);
230  // Solve Newton system
231  this->linear_solve(s,dl,d,du,r);
232  // Perform line search
233  tmp = rnorm;
234  alpha = 1.0;
235  utmp.assign(up->begin(),up->end());
236  this->update(utmp,s,-alpha);
237  this->compute_residual(r,utmp,*zp);
238  rnorm = this->compute_norm(r);
239  while ( rnorm > (1.0-1.e-4*alpha)*tmp && alpha > std::sqrt(ROL::ROL_EPSILON) ) {
240  alpha /= 2.0;
241  utmp.assign(up->begin(),up->end());
242  this->update(utmp,s,-alpha);
243  this->compute_residual(r,utmp,*zp);
244  rnorm = this->compute_norm(r);
245  }
246  // Update iterate
247  up->assign(utmp.begin(),utmp.end());
248  if ( rnorm < rtol ) {
249  break;
250  }
251  }
252  }
253 
254 
255 
257  const ROL::Vector<Real> &z, Real &tol) {
258  Teuchos::RCP<std::vector<Real> > jvp =
259  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(jv)).getVector());
260  Teuchos::RCP<const std::vector<Real> > vp =
261  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(v))).getVector();
262  Teuchos::RCP<const std::vector<Real> > up =
263  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(u))).getVector();
264  Teuchos::RCP<const std::vector<Real> > zp =
265  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(z))).getVector();
266  // Fill jvp
267  for (int i = 0; i < this->nx_; i++) {
268  (*jvp)[i] = this->nu_/this->dx_*2.0*(*vp)[i];
269  if ( i > 0 ) {
270  (*jvp)[i] += -this->nu_/this->dx_*(*vp)[i-1]
271  -(*up)[i-1]/6.0*(*vp)[i]
272  -((*up)[i]+2.0*(*up)[i-1])/6.0*(*vp)[i-1];
273  }
274  if ( i < this->nx_-1 ) {
275  (*jvp)[i] += -this->nu_/this->dx_*(*vp)[i+1]
276  +(*up)[i+1]/6.0*(*vp)[i]
277  +((*up)[i]+2.0*(*up)[i+1])/6.0*(*vp)[i+1];
278  }
279  }
280  (*jvp)[0] -= this->u0_/6.0*(*vp)[0];
281  (*jvp)[this->nx_-1] += this->u1_/6.0*(*vp)[this->nx_-1];
282  }
283 
285  const ROL::Vector<Real> &z, Real &tol) {
286  Teuchos::RCP<std::vector<Real> > jvp =
287  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(jv)).getVector());
288  Teuchos::RCP<const std::vector<Real> > vp =
289  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(v))).getVector();
290  Teuchos::RCP<const std::vector<Real> > up =
291  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(u))).getVector();
292  Teuchos::RCP<const std::vector<Real> > zp =
293  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(z))).getVector();
294  for (int i=0; i<this->nx_; i++) {
295  // Contribution from control
296  (*jvp)[i] = -this->dx_/6.0*((*vp)[i]+4.0*(*vp)[i+1]+(*vp)[i+2]);
297  }
298  }
299 
301  const ROL::Vector<Real> &z, Real &tol) {
302  Teuchos::RCP<std::vector<Real> > ijvp =
303  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(ijv)).getVector());
304  Teuchos::RCP<const std::vector<Real> > vp =
305  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(v))).getVector();
306  Teuchos::RCP<const std::vector<Real> > up =
307  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(u))).getVector();
308  Teuchos::RCP<const std::vector<Real> > zp =
309  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(z))).getVector();
310  // Get PDE Jacobian
311  std::vector<Real> d(this->nx_,0.0);
312  std::vector<Real> dl(this->nx_-1,0.0);
313  std::vector<Real> du(this->nx_-1,0.0);
314  this->compute_pde_jacobian(dl,d,du,*up);
315  // Solve solve state sensitivity system at current time step
316  this->linear_solve(*ijvp,dl,d,du,*vp);
317  }
318 
320  const ROL::Vector<Real> &z, Real &tol) {
321  Teuchos::RCP<std::vector<Real> > jvp =
322  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(ajv)).getVector());
323  Teuchos::RCP<const std::vector<Real> > vp =
324  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(v))).getVector();
325  Teuchos::RCP<const std::vector<Real> > up =
326  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(u))).getVector();
327  Teuchos::RCP<const std::vector<Real> > zp =
328  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(z))).getVector();
329  // Fill jvp
330  for (int i = 0; i < this->nx_; i++) {
331  (*jvp)[i] = this->nu_/this->dx_*2.0*(*vp)[i];
332  if ( i > 0 ) {
333  (*jvp)[i] += -this->nu_/this->dx_*(*vp)[i-1]
334  -(*up)[i-1]/6.0*(*vp)[i]
335  +((*up)[i-1]+2.0*(*up)[i])/6.0*(*vp)[i-1];
336  }
337  if ( i < this->nx_-1 ) {
338  (*jvp)[i] += -this->nu_/this->dx_*(*vp)[i+1]
339  +(*up)[i+1]/6.0*(*vp)[i]
340  -((*up)[i+1]+2.0*(*up)[i])/6.0*(*vp)[i+1];
341  }
342  }
343  (*jvp)[0] -= this->u0_/6.0*(*vp)[0];
344  (*jvp)[this->nx_-1] += this->u1_/6.0*(*vp)[this->nx_-1];
345  }
346 
348  const ROL::Vector<Real> &z, Real &tol) {
349  Teuchos::RCP<std::vector<Real> > jvp =
350  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(jv)).getVector());
351  Teuchos::RCP<const std::vector<Real> > vp =
352  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(v))).getVector();
353  Teuchos::RCP<const std::vector<Real> > up =
354  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(u))).getVector();
355  Teuchos::RCP<const std::vector<Real> > zp =
356  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(z))).getVector();
357  for (int i=0; i<this->nx_+2; i++) {
358  if ( i == 0 ) {
359  (*jvp)[i] = -this->dx_/6.0*(*vp)[i];
360  }
361  else if ( i == 1 ) {
362  (*jvp)[i] = -this->dx_/6.0*(4.0*(*vp)[i-1]+(*vp)[i]);
363  }
364  else if ( i == this->nx_ ) {
365  (*jvp)[i] = -this->dx_/6.0*(4.0*(*vp)[i-1]+(*vp)[i-2]);
366  }
367  else if ( i == this->nx_+1 ) {
368  (*jvp)[i] = -this->dx_/6.0*(*vp)[i-2];
369  }
370  else {
371  (*jvp)[i] = -this->dx_/6.0*((*vp)[i-2]+4.0*(*vp)[i-1]+(*vp)[i]);
372  }
373  }
374  }
375 
377  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
378  Teuchos::RCP<std::vector<Real> > iajvp =
379  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(iajv)).getVector());
380  Teuchos::RCP<const std::vector<Real> > vp =
381  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(v))).getVector();
382  Teuchos::RCP<const std::vector<Real> > up =
383  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(u))).getVector();
384  // Get PDE Jacobian
385  std::vector<Real> d(this->nx_,0.0);
386  std::vector<Real> du(this->nx_-1,0.0);
387  std::vector<Real> dl(this->nx_-1,0.0);
388  this->compute_pde_jacobian(dl,d,du,*up);
389  // Solve solve adjoint system at current time step
390  this->linear_solve(*iajvp,dl,d,du,*vp,true);
391  }
392 
394  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
395  Teuchos::RCP<std::vector<Real> > ahwvp =
396  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(ahwv)).getVector());
397  Teuchos::RCP<const std::vector<Real> > wp =
398  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(w))).getVector();
399  Teuchos::RCP<const std::vector<Real> > vp =
400  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(v))).getVector();
401  Teuchos::RCP<const std::vector<Real> > up =
402  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(u))).getVector();
403  Teuchos::RCP<const std::vector<Real> > zp =
404  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(z))).getVector();
405  for (int i=0; i<this->nx_; i++) {
406  // Contribution from nonlinear term
407  (*ahwvp)[i] = 0.0;
408  if (i<this->nx_-1){
409  (*ahwvp)[i] += ((*wp)[i]*(*vp)[i+1] - (*wp)[i+1]*(2.0*(*vp)[i]+(*vp)[i+1]))/6.0;
410  }
411  if (i>0) {
412  (*ahwvp)[i] += ((*wp)[i-1]*((*vp)[i-1]+2.0*(*vp)[i]) - (*wp)[i]*(*vp)[i-1])/6.0;
413  }
414  }
415  }
416 
418  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
419  ahwv.zero();
420  }
422  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
423  ahwv.zero();
424  }
426  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
427  ahwv.zero();
428  }
429 
430  //void solveAugmentedSystem(ROL::Vector<Real> &v1, ROL::Vector<Real> &v2, const ROL::Vector<Real> &b1,
431  // const ROL::Vector<Real> &b2, const ROL::Vector<Real> &x, Real &tol) {}
432 };
433 
434 template<class Real>
435 class Objective_BurgersControl : public ROL::Objective_SimOpt<Real> {
436 private:
437  Real alpha_; // Penalty Parameter
438 
439  int nx_; // Number of interior nodes
440  Real dx_; // Mesh spacing (i.e. 1/(nx+1))
441 
442 /***************************************************************/
443 /********** BEGIN PRIVATE MEMBER FUNCTION DECLARATION **********/
444 /***************************************************************/
445  Real evaluate_target(Real x) {
446  Real val = 0.0;
447  int example = 2;
448  switch (example) {
449  case 1: val = ((x<0.5) ? 1.0 : 0.0); break;
450  case 2: val = 1.0; break;
451  case 3: val = std::abs(std::sin(8.0*M_PI*x)); break;
452  case 4: val = std::exp(-0.5*(x-0.5)*(x-0.5)); break;
453  }
454  return val;
455  }
456 
457  Real dot(const std::vector<Real> &x, const std::vector<Real> &y) {
458  Real ip = 0.0;
459  Real c = (((int)x.size()==this->nx_) ? 4.0 : 2.0);
460  for (unsigned i=0; i<x.size(); i++) {
461  if ( i == 0 ) {
462  ip += this->dx_/6.0*(c*x[i] + x[i+1])*y[i];
463  }
464  else if ( i == x.size()-1 ) {
465  ip += this->dx_/6.0*(x[i-1] + c*x[i])*y[i];
466  }
467  else {
468  ip += this->dx_/6.0*(x[i-1] + 4.0*x[i] + x[i+1])*y[i];
469  }
470  }
471  return ip;
472  }
473 
474  void apply_mass(std::vector<Real> &Mu, const std::vector<Real> &u ) {
475  Mu.resize(u.size(),0.0);
476  Real c = (((int)u.size()==this->nx_) ? 4.0 : 2.0);
477  for (unsigned i=0; i<u.size(); i++) {
478  if ( i == 0 ) {
479  Mu[i] = this->dx_/6.0*(c*u[i] + u[i+1]);
480  }
481  else if ( i == u.size()-1 ) {
482  Mu[i] = this->dx_/6.0*(u[i-1] + c*u[i]);
483  }
484  else {
485  Mu[i] = this->dx_/6.0*(u[i-1] + 4.0*u[i] + u[i+1]);
486  }
487  }
488  }
489 /*************************************************************/
490 /********** END PRIVATE MEMBER FUNCTION DECLARATION **********/
491 /*************************************************************/
492 
493 public:
494 
495  Objective_BurgersControl(Real alpha = 1.e-4, int nx = 128) : alpha_(alpha), nx_(nx) {
496  dx_ = 1.0/((Real)nx+1.0);
497  }
498 
499  Real value( const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
500  Teuchos::RCP<const std::vector<Real> > up =
501  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(u))).getVector();
502  Teuchos::RCP<const std::vector<Real> > zp =
503  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(z))).getVector();
504  // COMPUTE RESIDUAL
505  Real res1 = 0.0, res2 = 0.0, res3 = 0.0;
506  Real valu = 0.0, valz = this->dot(*zp,*zp);
507  for (int i=0; i<this->nx_; i++) {
508  if ( i == 0 ) {
509  res1 = (*up)[i]-evaluate_target((Real)(i+1)*this->dx_);
510  res2 = (*up)[i+1]-evaluate_target((Real)(i+2)*this->dx_);
511  valu += this->dx_/6.0*(4.0*res1 + res2)*res1;
512  }
513  else if ( i == this->nx_-1 ) {
514  res1 = (*up)[i-1]-evaluate_target((Real)i*this->dx_);
515  res2 = (*up)[i]-evaluate_target((Real)(i+1)*this->dx_);
516  valu += this->dx_/6.0*(res1 + 4.0*res2)*res2;
517  }
518  else {
519  res1 = (*up)[i-1]-evaluate_target((Real)i*this->dx_);
520  res2 = (*up)[i]-evaluate_target((Real)(i+1)*this->dx_);
521  res3 = (*up)[i+1]-evaluate_target((Real)(i+2)*this->dx_);
522  valu += this->dx_/6.0*(res1 + 4.0*res2 + res3)*res2;
523  }
524  }
525  return 0.5*(valu + this->alpha_*valz);
526  }
527 
528  void gradient_1( ROL::Vector<Real> &g, const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
529  // Unwrap g
530  Teuchos::RCP<std::vector<Real> > gup = Teuchos::rcp_const_cast<std::vector<Real> >(
531  (Teuchos::dyn_cast<const ROL::StdVector<Real> >(g)).getVector());
532  // Unwrap x
533  Teuchos::RCP<const std::vector<Real> > up =
534  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(u))).getVector();
535  Teuchos::RCP<const std::vector<Real> > zp =
536  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(z))).getVector();
537  // COMPUTE GRADIENT WRT U
538  std::vector<Real> diff(this->nx_,0.0);
539  for (int i=0; i<this->nx_; i++) {
540  diff[i] = ((*up)[i]-this->evaluate_target((Real)(i+1)*this->dx_));
541  }
542  this->apply_mass(*gup,diff);
543  }
544 
545  void gradient_2( ROL::Vector<Real> &g, const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
546  // Unwrap g
547  Teuchos::RCP<std::vector<Real> > gzp = Teuchos::rcp_const_cast<std::vector<Real> >(
548  (Teuchos::dyn_cast<const ROL::StdVector<Real> >(g)).getVector());
549  // Unwrap x
550  Teuchos::RCP<const std::vector<Real> > up =
551  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(u))).getVector();
552  Teuchos::RCP<const std::vector<Real> > zp =
553  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(z))).getVector();
554  // COMPUTE GRADIENT WRT Z
555  for (int i=0; i<this->nx_+2; i++) {
556  if (i==0) {
557  (*gzp)[i] = this->alpha_*this->dx_/6.0*(2.0*(*zp)[i]+(*zp)[i+1]);
558  }
559  else if (i==this->nx_+1) {
560  (*gzp)[i] = this->alpha_*this->dx_/6.0*(2.0*(*zp)[i]+(*zp)[i-1]);
561  }
562  else {
563  (*gzp)[i] = this->alpha_*this->dx_/6.0*((*zp)[i-1]+4.0*(*zp)[i]+(*zp)[i+1]);
564  }
565  }
566  }
567 
569  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
570  Teuchos::RCP<std::vector<Real> > hvup = Teuchos::rcp_const_cast<std::vector<Real> >(
571  (Teuchos::dyn_cast<const ROL::StdVector<Real> >(hv)).getVector());
572  // Unwrap v
573  Teuchos::RCP<const std::vector<Real> > vup =
574  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(v))).getVector();
575  // COMPUTE GRADIENT WRT U
576  this->apply_mass(*hvup,*vup);
577  }
578 
580  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
581  hv.zero();
582  }
583 
585  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
586  hv.zero();
587  }
588 
590  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
591  Teuchos::RCP<std::vector<Real> > hvzp = Teuchos::rcp_const_cast<std::vector<Real> >(
592  (Teuchos::dyn_cast<const ROL::StdVector<Real> >(hv)).getVector());
593  // Unwrap v
594  Teuchos::RCP<const std::vector<Real> > vzp =
595  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(v))).getVector();
596  // COMPUTE GRADIENT WRT Z
597  for (int i=0; i<this->nx_+2; i++) {
598  if (i==0) {
599  (*hvzp)[i] = this->alpha_*this->dx_/6.0*(2.0*(*vzp)[i]+(*vzp)[i+1]);
600  }
601  else if (i==this->nx_+1) {
602  (*hvzp)[i] = this->alpha_*this->dx_/6.0*(2.0*(*vzp)[i]+(*vzp)[i-1]);
603  }
604  else {
605  (*hvzp)[i] = this->alpha_*this->dx_/6.0*((*vzp)[i-1]+4.0*(*vzp)[i]+(*vzp)[i+1]);
606  }
607  }
608  }
609 };
610 
611 typedef double RealT;
612 
613 int main(int argc, char *argv[]) {
614 
615  Teuchos::GlobalMPISession mpiSession(&argc, &argv);
616 
617  // This little trick lets us print to std::cout only if a (dummy) command-line argument is provided.
618  int iprint = argc - 1;
619  Teuchos::RCP<std::ostream> outStream;
620  Teuchos::oblackholestream bhs; // outputs nothing
621  if (iprint > 0)
622  outStream = Teuchos::rcp(&std::cout, false);
623  else
624  outStream = Teuchos::rcp(&bhs, false);
625 
626  int errorFlag = 0;
627 
628  // *** Example body.
629 
630  try {
631  // Initialize full objective function.
632  int nx = 256; // Set spatial discretization.
633  RealT alpha = 1.e-3; // Set penalty parameter.
634  RealT nu = 1.e-2; // Set viscosity parameter.
635  Objective_BurgersControl<RealT> obj(alpha,nx);
636  // Initialize equality constraints
638  // Initialize iteration vectors.
639  Teuchos::RCP<std::vector<RealT> > z_rcp = Teuchos::rcp( new std::vector<RealT> (nx+2, 1.0) );
640  Teuchos::RCP<std::vector<RealT> > yz_rcp = Teuchos::rcp( new std::vector<RealT> (nx+2, 1.0) );
641  for (int i=0; i<nx+2; i++) {
642  (*z_rcp)[i] = (RealT)rand()/(RealT)RAND_MAX;
643  (*yz_rcp)[i] = (RealT)rand()/(RealT)RAND_MAX;
644  }
645  ROL::StdVector<RealT> z(z_rcp);
646  ROL::StdVector<RealT> yz(yz_rcp);
647  Teuchos::RCP<ROL::Vector<RealT> > zp = Teuchos::rcp(&z,false);
648  Teuchos::RCP<ROL::Vector<RealT> > yzp = Teuchos::rcp(&yz,false);
649 
650  Teuchos::RCP<std::vector<RealT> > u_rcp = Teuchos::rcp( new std::vector<RealT> (nx, 1.0) );
651  Teuchos::RCP<std::vector<RealT> > yu_rcp = Teuchos::rcp( new std::vector<RealT> (nx, 1.0) );
652  for (int i=0; i<nx; i++) {
653  (*u_rcp)[i] = (RealT)rand()/(RealT)RAND_MAX;
654  (*yu_rcp)[i] = (RealT)rand()/(RealT)RAND_MAX;
655  }
656  ROL::StdVector<RealT> u(u_rcp);
657  ROL::StdVector<RealT> yu(yu_rcp);
658  Teuchos::RCP<ROL::Vector<RealT> > up = Teuchos::rcp(&u,false);
659  Teuchos::RCP<ROL::Vector<RealT> > yup = Teuchos::rcp(&yu,false);
660 
661  Teuchos::RCP<std::vector<RealT> > jv_rcp = Teuchos::rcp( new std::vector<RealT> (nx, 1.0) );
662  ROL::StdVector<RealT> jv(jv_rcp);
663  Teuchos::RCP<ROL::Vector<RealT> > jvp = Teuchos::rcp(&jv,false);
664 
665  ROL::Vector_SimOpt<RealT> x(up,zp);
666  ROL::Vector_SimOpt<RealT> y(yup,yzp);
667  // Check derivatives.
668  obj.checkGradient(x,x,y,true,*outStream);
669  obj.checkHessVec(x,x,y,true,*outStream);
670  con.checkApplyJacobian(x,y,jv,true,*outStream);
671  con.checkApplyAdjointHessian(x,yu,y,x,true,*outStream);
672 
673  // Initialize reduced objective function.
674  Teuchos::RCP<std::vector<RealT> > p_rcp = Teuchos::rcp( new std::vector<RealT> (nx, 1.0) );
675  ROL::StdVector<RealT> p(p_rcp);
676  Teuchos::RCP<ROL::Vector<RealT> > pp = Teuchos::rcp(&p,false);
677  Teuchos::RCP<ROL::Objective_SimOpt<RealT> > pobj = Teuchos::rcp(&obj,false);
678  Teuchos::RCP<ROL::EqualityConstraint_SimOpt<RealT> > pcon = Teuchos::rcp(&con,false);
679  ROL::Reduced_Objective_SimOpt<RealT> robj(pobj,pcon,up,pp);
680  // Check derivatives.
681  robj.checkGradient(z,z,yz,true,*outStream);
682  robj.checkHessVec(z,z,yz,true,*outStream);
683  // Optimization
684  std::string filename = "input.xml";
685  Teuchos::RCP<Teuchos::ParameterList> parlist_tr = Teuchos::rcp( new Teuchos::ParameterList() );
686  Teuchos::updateParametersFromXmlFile( filename, Teuchos::Ptr<Teuchos::ParameterList>(&*parlist_tr) );
687  // Define status test.
688  RealT gtol = 1e-14; // norm of gradient tolerance
689  RealT stol = 1e-16; // norm of step tolerance
690  int maxit = 1000; // maximum number of iterations
691  ROL::StatusTest<RealT> status(gtol, stol, maxit);
692  // Define step.
693  ROL::TrustRegionStep<RealT> stepr(*parlist_tr);
694  // Define algorithm.
695  ROL::DefaultAlgorithm<RealT> algo(stepr,status,false);
696  // Run Algorithm
697  //z.zero();
698  z.scale(50.0);
699  algo.run(z,robj,true,*outStream);
700  }
701  catch (std::logic_error err) {
702  *outStream << err.what() << "\n";
703  errorFlag = -1000;
704  }; // end try
705 
706  if (errorFlag != 0)
707  std::cout << "End Result: TEST FAILED\n";
708  else
709  std::cout << "End Result: TEST PASSED\n";
710 
711  return 0;
712 
713 }
714 
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 scale(const Real alpha)
Compute where .
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)
Defines the linear algebra or vector space interface for simulation-based optimization.
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:155
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:72
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.
virtual std::vector< std::vector< Real > > checkGradient(const Vector< Real > &x, const Vector< Real > &d, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
Finite-difference gradient check.
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.
virtual std::vector< std::vector< Real > > checkApplyJacobian(const Vector< Real > &x, 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)
Finite-difference check for the constraint Jacobian application.
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 ...
Provides the interface to evaluate simulation-based reduced 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 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 .
virtual std::vector< std::string > run(Vector< Real > &x, Objective< Real > &obj, bool print=false, std::ostream &outStream=std::cout)
Run algorithm on unconstrained problems (Type-U). This is the primary Type-U interface.
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)
Provides an interface to check status of optimization algorithms.
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.
int main(int argc, char *argv[])
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 std::vector< std::vector< Real > > checkHessVec(const Vector< Real > &x, const Vector< Real > &v, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
Finite-difference Hessian-applied-to-vector check.
double RealT
virtual Real norm() const =0
Returns where .
virtual std::vector< std::vector< Real > > checkApplyAdjointHessian(const Vector< Real > &x, const Vector< Real > &u, const Vector< Real > &v, const Vector< Real > &hv, const std::vector< Real > &step, const bool printToScreen=true, std::ostream &outStream=std::cout, const int order=1)
Finite-difference check for the application of the adjoint of constraint Hessian. ...
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 ...
Provides the interface to compute optimization steps with trust regions.
static const double ROL_EPSILON
Platform-dependent machine epsilon.
Definition: ROL_Types.hpp:115