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