ROL
example_07.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_Types.hpp"
50 #include "ROL_Vector.hpp"
51 #include "ROL_BoundConstraint.hpp"
55 
56 #include "Teuchos_LAPACK.hpp"
57 
58 template<class Real>
59 class L2VectorPrimal;
60 
61 template<class Real>
62 class L2VectorDual;
63 
64 template<class Real>
65 class H1VectorPrimal;
66 
67 template<class Real>
68 class H1VectorDual;
69 
70 template<class Real>
71 class BurgersFEM {
72 private:
73  int nx_;
74  Real dx_;
75  Real nu_;
76  Real nl_;
77  Real u0_;
78  Real u1_;
79  Real f_;
80  Real cH1_;
81  Real cL2_;
82 
83 private:
84  void update(std::vector<Real> &u, const std::vector<Real> &s, const Real alpha=1.0) const {
85  for (unsigned i=0; i<u.size(); i++) {
86  u[i] += alpha*s[i];
87  }
88  }
89 
90  void axpy(std::vector<Real> &out, const Real a, const std::vector<Real> &x, const std::vector<Real> &y) const {
91  for (unsigned i=0; i < x.size(); i++) {
92  out[i] = a*x[i] + y[i];
93  }
94  }
95 
96  void scale(std::vector<Real> &u, const Real alpha=0.0) const {
97  for (unsigned i=0; i<u.size(); i++) {
98  u[i] *= alpha;
99  }
100  }
101 
102  void linear_solve(std::vector<Real> &u, std::vector<Real> &dl, std::vector<Real> &d, std::vector<Real> &du,
103  const std::vector<Real> &r, const bool transpose = false) const {
104  if ( r.size() == 1 ) {
105  u.resize(1,r[0]/d[0]);
106  }
107  else if ( r.size() == 2 ) {
108  u.resize(2,0.0);
109  Real det = d[0]*d[1] - dl[0]*du[0];
110  u[0] = (d[1]*r[0] - du[0]*r[1])/det;
111  u[1] = (d[0]*r[1] - dl[0]*r[0])/det;
112  }
113  else {
114  u.assign(r.begin(),r.end());
115  // Perform LDL factorization
116  Teuchos::LAPACK<int,Real> lp;
117  std::vector<Real> du2(r.size()-2,0.0);
118  std::vector<int> ipiv(r.size(),0);
119  int info;
120  int dim = r.size();
121  int ldb = r.size();
122  int nhrs = 1;
123  lp.GTTRF(dim,&dl[0],&d[0],&du[0],&du2[0],&ipiv[0],&info);
124  char trans = 'N';
125  if ( transpose ) {
126  trans = 'T';
127  }
128  lp.GTTRS(trans,dim,nhrs,&dl[0],&d[0],&du[0],&du2[0],&ipiv[0],&u[0],ldb,&info);
129  }
130  }
131 
132 public:
133  BurgersFEM(int nx = 128, Real nl = 1.0, Real cH1 = 1.0, Real cL2 = 1.0)
134  : nx_(nx), dx_(1.0/((Real)nx+1.0)), nl_(nl), cH1_(cH1), cL2_(cL2) {}
135 
136  void set_problem_data(const Real nu, const Real f, const Real u0, const Real u1) {
137  nu_ = std::pow(10.0,nu-2.0);
138  f_ = 0.01*f;
139  u0_ = 1.0+0.001*u0;
140  u1_ = 0.001*u1;
141  }
142 
143  Real get_viscosity(void) const {
144  return nu_;
145  }
146 
147  int num_dof(void) const {
148  return nx_;
149  }
150 
151  Real mesh_spacing(void) const {
152  return dx_;
153  }
154 
155  /***************************************************************************/
156  /*********************** L2 INNER PRODUCT FUNCTIONS ************************/
157  /***************************************************************************/
158  // Compute L2 inner product
159  Real compute_L2_dot(const std::vector<Real> &x, const std::vector<Real> &y) const {
160  Real ip = 0.0;
161  Real c = (((int)x.size()==nx_) ? 4.0 : 2.0);
162  for (unsigned i=0; i<x.size(); i++) {
163  if ( i == 0 ) {
164  ip += dx_/6.0*(c*x[i] + x[i+1])*y[i];
165  }
166  else if ( i == x.size()-1 ) {
167  ip += dx_/6.0*(x[i-1] + c*x[i])*y[i];
168  }
169  else {
170  ip += dx_/6.0*(x[i-1] + 4.0*x[i] + x[i+1])*y[i];
171  }
172  }
173  return ip;
174  }
175 
176  // compute L2 norm
177  Real compute_L2_norm(const std::vector<Real> &r) const {
178  return std::sqrt(compute_L2_dot(r,r));
179  }
180 
181  // Apply L2 Reisz operator
182  void apply_mass(std::vector<Real> &Mu, const std::vector<Real> &u ) const {
183  Mu.resize(u.size(),0.0);
184  Real c = (((int)u.size()==nx_) ? 4.0 : 2.0);
185  for (unsigned i=0; i<u.size(); i++) {
186  if ( i == 0 ) {
187  Mu[i] = dx_/6.0*(c*u[i] + u[i+1]);
188  }
189  else if ( i == u.size()-1 ) {
190  Mu[i] = dx_/6.0*(u[i-1] + c*u[i]);
191  }
192  else {
193  Mu[i] = dx_/6.0*(u[i-1] + 4.0*u[i] + u[i+1]);
194  }
195  }
196  }
197 
198  // Apply L2 inverse Reisz operator
199  void apply_inverse_mass(std::vector<Real> &Mu, const std::vector<Real> &u) const {
200  unsigned nx = u.size();
201  // Build mass matrix
202  std::vector<Real> dl(nx-1,dx_/6.0);
203  std::vector<Real> d(nx,2.0*dx_/3.0);
204  std::vector<Real> du(nx-1,dx_/6.0);
205  if ( (int)nx != nx_ ) {
206  d[ 0] = dx_/3.0;
207  d[nx-1] = dx_/3.0;
208  }
209  linear_solve(Mu,dl,d,du,u);
210  }
211 
212  void test_inverse_mass(std::ostream &outStream = std::cout) {
213  // State Mass Matrix
214  std::vector<Real> u(nx_,0.0), Mu(nx_,0.0), iMMu(nx_,0.0), diff(nx_,0.0);
215  for (int i = 0; i < nx_; i++) {
216  u[i] = 2.0*(Real)rand()/(Real)RAND_MAX - 1.0;
217  }
218  apply_mass(Mu,u);
219  apply_inverse_mass(iMMu,Mu);
220  axpy(diff,-1.0,iMMu,u);
221  Real error = compute_L2_norm(diff);
222  Real normu = compute_L2_norm(u);
223  outStream << "Test Inverse State Mass Matrix\n";
224  outStream << " ||u - inv(M)Mu|| = " << error << "\n";
225  outStream << " ||u|| = " << normu << "\n";
226  outStream << " Relative Error = " << error/normu << "\n";
227  outStream << "\n";
228  // Control Mass Matrix
229  u.resize(nx_+2,0.0); Mu.resize(nx_+2,0.0); iMMu.resize(nx_+2,0.0); diff.resize(nx_+2,0.0);
230  for (int i = 0; i < nx_+2; i++) {
231  u[i] = 2.0*(Real)rand()/(Real)RAND_MAX - 1.0;
232  }
233  apply_mass(Mu,u);
234  apply_inverse_mass(iMMu,Mu);
235  axpy(diff,-1.0,iMMu,u);
236  error = compute_L2_norm(diff);
237  normu = compute_L2_norm(u);
238  outStream << "Test Inverse Control Mass Matrix\n";
239  outStream << " ||z - inv(M)Mz|| = " << error << "\n";
240  outStream << " ||z|| = " << normu << "\n";
241  outStream << " Relative Error = " << error/normu << "\n";
242  outStream << "\n";
243  }
244 
245  /***************************************************************************/
246  /*********************** H1 INNER PRODUCT FUNCTIONS ************************/
247  /***************************************************************************/
248  // Compute H1 inner product
249  Real compute_H1_dot(const std::vector<Real> &x, const std::vector<Real> &y) const {
250  Real ip = 0.0;
251  for (int i=0; i<nx_; i++) {
252  if ( i == 0 ) {
253  ip += cL2_*dx_/6.0*(4.0*x[i] + x[i+1])*y[i]; // Mass term
254  ip += cH1_*(2.0*x[i] - x[i+1])/dx_*y[i]; // Stiffness term
255  }
256  else if ( i == nx_-1 ) {
257  ip += cL2_*dx_/6.0*(x[i-1] + 4.0*x[i])*y[i]; // Mass term
258  ip += cH1_*(2.0*x[i] - x[i-1])/dx_*y[i]; // Stiffness term
259  }
260  else {
261  ip += cL2_*dx_/6.0*(x[i-1] + 4.0*x[i] + x[i+1])*y[i]; // Mass term
262  ip += cH1_*(2.0*x[i] - x[i-1] - x[i+1])/dx_*y[i]; // Stiffness term
263  }
264  }
265  return ip;
266  }
267 
268  // compute H1 norm
269  Real compute_H1_norm(const std::vector<Real> &r) const {
270  return std::sqrt(compute_H1_dot(r,r));
271  }
272 
273  // Apply H2 Reisz operator
274  void apply_H1(std::vector<Real> &Mu, const std::vector<Real> &u ) const {
275  Mu.resize(nx_,0.0);
276  for (int i=0; i<nx_; i++) {
277  if ( i == 0 ) {
278  Mu[i] = cL2_*dx_/6.0*(4.0*u[i] + u[i+1])
279  + cH1_*(2.0*u[i] - u[i+1])/dx_;
280  }
281  else if ( i == nx_-1 ) {
282  Mu[i] = cL2_*dx_/6.0*(u[i-1] + 4.0*u[i])
283  + cH1_*(2.0*u[i] - u[i-1])/dx_;
284  }
285  else {
286  Mu[i] = cL2_*dx_/6.0*(u[i-1] + 4.0*u[i] + u[i+1])
287  + cH1_*(2.0*u[i] - u[i-1] - u[i+1])/dx_;
288  }
289  }
290  }
291 
292  // Apply H1 inverse Reisz operator
293  void apply_inverse_H1(std::vector<Real> &Mu, const std::vector<Real> &u) const {
294  // Build mass matrix
295  std::vector<Real> dl(nx_-1,cL2_*dx_/6.0 - cH1_/dx_);
296  std::vector<Real> d(nx_,2.0*(cL2_*dx_/3.0 + cH1_/dx_));
297  std::vector<Real> du(nx_-1,cL2_*dx_/6.0 - cH1_/dx_);
298  linear_solve(Mu,dl,d,du,u);
299  }
300 
301  void test_inverse_H1(std::ostream &outStream = std::cout) {
302  std::vector<Real> u(nx_,0.0), Mu(nx_,0.0), iMMu(nx_,0.0), diff(nx_,0.0);
303  for (int i = 0; i < nx_; i++) {
304  u[i] = 2.0*(Real)rand()/(Real)RAND_MAX - 1.0;
305  }
306  apply_H1(Mu,u);
307  apply_inverse_H1(iMMu,Mu);
308  axpy(diff,-1.0,iMMu,u);
309  Real error = compute_H1_norm(diff);
310  Real normu = compute_H1_norm(u);
311  outStream << "Test Inverse State H1 Matrix\n";
312  outStream << " ||u - inv(M)Mu|| = " << error << "\n";
313  outStream << " ||u|| = " << normu << "\n";
314  outStream << " Relative Error = " << error/normu << "\n";
315  outStream << "\n";
316  }
317 
318  /***************************************************************************/
319  /*********************** PDE RESIDUAL AND SOLVE ****************************/
320  /***************************************************************************/
321  // Compute current PDE residual
322  void compute_residual(std::vector<Real> &r, const std::vector<Real> &u,
323  const std::vector<Real> &z) const {
324  r.clear();
325  r.resize(nx_,0.0);
326  for (int i=0; i<nx_; i++) {
327  // Contribution from stiffness term
328  if ( i==0 ) {
329  r[i] = nu_/dx_*(2.0*u[i]-u[i+1]);
330  }
331  else if (i==nx_-1) {
332  r[i] = nu_/dx_*(2.0*u[i]-u[i-1]);
333  }
334  else {
335  r[i] = nu_/dx_*(2.0*u[i]-u[i-1]-u[i+1]);
336  }
337  // Contribution from nonlinear term
338  if (i<nx_-1){
339  r[i] += nl_*u[i+1]*(u[i]+u[i+1])/6.0;
340  }
341  if (i>0) {
342  r[i] -= nl_*u[i-1]*(u[i-1]+u[i])/6.0;
343  }
344  // Contribution from control
345  r[i] -= dx_/6.0*(z[i]+4.0*z[i+1]+z[i+2]);
346  // Contribution from right-hand side
347  r[i] -= dx_*f_;
348  }
349  // Contribution from Dirichlet boundary terms
350  r[0] -= nl_*(u0_*u[ 0]/6.0 + u0_*u0_/6.0) + nu_*u0_/dx_;
351  r[nx_-1] += nl_*(u1_*u[nx_-1]/6.0 + u1_*u1_/6.0) - nu_*u1_/dx_;
352  }
353 
354  // Solve PDE
355  void solve(std::vector<Real> &u, const std::vector<Real> &z) const {
356  u.assign(nx_,1.0);
357  // Compute residual and residual norm
358  std::vector<Real> r(u.size(),0.0);
359  compute_residual(r,u,z);
360  Real rnorm = compute_L2_norm(r);
361  // Define tolerances
362  Real rtol = 1.e2*ROL::ROL_EPSILON;
363  Real maxit = 500;
364  // Initialize Jacobian storage
365  std::vector<Real> d(nx_,0.0);
366  std::vector<Real> dl(nx_-1,0.0);
367  std::vector<Real> du(nx_-1,0.0);
368  // Iterate Newton's method
369  Real alpha = 1.0, tmp = 0.0;
370  std::vector<Real> s(nx_,0.0);
371  std::vector<Real> utmp(nx_,0.0);
372  for (int i=0; i<maxit; i++) {
373  //std::cout << i << " " << rnorm << "\n";
374  // Get Jacobian
375  compute_pde_jacobian(dl,d,du,u);
376  // Solve Newton system
377  linear_solve(s,dl,d,du,r);
378  // Perform line search
379  tmp = rnorm;
380  alpha = 1.0;
381  utmp.assign(u.begin(),u.end());
382  update(utmp,s,-alpha);
383  compute_residual(r,utmp,z);
384  rnorm = compute_L2_norm(r);
385  while ( rnorm > (1.0-1.e-4*alpha)*tmp && alpha > std::sqrt(ROL::ROL_EPSILON) ) {
386  alpha /= 2.0;
387  utmp.assign(u.begin(),u.end());
388  update(utmp,s,-alpha);
389  compute_residual(r,utmp,z);
390  rnorm = compute_L2_norm(r);
391  }
392  // Update iterate
393  u.assign(utmp.begin(),utmp.end());
394  if ( rnorm < rtol ) {
395  break;
396  }
397  }
398  }
399 
400  /***************************************************************************/
401  /*********************** PDE JACOBIAN FUNCTIONS ****************************/
402  /***************************************************************************/
403  // Build PDE Jacobian trigiagonal matrix
404  void compute_pde_jacobian(std::vector<Real> &dl, // Lower diagonal
405  std::vector<Real> &d, // Diagonal
406  std::vector<Real> &du, // Upper diagonal
407  const std::vector<Real> &u) const { // State variable
408  // Get Diagonal and Off-Diagonal Entries of linear PDE Jacobian
409  d.clear();
410  d.resize(nx_,nu_*2.0/dx_);
411  dl.clear();
412  dl.resize(nx_-1,-nu_/dx_);
413  du.clear();
414  du.resize(nx_-1,-nu_/dx_);
415  // Contribution from nonlinearity
416  for (int i=0; i<nx_; i++) {
417  if (i<nx_-1) {
418  dl[i] += nl_*(-2.0*u[i]-u[i+1])/6.0;
419  d[i] += nl_*u[i+1]/6.0;
420  }
421  if (i>0) {
422  d[i] -= nl_*u[i-1]/6.0;
423  du[i-1] += nl_*(u[i-1]+2.0*u[i])/6.0;
424  }
425  }
426  // Contribution from Dirichlet boundary conditions
427  d[ 0] -= nl_*u0_/6.0;
428  d[nx_-1] += nl_*u1_/6.0;
429  }
430 
431  // Apply PDE Jacobian to a vector
432  void apply_pde_jacobian(std::vector<Real> &jv,
433  const std::vector<Real> &v,
434  const std::vector<Real> &u,
435  const std::vector<Real> &z) const {
436  // Fill jv
437  for (int i = 0; i < nx_; i++) {
438  jv[i] = nu_/dx_*2.0*v[i];
439  if ( i > 0 ) {
440  jv[i] += -nu_/dx_*v[i-1]-nl_*(u[i-1]/6.0*v[i]+(u[i]+2.0*u[i-1])/6.0*v[i-1]);
441  }
442  if ( i < nx_-1 ) {
443  jv[i] += -nu_/dx_*v[i+1]+nl_*(u[i+1]/6.0*v[i]+(u[i]+2.0*u[i+1])/6.0*v[i+1]);
444  }
445  }
446  jv[ 0] -= nl_*u0_/6.0*v[0];
447  jv[nx_-1] += nl_*u1_/6.0*v[nx_-1];
448  }
449 
450  // Apply inverse PDE Jacobian to a vector
451  void apply_inverse_pde_jacobian(std::vector<Real> &ijv,
452  const std::vector<Real> &v,
453  const std::vector<Real> &u,
454  const std::vector<Real> &z) const {
455  // Get PDE Jacobian
456  std::vector<Real> d(nx_,0.0);
457  std::vector<Real> dl(nx_-1,0.0);
458  std::vector<Real> du(nx_-1,0.0);
459  compute_pde_jacobian(dl,d,du,u);
460  // Solve solve state sensitivity system at current time step
461  linear_solve(ijv,dl,d,du,v);
462  }
463 
464  // Apply adjoint PDE Jacobian to a vector
465  void apply_adjoint_pde_jacobian(std::vector<Real> &ajv,
466  const std::vector<Real> &v,
467  const std::vector<Real> &u,
468  const std::vector<Real> &z) const {
469  // Fill jvp
470  for (int i = 0; i < nx_; i++) {
471  ajv[i] = nu_/dx_*2.0*v[i];
472  if ( i > 0 ) {
473  ajv[i] += -nu_/dx_*v[i-1]-nl_*(u[i-1]/6.0*v[i]
474  -(u[i-1]+2.0*u[i])/6.0*v[i-1]);
475  }
476  if ( i < nx_-1 ) {
477  ajv[i] += -nu_/dx_*v[i+1]+nl_*(u[i+1]/6.0*v[i]
478  -(u[i+1]+2.0*u[i])/6.0*v[i+1]);
479  }
480  }
481  ajv[ 0] -= nl_*u0_/6.0*v[0];
482  ajv[nx_-1] += nl_*u1_/6.0*v[nx_-1];
483  }
484 
485  // Apply inverse adjoint PDE Jacobian to a vector
486  void apply_inverse_adjoint_pde_jacobian(std::vector<Real> &iajv,
487  const std::vector<Real> &v,
488  const std::vector<Real> &u,
489  const std::vector<Real> &z) const {
490  // Get PDE Jacobian
491  std::vector<Real> d(nx_,0.0);
492  std::vector<Real> du(nx_-1,0.0);
493  std::vector<Real> dl(nx_-1,0.0);
494  compute_pde_jacobian(dl,d,du,u);
495  // Solve solve adjoint system at current time step
496  linear_solve(iajv,dl,d,du,v,true);
497  }
498 
499  /***************************************************************************/
500  /*********************** CONTROL JACOBIAN FUNCTIONS ************************/
501  /***************************************************************************/
502  // Apply control Jacobian to a vector
503  void apply_control_jacobian(std::vector<Real> &jv,
504  const std::vector<Real> &v,
505  const std::vector<Real> &u,
506  const std::vector<Real> &z) const {
507  for (int i=0; i<nx_; i++) {
508  // Contribution from control
509  jv[i] = -dx_/6.0*(v[i]+4.0*v[i+1]+v[i+2]);
510  }
511  }
512 
513  // Apply adjoint control Jacobian to a vector
514  void apply_adjoint_control_jacobian(std::vector<Real> &jv,
515  const std::vector<Real> &v,
516  const std::vector<Real> &u,
517  const std::vector<Real> &z) const {
518  for (int i=0; i<nx_+2; i++) {
519  if ( i == 0 ) {
520  jv[i] = -dx_/6.0*v[i];
521  }
522  else if ( i == 1 ) {
523  jv[i] = -dx_/6.0*(4.0*v[i-1]+v[i]);
524  }
525  else if ( i == nx_ ) {
526  jv[i] = -dx_/6.0*(4.0*v[i-1]+v[i-2]);
527  }
528  else if ( i == nx_+1 ) {
529  jv[i] = -dx_/6.0*v[i-2];
530  }
531  else {
532  jv[i] = -dx_/6.0*(v[i-2]+4.0*v[i-1]+v[i]);
533  }
534  }
535  }
536 
537  /***************************************************************************/
538  /*********************** AJDOINT HESSIANS **********************************/
539  /***************************************************************************/
540  void apply_adjoint_pde_hessian(std::vector<Real> &ahwv,
541  const std::vector<Real> &w,
542  const std::vector<Real> &v,
543  const std::vector<Real> &u,
544  const std::vector<Real> &z) const {
545  for (int i=0; i<nx_; i++) {
546  // Contribution from nonlinear term
547  ahwv[i] = 0.0;
548  if (i<nx_-1){
549  ahwv[i] += (w[i]*v[i+1] - w[i+1]*(2.0*v[i]+v[i+1]))/6.0;
550  }
551  if (i>0) {
552  ahwv[i] += (w[i-1]*(v[i-1]+2.0*v[i]) - w[i]*v[i-1])/6.0;
553  }
554  }
555  //ahwv.assign(u.size(),0.0);
556  }
557  void apply_adjoint_pde_control_hessian(std::vector<Real> &ahwv,
558  const std::vector<Real> &w,
559  const std::vector<Real> &v,
560  const std::vector<Real> &u,
561  const std::vector<Real> &z) {
562  ahwv.assign(u.size(),0.0);
563  }
564  void apply_adjoint_control_pde_hessian(std::vector<Real> &ahwv,
565  const std::vector<Real> &w,
566  const std::vector<Real> &v,
567  const std::vector<Real> &u,
568  const std::vector<Real> &z) {
569  ahwv.assign(z.size(),0.0);
570  }
571  void apply_adjoint_control_hessian(std::vector<Real> &ahwv,
572  const std::vector<Real> &w,
573  const std::vector<Real> &v,
574  const std::vector<Real> &u,
575  const std::vector<Real> &z) {
576  ahwv.assign(z.size(),0.0);
577  }
578 };
579 
580 template<class Real>
581 class L2VectorPrimal : public ROL::Vector<Real> {
582 private:
583  Teuchos::RCP<std::vector<Real> > vec_;
584  Teuchos::RCP<BurgersFEM<Real> > fem_;
585 
586  mutable Teuchos::RCP<L2VectorDual<Real> > dual_vec_;
587 
588 public:
589  L2VectorPrimal(const Teuchos::RCP<std::vector<Real> > & vec,
590  const Teuchos::RCP<BurgersFEM<Real> > &fem)
591  : vec_(vec), fem_(fem), dual_vec_(Teuchos::null) {}
592 
593  void set( const ROL::Vector<Real> &x ) {
594  const L2VectorPrimal &ex = Teuchos::dyn_cast<const L2VectorPrimal>(x);
595  const std::vector<Real>& xval = *ex.getVector();
596  std::copy(xval.begin(),xval.end(),vec_->begin());
597  }
598 
599  void plus( const ROL::Vector<Real> &x ) {
600  const L2VectorPrimal &ex = Teuchos::dyn_cast<const L2VectorPrimal>(x);
601  const std::vector<Real>& xval = *ex.getVector();
602  unsigned dimension = vec_->size();
603  for (unsigned i=0; i<dimension; i++) {
604  (*vec_)[i] += xval[i];
605  }
606  }
607 
608  void scale( const Real alpha ) {
609  unsigned dimension = vec_->size();
610  for (unsigned i=0; i<dimension; i++) {
611  (*vec_)[i] *= alpha;
612  }
613  }
614 
615  Real dot( const ROL::Vector<Real> &x ) const {
616  const L2VectorPrimal & ex = Teuchos::dyn_cast<const L2VectorPrimal>(x);
617  const std::vector<Real>& xval = *ex.getVector();
618  return fem_->compute_L2_dot(xval,*vec_);
619  }
620 
621  Real norm() const {
622  Real val = 0;
623  val = std::sqrt( dot(*this) );
624  return val;
625  }
626 
627  Teuchos::RCP<ROL::Vector<Real> > clone() const {
628  return Teuchos::rcp( new L2VectorPrimal( Teuchos::rcp(new std::vector<Real>(vec_->size(),0.0)),fem_));
629  }
630 
631  Teuchos::RCP<const std::vector<Real> > getVector() const {
632  return vec_;
633  }
634 
635  Teuchos::RCP<std::vector<Real> > getVector() {
636  return vec_;
637  }
638 
639  Teuchos::RCP<ROL::Vector<Real> > basis( const int i ) const {
640  Teuchos::RCP<L2VectorPrimal> e
641  = Teuchos::rcp( new L2VectorPrimal( Teuchos::rcp(new std::vector<Real>(vec_->size(),0.0)),fem_));
642  (*e->getVector())[i] = 1.0;
643  return e;
644  }
645 
646  int dimension() const {
647  return vec_->size();
648  }
649 
650  const ROL::Vector<Real>& dual() const {
651  dual_vec_ = Teuchos::rcp(new L2VectorDual<Real>(
652  Teuchos::rcp(new std::vector<Real>(*vec_)),fem_));
653 
654  fem_->apply_mass(*(Teuchos::rcp_const_cast<std::vector<Real> >(dual_vec_->getVector())),*vec_);
655  return *dual_vec_;
656  }
657 
658 };
659 
660 template<class Real>
661 class L2VectorDual : public ROL::Vector<Real> {
662 private:
663  Teuchos::RCP<std::vector<Real> > vec_;
664  Teuchos::RCP<BurgersFEM<Real> > fem_;
665 
666  mutable Teuchos::RCP<L2VectorPrimal<Real> > dual_vec_;
667 
668 public:
669  L2VectorDual(const Teuchos::RCP<std::vector<Real> > & vec,
670  const Teuchos::RCP<BurgersFEM<Real> > &fem)
671  : vec_(vec), fem_(fem), dual_vec_(Teuchos::null) {}
672 
673  void set( const ROL::Vector<Real> &x ) {
674  const L2VectorDual &ex = Teuchos::dyn_cast<const L2VectorDual>(x);
675  const std::vector<Real>& xval = *ex.getVector();
676  std::copy(xval.begin(),xval.end(),vec_->begin());
677  }
678 
679  void plus( const ROL::Vector<Real> &x ) {
680  const L2VectorDual &ex = Teuchos::dyn_cast<const L2VectorDual>(x);
681  const std::vector<Real>& xval = *ex.getVector();
682  unsigned dimension = vec_->size();
683  for (unsigned i=0; i<dimension; i++) {
684  (*vec_)[i] += xval[i];
685  }
686  }
687 
688  void scale( const Real alpha ) {
689  unsigned dimension = vec_->size();
690  for (unsigned i=0; i<dimension; i++) {
691  (*vec_)[i] *= alpha;
692  }
693  }
694 
695  Real dot( const ROL::Vector<Real> &x ) const {
696  const L2VectorDual & ex = Teuchos::dyn_cast<const L2VectorDual>(x);
697  const std::vector<Real>& xval = *ex.getVector();
698  unsigned dimension = vec_->size();
699  std::vector<Real> Mx(dimension,0.0);
700  fem_->apply_inverse_mass(Mx,xval);
701  Real val = 0.0;
702  for (unsigned i = 0; i < dimension; i++) {
703  val += Mx[i]*(*vec_)[i];
704  }
705  return val;
706  }
707 
708  Real norm() const {
709  Real val = 0;
710  val = std::sqrt( dot(*this) );
711  return val;
712  }
713 
714  Teuchos::RCP<ROL::Vector<Real> > clone() const {
715  return Teuchos::rcp( new L2VectorDual( Teuchos::rcp(new std::vector<Real>(vec_->size(),0.0)),fem_));
716  }
717 
718  Teuchos::RCP<const std::vector<Real> > getVector() const {
719  return vec_;
720  }
721 
722  Teuchos::RCP<std::vector<Real> > getVector() {
723  return vec_;
724  }
725 
726  Teuchos::RCP<ROL::Vector<Real> > basis( const int i ) const {
727  Teuchos::RCP<L2VectorDual> e
728  = Teuchos::rcp( new L2VectorDual( Teuchos::rcp(new std::vector<Real>(vec_->size(),0.0)),fem_));
729  (*e->getVector())[i] = 1.0;
730  return e;
731  }
732 
733  int dimension() const {
734  return vec_->size();
735  }
736 
737  const ROL::Vector<Real>& dual() const {
738  dual_vec_ = Teuchos::rcp(new L2VectorPrimal<Real>(
739  Teuchos::rcp(new std::vector<Real>(*vec_)),fem_));
740 
741  fem_->apply_inverse_mass(*(Teuchos::rcp_const_cast<std::vector<Real> >(dual_vec_->getVector())),*vec_);
742  return *dual_vec_;
743  }
744 
745 };
746 
747 template<class Real>
748 class H1VectorPrimal : public ROL::Vector<Real> {
749 private:
750  Teuchos::RCP<std::vector<Real> > vec_;
751  Teuchos::RCP<BurgersFEM<Real> > fem_;
752 
753  mutable Teuchos::RCP<H1VectorDual<Real> > dual_vec_;
754 
755 public:
756  H1VectorPrimal(const Teuchos::RCP<std::vector<Real> > & vec,
757  const Teuchos::RCP<BurgersFEM<Real> > &fem)
758  : vec_(vec), fem_(fem), dual_vec_(Teuchos::null) {}
759 
760  void set( const ROL::Vector<Real> &x ) {
761  const H1VectorPrimal &ex = Teuchos::dyn_cast<const H1VectorPrimal>(x);
762  const std::vector<Real>& xval = *ex.getVector();
763  std::copy(xval.begin(),xval.end(),vec_->begin());
764  }
765 
766  void plus( const ROL::Vector<Real> &x ) {
767  const H1VectorPrimal &ex = Teuchos::dyn_cast<const H1VectorPrimal>(x);
768  const std::vector<Real>& xval = *ex.getVector();
769  unsigned dimension = vec_->size();
770  for (unsigned i=0; i<dimension; i++) {
771  (*vec_)[i] += xval[i];
772  }
773  }
774 
775  void scale( const Real alpha ) {
776  unsigned dimension = vec_->size();
777  for (unsigned i=0; i<dimension; i++) {
778  (*vec_)[i] *= alpha;
779  }
780  }
781 
782  Real dot( const ROL::Vector<Real> &x ) const {
783  const H1VectorPrimal & ex = Teuchos::dyn_cast<const H1VectorPrimal>(x);
784  const std::vector<Real>& xval = *ex.getVector();
785  return fem_->compute_H1_dot(xval,*vec_);
786  }
787 
788  Real norm() const {
789  Real val = 0;
790  val = std::sqrt( dot(*this) );
791  return val;
792  }
793 
794  Teuchos::RCP<ROL::Vector<Real> > clone() const {
795  return Teuchos::rcp( new H1VectorPrimal( Teuchos::rcp(new std::vector<Real>(vec_->size(),0.0)),fem_));
796  }
797 
798  Teuchos::RCP<const std::vector<Real> > getVector() const {
799  return vec_;
800  }
801 
802  Teuchos::RCP<std::vector<Real> > getVector() {
803  return vec_;
804  }
805 
806  Teuchos::RCP<ROL::Vector<Real> > basis( const int i ) const {
807  Teuchos::RCP<H1VectorPrimal> e
808  = Teuchos::rcp( new H1VectorPrimal( Teuchos::rcp(new std::vector<Real>(vec_->size(),0.0)),fem_));
809  (*e->getVector())[i] = 1.0;
810  return e;
811  }
812 
813  int dimension() const {
814  return vec_->size();
815  }
816 
817  const ROL::Vector<Real>& dual() const {
818  dual_vec_ = Teuchos::rcp(new H1VectorDual<Real>(
819  Teuchos::rcp(new std::vector<Real>(*vec_)),fem_));
820 
821  fem_->apply_H1(*(Teuchos::rcp_const_cast<std::vector<Real> >(dual_vec_->getVector())),*vec_);
822  return *dual_vec_;
823  }
824 
825 };
826 
827 template<class Real>
828 class H1VectorDual : public ROL::Vector<Real> {
829 private:
830  Teuchos::RCP<std::vector<Real> > vec_;
831  Teuchos::RCP<BurgersFEM<Real> > fem_;
832 
833  mutable Teuchos::RCP<H1VectorPrimal<Real> > dual_vec_;
834 
835 public:
836  H1VectorDual(const Teuchos::RCP<std::vector<Real> > & vec,
837  const Teuchos::RCP<BurgersFEM<Real> > &fem)
838  : vec_(vec), fem_(fem), dual_vec_(Teuchos::null) {}
839 
840  void set( const ROL::Vector<Real> &x ) {
841  const H1VectorDual &ex = Teuchos::dyn_cast<const H1VectorDual>(x);
842  const std::vector<Real>& xval = *ex.getVector();
843  std::copy(xval.begin(),xval.end(),vec_->begin());
844  }
845 
846  void plus( const ROL::Vector<Real> &x ) {
847  const H1VectorDual &ex = Teuchos::dyn_cast<const H1VectorDual>(x);
848  const std::vector<Real>& xval = *ex.getVector();
849  unsigned dimension = vec_->size();
850  for (unsigned i=0; i<dimension; i++) {
851  (*vec_)[i] += xval[i];
852  }
853  }
854 
855  void scale( const Real alpha ) {
856  unsigned dimension = vec_->size();
857  for (unsigned i=0; i<dimension; i++) {
858  (*vec_)[i] *= alpha;
859  }
860  }
861 
862  Real dot( const ROL::Vector<Real> &x ) const {
863  const H1VectorDual & ex = Teuchos::dyn_cast<const H1VectorDual>(x);
864  const std::vector<Real>& xval = *ex.getVector();
865  unsigned dimension = vec_->size();
866  std::vector<Real> Mx(dimension,0.0);
867  fem_->apply_inverse_H1(Mx,xval);
868  Real val = 0.0;
869  for (unsigned i = 0; i < dimension; i++) {
870  val += Mx[i]*(*vec_)[i];
871  }
872  return val;
873  }
874 
875  Real norm() const {
876  Real val = 0;
877  val = std::sqrt( dot(*this) );
878  return val;
879  }
880 
881  Teuchos::RCP<ROL::Vector<Real> > clone() const {
882  return Teuchos::rcp( new H1VectorDual( Teuchos::rcp(new std::vector<Real>(vec_->size(),0.0)),fem_));
883  }
884 
885  Teuchos::RCP<const std::vector<Real> > getVector() const {
886  return vec_;
887  }
888 
889  Teuchos::RCP<std::vector<Real> > getVector() {
890  return vec_;
891  }
892 
893  Teuchos::RCP<ROL::Vector<Real> > basis( const int i ) const {
894  Teuchos::RCP<H1VectorDual> e
895  = Teuchos::rcp( new H1VectorDual( Teuchos::rcp(new std::vector<Real>(vec_->size(),0.0)),fem_));
896  (*e->getVector())[i] = 1.0;
897  return e;
898  }
899 
900  int dimension() const {
901  return vec_->size();
902  }
903 
904  const ROL::Vector<Real>& dual() const {
905  dual_vec_ = Teuchos::rcp(new H1VectorPrimal<Real>(
906  Teuchos::rcp(new std::vector<Real>(*vec_)),fem_));
907 
908  fem_->apply_inverse_H1(*(Teuchos::rcp_const_cast<std::vector<Real> >(dual_vec_->getVector())),*vec_);
909  return *dual_vec_;
910  }
911 
912 };
913 
914 template<class Real>
916 private:
917 
920 
923 
926 
927  typedef typename std::vector<Real>::size_type uint;
928 
929  Teuchos::RCP<BurgersFEM<Real> > fem_;
930  bool useHessian_;
931 
932 public:
933  EqualityConstraint_BurgersControl(Teuchos::RCP<BurgersFEM<Real> > &fem, bool useHessian = true)
934  : fem_(fem), useHessian_(useHessian) {}
935 
937  const ROL::Vector<Real> &z, Real &tol) {
938  Teuchos::RCP<std::vector<Real> > cp =
939  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<PrimalConstraintVector>(c)).getVector());
940  Teuchos::RCP<const std::vector<Real> > up =
941  (Teuchos::dyn_cast<PrimalStateVector>(const_cast<ROL::Vector<Real> &>(u))).getVector();
942  Teuchos::RCP<const std::vector<Real> > zp =
943  (Teuchos::dyn_cast<PrimalControlVector>(const_cast<ROL::Vector<Real> &>(z))).getVector();
944 
945  const std::vector<Real> param
947  fem_->set_problem_data(param[0],param[1],param[2],param[3]);
948 
949  fem_->compute_residual(*cp,*up,*zp);
950  }
951 
952  void solve(ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
953  Teuchos::RCP<std::vector<Real> > up =
954  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<PrimalStateVector>(u)).getVector());
955  up->assign(up->size(),z.norm()/up->size());
956  Teuchos::RCP<const std::vector<Real> > zp =
957  (Teuchos::dyn_cast<PrimalControlVector>(const_cast<ROL::Vector<Real> &>(z))).getVector();
958 
959  const std::vector<Real> param
961  fem_->set_problem_data(param[0],param[1],param[2],param[3]);
962 
963  fem_->solve(*up,*zp);
964  }
965 
967  const ROL::Vector<Real> &z, Real &tol) {
968  Teuchos::RCP<std::vector<Real> > jvp =
969  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<PrimalConstraintVector>(jv)).getVector());
970  Teuchos::RCP<const std::vector<Real> > vp =
971  (Teuchos::dyn_cast<PrimalStateVector>(const_cast<ROL::Vector<Real> &>(v))).getVector();
972  Teuchos::RCP<const std::vector<Real> > up =
973  (Teuchos::dyn_cast<PrimalStateVector>(const_cast<ROL::Vector<Real> &>(u))).getVector();
974  Teuchos::RCP<const std::vector<Real> > zp =
975  (Teuchos::dyn_cast<PrimalControlVector>(const_cast<ROL::Vector<Real> &>(z))).getVector();
976 
977  const std::vector<Real> param
979  fem_->set_problem_data(param[0],param[1],param[2],param[3]);
980 
981  fem_->apply_pde_jacobian(*jvp,*vp,*up,*zp);
982  }
983 
985  const ROL::Vector<Real> &z, Real &tol) {
986  Teuchos::RCP<std::vector<Real> > jvp =
987  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<PrimalConstraintVector>(jv)).getVector());
988  Teuchos::RCP<const std::vector<Real> > vp =
989  (Teuchos::dyn_cast<PrimalControlVector>(const_cast<ROL::Vector<Real> &>(v))).getVector();
990  Teuchos::RCP<const std::vector<Real> > up =
991  (Teuchos::dyn_cast<PrimalStateVector>(const_cast<ROL::Vector<Real> &>(u))).getVector();
992  Teuchos::RCP<const std::vector<Real> > zp =
993  (Teuchos::dyn_cast<PrimalControlVector>(const_cast<ROL::Vector<Real> &>(z))).getVector();
994 
995  const std::vector<Real> param
997  fem_->set_problem_data(param[0],param[1],param[2],param[3]);
998 
999  fem_->apply_control_jacobian(*jvp,*vp,*up,*zp);
1000  }
1001 
1003  const ROL::Vector<Real> &z, Real &tol) {
1004  Teuchos::RCP<std::vector<Real> > ijvp =
1005  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<PrimalStateVector>(ijv)).getVector());
1006  Teuchos::RCP<const std::vector<Real> > vp =
1007  (Teuchos::dyn_cast<PrimalConstraintVector>(const_cast<ROL::Vector<Real> &>(v))).getVector();
1008  Teuchos::RCP<const std::vector<Real> > up =
1009  (Teuchos::dyn_cast<PrimalStateVector>(const_cast<ROL::Vector<Real> &>(u))).getVector();
1010  Teuchos::RCP<const std::vector<Real> > zp =
1011  (Teuchos::dyn_cast<PrimalControlVector>(const_cast<ROL::Vector<Real> &>(z))).getVector();
1012 
1013  const std::vector<Real> param
1015  fem_->set_problem_data(param[0],param[1],param[2],param[3]);
1016 
1017  fem_->apply_inverse_pde_jacobian(*ijvp,*vp,*up,*zp);
1018  }
1019 
1021  const ROL::Vector<Real> &z, Real &tol) {
1022  Teuchos::RCP<std::vector<Real> > jvp =
1023  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<DualStateVector>(ajv)).getVector());
1024  Teuchos::RCP<const std::vector<Real> > vp =
1025  (Teuchos::dyn_cast<DualConstraintVector>(const_cast<ROL::Vector<Real> &>(v))).getVector();
1026  Teuchos::RCP<const std::vector<Real> > up =
1027  (Teuchos::dyn_cast<PrimalStateVector>(const_cast<ROL::Vector<Real> &>(u))).getVector();
1028  Teuchos::RCP<const std::vector<Real> > zp =
1029  (Teuchos::dyn_cast<PrimalControlVector>(const_cast<ROL::Vector<Real> &>(z))).getVector();
1030 
1031  const std::vector<Real> param
1033  fem_->set_problem_data(param[0],param[1],param[2],param[3]);
1034 
1035  fem_->apply_adjoint_pde_jacobian(*jvp,*vp,*up,*zp);
1036  }
1037 
1039  const ROL::Vector<Real> &z, Real &tol) {
1040  Teuchos::RCP<std::vector<Real> > jvp =
1041  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<DualControlVector>(jv)).getVector());
1042  Teuchos::RCP<const std::vector<Real> > vp =
1043  (Teuchos::dyn_cast<DualConstraintVector>(const_cast<ROL::Vector<Real> &>(v))).getVector();
1044  Teuchos::RCP<const std::vector<Real> > up =
1045  (Teuchos::dyn_cast<PrimalStateVector>(const_cast<ROL::Vector<Real> &>(u))).getVector();
1046  Teuchos::RCP<const std::vector<Real> > zp =
1047  (Teuchos::dyn_cast<PrimalControlVector>(const_cast<ROL::Vector<Real> &>(z))).getVector();
1048 
1049  const std::vector<Real> param
1051  fem_->set_problem_data(param[0],param[1],param[2],param[3]);
1052 
1053  fem_->apply_adjoint_control_jacobian(*jvp,*vp,*up,*zp);
1054  }
1055 
1057  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
1058  Teuchos::RCP<std::vector<Real> > iajvp =
1059  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<DualConstraintVector>(iajv)).getVector());
1060  Teuchos::RCP<const std::vector<Real> > vp =
1061  (Teuchos::dyn_cast<DualStateVector>(const_cast<ROL::Vector<Real> &>(v))).getVector();
1062  Teuchos::RCP<const std::vector<Real> > up =
1063  (Teuchos::dyn_cast<PrimalStateVector>(const_cast<ROL::Vector<Real> &>(u))).getVector();
1064  Teuchos::RCP<const std::vector<Real> > zp =
1065  (Teuchos::dyn_cast<PrimalControlVector>(const_cast<ROL::Vector<Real> &>(z))).getVector();
1066 
1067  const std::vector<Real> param
1069  fem_->set_problem_data(param[0],param[1],param[2],param[3]);
1070 
1071  fem_->apply_inverse_adjoint_pde_jacobian(*iajvp,*vp,*up,*zp);
1072  }
1073 
1075  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
1076  if ( useHessian_ ) {
1077  Teuchos::RCP<std::vector<Real> > ahwvp =
1078  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<DualStateVector>(ahwv)).getVector());
1079  Teuchos::RCP<const std::vector<Real> > wp =
1080  (Teuchos::dyn_cast<DualConstraintVector>(const_cast<ROL::Vector<Real> &>(w))).getVector();
1081  Teuchos::RCP<const std::vector<Real> > vp =
1082  (Teuchos::dyn_cast<PrimalStateVector>(const_cast<ROL::Vector<Real> &>(v))).getVector();
1083  Teuchos::RCP<const std::vector<Real> > up =
1084  (Teuchos::dyn_cast<PrimalStateVector>(const_cast<ROL::Vector<Real> &>(u))).getVector();
1085  Teuchos::RCP<const std::vector<Real> > zp =
1086  (Teuchos::dyn_cast<PrimalControlVector>(const_cast<ROL::Vector<Real> &>(z))).getVector();
1087 
1088  const std::vector<Real> param
1090  fem_->set_problem_data(param[0],param[1],param[2],param[3]);
1091 
1092  fem_->apply_adjoint_pde_hessian(*ahwvp,*wp,*vp,*up,*zp);
1093  }
1094  else {
1095  ahwv.zero();
1096  }
1097  }
1098 
1100  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
1101  if ( useHessian_ ) {
1102  Teuchos::RCP<std::vector<Real> > ahwvp =
1103  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<DualControlVector>(ahwv)).getVector());
1104  Teuchos::RCP<const std::vector<Real> > wp =
1105  (Teuchos::dyn_cast<DualConstraintVector>(const_cast<ROL::Vector<Real> &>(w))).getVector();
1106  Teuchos::RCP<const std::vector<Real> > vp =
1107  (Teuchos::dyn_cast<PrimalStateVector>(const_cast<ROL::Vector<Real> &>(v))).getVector();
1108  Teuchos::RCP<const std::vector<Real> > up =
1109  (Teuchos::dyn_cast<PrimalStateVector>(const_cast<ROL::Vector<Real> &>(u))).getVector();
1110  Teuchos::RCP<const std::vector<Real> > zp =
1111  (Teuchos::dyn_cast<PrimalControlVector>(const_cast<ROL::Vector<Real> &>(z))).getVector();
1112 
1113  const std::vector<Real> param
1115  fem_->set_problem_data(param[0],param[1],param[2],param[3]);
1116 
1117  fem_->apply_adjoint_control_pde_hessian(*ahwvp,*wp,*vp,*up,*zp);
1118  }
1119  else {
1120  ahwv.zero();
1121  }
1122  }
1124  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
1125  if ( useHessian_ ) {
1126  Teuchos::RCP<std::vector<Real> > ahwvp =
1127  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<DualStateVector>(ahwv)).getVector());
1128  Teuchos::RCP<const std::vector<Real> > wp =
1129  (Teuchos::dyn_cast<DualConstraintVector>(const_cast<ROL::Vector<Real> &>(w))).getVector();
1130  Teuchos::RCP<const std::vector<Real> > vp =
1131  (Teuchos::dyn_cast<PrimalControlVector>(const_cast<ROL::Vector<Real> &>(v))).getVector();
1132  Teuchos::RCP<const std::vector<Real> > up =
1133  (Teuchos::dyn_cast<PrimalStateVector>(const_cast<ROL::Vector<Real> &>(u))).getVector();
1134  Teuchos::RCP<const std::vector<Real> > zp =
1135  (Teuchos::dyn_cast<PrimalControlVector>(const_cast<ROL::Vector<Real> &>(z))).getVector();
1136 
1137  const std::vector<Real> param
1139  fem_->set_problem_data(param[0],param[1],param[2],param[3]);
1140 
1141  fem_->apply_adjoint_pde_control_hessian(*ahwvp,*wp,*vp,*up,*zp);
1142  }
1143  else {
1144  ahwv.zero();
1145  }
1146  }
1148  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol) {
1149  if ( useHessian_ ) {
1150  Teuchos::RCP<std::vector<Real> > ahwvp =
1151  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<DualControlVector>(ahwv)).getVector());
1152  Teuchos::RCP<const std::vector<Real> > wp =
1153  (Teuchos::dyn_cast<DualConstraintVector>(const_cast<ROL::Vector<Real> &>(w))).getVector();
1154  Teuchos::RCP<const std::vector<Real> > vp =
1155  (Teuchos::dyn_cast<PrimalControlVector>(const_cast<ROL::Vector<Real> &>(v))).getVector();
1156  Teuchos::RCP<const std::vector<Real> > up =
1157  (Teuchos::dyn_cast<PrimalStateVector>(const_cast<ROL::Vector<Real> &>(u))).getVector();
1158  Teuchos::RCP<const std::vector<Real> > zp =
1159  (Teuchos::dyn_cast<PrimalControlVector>(const_cast<ROL::Vector<Real> &>(z))).getVector();
1160 
1161  const std::vector<Real> param
1163  fem_->set_problem_data(param[0],param[1],param[2],param[3]);
1164 
1165  fem_->apply_adjoint_control_hessian(*ahwvp,*wp,*vp,*up,*zp);
1166  }
1167  else {
1168  ahwv.zero();
1169  }
1170  }
1171 };
1172 
1173 template<class Real>
1175 private:
1176 
1179 
1182 
1183  Teuchos::RCP<BurgersFEM<Real> > fem_;
1184 
1185  Real x_;
1186  std::vector<int> indices_;
1187 
1188 public:
1189  Objective_BurgersControl(const Teuchos::RCP<BurgersFEM<Real> > &fem,
1190  Real x = 0.0) : fem_(fem), x_(x) {
1191  for (int i = 1; i < fem_->num_dof()+1; i++) {
1192  if ( (Real)i*(fem_->mesh_spacing()) >= x_ ) {
1193  indices_.push_back(i-1);
1194  }
1195  }
1196  }
1197 
1198  Real value( const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
1199  Teuchos::RCP<const std::vector<Real> > up =
1200  (Teuchos::dyn_cast<PrimalStateVector>(const_cast<ROL::Vector<Real> &>(u))).getVector();
1201 
1202 // const std::vector<Real> param
1203 // = ROL::ParametrizedObjective_SimOpt<Real>::getParameter();
1204 // fem_->set_problem_data(param[0],param[1],param[2],param[3]);
1205 // Real nu = fem_->get_viscosity();
1206 //
1207 // return 0.5*nu*fem_->compute_H1_dot(*up,*up);
1208 
1209  Real val = 0.5*((((Real)indices_[0]+1.)*(fem_->mesh_spacing())-x_)
1210  *(x_+(2.-((Real)indices_[0]+1.))*(fem_->mesh_spacing()))/(fem_->mesh_spacing())
1211  +(fem_->mesh_spacing())) * (*up)[indices_[0]];
1212  for (uint i = 1; i < indices_.size(); i++) {
1213  val += (fem_->mesh_spacing())*(*up)[indices_[i]];
1214  }
1215  return -val;
1216  }
1217 
1218  void gradient_1( ROL::Vector<Real> &g, const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
1219  Teuchos::RCP<std::vector<Real> > gp =
1220  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<DualStateVector>(g)).getVector());
1221  Teuchos::RCP<const std::vector<Real> > up =
1222  (Teuchos::dyn_cast<PrimalStateVector>(const_cast<ROL::Vector<Real> &>(u))).getVector();
1223 
1224 // const std::vector<Real> param
1225 // = ROL::ParametrizedObjective_SimOpt<Real>::getParameter();
1226 // fem_->set_problem_data(param[0],param[1],param[2],param[3]);
1227 // Real nu = fem_->get_viscosity();
1228 //
1229 // fem_->apply_H1(*gp,*up);
1230 // g.scale(nu);
1231 
1232  g.zero();
1233  (*gp)[indices_[0]] = -0.5*((((Real)indices_[0]+1.)*(fem_->mesh_spacing())-x_)
1234  *(x_+(2.-((Real)indices_[0]+1.))*(fem_->mesh_spacing()))/(fem_->mesh_spacing())
1235  +(fem_->mesh_spacing()));
1236 
1237 
1238  for (uint i = 1; i < indices_.size(); i++) {
1239  (*gp)[indices_[i]] = -(fem_->mesh_spacing());
1240  }
1241  }
1242 
1243  void gradient_2( ROL::Vector<Real> &g, const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
1244  g.zero();
1245  }
1246 
1248  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
1249 // Teuchos::RCP<std::vector<Real> > hvp =
1250 // Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<DualStateVector>(hv)).getVector());
1251 // Teuchos::RCP<const std::vector<Real> > vp =
1252 // (Teuchos::dyn_cast<PrimalStateVector>(const_cast<ROL::Vector<Real> &>(v))).getVector();
1253 //
1254 // const std::vector<Real> param
1255 // = ROL::ParametrizedObjective_SimOpt<Real>::getParameter();
1256 // fem_->set_problem_data(param[0],param[1],param[2],param[3]);
1257 // Real nu = fem_->get_viscosity();
1258 //
1259 // fem_->apply_H1(*hvp,*vp);
1260 // hv.scale(nu);
1261 
1262  hv.zero();
1263  }
1264 
1266  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
1267  hv.zero();
1268  }
1269 
1271  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
1272  hv.zero();
1273  }
1274 
1276  const ROL::Vector<Real> &u, const ROL::Vector<Real> &z, Real &tol ) {
1277  hv.zero();
1278  }
1279 };
1280 
1281 template<class Real>
1282 class L2BoundConstraint : public ROL::BoundConstraint<Real> {
1283 private:
1284  int dim_;
1285  std::vector<Real> x_lo_;
1286  std::vector<Real> x_up_;
1287  Real min_diff_;
1288  Real scale_;
1289  Teuchos::RCP<BurgersFEM<Real> > fem_;
1290 
1291  void cast_vector(Teuchos::RCP<std::vector<Real> > &xvec,
1292  ROL::Vector<Real> &x) const {
1293  try {
1294  xvec = Teuchos::rcp_const_cast<std::vector<Real> >(
1295  (Teuchos::dyn_cast<L2VectorPrimal<Real> >(x)).getVector());
1296  }
1297  catch (std::exception &e) {
1298  xvec = Teuchos::rcp_const_cast<std::vector<Real> >(
1299  (Teuchos::dyn_cast<L2VectorDual<Real> >(x)).getVector());
1300  }
1301  }
1302 
1303  void cast_const_vector(Teuchos::RCP<const std::vector<Real> > &xvec,
1304  const ROL::Vector<Real> &x) const {
1305  try {
1306  xvec = (Teuchos::dyn_cast<L2VectorPrimal<Real> >(
1307  const_cast<ROL::Vector<Real> &>(x))).getVector();
1308  }
1309  catch (std::exception &e) {
1310  xvec = (Teuchos::dyn_cast<L2VectorDual<Real> >(
1311  const_cast<ROL::Vector<Real> &>(x))).getVector();
1312  }
1313  }
1314 
1315  void axpy(std::vector<Real> &out, const Real a,
1316  const std::vector<Real> &x, const std::vector<Real> &y) const{
1317  out.resize(dim_,0.0);
1318  for (unsigned i = 0; i < dim_; i++) {
1319  out[i] = a*x[i] + y[i];
1320  }
1321  }
1322 
1323  void projection(std::vector<Real> &x) {
1324  for ( int i = 0; i < dim_; i++ ) {
1325  x[i] = std::max(x_lo_[i],std::min(x_up_[i],x[i]));
1326  }
1327  }
1328 
1329 public:
1330  L2BoundConstraint(std::vector<Real> &l, std::vector<Real> &u,
1331  const Teuchos::RCP<BurgersFEM<Real> > &fem, Real scale = 1.0)
1332  : x_lo_(l), x_up_(u), scale_(scale), fem_(fem) {
1333  dim_ = x_lo_.size();
1334  for ( int i = 0; i < dim_; i++ ) {
1335  if ( i == 0 ) {
1336  min_diff_ = x_up_[i] - x_lo_[i];
1337  }
1338  else {
1339  min_diff_ = ( (min_diff_ < (x_up_[i] - x_lo_[i])) ? min_diff_ : (x_up_[i] - x_lo_[i]) );
1340  }
1341  }
1342  min_diff_ *= 0.5;
1343  }
1344 
1345  bool isFeasible( const ROL::Vector<Real> &x ) {
1346  Teuchos::RCP<const std::vector<Real> > ex; cast_const_vector(ex,x);
1347  bool val = true;
1348  int cnt = 1;
1349  for ( int i = 0; i < dim_; i++ ) {
1350  if ( (*ex)[i] >= x_lo_[i] && (*ex)[i] <= x_up_[i] ) { cnt *= 1; }
1351  else { cnt *= 0; }
1352  }
1353  if ( cnt == 0 ) { val = false; }
1354  return val;
1355  }
1356 
1358  Teuchos::RCP<std::vector<Real> > ex; cast_vector(ex,x);
1359  projection(*ex);
1360  }
1361 
1363  Teuchos::RCP<const std::vector<Real> > ex; cast_const_vector(ex,x);
1364  Teuchos::RCP<std::vector<Real> > ev; cast_vector(ev,v);
1365  Real epsn = std::min(scale_*eps,min_diff_);
1366  for ( int i = 0; i < dim_; i++ ) {
1367  if ( ((*ex)[i] <= x_lo_[i]+epsn) ) {
1368  (*ev)[i] = 0.0;
1369  }
1370  }
1371  }
1372 
1374  Teuchos::RCP<const std::vector<Real> > ex; cast_const_vector(ex,x);
1375  Teuchos::RCP<std::vector<Real> > ev; cast_vector(ev,v);
1376  Real epsn = std::min(scale_*eps,min_diff_);
1377  for ( int i = 0; i < dim_; i++ ) {
1378  if ( ((*ex)[i] >= x_up_[i]-epsn) ) {
1379  (*ev)[i] = 0.0;
1380  }
1381  }
1382  }
1383 
1384  void pruneActive(ROL::Vector<Real> &v, const ROL::Vector<Real> &x, Real eps) {
1385  Teuchos::RCP<const std::vector<Real> > ex; cast_const_vector(ex,x);
1386  Teuchos::RCP<std::vector<Real> > ev; cast_vector(ev,v);
1387  Real epsn = std::min(scale_*eps,min_diff_);
1388  for ( int i = 0; i < dim_; i++ ) {
1389  if ( ((*ex)[i] <= x_lo_[i]+epsn) ||
1390  ((*ex)[i] >= x_up_[i]-epsn) ) {
1391  (*ev)[i] = 0.0;
1392  }
1393  }
1394  }
1395 
1397  Teuchos::RCP<const std::vector<Real> > ex; cast_const_vector(ex,x);
1398  Teuchos::RCP<const std::vector<Real> > eg; cast_const_vector(eg,g);
1399  Teuchos::RCP<std::vector<Real> > ev; cast_vector(ev,v);
1400  Real epsn = std::min(scale_*eps,min_diff_);
1401  for ( int i = 0; i < dim_; i++ ) {
1402  if ( ((*ex)[i] <= x_lo_[i]+epsn && (*eg)[i] > 0.0) ) {
1403  (*ev)[i] = 0.0;
1404  }
1405  }
1406  }
1407 
1409  Teuchos::RCP<const std::vector<Real> > ex; cast_const_vector(ex,x);
1410  Teuchos::RCP<const std::vector<Real> > eg; cast_const_vector(eg,g);
1411  Teuchos::RCP<std::vector<Real> > ev; cast_vector(ev,v);
1412  Real epsn = std::min(scale_*eps,min_diff_);
1413  for ( int i = 0; i < dim_; i++ ) {
1414  if ( ((*ex)[i] >= x_up_[i]-epsn && (*eg)[i] < 0.0) ) {
1415  (*ev)[i] = 0.0;
1416  }
1417  }
1418  }
1419 
1420  void pruneActive(ROL::Vector<Real> &v, const ROL::Vector<Real> &g, const ROL::Vector<Real> &x, Real eps) {
1421  Teuchos::RCP<const std::vector<Real> > ex; cast_const_vector(ex,x);
1422  Teuchos::RCP<const std::vector<Real> > eg; cast_const_vector(eg,g);
1423  Teuchos::RCP<std::vector<Real> > ev; cast_vector(ev,v);
1424  Real epsn = std::min(scale_*eps,min_diff_);
1425  for ( int i = 0; i < dim_; i++ ) {
1426  if ( ((*ex)[i] <= x_lo_[i]+epsn && (*eg)[i] > 0.0) ||
1427  ((*ex)[i] >= x_up_[i]-epsn && (*eg)[i] < 0.0) ) {
1428  (*ev)[i] = 0.0;
1429  }
1430  }
1431  }
1432 
1434  Teuchos::RCP<std::vector<Real> > us = Teuchos::rcp( new std::vector<Real>(dim_,0.0) );
1435  us->assign(x_up_.begin(),x_up_.end());
1436  Teuchos::RCP<ROL::Vector<Real> > up = Teuchos::rcp( new L2VectorPrimal<Real>(us,fem_) );
1437  u.set(*up);
1438  }
1439 
1441  Teuchos::RCP<std::vector<Real> > ls = Teuchos::rcp( new std::vector<Real>(dim_,0.0) );
1442  ls->assign(x_lo_.begin(),x_lo_.end());
1443  Teuchos::RCP<ROL::Vector<Real> > lp = Teuchos::rcp( new L2VectorPrimal<Real>(ls,fem_) );
1444  l.set(*lp);
1445  }
1446 };
1447 
1448 template<class Real>
1449 class H1BoundConstraint : public ROL::BoundConstraint<Real> {
1450 private:
1451  int dim_;
1452  std::vector<Real> x_lo_;
1453  std::vector<Real> x_up_;
1454  Real min_diff_;
1455  Real scale_;
1456  Teuchos::RCP<BurgersFEM<Real> > fem_;
1457 
1458  void cast_vector(Teuchos::RCP<std::vector<Real> > &xvec,
1459  ROL::Vector<Real> &x) const {
1460  try {
1461  xvec = Teuchos::rcp_const_cast<std::vector<Real> >(
1462  (Teuchos::dyn_cast<H1VectorPrimal<Real> >(x)).getVector());
1463  }
1464  catch (std::exception &e) {
1465  xvec = Teuchos::rcp_const_cast<std::vector<Real> >(
1466  (Teuchos::dyn_cast<H1VectorDual<Real> >(x)).getVector());
1467  }
1468  }
1469 
1470  void cast_const_vector(Teuchos::RCP<const std::vector<Real> > &xvec,
1471  const ROL::Vector<Real> &x) const {
1472  try {
1473  xvec = (Teuchos::dyn_cast<H1VectorPrimal<Real> >(
1474  const_cast<ROL::Vector<Real> &>(x))).getVector();
1475  }
1476  catch (std::exception &e) {
1477  xvec = (Teuchos::dyn_cast<H1VectorDual<Real> >(
1478  const_cast<ROL::Vector<Real> &>(x))).getVector();
1479  }
1480  }
1481 
1482  void axpy(std::vector<Real> &out, const Real a,
1483  const std::vector<Real> &x, const std::vector<Real> &y) const{
1484  out.resize(dim_,0.0);
1485  for (unsigned i = 0; i < dim_; i++) {
1486  out[i] = a*x[i] + y[i];
1487  }
1488  }
1489 
1490  void projection(std::vector<Real> &x) {
1491  for ( int i = 0; i < dim_; i++ ) {
1492  x[i] = std::max(x_lo_[i],std::min(x_up_[i],x[i]));
1493  }
1494  }
1495 
1496 public:
1497  H1BoundConstraint(std::vector<Real> &l, std::vector<Real> &u,
1498  const Teuchos::RCP<BurgersFEM<Real> > &fem, Real scale = 1.0)
1499  : x_lo_(l), x_up_(u), scale_(scale), fem_(fem) {
1500  dim_ = x_lo_.size();
1501  for ( int i = 0; i < dim_; i++ ) {
1502  if ( i == 0 ) {
1503  min_diff_ = x_up_[i] - x_lo_[i];
1504  }
1505  else {
1506  min_diff_ = ( (min_diff_ < (x_up_[i] - x_lo_[i])) ? min_diff_ : (x_up_[i] - x_lo_[i]) );
1507  }
1508  }
1509  min_diff_ *= 0.5;
1510  }
1511 
1512  bool isFeasible( const ROL::Vector<Real> &x ) {
1513  Teuchos::RCP<const std::vector<Real> > ex; cast_const_vector(ex,x);
1514  bool val = true;
1515  int cnt = 1;
1516  for ( int i = 0; i < dim_; i++ ) {
1517  if ( (*ex)[i] >= x_lo_[i] && (*ex)[i] <= x_up_[i] ) { cnt *= 1; }
1518  else { cnt *= 0; }
1519  }
1520  if ( cnt == 0 ) { val = false; }
1521  return val;
1522  }
1523 
1525  Teuchos::RCP<std::vector<Real> > ex; cast_vector(ex,x);
1526  projection(*ex);
1527  }
1528 
1530  Teuchos::RCP<const std::vector<Real> > ex; cast_const_vector(ex,x);
1531  Teuchos::RCP<std::vector<Real> > ev; cast_vector(ev,v);
1532  Real epsn = std::min(scale_*eps,min_diff_);
1533  for ( int i = 0; i < dim_; i++ ) {
1534  if ( ((*ex)[i] <= x_lo_[i]+epsn) ) {
1535  (*ev)[i] = 0.0;
1536  }
1537  }
1538  }
1539 
1541  Teuchos::RCP<const std::vector<Real> > ex; cast_const_vector(ex,x);
1542  Teuchos::RCP<std::vector<Real> > ev; cast_vector(ev,v);
1543  Real epsn = std::min(scale_*eps,min_diff_);
1544  for ( int i = 0; i < dim_; i++ ) {
1545  if ( ((*ex)[i] >= x_up_[i]-epsn) ) {
1546  (*ev)[i] = 0.0;
1547  }
1548  }
1549  }
1550 
1551  void pruneActive(ROL::Vector<Real> &v, const ROL::Vector<Real> &x, Real eps) {
1552  Teuchos::RCP<const std::vector<Real> > ex; cast_const_vector(ex,x);
1553  Teuchos::RCP<std::vector<Real> > ev; cast_vector(ev,v);
1554  Real epsn = std::min(scale_*eps,min_diff_);
1555  for ( int i = 0; i < dim_; i++ ) {
1556  if ( ((*ex)[i] <= x_lo_[i]+epsn) ||
1557  ((*ex)[i] >= x_up_[i]-epsn) ) {
1558  (*ev)[i] = 0.0;
1559  }
1560  }
1561  }
1562 
1564  Teuchos::RCP<const std::vector<Real> > ex; cast_const_vector(ex,x);
1565  Teuchos::RCP<const std::vector<Real> > eg; cast_const_vector(eg,g);
1566  Teuchos::RCP<std::vector<Real> > ev; cast_vector(ev,v);
1567  Real epsn = std::min(scale_*eps,min_diff_);
1568  for ( int i = 0; i < dim_; i++ ) {
1569  if ( ((*ex)[i] <= x_lo_[i]+epsn && (*eg)[i] > 0.0) ) {
1570  (*ev)[i] = 0.0;
1571  }
1572  }
1573  }
1574 
1576  Teuchos::RCP<const std::vector<Real> > ex; cast_const_vector(ex,x);
1577  Teuchos::RCP<const std::vector<Real> > eg; cast_const_vector(eg,g);
1578  Teuchos::RCP<std::vector<Real> > ev; cast_vector(ev,v);
1579  Real epsn = std::min(scale_*eps,min_diff_);
1580  for ( int i = 0; i < dim_; i++ ) {
1581  if ( ((*ex)[i] >= x_up_[i]-epsn && (*eg)[i] < 0.0) ) {
1582  (*ev)[i] = 0.0;
1583  }
1584  }
1585  }
1586 
1587  void pruneActive(ROL::Vector<Real> &v, const ROL::Vector<Real> &g, const ROL::Vector<Real> &x, Real eps) {
1588  Teuchos::RCP<const std::vector<Real> > ex; cast_const_vector(ex,x);
1589  Teuchos::RCP<const std::vector<Real> > eg; cast_const_vector(eg,g);
1590  Teuchos::RCP<std::vector<Real> > ev; cast_vector(ev,v);
1591  Real epsn = std::min(scale_*eps,min_diff_);
1592  for ( int i = 0; i < dim_; i++ ) {
1593  if ( ((*ex)[i] <= x_lo_[i]+epsn && (*eg)[i] > 0.0) ||
1594  ((*ex)[i] >= x_up_[i]-epsn && (*eg)[i] < 0.0) ) {
1595  (*ev)[i] = 0.0;
1596  }
1597  }
1598  }
1599 
1601  Teuchos::RCP<std::vector<Real> > us = Teuchos::rcp( new std::vector<Real>(dim_,0.0) );
1602  us->assign(x_up_.begin(),x_up_.end());
1603  Teuchos::RCP<ROL::Vector<Real> > up = Teuchos::rcp( new H1VectorPrimal<Real>(us,fem_) );
1604  u.set(*up);
1605  }
1606 
1608  Teuchos::RCP<std::vector<Real> > ls = Teuchos::rcp( new std::vector<Real>(dim_,0.0) );
1609  ls->assign(x_lo_.begin(),x_lo_.end());
1610  Teuchos::RCP<ROL::Vector<Real> > lp = Teuchos::rcp( new H1VectorPrimal<Real>(ls,fem_) );
1611  l.set(*lp);
1612  }
1613 };
1614 
1615 template<class Real, class Ordinal>
1616 class L2VectorBatchManager : public ROL::TeuchosBatchManager<Real,Ordinal> {
1617 private:
1618  void cast_vector(Teuchos::RCP<std::vector<Real> > &xvec,
1619  ROL::Vector<Real> &x) const {
1620  try {
1621  xvec = Teuchos::rcp_const_cast<std::vector<Real> >(
1622  (Teuchos::dyn_cast<L2VectorPrimal<Real> >(x)).getVector());
1623  }
1624  catch (std::exception &e) {
1625  xvec = Teuchos::rcp_const_cast<std::vector<Real> >(
1626  (Teuchos::dyn_cast<L2VectorDual<Real> >(x)).getVector());
1627  }
1628  }
1629 
1630 public:
1631  L2VectorBatchManager(const Teuchos::RCP<const Teuchos::Comm<Ordinal> > &comm)
1632  : ROL::TeuchosBatchManager<Real,Ordinal>(comm) {}
1634  Teuchos::RCP<std::vector<Real> > input_ptr;
1635  cast_vector(input_ptr,input);
1636  int dim_i = input_ptr->size();
1637  Teuchos::RCP<std::vector<Real> > output_ptr;
1638  cast_vector(output_ptr,output);
1639  int dim_o = output_ptr->size();
1640  if ( dim_i != dim_o ) {
1641  std::cout << "L2VectorBatchManager: DIMENSION MISMATCH ON RANK "
1643  }
1644  else {
1645  ROL::TeuchosBatchManager<Real,Ordinal>::sumAll(&(*input_ptr)[0],&(*output_ptr)[0],dim_i);
1646  }
1647  }
1648 };
1649 
1650 template<class Real, class Ordinal>
1651 class H1VectorBatchManager : public ROL::TeuchosBatchManager<Real,Ordinal> {
1652 private:
1653  void cast_vector(Teuchos::RCP<std::vector<Real> > &xvec,
1654  ROL::Vector<Real> &x) const {
1655  try {
1656  xvec = Teuchos::rcp_const_cast<std::vector<Real> >(
1657  (Teuchos::dyn_cast<H1VectorPrimal<Real> >(x)).getVector());
1658  }
1659  catch (std::exception &e) {
1660  xvec = Teuchos::rcp_const_cast<std::vector<Real> >(
1661  (Teuchos::dyn_cast<H1VectorDual<Real> >(x)).getVector());
1662  }
1663  }
1664 
1665 public:
1666  H1VectorBatchManager(const Teuchos::RCP<const Teuchos::Comm<Ordinal> > &comm)
1667  : ROL::TeuchosBatchManager<Real,Ordinal>(comm) {}
1669  Teuchos::RCP<std::vector<Real> > input_ptr;
1670  cast_vector(input_ptr,input);
1671  int dim_i = input_ptr->size();
1672  Teuchos::RCP<std::vector<Real> > output_ptr;
1673  cast_vector(output_ptr,output);
1674  int dim_o = output_ptr->size();
1675  if ( dim_i != dim_o ) {
1676  std::cout << "H1VectorBatchManager: DIMENSION MISMATCH ON RANK "
1678  }
1679  else {
1680  ROL::TeuchosBatchManager<Real,Ordinal>::sumAll(&(*input_ptr)[0],&(*output_ptr)[0],dim_i);
1681  }
1682  }
1683 };
1684 
1685 template<class Real>
1686 Real random(const Teuchos::RCP<const Teuchos::Comm<int> > &comm) {
1687  Real val = 0.0;
1688  if ( Teuchos::rank<int>(*comm)==0 ) {
1689  val = (Real)rand()/(Real)RAND_MAX;
1690  }
1691  Teuchos::broadcast<int,Real>(*comm,0,1,&val);
1692  return val;
1693 }
BurgersFEM(int nx=128, Real nl=1.0, Real cH1=1.0, Real cL2=1.0)
Definition: example_07.hpp:133
L2VectorBatchManager(const Teuchos::RCP< const Teuchos::Comm< Ordinal > > &comm)
void pruneUpperActive(ROL::Vector< Real > &v, const ROL::Vector< Real > &x, Real eps)
Set variables to zero if they correspond to the upper -active set.
void cast_vector(Teuchos::RCP< std::vector< Real > > &xvec, ROL::Vector< Real > &x) const
L2VectorPrimal< Real > PrimalControlVector
Definition: example_07.hpp:921
Real compute_H1_dot(const std::vector< Real > &x, const std::vector< Real > &y) const
Definition: example_07.hpp:249
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 ...
Real norm() const
Returns where .
Definition: example_07.hpp:788
void cast_const_vector(Teuchos::RCP< const std::vector< Real > > &xvec, const ROL::Vector< Real > &x) const
void apply_mass(std::vector< Real > &Mu, const std::vector< Real > &u) const
Definition: example_07.hpp:182
Teuchos::RCP< const std::vector< Real > > getVector() const
Definition: example_04.hpp:709
bool isFeasible(const ROL::Vector< Real > &x)
Check if the vector, v, is feasible.
void axpy(std::vector< Real > &out, const Real a, const std::vector< Real > &x, const std::vector< Real > &y) const
Definition: example_07.hpp:90
Teuchos::RCP< std::vector< Real > > getVector()
Definition: example_07.hpp:722
Teuchos::RCP< ROL::Vector< Real > > basis(const int i) const
Return i-th basis vector.
Definition: example_07.hpp:639
Objective_BurgersControl(const Teuchos::RCP< BurgersFEM< Real > > &fem, Real x=0.0)
Real dot(const ROL::Vector< Real > &x) const
Compute where .
Definition: example_07.hpp:695
void apply_adjoint_pde_control_hessian(std::vector< Real > &ahwv, const std::vector< Real > &w, const std::vector< Real > &v, const std::vector< Real > &u, const std::vector< Real > &z)
Definition: example_07.hpp:557
Teuchos::RCP< const std::vector< Real > > getVector() const
Definition: example_04.hpp:622
int dimension() const
Return dimension of the vector space.
Definition: example_07.hpp:646
Real get_viscosity(void) const
Definition: example_07.hpp:143
void plus(const ROL::Vector< Real > &x)
Compute , where .
Definition: example_07.hpp:766
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 .
void apply_adjoint_pde_hessian(std::vector< Real > &ahwv, const std::vector< Real > &w, const std::vector< Real > &v, const std::vector< Real > &u, const std::vector< Real > &z) const
Definition: example_07.hpp:540
L2VectorDual(const Teuchos::RCP< std::vector< Real > > &vec, const Teuchos::RCP< BurgersFEM< Real > > &fem)
Definition: example_07.hpp:669
Teuchos::RCP< ROL::Vector< Real > > basis(const int i) const
Return i-th basis vector.
Definition: example_07.hpp:893
Real norm() const
Returns where .
Definition: example_07.hpp:621
void pruneActive(ROL::Vector< Real > &v, const ROL::Vector< Real > &g, const ROL::Vector< Real > &x, Real eps)
Set variables to zero if they correspond to the -binding set.
Contains definitions of custom data types in ROL.
Real dot(const ROL::Vector< Real > &x) const
Compute where .
Definition: example_07.hpp:862
void plus(const ROL::Vector< Real > &x)
Compute , where .
Definition: example_07.hpp:846
void pruneLowerActive(ROL::Vector< Real > &v, const ROL::Vector< Real > &x, Real eps)
Set variables to zero if they correspond to the lower -active set.
Teuchos::RCP< const std::vector< Real > > getVector() const
Definition: example_04.hpp:789
Real compute_H1_norm(const std::vector< Real > &r) const
Definition: example_07.hpp:269
void solve(std::vector< Real > &u, const std::vector< Real > &z) const
Definition: example_07.hpp:355
Real norm() const
Returns where .
Definition: example_07.hpp:708
void pruneUpperActive(ROL::Vector< Real > &v, const ROL::Vector< Real > &g, const ROL::Vector< Real > &x, Real eps)
Set variables to zero if they correspond to the upper -binding set.
std::vector< Real >::size_type uint
Definition: example_07.hpp:927
Teuchos::RCP< std::vector< Real > > getVector()
Definition: example_07.hpp:889
Real random(const Teuchos::RCP< const Teuchos::Comm< int > > &comm)
void project(ROL::Vector< Real > &x)
Project optimization variables onto the bounds.
Teuchos::RCP< std::vector< Real > > getVector()
Definition: example_07.hpp:802
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 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 ...
H1VectorDual< Real > DualStateVector
void scale(const Real alpha)
Compute where .
Definition: example_07.hpp:608
L2VectorDual< Real > DualControlVector
Definition: example_07.hpp:922
Real dot(const ROL::Vector< Real > &x) const
Compute where .
Definition: example_07.hpp:782
int num_dof(void) const
Definition: example_07.hpp:147
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 pruneActive(ROL::Vector< Real > &v, const ROL::Vector< Real > &x, Real eps)
Set variables to zero if they correspond to the -active set.
void plus(const ROL::Vector< Real > &x)
Compute , where .
Definition: example_07.hpp:599
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.
Real mesh_spacing(void) const
Definition: example_07.hpp:151
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 ...
bool isFeasible(const ROL::Vector< Real > &x)
Check if the vector, v, is feasible.
int dimension() const
Return dimension of the vector space.
Definition: example_07.hpp:813
void test_inverse_mass(std::ostream &outStream=std::cout)
Definition: example_07.hpp:212
void apply_pde_jacobian(std::vector< Real > &jv, const std::vector< Real > &v, const std::vector< Real > &u, const std::vector< Real > &z) const
Definition: example_07.hpp:432
H1BoundConstraint(std::vector< Real > &l, std::vector< Real > &u, const Teuchos::RCP< BurgersFEM< Real > > &fem, Real scale=1.0)
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 pruneActive(ROL::Vector< Real > &v, const ROL::Vector< Real > &x, Real eps)
Set variables to zero if they correspond to the -active set.
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...
void pruneLowerActive(ROL::Vector< Real > &v, const ROL::Vector< Real > &g, const ROL::Vector< Real > &x, Real eps)
Set variables to zero if they correspond to the lower -binding set.
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_07.hpp:966
Teuchos::RCP< ROL::Vector< Real > > clone() const
Clone to make a new (uninitialized) vector.
Definition: example_07.hpp:794
void scale(const Real alpha)
Compute where .
Definition: example_07.hpp:688
void compute_pde_jacobian(std::vector< Real > &dl, std::vector< Real > &d, std::vector< Real > &du, const std::vector< Real > &u) const
Definition: example_04.hpp:395
H1VectorPrimal< Real > PrimalStateVector
Definition: example_07.hpp:918
Teuchos::RCP< ROL::Vector< Real > > basis(const int i) const
Return i-th basis vector.
Definition: example_07.hpp:806
Real compute_L2_dot(const std::vector< Real > &x, const std::vector< Real > &y) const
Definition: example_07.hpp:159
void pruneLowerActive(ROL::Vector< Real > &v, const ROL::Vector< Real > &x, Real eps)
Set variables to zero if they correspond to the lower -active set.
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 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) const
Definition: example_07.hpp:102
void cast_vector(Teuchos::RCP< std::vector< Real > > &xvec, ROL::Vector< Real > &x) const
void apply_inverse_mass(std::vector< Real > &Mu, const std::vector< Real > &u) const
Definition: example_07.hpp:199
void axpy(std::vector< Real > &out, const Real a, const std::vector< Real > &x, const std::vector< Real > &y) const
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 setVectorToUpperBound(ROL::Vector< Real > &u)
Set the input vector to the upper bound.
H1VectorDual< Real > PrimalConstraintVector
Definition: example_07.hpp:924
void update(std::vector< Real > &u, const std::vector< Real > &s, const Real alpha=1.0) const
Definition: example_07.hpp:84
void pruneLowerActive(ROL::Vector< Real > &v, const ROL::Vector< Real > &g, const ROL::Vector< Real > &x, Real eps)
Set variables to zero if they correspond to the lower -binding set.
void scale(const Real alpha)
Compute where .
Definition: example_07.hpp:775
void projection(std::vector< Real > &x)
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 cast_vector(Teuchos::RCP< std::vector< Real > > &xvec, ROL::Vector< Real > &x) const
H1VectorPrimal< Real > DualConstraintVector
Definition: example_07.hpp:925
void setVectorToLowerBound(ROL::Vector< Real > &l)
Set the input vector to the lower bound.
void apply_adjoint_control_hessian(std::vector< Real > &ahwv, const std::vector< Real > &w, const std::vector< Real > &v, const std::vector< Real > &u, const std::vector< Real > &z)
Definition: example_07.hpp:571
void solve(ROL::Vector< Real > &u, const ROL::Vector< Real > &z, Real &tol)
Given , solve for .
Definition: example_07.hpp:952
Provides the interface to apply upper and lower bound constraints.
void sumAll(Real *input, Real *output, int dim)
void compute_residual(std::vector< Real > &r, const std::vector< Real > &u, const std::vector< Real > &z) const
Definition: example_07.hpp:322
void projection(std::vector< Real > &x)
void cast_const_vector(Teuchos::RCP< const std::vector< Real > > &xvec, const ROL::Vector< Real > &x) const
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_07.hpp:984
H1VectorDual< Real > DualStateVector
Definition: example_07.hpp:919
void set_problem_data(const Real nu, const Real f, const Real u0, const Real u1)
Definition: example_07.hpp:136
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 apply_control_jacobian(std::vector< Real > &jv, const std::vector< Real > &v, const std::vector< Real > &u, const std::vector< Real > &z) const
Definition: example_07.hpp:503
void test_inverse_H1(std::ostream &outStream=std::cout)
Definition: example_07.hpp:301
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_07.hpp:936
const ROL::Vector< Real > & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
Definition: example_07.hpp:817
H1VectorDual(const Teuchos::RCP< std::vector< Real > > &vec, const Teuchos::RCP< BurgersFEM< Real > > &fem)
Definition: example_07.hpp:836
void apply_adjoint_control_jacobian(std::vector< Real > &jv, const std::vector< Real > &v, const std::vector< Real > &u, const std::vector< Real > &z) const
Definition: example_07.hpp:514
Real dot(const ROL::Vector< Real > &x) const
Compute where .
Definition: example_07.hpp:615
L2VectorPrimal< Real > PrimalControlVector
void pruneUpperActive(ROL::Vector< Real > &v, const ROL::Vector< Real > &g, const ROL::Vector< Real > &x, Real eps)
Set variables to zero if they correspond to the upper -binding set.
void pruneActive(ROL::Vector< Real > &v, const ROL::Vector< Real > &g, const ROL::Vector< Real > &x, Real eps)
Set variables to zero if they correspond to the -binding set.
H1VectorPrimal< Real > PrimalStateVector
L2VectorPrimal(const Teuchos::RCP< std::vector< Real > > &vec, const Teuchos::RCP< BurgersFEM< Real > > &fem)
Definition: example_07.hpp:589
void apply_inverse_pde_jacobian(std::vector< Real > &ijv, const std::vector< Real > &v, const std::vector< Real > &u, const std::vector< Real > &z) const
Definition: example_07.hpp:451
void axpy(std::vector< Real > &out, const Real a, const std::vector< Real > &x, const std::vector< Real > &y) const
void cast_vector(Teuchos::RCP< std::vector< Real > > &xvec, ROL::Vector< Real > &x) const
EqualityConstraint_BurgersControl(Teuchos::RCP< BurgersFEM< Real > > &fem, bool useHessian=true)
Definition: example_07.hpp:933
Teuchos::RCP< ROL::Vector< Real > > clone() const
Clone to make a new (uninitialized) vector.
Definition: example_07.hpp:627
std::vector< int > indices_
void apply_H1(std::vector< Real > &Mu, const std::vector< Real > &u) const
Definition: example_07.hpp:274
void plus(const ROL::Vector< Real > &x)
Compute , where .
Definition: example_07.hpp:679
virtual void set(const Vector &x)
Set where .
Definition: ROL_Vector.hpp:196
void sumAll(ROL::Vector< Real > &input, ROL::Vector< Real > &output)
void apply_inverse_H1(std::vector< Real > &Mu, const std::vector< Real > &u) const
Definition: example_07.hpp:293
void project(ROL::Vector< Real > &x)
Project optimization variables onto the bounds.
virtual Real norm() const =0
Returns where .
void apply_adjoint_pde_jacobian(std::vector< Real > &ajv, const std::vector< Real > &v, const std::vector< Real > &u, const std::vector< Real > &z) const
Definition: example_07.hpp:465
void apply_inverse_adjoint_pde_jacobian(std::vector< Real > &iajv, const std::vector< Real > &v, const std::vector< Real > &u, const std::vector< Real > &z) const
Definition: example_07.hpp:486
void apply_adjoint_control_pde_hessian(std::vector< Real > &ahwv, const std::vector< Real > &w, const std::vector< Real > &v, const std::vector< Real > &u, const std::vector< Real > &z)
Definition: example_07.hpp:564
int dimension() const
Return dimension of the vector space.
Definition: example_07.hpp:900
L2VectorDual< Real > DualControlVector
const ROL::Vector< Real > & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
Definition: example_07.hpp:904
Teuchos::RCP< ROL::Vector< Real > > basis(const int i) const
Return i-th basis vector.
Definition: example_07.hpp:726
H1VectorPrimal(const Teuchos::RCP< std::vector< Real > > &vec, const Teuchos::RCP< BurgersFEM< Real > > &fem)
Definition: example_07.hpp:756
void sumAll(ROL::Vector< Real > &input, ROL::Vector< Real > &output)
Real norm() const
Returns where .
Definition: example_07.hpp:875
void scale(std::vector< Real > &u, const Real alpha=0.0) const
Definition: example_07.hpp:96
Teuchos::RCP< std::vector< Real > > getVector()
Definition: example_07.hpp:635
const ROL::Vector< Real > & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
Definition: example_07.hpp:737
Teuchos::RCP< ROL::Vector< Real > > clone() const
Clone to make a new (uninitialized) vector.
Definition: example_07.hpp:881
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 ...
L2BoundConstraint(std::vector< Real > &l, std::vector< Real > &u, const Teuchos::RCP< BurgersFEM< Real > > &fem, Real scale=1.0)
H1VectorBatchManager(const Teuchos::RCP< const Teuchos::Comm< Ordinal > > &comm)
void pruneUpperActive(ROL::Vector< Real > &v, const ROL::Vector< Real > &x, Real eps)
Set variables to zero if they correspond to the upper -active set.
Teuchos::RCP< const std::vector< Real > > getVector() const
Definition: example_04.hpp:876
void setVectorToLowerBound(ROL::Vector< Real > &l)
Set the input vector to the lower bound.
void scale(const Real alpha)
Compute where .
Definition: example_07.hpp:855
int dimension() const
Return dimension of the vector space.
Definition: example_07.hpp:733
static const double ROL_EPSILON
Platform-dependent machine epsilon.
Definition: ROL_Types.hpp:118
const ROL::Vector< Real > & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
Definition: example_07.hpp:650
Real compute_L2_norm(const std::vector< Real > &r) const
Definition: example_07.hpp:177
Teuchos::RCP< ROL::Vector< Real > > clone() const
Clone to make a new (uninitialized) vector.
Definition: example_07.hpp:714
void setVectorToUpperBound(ROL::Vector< Real > &u)
Set the input vector to the upper bound.