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