ROL
example_03.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 
52 #include "ROL_StdVector.hpp"
53 #include "ROL_Vector_SimOpt.hpp"
55 #include "ROL_Objective_SimOpt.hpp"
57 
58 #include "Teuchos_oblackholestream.hpp"
59 #include "Teuchos_GlobalMPISession.hpp"
60 #include "Teuchos_XMLParameterListHelpers.hpp"
61 #include "Teuchos_LAPACK.hpp"
62 
63 #include <iostream>
64 #include <algorithm>
65 #include <ctime>
66 
67 template<class Real>
69 private:
70  unsigned nx_;
71  unsigned nt_;
72 
73  Real dx_;
74  Real T_;
75  Real dt_;
76 
77  Real nu_;
78  Real u0_;
79  Real u1_;
80  Real f_;
81  std::vector<Real> u_init_;
82 
83 private:
84  Real compute_norm(const std::vector<Real> &r) {
85  return std::sqrt(dot(r,r));
86  }
87 
88  Real dot(const std::vector<Real> &x, const std::vector<Real> &y) {
89  Real ip = 0.0;
90  Real c = ((x.size()==nx_) ? 4.0 : 2.0);
91  for (unsigned i = 0; i < x.size(); i++) {
92  if ( i == 0 ) {
93  ip += dx_/6.0*(c*x[i] + x[i+1])*y[i];
94  }
95  else if ( i == x.size()-1 ) {
96  ip += dx_/6.0*(x[i-1] + c*x[i])*y[i];
97  }
98  else {
99  ip += dx_/6.0*(x[i-1] + 4.0*x[i] + x[i+1])*y[i];
100  }
101  }
102  return ip;
103  }
104 
106 
107  void update(std::vector<Real> &u, const std::vector<Real> &s, const Real alpha=1.0) {
108  for (unsigned i = 0; i < u.size(); i++) {
109  u[i] += alpha*s[i];
110  }
111  }
112 
113  void scale(std::vector<Real> &u, const Real alpha=0.0) {
114  for (unsigned i = 0; i < u.size(); i++) {
115  u[i] *= alpha;
116  }
117  }
118 
119  void compute_residual(std::vector<Real> &r, const std::vector<Real> &uold, const std::vector<Real> &zold,
120  const std::vector<Real> &unew, const std::vector<Real> &znew) {
121  r.clear();
122  r.resize(nx_,0.0);
123  for (unsigned n = 0; n < nx_; n++) {
124  // Contribution from mass & stiffness term at time step t and t-1
125  r[n] += (4.0*dx_/6.0 + 0.5*dt_*2.0*nu_/dx_)*unew[n];
126  r[n] += (-4.0*dx_/6.0 + 0.5*dt_*2.0*nu_/dx_)*uold[n];
127  if ( n > 0 ) {
128  r[n] += (dx_/6.0 - 0.5*dt_*nu_/dx_)*unew[n-1];
129  r[n] -= (dx_/6.0 + 0.5*dt_*nu_/dx_)*uold[n-1];
130  }
131  if ( n < nx_-1 ) {
132  r[n] += (dx_/6.0 - 0.5*dt_*nu_/dx_)*unew[n+1];
133  r[n] -= (dx_/6.0 + 0.5*dt_*nu_/dx_)*uold[n+1];
134  }
135  // Contribution from nonlinear term
136  if ( n > 0 ) {
137  r[n] -= 0.5*dt_*unew[n-1]*(unew[n-1]+unew[n])/6.0;
138  r[n] -= 0.5*dt_*uold[n-1]*(uold[n-1]+uold[n])/6.0;
139  }
140  if ( n < nx_-1 ){
141  r[n] += 0.5*dt_*unew[n+1]*(unew[n]+unew[n+1])/6.0;
142  r[n] += 0.5*dt_*uold[n+1]*(uold[n]+uold[n+1])/6.0;
143  }
144  // Contribution from control
145  r[n] -= 0.5*dt_*dx_/6.0*(znew[n]+4.0*znew[n+1]+znew[n+2]);
146  r[n] -= 0.5*dt_*dx_/6.0*(zold[n]+4.0*zold[n+1]+zold[n+2]);
147  // Contribution from right-hand side
148  r[n] -= dt_*dx_*f_;
149  }
150  // Contribution from Dirichlet boundary terms
151  r[0] -= dt_*(0.5*u0_*(unew[ 0] + uold[ 0])/6.0 + u0_*u0_/6.0 + nu_*u0_/dx_);
152  r[nx_-1] += dt_*(0.5*u1_*(unew[nx_-1] + uold[nx_-1])/6.0 + u1_*u1_/6.0 - nu_*u1_/dx_);
153  }
154 
155  void compute_pde_jacobian(std::vector<Real> &dl, std::vector<Real> &d, std::vector<Real> &du,
156  const std::vector<Real> &u) {
157  // Get Diagonal and Off-Diagonal Entries of linear PDE Jacobian
158  d.clear();
159  d.resize(nx_,4.0*dx_/6.0 + 0.5*dt_*nu_*2.0/dx_);
160  dl.clear();
161  dl.resize(nx_-1,dx_/6.0 - 0.5*dt_*nu_/dx_);
162  du.clear();
163  du.resize(nx_-1,dx_/6.0 - 0.5*dt_*nu_/dx_);
164  // Contribution from nonlinearity
165  for (unsigned n = 0; n < nx_; n++) {
166  if ( n < nx_-1 ) {
167  dl[n] += 0.5*dt_*(-2.0*u[n]-u[n+1])/6.0;
168  d[n] += 0.5*dt_*u[n+1]/6.0;
169  }
170  if ( n > 0 ) {
171  d[n] -= 0.5*dt_*u[n-1]/6.0;
172  du[n-1] += 0.5*dt_*(u[n-1]+2.0*u[n])/6.0;
173  }
174  }
175  // Contribution from Dirichlet boundary conditions
176  d[0] -= 0.5*dt_*u0_/6.0;
177  d[nx_-1] += 0.5*dt_*u1_/6.0;
178  }
179 
180  void apply_pde_jacobian_new(std::vector<Real> &jv, const std::vector<Real> &v, const std::vector<Real> &u,
181  bool adjoint = false) {
182  jv.clear();
183  jv.resize(nx_,0.0);
184  // Fill Jacobian times a vector
185  for (unsigned n = 0; n < nx_; n++) {
186  jv[n] = (4.0*dx_/6.0 + 0.5*dt_*nu_/dx_*2.0)*v[n]; // Mass & Stiffness
187  if ( n > 0 ) {
188  jv[n] += dx_/6.0*v[n-1]-0.5*dt_*nu_/dx_*v[n-1]; // Mass & Stiffness
189  if ( adjoint ) {
190  jv[n] -= 0.5*dt_*(u[n-1]/6.0*v[n]-(u[n-1]+2.0*u[n])/6.0*v[n-1]); // Nonlinearity
191  }
192  else {
193  jv[n] -= 0.5*dt_*(u[n-1]/6.0*v[n]+(u[n]+2.0*u[n-1])/6.0*v[n-1]); // Nonlinearity
194  }
195  }
196  if ( n < nx_-1 ) {
197  jv[n] += dx_/6.0*v[n+1]-0.5*dt_*nu_/dx_*v[n+1]; // Mass & Stiffness
198  if ( adjoint ) {
199  jv[n] += 0.5*dt_*(u[n+1]/6.0*v[n]-(u[n+1]+2.0*u[n])/6.0*v[n+1]); // Nonlinearity
200  }
201  else {
202  jv[n] += 0.5*dt_*(u[n+1]/6.0*v[n]+(u[n]+2.0*u[n+1])/6.0*v[n+1]); // Nonlinearity
203  }
204  }
205  }
206  jv[0] -= 0.5*dt_*u0_/6.0*v[0]; // Nonlinearity
207  jv[nx_-1] += 0.5*dt_*u1_/6.0*v[nx_-1]; // Nonlinearity
208  }
209 
210  void apply_pde_jacobian_old(std::vector<Real> &jv, const std::vector<Real> &v, const std::vector<Real> &u,
211  bool adjoint = false) {
212  jv.clear();
213  jv.resize(nx_,0.0);
214  // Fill Jacobian times a vector
215  for (unsigned n = 0; n < nx_; n++) {
216  jv[n] = (-4.0*dx_/6.0 + 0.5*dt_*nu_/dx_*2.0)*v[n]; // Mass & Stiffness
217  if ( n > 0 ) {
218  jv[n] += -dx_/6.0*v[n-1]-0.5*dt_*nu_/dx_*v[n-1]; // Mass & Stiffness
219  if ( adjoint ) {
220  jv[n] -= 0.5*dt_*(u[n-1]/6.0*v[n]-(u[n-1]+2.0*u[n])/6.0*v[n-1]); // Nonlinearity
221  }
222  else {
223  jv[n] -= 0.5*dt_*(u[n-1]/6.0*v[n]+(u[n]+2.0*u[n-1])/6.0*v[n-1]); // Nonlinearity
224  }
225  }
226  if ( n < nx_-1 ) {
227  jv[n] += -dx_/6.0*v[n+1]-0.5*dt_*nu_/dx_*v[n+1]; // Mass & Stiffness
228  if ( adjoint ) {
229  jv[n] += 0.5*dt_*(u[n+1]/6.0*v[n]-(u[n+1]+2.0*u[n])/6.0*v[n+1]); // Nonlinearity
230  }
231  else {
232  jv[n] += 0.5*dt_*(u[n+1]/6.0*v[n]+(u[n]+2.0*u[n+1])/6.0*v[n+1]); // Nonlinearity
233  }
234  }
235  }
236  jv[0] -= 0.5*dt_*u0_/6.0*v[0]; // Nonlinearity
237  jv[nx_-1] += 0.5*dt_*u1_/6.0*v[nx_-1]; // Nonlinearity
238  }
239 
240  void apply_pde_jacobian(std::vector<Real> &jv, const std::vector<Real> &vold, const std::vector<Real> &uold,
241  const std::vector<Real> &vnew, const std::vector<Real> unew, bool adjoint = false) {
242  jv.clear();
243  jv.resize(nx_,0.0);
244  // Fill Jacobian times a vector
245  for (unsigned n = 0; n < nx_; n++) {
246  jv[n] += (4.0*dx_/6.0+0.5*dt_*nu_/dx_*2.0)*vnew[n]; // Mass & Stiffness
247  jv[n] -= (4.0*dx_/6.0-0.5*dt_*nu_/dx_*2.0)*vold[n]; // Mass & Stiffness
248  if ( n > 0 ) {
249  jv[n] += dx_/6.0*vnew[n-1]-0.5*dt_*nu_/dx_*vnew[n-1]; // Mass & Stiffness
250  jv[n] -= dx_/6.0*vold[n-1]+0.5*dt_*nu_/dx_*vold[n-1]; // Mass & Stiffness
251  if ( adjoint ) {
252  jv[n] -= 0.5*dt_*(unew[n-1]/6.0*vnew[n]-(unew[n-1]+2.0*unew[n])/6.0*vnew[n-1]); // Nonlinearity
253  jv[n] -= 0.5*dt_*(uold[n-1]/6.0*vold[n]-(uold[n-1]+2.0*uold[n])/6.0*vold[n-1]); // Nonlinearity
254  }
255  else {
256  jv[n] -= 0.5*dt_*(unew[n-1]/6.0*vnew[n]+(unew[n]+2.0*unew[n-1])/6.0*vnew[n-1]); // Nonlinearity
257  jv[n] -= 0.5*dt_*(uold[n-1]/6.0*vold[n]+(uold[n]+2.0*uold[n-1])/6.0*vold[n-1]); // Nonlinearity
258  }
259  }
260  if ( n < nx_-1 ) {
261  jv[n] += dx_/6.0*vnew[n+1]-0.5*dt_*nu_/dx_*vnew[n+1]; // Mass & Stiffness
262  jv[n] -= dx_/6.0*vold[n+1]+0.5*dt_*nu_/dx_*vold[n+1]; // Mass & Stiffness
263  if ( adjoint ) {
264  jv[n] += 0.5*dt_*(unew[n+1]/6.0*vnew[n]-(unew[n+1]+2.0*unew[n])/6.0*vnew[n+1]); // Nonlinearity
265  jv[n] += 0.5*dt_*(uold[n+1]/6.0*vold[n]-(uold[n+1]+2.0*uold[n])/6.0*vold[n+1]); // Nonlinearity
266  }
267  else {
268  jv[n] += 0.5*dt_*(unew[n+1]/6.0*vnew[n]+(unew[n]+2.0*unew[n+1])/6.0*vnew[n+1]); // Nonlinearity
269  jv[n] += 0.5*dt_*(uold[n+1]/6.0*vold[n]+(uold[n]+2.0*uold[n+1])/6.0*vold[n+1]); // Nonlinearity
270  }
271  }
272  }
273  jv[0] -= 0.5*dt_*u0_/6.0*vnew[0]; // Nonlinearity
274  jv[0] -= 0.5*dt_*u0_/6.0*vold[0]; // Nonlinearity
275  jv[nx_-1] += 0.5*dt_*u1_/6.0*vnew[nx_-1]; // Nonlinearity
276  jv[nx_-1] += 0.5*dt_*u1_/6.0*vold[nx_-1]; // Nonlinearity
277  }
278 
279  void apply_pde_hessian(std::vector<Real> &hv, const std::vector<Real> &wold, const std::vector<Real> &vold,
280  const std::vector<Real> &wnew, const std::vector<Real> &vnew) {
281  hv.clear();
282  hv.resize(nx_,0.0);
283  for (unsigned n = 0; n < nx_; n++) {
284  if ( n > 0 ) {
285  hv[n] += 0.5*dt_*((wnew[n-1]*(vnew[n-1]+2.0*vnew[n]) - wnew[n]*vnew[n-1])/6.0);
286  hv[n] += 0.5*dt_*((wold[n-1]*(vold[n-1]+2.0*vold[n]) - wold[n]*vold[n-1])/6.0);
287  }
288  if ( n < nx_-1 ){
289  hv[n] += 0.5*dt_*((wnew[n]*vnew[n+1] - wnew[n+1]*(2.0*vnew[n]+vnew[n+1]))/6.0);
290  hv[n] += 0.5*dt_*((wold[n]*vold[n+1] - wold[n+1]*(2.0*vold[n]+vold[n+1]))/6.0);
291  }
292  }
293  }
294 
295  void apply_control_jacobian(std::vector<Real> &jv, const std::vector<Real> &v, bool adjoint = false) {
296  jv.clear();
297  unsigned dim = ((adjoint == true) ? nx_+2 : nx_);
298  jv.resize(dim,0.0);
299  for (unsigned n = 0; n < dim; n++) {
300  if ( adjoint ) {
301  if ( n == 0 ) {
302  jv[n] = -0.5*dt_*dx_/6.0*v[n];
303  }
304  else if ( n == 1 ) {
305  jv[n] = -0.5*dt_*dx_/6.0*(4.0*v[n-1]+v[n]);
306  }
307  else if ( n == nx_ ) {
308  jv[n] = -0.5*dt_*dx_/6.0*(4.0*v[n-1]+v[n-2]);
309  }
310  else if ( n == nx_+1 ) {
311  jv[n] = -0.5*dt_*dx_/6.0*v[n-2];
312  }
313  else {
314  jv[n] = -0.5*dt_*dx_/6.0*(v[n-2]+4.0*v[n-1]+v[n]);
315  }
316  }
317  else {
318  jv[n] -= 0.5*dt_*dx_/6.0*(v[n]+4.0*v[n+1]+v[n+2]);
319  }
320  }
321  }
322 
323  void run_newton(std::vector<Real> &u, const std::vector<Real> &znew,
324  const std::vector<Real> &uold, const std::vector<Real> &zold) {
325  u.clear();
326  u.resize(nx_,0.0);
327  // Compute residual and residual norm
328  std::vector<Real> r(nx_,0.0);
329  compute_residual(r,uold,zold,u,znew);
330  Real rnorm = compute_norm(r);
331  // Define tolerances
332  Real rtol = 1.e2*ROL::ROL_EPSILON<Real>();
333  unsigned maxit = 500;
334  // Initialize Jacobian storage
335  std::vector<Real> d(nx_,0.0);
336  std::vector<Real> dl(nx_-1,0.0);
337  std::vector<Real> du(nx_-1,0.0);
338  // Iterate Newton's method
339  Real alpha = 1.0, tmp = 0.0;
340  std::vector<Real> s(nx_,0.0);
341  std::vector<Real> utmp(nx_,0.0);
342  for (unsigned i = 0; i < maxit; i++) {
343  //std::cout << i << " " << rnorm << "\n";
344  // Get Jacobian
345  compute_pde_jacobian(dl,d,du,u);
346  // Solve Newton system
347  linear_solve(s,dl,d,du,r);
348  // Perform line search
349  tmp = rnorm;
350  alpha = 1.0;
351  utmp.assign(u.begin(),u.end());
352  update(utmp,s,-alpha);
353  compute_residual(r,uold,zold,utmp,znew);
354  rnorm = compute_norm(r);
355  while ( rnorm > (1.0-1.e-4*alpha)*tmp && alpha > std::sqrt(ROL::ROL_EPSILON<Real>()) ) {
356  alpha *= 0.5;
357  utmp.assign(u.begin(),u.end());
358  update(utmp,s,-alpha);
359  compute_residual(r,uold,zold,utmp,znew);
360  rnorm = compute_norm(r);
361  }
362  // Update iterate
363  u.assign(utmp.begin(),utmp.end());
364  if ( rnorm < rtol ) {
365  break;
366  }
367  }
368  }
369 
370  void linear_solve(std::vector<Real> &u,
371  const std::vector<Real> &dl, const std::vector<Real> &d, const std::vector<Real> &du,
372  const std::vector<Real> &r, const bool transpose = false) {
373  bool useLAPACK = false;
374  if ( useLAPACK ) { // DIRECT SOLVE: USE LAPACK
375  u.assign(r.begin(),r.end());
376  // Store matrix diagonal & off-diagonals.
377  std::vector<Real> Dl(dl);
378  std::vector<Real> Du(du);
379  std::vector<Real> D(d);
380  // Perform LDL factorization
381  Teuchos::LAPACK<int,Real> lp;
382  std::vector<Real> Du2(nx_-2,0.0);
383  std::vector<int> ipiv(nx_,0);
384  int info;
385  int ldb = nx_;
386  int nhrs = 1;
387  lp.GTTRF(nx_,&Dl[0],&D[0],&Du[0],&Du2[0],&ipiv[0],&info);
388  char trans = ((transpose == true) ? 'T' : 'N');
389  lp.GTTRS(trans,nx_,nhrs,&Dl[0],&D[0],&Du[0],&Du2[0],&ipiv[0],&u[0],ldb,&info);
390  }
391  else { // ITERATIVE SOLVE: USE GAUSS-SEIDEL
392  u.clear();
393  u.resize(nx_,0.0);
394  unsigned maxit = 100;
395  Real rtol = std::min(1.e-12,1.e-4*std::sqrt(dot(r,r)));
396  Real resid = 0.0;
397  Real rnorm = 10.0*rtol;
398  for (unsigned i = 0; i < maxit; i++) {
399  for (unsigned n = 0; n < nx_; n++) {
400  u[n] = r[n]/d[n];
401  if ( n < nx_-1 ) {
402  u[n] -= ((transpose == false) ? du[n] : dl[n])*u[n+1]/d[n];
403  }
404  if ( n > 0 ) {
405  u[n] -= ((transpose == false) ? dl[n-1] : du[n-1])*u[n-1]/d[n];
406  }
407  }
408  // Compute Residual
409  rnorm = 0.0;
410  for (unsigned n = 0; n < nx_; n++) {
411  resid = r[n] - d[n]*u[n];
412  if ( n < nx_-1 ) {
413  resid -= ((transpose == false) ? du[n] : dl[n])*u[n+1];
414  }
415  if ( n > 0 ) {
416  resid -= ((transpose == false) ? dl[n-1] : du[n-1])*u[n-1];
417  }
418  rnorm += resid*resid;
419  }
420  rnorm = std::sqrt(rnorm);
421  if ( rnorm < rtol ) {
422  //std::cout << "\ni = " << i+1 << " rnorm = " << rnorm << "\n";
423  break;
424  }
425  }
426  }
427  }
428 
429 public:
430 
431  EqualityConstraint_BurgersControl(int nx = 128, int nt = 100, Real T = 1,
432  Real nu = 1.e-2, Real u0 = 0.0, Real u1 = 0.0, Real f = 0.0)
433  : nx_((unsigned)nx), nt_((unsigned)nt), T_(T), nu_(nu), u0_(u0), u1_(u1), f_(f) {
434  dx_ = 1.0/((Real)nx+1.0);
435  dt_ = T/((Real)nt);
436  u_init_.clear();
437  u_init_.resize(nx_,0.0);
438  Real x = 0.0;
439  for (unsigned n = 0; n < nx_; n++) {
440  x = (Real)(n+1)*dx_;
441  u_init_[n] = ((x <= 0.5) ? 1.0 : 0.0);
442  }
443  }
444 
445  void value(ROL::Vector<Real> &c, const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
446  Teuchos::RCP<std::vector<Real> > cp =
447  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(c)).getVector());
448  Teuchos::RCP<const std::vector<Real> > up =
449  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(u))).getVector();
450  Teuchos::RCP<const std::vector<Real> > zp =
451  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(z))).getVector();
452  // Initialize storage
453  std::vector<Real> C(nx_,0.0);
454  std::vector<Real> uold(u_init_);
455  std::vector<Real> unew(u_init_);
456  std::vector<Real> zold(nx_+2,0.0);
457  std::vector<Real> znew(nx_+2,0.0);
458  // Copy initial control
459  for (unsigned n = 0; n < nx_+2; n++) {
460  zold[n] = (*zp)[n];
461  }
462  // Evaluate residual
463  for (unsigned t = 0; t < nt_; t++) {
464  // Copy current state at time step t
465  for (unsigned n = 0; n < nx_; n++) {
466  unew[n] = (*up)[t*nx_+n];
467  }
468  // Copy current control at time step t
469  for (unsigned n = 0; n < nx_+2; n++) {
470  znew[n] = (*zp)[(t+1)*(nx_+2)+n];
471  }
472  // Compute residual at time step t
473  compute_residual(C,uold,zold,unew,znew);
474  // Store residual at time step t
475  for (unsigned n = 0; n < nx_; n++) {
476  (*cp)[t*nx_+n] = C[n];
477  }
478  // Copy previous state and control at time step t+1
479  uold.assign(unew.begin(),unew.end());
480  zold.assign(znew.begin(),znew.end());
481  }
482  }
483 
484  void solve(ROL::Vector<Real> &c, ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
485  Teuchos::RCP<std::vector<Real> > up =
486  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(u)).getVector());
487  up->assign(up->size(),z.norm()/up->size());
488  Teuchos::RCP<const std::vector<Real> > zp =
489  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(z))).getVector();
490  // Initialize storage
491  std::vector<Real> uold(u_init_);
492  std::vector<Real> unew(u_init_);
493  std::vector<Real> zold(nx_+2,0.0);
494  std::vector<Real> znew(nx_+2,0.0);
495  // Copy initial control
496  for (unsigned n = 0; n < nx_+2; n++) {
497  zold[n] = (*zp)[n];
498  }
499  // Time step using the trapezoidal rule
500  for (unsigned t = 0; t < nt_; t++) {
501  // Copy current control at time step t
502  for (unsigned n = 0; n < nx_+2; n++) {
503  znew[n] = (*zp)[(t+1)*(nx_+2)+n];
504  }
505  // Solve nonlinear system at time step t
506  run_newton(unew,znew,uold,zold);
507  // store state at time step t
508  for (unsigned n = 0; n < nx_; n++) {
509  (*up)[t*nx_+n] = unew[n];
510  }
511  // Copy previous state and control at time step t+1
512  uold.assign(unew.begin(),unew.end());
513  zold.assign(znew.begin(),znew.end());
514  }
515  value(c,u,z,tol);
516  }
517 
519  const ROL::Vector<Real> &z, Real &tol) {
520  Teuchos::RCP<std::vector<Real> > jvp =
521  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(jv)).getVector());
522  Teuchos::RCP<const std::vector<Real> > vp =
523  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(v))).getVector();
524  Teuchos::RCP<const std::vector<Real> > up =
525  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(u))).getVector();
526  std::vector<Real> J(nx_,0.0);
527  std::vector<Real> vold(nx_,0.0);
528  std::vector<Real> vnew(nx_,0.0);
529  std::vector<Real> uold(nx_,0.0);
530  std::vector<Real> unew(nx_,0.0);
531  for (unsigned t = 0; t < nt_; t++) {
532  for (unsigned n = 0; n < nx_; n++) {
533  unew[n] = (*up)[t*nx_+n];
534  vnew[n] = (*vp)[t*nx_+n];
535  }
536  apply_pde_jacobian(J,vold,uold,vnew,unew);
537  for (unsigned n = 0; n < nx_; n++) {
538  (*jvp)[t*nx_+n] = J[n];
539  }
540  vold.assign(vnew.begin(),vnew.end());
541  uold.assign(unew.begin(),unew.end());
542  }
543  }
544 
546  const ROL::Vector<Real> &z, Real &tol) {
547  Teuchos::RCP<std::vector<Real> > jvp =
548  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(jv)).getVector());
549  Teuchos::RCP<const std::vector<Real> > vp =
550  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(v))).getVector();
551  Teuchos::RCP<const std::vector<Real> > zp =
552  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(z))).getVector();
553  std::vector<Real> vnew(nx_+2,0.0);
554  std::vector<Real> vold(nx_+2,0.0);
555  std::vector<Real> jnew(nx_,0.0);
556  std::vector<Real> jold(nx_,0.0);
557  for (unsigned n = 0; n < nx_+2; n++) {
558  vold[n] = (*vp)[n];
559  }
560  apply_control_jacobian(jold,vold);
561  for (unsigned t = 0; t < nt_; t++) {
562  for (unsigned n = 0; n < nx_+2; n++) {
563  vnew[n] = (*vp)[(t+1)*(nx_+2)+n];
564  }
565  apply_control_jacobian(jnew,vnew);
566  for (unsigned n = 0; n < nx_; n++) {
567  // Contribution from control
568  (*jvp)[t*nx_+n] = jnew[n] + jold[n];
569  }
570  jold.assign(jnew.begin(),jnew.end());
571  }
572  }
573 
575  const ROL::Vector<Real> &z, Real &tol) {
576  Teuchos::RCP<std::vector<Real> > ijvp =
577  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(ijv)).getVector());
578  Teuchos::RCP<const std::vector<Real> > vp =
579  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(v))).getVector();
580  Teuchos::RCP<const std::vector<Real> > up =
581  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(u))).getVector();
582  std::vector<Real> J(nx_,0.0);
583  std::vector<Real> r(nx_,0.0);
584  std::vector<Real> s(nx_,0.0);
585  std::vector<Real> uold(nx_,0.0);
586  std::vector<Real> unew(nx_,0.0);
587  std::vector<Real> d(nx_,0.0);
588  std::vector<Real> dl(nx_-1,0.0);
589  std::vector<Real> du(nx_-1,0.0);
590  for (unsigned t = 0; t < nt_; t++) {
591  // Build rhs.
592  for (unsigned n = 0; n < nx_; n++) {
593  r[n] = (*vp)[t*nx_+n];
594  unew[n] = (*up)[t*nx_+n];
595  }
596  apply_pde_jacobian_old(J,s,uold);
597  update(r,J,-1.0);
598  // Build Jacobian.
599  compute_pde_jacobian(dl,d,du,unew);
600  // Solve linear system.
601  linear_solve(s,dl,d,du,r);
602  // Copy solution.
603  for (unsigned n = 0; n < nx_; n++) {
604  (*ijvp)[t*nx_+n] = s[n];
605  }
606  uold.assign(unew.begin(),unew.end());
607  }
608  }
609 
611  const ROL::Vector<Real> &z, Real &tol) {
612  Teuchos::RCP<std::vector<Real> > jvp =
613  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(ajv)).getVector());
614  Teuchos::RCP<const std::vector<Real> > vp =
615  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(v))).getVector();
616  Teuchos::RCP<const std::vector<Real> > up =
617  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(u))).getVector();
618  std::vector<Real> J(nx_,0.0);
619  std::vector<Real> vold(nx_,0.0);
620  std::vector<Real> vnew(nx_,0.0);
621  std::vector<Real> unew(nx_,0.0);
622  for (unsigned t = nt_; t > 0; t--) {
623  for (unsigned n = 0; n < nx_; n++) {
624  unew[n] = (*up)[(t-1)*nx_+n];
625  vnew[n] = (*vp)[(t-1)*nx_+n];
626  }
627  apply_pde_jacobian(J,vold,unew,vnew,unew,true);
628  for (unsigned n = 0; n < nx_; n++) {
629  (*jvp)[(t-1)*nx_+n] = J[n];
630  }
631  vold.assign(vnew.begin(),vnew.end());
632  }
633  }
634 
636  const ROL::Vector<Real> &z, Real &tol) {
637  Teuchos::RCP<std::vector<Real> > jvp =
638  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(jv)).getVector());
639  Teuchos::RCP<const std::vector<Real> > vp =
640  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(v))).getVector();
641  Teuchos::RCP<const std::vector<Real> > zp =
642  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(z))).getVector();
643  std::vector<Real> vnew(nx_,0.0);
644  std::vector<Real> vold(nx_,0.0);
645  std::vector<Real> jnew(nx_+2,0.0);
646  std::vector<Real> jold(nx_+2,0.0);
647  for (unsigned t = nt_+1; t > 0; t--) {
648  for (unsigned n = 0; n < nx_; n++) {
649  if ( t > 1 ) {
650  vnew[n] = (*vp)[(t-2)*nx_+n];
651  }
652  else {
653  vnew[n] = 0.0;
654  }
655  }
656  apply_control_jacobian(jnew,vnew,true);
657  for (unsigned n = 0; n < nx_+2; n++) {
658  // Contribution from control
659  (*jvp)[(t-1)*(nx_+2)+n] = jnew[n] + jold[n];
660  }
661  jold.assign(jnew.begin(),jnew.end());
662  }
663  }
664 
666  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
667  Teuchos::RCP<std::vector<Real> > ijvp =
668  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(ijv)).getVector());
669  Teuchos::RCP<const std::vector<Real> > vp =
670  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(v))).getVector();
671  Teuchos::RCP<const std::vector<Real> > up =
672  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(u))).getVector();
673  std::vector<Real> J(nx_,0.0);
674  std::vector<Real> r(nx_,0.0);
675  std::vector<Real> s(nx_,0.0);
676  std::vector<Real> unew(nx_,0.0);
677  std::vector<Real> d(nx_,0.0);
678  std::vector<Real> dl(nx_-1,0.0);
679  std::vector<Real> du(nx_-1,0.0);
680  for (unsigned t = nt_; t > 0; t--) {
681  // Build rhs.
682  for (unsigned n = 0; n < nx_; n++) {
683  r[n] = (*vp)[(t-1)*nx_+n];
684  unew[n] = (*up)[(t-1)*nx_+n];
685  }
686  apply_pde_jacobian_old(J,s,unew,true);
687  update(r,J,-1.0);
688  // Build Jacobian.
689  compute_pde_jacobian(dl,d,du,unew);
690  // Solve linear system.
691  linear_solve(s,dl,d,du,r,true);
692  // Copy solution.
693  for (unsigned n = 0; n < nx_; n++) {
694  (*ijvp)[(t-1)*nx_+n] = s[n];
695  }
696  }
697  }
698 
700  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
701  Teuchos::RCP<std::vector<Real> > hwvp =
702  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(hwv)).getVector());
703  Teuchos::RCP<const std::vector<Real> > wp =
704  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(w))).getVector();
705  Teuchos::RCP<const std::vector<Real> > vp =
706  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(v))).getVector();
707  std::vector<Real> snew(nx_,0.0);
708  std::vector<Real> wnew(nx_,0.0);
709  std::vector<Real> wold(nx_,0.0);
710  std::vector<Real> vnew(nx_,0.0);
711  for (unsigned t = nt_; t > 0; t--) {
712  for (unsigned n = 0; n < nx_; n++) {
713  vnew[n] = (*vp)[(t-1)*nx_+n];
714  wnew[n] = (*wp)[(t-1)*nx_+n];
715  }
716  apply_pde_hessian(snew,wold,vnew,wnew,vnew);
717  for (unsigned n = 0; n < nx_; n++) {
718  (*hwvp)[(t-1)*nx_+n] = snew[n];
719  }
720  wold.assign(wnew.begin(),wnew.end());
721  }
722  }
723 
725  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
726  ahwv.zero();
727  }
729  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
730  ahwv.zero();
731  }
733  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
734  ahwv.zero();
735  }
736 };
737 
738 template<class Real>
739 class Objective_BurgersControl : public ROL::Objective_SimOpt<Real> {
740 private:
741  Real alpha_; // Penalty Parameter
742 
743  unsigned nx_; // Number of interior nodes
744  Real dx_; // Mesh spacing (i.e. 1/(nx+1))
745  unsigned nt_; // Number of time steps
746  Real dt_; // Time step size
747  Real T_; // Final time
748 
749 /***************************************************************/
750 /********** BEGIN PRIVATE MEMBER FUNCTION DECLARATION **********/
751 /***************************************************************/
752  Real evaluate_target(Real x) {
753  Real val = 0.0;
754  int example = 1;
755  switch (example) {
756  case 1: val = ((x<=0.5) ? 1.0 : 0.0); break;
757  case 2: val = 1.0; break;
758  case 3: val = std::abs(std::sin(8.0*M_PI*x)); break;
759  case 4: val = std::exp(-0.5*(x-0.5)*(x-0.5)); break;
760  }
761  return val;
762  }
763 
764  Real dot(const std::vector<Real> &x, const std::vector<Real> &y) {
765  Real ip = 0.0;
766  Real c = ((x.size()==nx_) ? 4.0 : 2.0);
767  for (unsigned i=0; i<x.size(); i++) {
768  if ( i == 0 ) {
769  ip += dx_/6.0*(c*x[i] + x[i+1])*y[i];
770  }
771  else if ( i == x.size()-1 ) {
772  ip += dx_/6.0*(x[i-1] + c*x[i])*y[i];
773  }
774  else {
775  ip += dx_/6.0*(x[i-1] + 4.0*x[i] + x[i+1])*y[i];
776  }
777  }
778  return ip;
779  }
780 
781  void apply_mass(std::vector<Real> &Mu, const std::vector<Real> &u ) {
782  Mu.resize(u.size(),0.0);
783  Real c = ((u.size()==nx_) ? 4.0 : 2.0);
784  for (unsigned i=0; i<u.size(); i++) {
785  if ( i == 0 ) {
786  Mu[i] = dx_/6.0*(c*u[i] + u[i+1]);
787  }
788  else if ( i == u.size()-1 ) {
789  Mu[i] = dx_/6.0*(u[i-1] + c*u[i]);
790  }
791  else {
792  Mu[i] = dx_/6.0*(u[i-1] + 4.0*u[i] + u[i+1]);
793  }
794  }
795  }
796 /*************************************************************/
797 /********** END PRIVATE MEMBER FUNCTION DECLARATION **********/
798 /*************************************************************/
799 
800 public:
801 
802  Objective_BurgersControl(Real alpha = 1.e-4, int nx = 128, int nt = 100, Real T = 1.0)
803  : alpha_(alpha), nx_((unsigned)nx), nt_((unsigned)nt), T_(T) {
804  dx_ = 1.0/((Real)nx+1.0);
805  dt_ = T/((Real)nt);
806  }
807 
808  Real value( const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
809  Teuchos::RCP<const std::vector<Real> > up =
810  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(u))).getVector();
811  Teuchos::RCP<const std::vector<Real> > zp =
812  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(z))).getVector();
813  // COMPUTE RESIDUAL
814  std::vector<Real> U(nx_,0.0);
815  std::vector<Real> Z(nx_+2,0.0);
816  for (unsigned n = 0; n < nx_+2; n++) {
817  Z[n] = (*zp)[n];
818  }
819  Real ss = 0.5*dt_;
820  Real val = 0.5*ss*alpha_*dot(Z,Z);
821  for (unsigned t = 0; t < nt_; t++) {
822  ss = ((t < nt_-1) ? dt_ : 0.5*dt_);
823  for (unsigned n = 0; n < nx_; n++) {
824  U[n] = (*up)[t*nx_+n]-evaluate_target((Real)(n+1)*dx_);
825  }
826  val += 0.5*ss*dot(U,U);
827  for (unsigned n = 0; n < nx_+2; n++) {
828  Z[n] = (*zp)[(t+1)*(nx_+2)+n];
829  }
830  val += 0.5*ss*alpha_*dot(Z,Z);
831  }
832  return val;
833  }
834 
835  void gradient_1( ROL::Vector<Real> &g, const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
836  Teuchos::RCP<std::vector<Real> > gp = Teuchos::rcp_const_cast<std::vector<Real> >(
837  (Teuchos::dyn_cast<const ROL::StdVector<Real> >(g)).getVector());
838  Teuchos::RCP<const std::vector<Real> > up =
839  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(u))).getVector();
840  // COMPUTE GRADIENT WRT U
841  std::vector<Real> U(nx_,0.0);
842  std::vector<Real> M(nx_,0.0);
843  Real ss = dt_;
844  for (unsigned t = 0; t < nt_; t++) {
845  ss = ((t < nt_-1) ? dt_ : 0.5*dt_);
846  for (unsigned n = 0; n < nx_; n++) {
847  U[n] = (*up)[t*nx_+n]-evaluate_target((Real)(n+1)*dx_);
848  }
849  apply_mass(M,U);
850  for (unsigned n = 0; n < nx_; n++) {
851  (*gp)[t*nx_+n] = ss*M[n];
852  }
853  }
854  }
855 
856  void gradient_2( ROL::Vector<Real> &g, const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
857  Teuchos::RCP<std::vector<Real> > gp = Teuchos::rcp_const_cast<std::vector<Real> >(
858  (Teuchos::dyn_cast<const ROL::StdVector<Real> >(g)).getVector());
859  Teuchos::RCP<const std::vector<Real> > zp =
860  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(z))).getVector();
861  // COMPUTE GRADIENT WRT Z
862  std::vector<Real> Z(nx_+2,0.0);
863  std::vector<Real> M(nx_+2,0.0);
864  Real ss = dt_;
865  for (unsigned t = 0; t < nt_+1; t++) {
866  ss = ((t < nt_ && t > 0) ? dt_ : 0.5*dt_);
867  for (unsigned n = 0; n < nx_+2; n++) {
868  Z[n] = (*zp)[t*(nx_+2)+n];
869  }
870  apply_mass(M,Z);
871  for (unsigned n = 0; n < nx_+2; n++) {
872  (*gp)[t*(nx_+2)+n] = ss*alpha_*M[n];
873  }
874  }
875  }
876 
878  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
879  Teuchos::RCP<std::vector<Real> > hvp = Teuchos::rcp_const_cast<std::vector<Real> >(
880  (Teuchos::dyn_cast<const ROL::StdVector<Real> >(hv)).getVector());
881  Teuchos::RCP<const std::vector<Real> > vp =
882  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(v))).getVector();
883  // COMPUTE GRADIENT WRT U
884  std::vector<Real> V(nx_,0.0);
885  std::vector<Real> M(nx_,0.0);
886  Real ss = 0.5*dt_;
887  for (unsigned t = 0; t < nt_; t++) {
888  ss = ((t < nt_-1) ? dt_ : 0.5*dt_);
889  for (unsigned n = 0; n < nx_; n++) {
890  V[n] = (*vp)[t*nx_+n];
891  }
892  apply_mass(M,V);
893  for (unsigned n = 0; n < nx_; n++) {
894  (*hvp)[t*nx_+n] = ss*M[n];
895  }
896  }
897  }
898 
900  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
901  hv.zero();
902  }
903 
905  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
906  hv.zero();
907  }
908 
910  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
911  Teuchos::RCP<std::vector<Real> > hvp = Teuchos::rcp_const_cast<std::vector<Real> >(
912  (Teuchos::dyn_cast<const ROL::StdVector<Real> >(hv)).getVector());
913  Teuchos::RCP<const std::vector<Real> > vp =
914  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(v))).getVector();
915  // COMPUTE GRADIENT WRT Z
916  std::vector<Real> V(nx_+2,0.0);
917  std::vector<Real> M(nx_+2,0.0);
918  Real ss = 0.0;
919  for (unsigned t = 0; t < nt_+1; t++) {
920  ss = ((t < nt_ && t > 0) ? dt_ : 0.5*dt_);
921  for (unsigned n = 0; n < nx_+2; n++) {
922  V[n] = (*vp)[t*(nx_+2)+n];
923  }
924  apply_mass(M,V);
925  for (unsigned n = 0; n < nx_+2; n++) {
926  (*hvp)[t*(nx_+2)+n] = ss*alpha_*M[n];
927  }
928  }
929  }
930 };
Provides the interface to evaluate simulation-based objective functions.
void run_newton(std::vector< Real > &u, const std::vector< Real > &znew, const std::vector< Real > &uold, const std::vector< Real > &zold)
Definition: example_03.hpp:323
void compute_pde_jacobian(std::vector< Real > &dl, std::vector< Real > &d, std::vector< Real > &du, const std::vector< Real > &u)
Definition: example_03.hpp:155
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...
Definition: example_03.hpp:728
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)
Real evaluate_target(Real x)
Definition: example_03.hpp:752
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 ...
Definition: example_03.hpp:732
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 .
Definition: example_03.hpp:574
Real dot(const std::vector< Real > &x, const std::vector< Real > &y)
Contains definitions of custom data types in ROL.
void compute_residual(std::vector< Real > &r, const std::vector< Real > &uold, const std::vector< Real > &zold, const std::vector< Real > &unew, const std::vector< Real > &znew)
Definition: example_03.hpp:119
Real compute_norm(const std::vector< Real > &r)
Definition: example_03.hpp:84
EqualityConstraint_BurgersControl(int nx=128, int nt=100, Real T=1, Real nu=1.e-2, Real u0=0.0, Real u1=0.0, Real f=0.0)
Definition: example_03.hpp:431
void solve(ROL::Vector< Real > &c, ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Given , solve for .
Definition: example_03.hpp:484
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 apply_control_jacobian(std::vector< Real > &jv, const std::vector< Real > &v, bool adjoint=false)
Definition: example_03.hpp:295
void compute_residual(std::vector< Real > &r, const std::vector< Real > &u, const std::vector< Real > &z)
void apply_pde_hessian(std::vector< Real > &hv, const std::vector< Real > &wold, const std::vector< Real > &vold, const std::vector< Real > &wnew, const std::vector< Real > &vnew)
Definition: example_03.hpp:279
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)
Definition: example_03.hpp:909
Real value(const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Compute value.
Definition: example_03.hpp:808
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.
Definition: example_03.hpp:835
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...
Definition: example_03.hpp:610
void apply_pde_jacobian_new(std::vector< Real > &jv, const std::vector< Real > &v, const std::vector< Real > &u, bool adjoint=false)
Definition: example_03.hpp:180
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...
Definition: example_03.hpp:635
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 .
Definition: example_03.hpp:518
Real dot(const std::vector< Real > &x, const std::vector< Real > &y)
Definition: example_03.hpp:764
void hessVec_21(ROL::Vector< Real > &hv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Definition: example_03.hpp:904
void hessVec_12(ROL::Vector< Real > &hv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Definition: example_03.hpp:899
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.
Definition: example_03.hpp:856
void scale(std::vector< Real > &u, const Real alpha=0.0)
Definition: example_03.hpp:113
Objective_BurgersControl(Real alpha=1.e-4, int nx=128, int nt=100, Real T=1.0)
Definition: example_03.hpp:802
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 .
Definition: example_03.hpp:545
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.
Definition: example_03.hpp:877
void value(ROL::Vector< Real > &c, const ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Evaluate the constraint operator at .
Definition: example_03.hpp:445
void applyInverseAdjointJacobian_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 of the adjoint of the partial constraint Jacobian at , , to the vector ...
Definition: example_03.hpp:665
void apply_mass(std::vector< Real > &Mu, const std::vector< Real > &u)
Definition: example_03.hpp:781
void linear_solve(std::vector< Real > &u, const std::vector< Real > &dl, const std::vector< Real > &d, const std::vector< Real > &du, const std::vector< Real > &r, const bool transpose=false)
Definition: example_03.hpp:370
void apply_pde_jacobian(std::vector< Real > &jv, const std::vector< Real > &vold, const std::vector< Real > &uold, const std::vector< Real > &vnew, const std::vector< Real > unew, bool adjoint=false)
Definition: example_03.hpp:240
virtual Real norm() const =0
Returns where .
void update(std::vector< Real > &u, const std::vector< Real > &s, const Real alpha=1.0)
Definition: example_03.hpp:107
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...
Definition: example_03.hpp:724
void apply_pde_jacobian_old(std::vector< Real > &jv, const std::vector< Real > &v, const std::vector< Real > &u, bool adjoint=false)
Definition: example_03.hpp:210
void applyAdjointHessian_11(ROL::Vector< Real > &hwv, 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 ...
Definition: example_03.hpp:699