ROL
burgers-control/example_01.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ************************************************************************
3 //
4 // Rapid Optimization Library (ROL) Package
5 // Copyright (2014) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact lead developers:
38 // Drew Kouri (dpkouri@sandia.gov) and
39 // Denis Ridzal (dridzal@sandia.gov)
40 //
41 // ************************************************************************
42 // @HEADER
43 
49 #include "ROL_Algorithm.hpp"
51 #include "ROL_TrustRegionStep.hpp"
52 #include "ROL_Types.hpp"
53 #include "ROL_Stream.hpp"
54 #include "Teuchos_GlobalMPISession.hpp"
55 #include "Teuchos_XMLParameterListHelpers.hpp"
56 #include "Teuchos_LAPACK.hpp"
57 
58 #include <iostream>
59 #include <fstream>
60 #include <algorithm>
61 
62 #include "ROL_StdVector.hpp"
63 #include "ROL_Objective.hpp"
64 #include "ROL_Bounds.hpp"
65 
66 template<class Real>
68 private:
69  Real alpha_; // Penalty Parameter
70 
71  int nx_; // Number of interior nodes
72  Real dx_; // Mesh spacing (i.e. 1/(nx+1))
73 
74 /***************************************************************/
75 /********** BEGIN PRIVATE MEMBER FUNCTION DECLARATION **********/
76 /***************************************************************/
77  Real evaluate_target(Real x) {
78  Real val = 0.0;
79  int example = 2;
80  switch (example) {
81  case 1: val = ((x<0.5) ? 1.0 : 0.0); break;
82  case 2: val = 1.0; break;
83  case 3: val = std::abs(std::sin(8.0*M_PI*x)); break;
84  case 4: val = std::exp(-0.5*(x-0.5)*(x-0.5)); break;
85  }
86  return val;
87  }
88 
89  Real compute_norm(const std::vector<Real> &r) {
90  return std::sqrt(this->dot(r,r));
91  }
92 
93  Real dot(const std::vector<Real> &x, const std::vector<Real> &y) {
94  Real ip = 0.0;
95  Real c = (((int)x.size()==this->nx_) ? 4.0 : 2.0);
96  for (unsigned i=0; i<x.size(); i++) {
97  if ( i == 0 ) {
98  ip += this->dx_/6.0*(c*x[i] + x[i+1])*y[i];
99  }
100  else if ( i == x.size()-1 ) {
101  ip += this->dx_/6.0*(x[i-1] + c*x[i])*y[i];
102  }
103  else {
104  ip += this->dx_/6.0*(x[i-1] + 4.0*x[i] + x[i+1])*y[i];
105  }
106  }
107  return ip;
108  }
109 
111 
112  void update(std::vector<Real> &u, const std::vector<Real> &s, const Real alpha=1.0) {
113  for (unsigned i=0; i<u.size(); i++) {
114  u[i] += alpha*s[i];
115  }
116  }
117 
118  void scale(std::vector<Real> &u, const Real alpha=0.0) {
119  for (unsigned i=0; i<u.size(); i++) {
120  u[i] *= alpha;
121  }
122  }
123 
124  void apply_mass(std::vector<Real> &Mu, const std::vector<Real> &u ) {
125  Mu.resize(u.size(),0.0);
126  Real c = (((int)u.size()==this->nx_) ? 4.0 : 2.0);
127  for (unsigned i=0; i<u.size(); i++) {
128  if ( i == 0 ) {
129  Mu[i] = this->dx_/6.0*(c*u[i] + u[i+1]);
130  }
131  else if ( i == u.size()-1 ) {
132  Mu[i] = this->dx_/6.0*(u[i-1] + c*u[i]);
133  }
134  else {
135  Mu[i] = this->dx_/6.0*(u[i-1] + 4.0*u[i] + u[i+1]);
136  }
137  }
138  }
139 
140  void compute_residual(std::vector<Real> &r, const std::vector<Real> &u,
141  const std::vector<Real> &z, const std::vector<Real> &param) {
142  r.clear();
143  r.resize(this->nx_,0.0);
144  Real nu = std::pow(10.0,param[0]-2.0);
145  Real f = param[1]/100.0;
146  Real u0 = 1.0+param[2]/1000.0;
147  Real u1 = param[3]/1000.0;
148  for (int i=0; i<this->nx_; i++) {
149  // Contribution from stiffness term
150  if ( i==0 ) {
151  r[i] = nu/this->dx_*(2.0*u[i]-u[i+1]);
152  }
153  else if (i==this->nx_-1) {
154  r[i] = nu/this->dx_*(2.0*u[i]-u[i-1]);
155  }
156  else {
157  r[i] = nu/this->dx_*(2.0*u[i]-u[i-1]-u[i+1]);
158  }
159  // Contribution from nonlinear term
160  if (i<this->nx_-1){
161  r[i] += u[i+1]*(u[i]+u[i+1])/6.0;
162  }
163  if (i>0) {
164  r[i] -= u[i-1]*(u[i-1]+u[i])/6.0;
165  }
166  // Contribution from control
167  r[i] -= this->dx_/6.0*(z[i]+4.0*z[i+1]+z[i+2]);
168  // Contribution from right-hand side
169  r[i] -= this->dx_*f;
170  }
171  // Contribution from Dirichlet boundary terms
172  r[0] -= u0*u[ 0]/6.0 + u0*u0/6.0 + nu*u0/this->dx_;
173  r[this->nx_-1] += u1*u[this->nx_-1]/6.0 + u1*u1/6.0 - nu*u1/this->dx_;
174  }
175 
176  void compute_pde_jacobian(std::vector<Real> &dl, std::vector<Real> &d, std::vector<Real> &du,
177  const std::vector<Real> &u, const std::vector<Real> &param) {
178  Real nu = std::pow(10.0,param[0]-2.0);
179  Real u0 = 1.0+param[2]/1000.0;
180  Real u1 = param[3]/1000.0;
181  // Get Diagonal and Off-Diagonal Entries of linear PDE Jacobian
182  d.clear();
183  d.resize(this->nx_,nu*2.0/this->dx_);
184  dl.clear();
185  dl.resize(this->nx_-1,-nu/this->dx_);
186  du.clear();
187  du.resize(this->nx_-1,-nu/this->dx_);
188  // Contribution from nonlinearity
189  for (int i=0; i<this->nx_; i++) {
190  if (i<this->nx_-1) {
191  dl[i] += (-2.0*u[i]-u[i+1])/6.0;
192  d[i] += u[i+1]/6.0;
193  }
194  if (i>0) {
195  d[i] += -u[i-1]/6.0;
196  du[i-1] += (u[i-1]+2.0*u[i])/6.0;
197  }
198  }
199  // Contribution from Dirichlet boundary conditions
200  d[0] -= u0/6.0;
201  d[this->nx_-1] += u1/6.0;
202  }
203 
204  void add_pde_hessian(std::vector<Real> &r, const std::vector<Real> &u, const std::vector<Real> &p,
205  const std::vector<Real> &s, Real alpha = 1.0) {
206  for (int i=0; i<this->nx_; i++) {
207  // Contribution from nonlinear term
208  if (i<this->nx_-1){
209  r[i] += alpha*(p[i]*s[i+1] - p[i+1]*(2.0*s[i]+s[i+1]))/6.0;
210  }
211  if (i>0) {
212  r[i] += alpha*(p[i-1]*(s[i-1]+2.0*s[i]) - p[i]*s[i-1])/6.0;
213  }
214  }
215  }
216 
217  void linear_solve(std::vector<Real> &u, std::vector<Real> &dl, std::vector<Real> &d, std::vector<Real> &du,
218  const std::vector<Real> &r, const bool transpose = false) {
219  u.assign(r.begin(),r.end());
220  // Perform LDL factorization
221  Teuchos::LAPACK<int,Real> lp;
222  std::vector<Real> du2(this->nx_-2,0.0);
223  std::vector<int> ipiv(this->nx_,0);
224  int info;
225  int ldb = this->nx_;
226  int nhrs = 1;
227  lp.GTTRF(this->nx_,&dl[0],&d[0],&du[0],&du2[0],&ipiv[0],&info);
228  char trans = 'N';
229  if ( transpose ) {
230  trans = 'T';
231  }
232  lp.GTTRS(trans,this->nx_,nhrs,&dl[0],&d[0],&du[0],&du2[0],&ipiv[0],&u[0],ldb,&info);
233  }
234 
235  void run_newton(std::vector<Real> &u, const std::vector<Real> &z, const std::vector<Real> &param) {
236  // Compute residual and residual norm
237  std::vector<Real> r(u.size(),0.0);
238  this->compute_residual(r,u,z,param);
239  Real rnorm = this->compute_norm(r);
240  // Define tolerances
241  Real tol = 1.e2*ROL::ROL_EPSILON<Real>();
242  Real maxit = 500;
243  // Initialize Jacobian storage
244  std::vector<Real> d(this->nx_,0.0);
245  std::vector<Real> dl(this->nx_-1,0.0);
246  std::vector<Real> du(this->nx_-1,0.0);
247  // Iterate Newton's method
248  Real alpha = 1.0, tmp = 0.0;
249  std::vector<Real> s(this->nx_,0.0);
250  std::vector<Real> utmp(this->nx_,0.0);
251  for (int i=0; i<maxit; i++) {
252  //std::cout << i << " " << rnorm << "\n";
253  // Get Jacobian
254  this->compute_pde_jacobian(dl,d,du,u,param);
255  // Solve Newton system
256  this->linear_solve(s,dl,d,du,r);
257  // Perform line search
258  tmp = rnorm;
259  alpha = 1.0;
260  utmp.assign(u.begin(),u.end());
261  this->update(utmp,s,-alpha);
262  this->compute_residual(r,utmp,z,param);
263  rnorm = this->compute_norm(r);
264  while ( rnorm > (1.0-1.e-4*alpha)*tmp && alpha > std::sqrt(ROL::ROL_EPSILON<Real>()) ) {
265  alpha /= 2.0;
266  utmp.assign(u.begin(),u.end());
267  this->update(utmp,s,-alpha);
268  this->compute_residual(r,utmp,z,param);
269  rnorm = this->compute_norm(r);
270  }
271  // Update iterate
272  u.assign(utmp.begin(),utmp.end());
273  if ( rnorm < tol ) {
274  break;
275  }
276  }
277  }
278 /*************************************************************/
279 /********** END PRIVATE MEMBER FUNCTION DECLARATION **********/
280 /*************************************************************/
281 
282 public:
283 
284  Objective_BurgersControl(Real alpha = 1.e-4, int nx = 128) : alpha_(alpha), nx_(nx) {
285  dx_ = 1.0/((Real)nx+1.0);
286  }
287 
288  void solve_state(std::vector<Real> &u, const std::vector<Real> &z, const std::vector<Real> &param) {
289  u.clear();
290  u.resize(this->nx_,1.0);
291  this->run_newton(u,z,param);
292  }
293 
294  void solve_adjoint(std::vector<Real> &p, const std::vector<Real> &u, const std::vector<Real> &param) {
295  // Initialize State Storage
296  p.clear();
297  p.resize(this->nx_);
298  // Get PDE Jacobian
299  std::vector<Real> d(this->nx_,0.0);
300  std::vector<Real> du(this->nx_-1,0.0);
301  std::vector<Real> dl(this->nx_-1,0.0);
302  this->compute_pde_jacobian(dl,d,du,u,param);
303  // Get Right Hand Side
304  std::vector<Real> r(this->nx_,0.0);
305  std::vector<Real> diff(this->nx_,0.0);
306  for (int i=0; i<this->nx_; i++) {
307  diff[i] = -(u[i]-this->evaluate_target((Real)(i+1)*this->dx_));
308  }
309  this->apply_mass(r,diff);
310  // Solve solve adjoint system at current time step
311  this->linear_solve(p,dl,d,du,r,true);
312  }
313 
314  void solve_state_sensitivity(std::vector<Real> &v, const std::vector<Real> &u,
315  const std::vector<Real> &z, const std::vector<Real> &param) {
316  // Initialize State Storage
317  v.clear();
318  v.resize(this->nx_);
319  // Get PDE Jacobian
320  std::vector<Real> d(this->nx_,0.0);
321  std::vector<Real> dl(this->nx_-1,0.0);
322  std::vector<Real> du(this->nx_-1,0.0);
323  this->compute_pde_jacobian(dl,d,du,u,param);
324  // Get Right Hand Side
325  std::vector<Real> r(this->nx_,0.0);
326  for (int i=0; i<this->nx_; i++) {
327  r[i] = this->dx_/6.0*(z[i]+4.0*z[i+1]+z[i+2]);
328  }
329  // Solve solve state sensitivity system at current time step
330  this->linear_solve(v,dl,d,du,r);
331  }
332 
333  void solve_adjoint_sensitivity(std::vector<Real> &q, const std::vector<Real> &u,
334  const std::vector<Real> &p, const std::vector<Real> &v,
335  const std::vector<Real> &z, const std::vector<Real> &param) {
336  // Initialize State Storage
337  q.clear();
338  q.resize(this->nx_);
339  // Get PDE Jacobian
340  std::vector<Real> d(this->nx_,0.0);
341  std::vector<Real> dl(this->nx_-1,0.0);
342  std::vector<Real> du(this->nx_-1,0.0);
343  this->compute_pde_jacobian(dl,d,du,u,param);
344  // Get Right Hand Side
345  std::vector<Real> r(this->nx_,0.0);
346  this->apply_mass(r,v);
347  this->scale(r,-1.0);
348  this->add_pde_hessian(r,u,p,v,-1.0);
349  // Solve Tridiagonal System Using LAPACK's SPD Tridiagonal Solver
350  this->linear_solve(q,dl,d,du,r,true);
351  }
352 
353  Real value( const ROL::Vector<Real> &z, Real &tol ) {
354  ROL::Ptr<const std::vector<Real> > zp =
355  dynamic_cast<const ROL::StdVector<Real>&>(z).getVector();
356  // SOLVE STATE EQUATION
357  std::vector<Real> param(4,0.0);
358  std::vector<Real> u;
359  this->solve_state(u,*zp,param);
360  // COMPUTE RESIDUAL
361  Real val = this->alpha_*0.5*this->dot(*zp,*zp);
362  Real res = 0.0, res1 = 0.0, res2 = 0.0, res3 = 0.0;
363  for (int i=0; i<this->nx_; i++) {
364  if ( i == 0 ) {
365  res1 = u[i]-evaluate_target((Real)(i+1)*this->dx_);
366  res2 = u[i+1]-evaluate_target((Real)(i+2)*this->dx_);
367  res = this->dx_/6.0*(4.0*res1 + res2)*res1;
368  }
369  else if ( i == this->nx_-1 ) {
370  res1 = u[i-1]-evaluate_target((Real)i*this->dx_);
371  res2 = u[i]-evaluate_target((Real)(i+1)*this->dx_);
372  res = this->dx_/6.0*(res1 + 4.0*res2)*res2;
373  }
374  else {
375  res1 = u[i-1]-evaluate_target((Real)i*this->dx_);
376  res2 = u[i]-evaluate_target((Real)(i+1)*this->dx_);
377  res3 = u[i+1]-evaluate_target((Real)(i+2)*this->dx_);
378  res = this->dx_/6.0*(res1 + 4.0*res2 + res3)*res2;
379  }
380  val += 0.5*res;
381  }
382  return val;
383  }
384 
385  void gradient( ROL::Vector<Real> &g, const ROL::Vector<Real> &z, Real &tol ) {
386  ROL::Ptr<const std::vector<Real> > zp =
387  dynamic_cast<const ROL::StdVector<Real>&>(z).getVector();
388  ROL::Ptr<std::vector<Real> > gp =
389  dynamic_cast<ROL::StdVector<Real>&>(g).getVector();
390  // SOLVE STATE EQUATION
391  std::vector<Real> param(4,0.0);
392  std::vector<Real> u;
393  this->solve_state(u,*zp,param);
394  // SOLVE ADJOINT EQUATION
395  std::vector<Real> p;
396  this->solve_adjoint(p,u,param);
397  // COMPUTE GRADIENT
398  for (int i=0; i<this->nx_+2; i++) {
399  if (i==0) {
400  (*gp)[i] = this->alpha_*this->dx_/6.0*(2.0*(*zp)[i]+(*zp)[i+1]);
401  (*gp)[i] -= this->dx_/6.0*(p[i]);
402  }
403  else if (i==this->nx_+1) {
404  (*gp)[i] = this->alpha_*this->dx_/6.0*(2.0*(*zp)[i]+(*zp)[i-1]);
405  (*gp)[i] -= this->dx_/6.0*(p[i-2]);
406  }
407  else {
408  (*gp)[i] = this->alpha_*this->dx_/6.0*((*zp)[i-1]+4.0*(*zp)[i]+(*zp)[i+1]);
409  if (i==1) {
410  (*gp)[i] -= this->dx_/6.0*(4.0*p[i-1]+p[i]);
411  }
412  else if (i==this->nx_) {
413  (*gp)[i] -= this->dx_/6.0*(4.0*p[i-1]+p[i-2]);
414  }
415  else {
416  (*gp)[i] -= this->dx_/6.0*(p[i-2]+4.0*p[i-1]+p[i]);
417  }
418  }
419  }
420  }
421 
422  void hessVec( ROL::Vector<Real> &hv, const ROL::Vector<Real> &v, const ROL::Vector<Real> &z, Real &tol ) {
423  ROL::Ptr<const std::vector<Real> > zp =
424  dynamic_cast<const ROL::StdVector<Real>&>(z).getVector();
425  ROL::Ptr<const std::vector<Real> > vp =
426  dynamic_cast<const ROL::StdVector<Real>&>(v).getVector();
427  ROL::Ptr<std::vector<Real> > hvp =
428  dynamic_cast<ROL::StdVector<Real>&>(hv).getVector();
429  // SOLVE STATE EQUATION
430  std::vector<Real> param(4,0.0);
431  std::vector<Real> u;
432  this->solve_state(u,*zp,param);
433  // SOLVE ADJOINT EQUATION
434  std::vector<Real> p;
435  this->solve_adjoint(p,u,param);
436  // SOLVE STATE SENSITIVITY EQUATION
437  std::vector<Real> s;
438  this->solve_state_sensitivity(s,u,*vp,param);
439  // SOLVE ADJOINT SENSITIVITY EQUATION
440  std::vector<Real> q;
441  this->solve_adjoint_sensitivity(q,u,p,s,*vp,param);
442  // COMPUTE HESSVEC
443  for (int i=0; i<this->nx_+2; i++) {
444  if (i==0) {
445  (*hvp)[i] = this->alpha_*this->dx_/6.0*(2.0*(*vp)[i]+(*vp)[i+1]);
446  (*hvp)[i] -= this->dx_/6.0*(q[i]);
447  }
448  else if (i==this->nx_+1) {
449  (*hvp)[i] = this->alpha_*this->dx_/6.0*(2.0*(*vp)[i]+(*vp)[i-1]);
450  (*hvp)[i] -= this->dx_/6.0*(q[i-2]);
451  }
452  else {
453  (*hvp)[i] = this->alpha_*this->dx_/6.0*((*vp)[i-1]+4.0*(*vp)[i]+(*vp)[i+1]);
454  if (i==1) {
455  (*hvp)[i] -= this->dx_/6.0*(4.0*q[i-1]+q[i]);
456  }
457  else if (i==this->nx_) {
458  (*hvp)[i] -= this->dx_/6.0*(4.0*q[i-1]+q[i-2]);
459  }
460  else {
461  (*hvp)[i] -= this->dx_/6.0*(q[i-2]+4.0*q[i-1]+q[i]);
462  }
463  }
464  }
465  }
466 };
467 
468 template<class Real>
470 private:
471  int dim_;
472  std::vector<Real> x_lo_;
473  std::vector<Real> x_up_;
474  Real min_diff_;
475 public:
477  x_lo_.resize(dim_,0.0);
478  x_up_.resize(dim_,1.0);
479  }
480  bool isFeasible( const ROL::Vector<Real> &x ) {
481  ROL::Ptr<const std::vector<Real> > ex =
482  dynamic_cast<const ROL::StdVector<Real>&>(x).getVector();
483  bool val = true;
484  int cnt = 1;
485  for ( int i = 0; i < this->dim_; i++ ) {
486  if ( (*ex)[i] >= this->x_lo_[i] && (*ex)[i] <= this->x_up_[i] ) { cnt *= 1; }
487  else { cnt *= 0; }
488  }
489  if ( cnt == 0 ) { val = false; }
490  return val;
491  }
493  ROL::Ptr<std::vector<Real> > ex =
494  dynamic_cast<ROL::StdVector<Real>&>(x).getVector();
495  for ( int i = 0; i < this->dim_; i++ ) {
496  (*ex)[i] = std::max(this->x_lo_[i],std::min(this->x_up_[i],(*ex)[i]));
497  }
498  }
499  void pruneLowerActive(ROL::Vector<Real> &v, const ROL::Vector<Real> &x, Real eps) {
500  ROL::Ptr<const std::vector<Real> > ex =
501  dynamic_cast<const ROL::StdVector<Real>&>(x).getVector();
502  ROL::Ptr<std::vector<Real> > ev =
503  dynamic_cast<ROL::StdVector<Real>&>(v).getVector();
504  Real epsn = std::min(eps,this->min_diff_);
505  for ( int i = 0; i < this->dim_; i++ ) {
506  if ( ((*ex)[i] <= this->x_lo_[i]+epsn) ) {
507  (*ev)[i] = 0.0;
508  }
509  }
510  }
511  void pruneUpperActive(ROL::Vector<Real> &v, const ROL::Vector<Real> &x, Real eps) {
512  ROL::Ptr<const std::vector<Real> > ex =
513  dynamic_cast<const ROL::StdVector<Real>&>(x).getVector();
514  ROL::Ptr<std::vector<Real> > ev =
515  dynamic_cast<ROL::StdVector<Real>&>(v).getVector();
516  Real epsn = std::min(eps,this->min_diff_);
517  for ( int i = 0; i < this->dim_; i++ ) {
518  if ( ((*ex)[i] >= this->x_up_[i]-epsn) ) {
519  (*ev)[i] = 0.0;
520  }
521  }
522  }
524  ROL::Ptr<const std::vector<Real> > ex =
525  dynamic_cast<ROL::StdVector<Real>&>(x).getVector();
526  ROL::Ptr<const std::vector<Real> > eg =
527  dynamic_cast<const ROL::StdVector<Real>&>(g).getVector();
528  ROL::Ptr<std::vector<Real> > ev =
529  dynamic_cast<ROL::StdVector<Real>&>(v).getVector();
530  Real epsn = std::min(eps,this->min_diff_);
531  for ( int i = 0; i < this->dim_; i++ ) {
532  if ( ((*ex)[i] <= this->x_lo_[i]+epsn && (*eg)[i] > 0.0) ){
533  (*ev)[i] = 0.0;
534  }
535  }
536  }
538  ROL::Ptr<const std::vector<Real> > ex =
539  dynamic_cast<const ROL::StdVector<Real>&>(x).getVector();
540  ROL::Ptr<const std::vector<Real> > eg =
541  dynamic_cast<const ROL::StdVector<Real>&>(g).getVector();
542  ROL::Ptr<std::vector<Real> > ev =
543  dynamic_cast<ROL::StdVector<Real>&>(v).getVector();
544  Real epsn = std::min(eps,this->min_diff_);
545  for ( int i = 0; i < this->dim_; i++ ) {
546  if ( ((*ex)[i] >= this->x_up_[i]-epsn && (*eg)[i] < 0.0) ) {
547  (*ev)[i] = 0.0;
548  }
549  }
550  }
552  ROL::Ptr<std::vector<Real> > us = ROL::makePtr<std::vector<Real>>(this->dim_,0.0);
553  us->assign(this->x_up_.begin(),this->x_up_.end());
554  ROL::Ptr<ROL::Vector<Real> > up = ROL::makePtr<ROL::StdVector<Real>>(us);
555  u.set(*up);
556  }
558  ROL::Ptr<std::vector<Real> > ls = ROL::makePtr<std::vector<Real>>(this->dim_,0.0);
559  ls->assign(this->x_lo_.begin(),this->x_lo_.end());
560  ROL::Ptr<ROL::Vector<Real> > lp = ROL::makePtr<ROL::StdVector<Real>>(ls);
561  l.set(*lp);
562  }
563 };
Provides the interface to evaluate 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 solve_adjoint_sensitivity(std::vector< Real > &q, const std::vector< Real > &u, const std::vector< Real > &p, const std::vector< Real > &v, const std::vector< Real > &z, const std::vector< Real > &param)
void solve_state(std::vector< Real > &u, const std::vector< Real > &z, const std::vector< Real > &param)
void run_newton(std::vector< Real > &u, const std::vector< Real > &z, const std::vector< Real > &param)
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.
Contains definitions of custom data types in ROL.
void solve_state_sensitivity(std::vector< Real > &v, const std::vector< Real > &u, const std::vector< Real > &z, const std::vector< Real > &param)
void compute_residual(std::vector< Real > &r, const std::vector< Real > &u, const std::vector< Real > &z, const std::vector< Real > &param)
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:80
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 -binding set.
Defines a no-output stream class ROL::NullStream and a function makeStreamPtr which either wraps a re...
void hessVec(ROL::Vector< Real > &hv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &z, Real &tol)
Apply Hessian approximation to vector.
void linear_solve(std::vector< Real > &u, std::vector< Real > &dl, std::vector< Real > &d, std::vector< Real > &du, const std::vector< Real > &r, const bool transpose=false)
Real dot(const std::vector< Real > &x, const std::vector< Real > &y)
Objective_BurgersControl(Real alpha=1.e-4, int nx=128)
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 compute_pde_jacobian(std::vector< Real > &dl, std::vector< Real > &d, std::vector< Real > &du, const std::vector< Real > &u, const std::vector< Real > &param)
void update(std::vector< Real > &u, const std::vector< Real > &s, const Real alpha=1.0)
void setVectorToLowerBound(ROL::Vector< Real > &l)
Provides the interface to apply upper and lower bound constraints.
void add_pde_hessian(std::vector< Real > &r, const std::vector< Real > &u, const std::vector< Real > &p, const std::vector< Real > &s, Real alpha=1.0)
void solve_adjoint(std::vector< Real > &p, const std::vector< Real > &u, const std::vector< Real > &param)
void project(ROL::Vector< Real > &x)
Project optimization variables onto the bounds.
Real value(const ROL::Vector< Real > &z, Real &tol)
Compute value.
void apply_mass(std::vector< Real > &Mu, const std::vector< Real > &u)
virtual void set(const Vector &x)
Set where .
Definition: ROL_Vector.hpp:209
void gradient(ROL::Vector< Real > &g, const ROL::Vector< Real > &z, Real &tol)
Compute gradient.
bool isFeasible(const ROL::Vector< Real > &x)
Check if the vector, v, is feasible.
constexpr auto dim
void setVectorToUpperBound(ROL::Vector< Real > &u)
Real compute_norm(const std::vector< Real > &r)
void scale(std::vector< Real > &u, const Real alpha=0.0)