ROL
poisson-inversion/example_02.cpp
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_StdVector.hpp"
50 #include "ROL_Objective.hpp"
52 #include "ROL_BoundConstraint.hpp"
53 #include "ROL_LineSearchStep.hpp"
54 #include "ROL_TrustRegionStep.hpp"
55 #include "ROL_Algorithm.hpp"
56 #include "ROL_Types.hpp"
57 #include "Teuchos_oblackholestream.hpp"
58 #include "Teuchos_GlobalMPISession.hpp"
59 
60 #include <iostream>
61 #include <algorithm>
62 
63 template<class Real>
65 private:
66  int nu_;
67  int nz_;
68 
69  Real hu_;
70  Real hz_;
71 
72  Real u0_;
73  Real u1_;
74 
75  Real alpha_;
76 
78  Teuchos::SerialDenseMatrix<int, Real> H_;
79 
80 public:
81 
82  /* CONSTRUCTOR */
83  Objective_PoissonInversion(int nz = 32, Real alpha = 1.e-4)
84  : nz_(nz), u0_(1.0), u1_(0.0), alpha_(alpha), useCorrection_(false) {
85  nu_ = nz_-1;
86  hu_ = 1.0/((Real)nu_+1.0);
87  hz_ = hu_;
88  }
89 
90  void apply_mass(std::vector<Real> &Mz, const std::vector<Real> &z ) {
91  Mz.resize(this->nu_,0.0);
92  for (int i=0; i<this->nu_; i++) {
93  if ( i == 0 ) {
94  Mz[i] = this->hu_/6.0*(2.0*z[i] + z[i+1]);
95  }
96  else if ( i == this->nu_-1 ) {
97  Mz[i] = this->hu_/6.0*(z[i-1] + 2.0*z[i]);
98  }
99  else {
100  Mz[i] = this->hu_/6.0*(z[i-1] + 4.0*z[i] + z[i+1]);
101  }
102  }
103  }
104 
105  Real evaluate_target(Real x) {
106  return (x <= 0.5) ? 1.0 : 0.0;
107  }
108 
109  void apply_linearized_control_operator( std::vector<Real> &Bd, const std::vector<Real> &z,
110  const std::vector<Real> &d, const std::vector<Real> &u,
111  bool addBC = true ) {
112  Bd.clear();
113  Bd.resize(this->nu_,0.0);
114  for (int i = 0; i < this->nu_; i++) {
115  if ( i == 0 ) {
116  Bd[i] = 1.0/this->hu_*( u[i]*d[i] + (u[i]-u[i+1])*d[i+1] );
117  }
118  else if ( i == this->nu_-1 ) {
119  Bd[i] = 1.0/this->hu_*( (u[i]-u[i-1])*d[i] + u[i]*d[i+1] );
120  }
121  else {
122  Bd[i] = 1.0/this->hu_*( (u[i]-u[i-1])*d[i] + (u[i]-u[i+1])*d[i+1] );
123  }
124  }
125  if ( addBC ) {
126  Bd[ 0] -= this->u0_*d[ 0]/this->hu_;
127  Bd[this->nu_-1] -= this->u1_*d[this-> nz_-1]/this->hu_;
128  }
129  }
130 
131  void apply_transposed_linearized_control_operator( std::vector<Real> &Bd, const std::vector<Real> &z,
132  const std::vector<Real> &d, const std::vector<Real> &u,
133  bool addBC = true ) {
134  Bd.clear();
135  Bd.resize(this->nz_,0.0);
136  for (int i = 0; i < this->nz_; i++) {
137  if ( i == 0 ) {
138  Bd[i] = 1.0/this->hu_*u[i]*d[i];
139  }
140  else if ( i == this->nz_-1 ) {
141  Bd[i] = 1.0/this->hu_*u[i-1]*d[i-1];
142  }
143  else {
144  Bd[i] = 1.0/this->hu_*( (u[i]-u[i-1])*(d[i]-d[i-1]) );
145  }
146  }
147  if ( addBC ) {
148  Bd[ 0] -= this->u0_*d[ 0]/this->hu_;
149  Bd[this->nz_-1] -= this->u1_*d[this-> nu_-1]/this->hu_;
150  }
151  }
152 
153  /* STATE AND ADJOINT EQUATION DEFINTIONS */
154  void solve_state_equation(std::vector<Real> &u, const std::vector<Real> &z) {
155  // Get Diagonal and Off-Diagonal Entries of PDE Jacobian
156  std::vector<Real> d(this->nu_,1.0);
157  std::vector<Real> o(this->nu_-1,1.0);
158  for ( int i = 0; i < this->nu_; i++ ) {
159  d[i] = (z[i] + z[i+1])/this->hu_;
160  if ( i < this->nu_-1 ) {
161  o[i] *= -z[i+1]/this->hu_;
162  }
163  }
164  // Set right hand side
165  u.clear();
166  u.resize(this->nu_,0.0);
167  u[ 0] = z[ 0]/this->hu_ * this->u0_;
168  u[this->nu_-1] = z[this->nz_-1]/this->hu_ * this->u1_;
169  // Solve Tridiagonal System Using LAPACK's SPD Tridiagonal Solver
170  Teuchos::LAPACK<int,Real> lp;
171  int info;
172  int ldb = this->nu_;
173  int nhrs = 1;
174  lp.PTTRF(this->nu_,&d[0],&o[0],&info);
175  lp.PTTRS(this->nu_,nhrs,&d[0],&o[0],&u[0],ldb,&info);
176  }
177 
178  void solve_adjoint_equation(std::vector<Real> &p, const std::vector<Real> &u, const std::vector<Real> &z) {
179  // Get Diagonal and Off-Diagonal Entries of PDE Jacobian
180  std::vector<Real> d(this->nu_,1.0);
181  std::vector<Real> o(this->nu_-1,1.0);
182  for ( int i = 0; i < this->nu_; i++ ) {
183  d[i] = (z[i] + z[i+1])/this->hu_;
184  if ( i < this->nu_-1 ) {
185  o[i] *= -z[i+1]/this->hu_;
186  }
187  }
188  // Set right hand side
189  std::vector<Real> r(this->nu_,0.0);
190  for (int i = 0; i < this->nu_; i++) {
191  r[i] = -(u[i]-this->evaluate_target((Real)(i+1)*this->hu_));
192  }
193  p.clear();
194  p.resize(this->nu_,0.0);
195  this->apply_mass(p,r);
196  // Solve Tridiagonal System Using LAPACK's SPD Tridiagonal Solver
197  Teuchos::LAPACK<int,Real> lp;
198  int info;
199  int ldb = this->nu_;
200  int nhrs = 1;
201  lp.PTTRF(this->nu_,&d[0],&o[0],&info);
202  lp.PTTRS(this->nu_,nhrs,&d[0],&o[0],&p[0],ldb,&info);
203  }
204 
205  void solve_state_sensitivity_equation(std::vector<Real> &w, const std::vector<Real> &v,
206  const std::vector<Real> &u, const std::vector<Real> &z) {
207  // Get Diagonal and Off-Diagonal Entries of PDE Jacobian
208  std::vector<Real> d(this->nu_,1.0);
209  std::vector<Real> o(this->nu_-1,1.0);
210  for ( int i = 0; i < this->nu_; i++ ) {
211  d[i] = (z[i] + z[i+1])/this->hu_;
212  if ( i < this->nu_-1 ) {
213  o[i] *= -z[i+1]/this->hu_;
214  }
215  }
216  // Set right hand side
217  w.clear();
218  w.resize(this->nu_,0.0);
219  this->apply_linearized_control_operator(w,z,v,u);
220  // Solve Tridiagonal System Using LAPACK's SPD Tridiagonal Solver
221  Teuchos::LAPACK<int,Real> lp;
222  int info;
223  int ldb = this->nu_;
224  int nhrs = 1;
225  lp.PTTRF(this->nu_,&d[0],&o[0],&info);
226  lp.PTTRS(this->nu_,nhrs,&d[0],&o[0],&w[0],ldb,&info);
227  }
228 
229  void solve_adjoint_sensitivity_equation(std::vector<Real> &q, const std::vector<Real> &w,
230  const std::vector<Real> &v, const std::vector<Real> &p,
231  const std::vector<Real> &u, const std::vector<Real> &z) {
232  // Get Diagonal and Off-Diagonal Entries of PDE Jacobian
233  std::vector<Real> d(this->nu_,1.0);
234  std::vector<Real> o(this->nu_-1,1.0);
235  for ( int i = 0; i < this->nu_; i++ ) {
236  d[i] = (z[i] + z[i+1])/this->hu_;
237  if ( i < this->nu_-1 ) {
238  o[i] *= -z[i+1]/this->hu_;
239  }
240  }
241  // Set right hand side
242  q.clear();
243  q.resize(this->nu_,0.0);
244  this->apply_mass(q,w);
245  std::vector<Real> res(this->nu_,0.0);
246  this->apply_linearized_control_operator(res,z,v,p,false);
247  for (int i = 0; i < this->nu_; i++) {
248  q[i] -= res[i];
249  }
250  // Solve Tridiagonal System Using LAPACK's SPD Tridiagonal Solver
251  Teuchos::LAPACK<int,Real> lp;
252  int info;
253  int ldb = this->nu_;
254  int nhrs = 1;
255  lp.PTTRF(this->nu_,&d[0],&o[0],&info);
256  lp.PTTRS(this->nu_,nhrs,&d[0],&o[0],&q[0],ldb,&info);
257  }
258 
259  void update(const ROL::Vector<Real> &z, bool flag, int iter) {
260  if ( flag && this->useCorrection_ ) {
261  Real tol = std::sqrt(ROL::ROL_EPSILON);
262  this->H_.shape(this->nz_,this->nz_);
263  Teuchos::RCP<ROL::Vector<Real> > e = z.clone();
264  Teuchos::RCP<ROL::Vector<Real> > h = z.clone();
265  for ( int i = 0; i < this->nz_; i++ ) {
266  e = z.basis(i);
267  this->hessVec_true(*h,*e,z,tol);
268  for ( int j = 0; j < this->nz_; j++ ) {
269  e = z.basis(j);
270  (this->H_)(j,i) = e->dot(*h);
271  }
272  }
273  std::vector<std::vector<Real> > eigenvals = ROL::computeEigenvalues<Real>(this->H_);
274  std::sort((eigenvals[0]).begin(), (eigenvals[0]).end());
275  Real inertia = (eigenvals[0])[0];
276  Real correction = 0.0;
277  if ( inertia <= 0.0 ) {
278  correction = (1.0+std::sqrt(ROL::ROL_EPSILON))*std::abs(inertia);
279  if ( inertia == 0.0 ) {
280  int cnt = 0;
281  while ( eigenvals[0][cnt] == 0.0 ) {
282  cnt++;
283  }
284  correction = std::sqrt(ROL::ROL_EPSILON)*eigenvals[0][cnt];
285  if ( cnt == this->nz_-1 ) {
286  correction = 1.0;
287  }
288  }
289  for ( int i = 0; i < this->nz_; i++ ) {
290  (this->H_)(i,i) += correction;
291  }
292  }
293  }
294  }
295 
296  /* OBJECTIVE FUNCTION DEFINITIONS */
297  Real value( const ROL::Vector<Real> &z, Real &tol ) {
298  Teuchos::RCP<const std::vector<Real> > zp =
299  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(z))).getVector();
300  // SOLVE STATE EQUATION
301  std::vector<Real> u(this->nu_,0.0);
302  this->solve_state_equation(u,*zp);
303  // EVALUATE OBJECTIVE
304  Real val = 0.0;
305  for (int i=0; i<this->nz_;i++) {
306  val += this->hz_*this->alpha_*0.5*(*zp)[i]*(*zp)[i];
307  }
308  Real res = 0.0;
309  Real res1 = 0.0;
310  Real res2 = 0.0;
311  Real res3 = 0.0;
312  for (int i=0; i<this->nu_; i++) {
313  if ( i == 0 ) {
314  res1 = u[i]-evaluate_target((Real)(i+1)*this->hu_);
315  res2 = u[i+1]-evaluate_target((Real)(i+2)*this->hu_);
316  res = this->hu_/6.0*(2.0*res1 + res2)*res1;
317  }
318  else if ( i == this->nu_-1 ) {
319  res1 = u[i-1]-evaluate_target((Real)i*this->hu_);
320  res2 = u[i]-evaluate_target((Real)(i+1)*this->hu_);
321  res = this->hu_/6.0*(res1 + 2.0*res2)*res2;
322  }
323  else {
324  res1 = u[i-1]-evaluate_target((Real)i*this->hu_);
325  res2 = u[i]-evaluate_target((Real)(i+1)*this->hu_);
326  res3 = u[i+1]-evaluate_target((Real)(i+2)*this->hu_);
327  res = this->hu_/6.0*(res1 + 4.0*res2 + res3)*res2;
328  }
329  val += 0.5*res;
330  }
331  return val;
332  }
333 
334  void gradient( ROL::Vector<Real> &g, const ROL::Vector<Real> &z, Real &tol ) {
335  Teuchos::RCP<const std::vector<Real> > zp =
336  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(z))).getVector();
337  Teuchos::RCP<std::vector<Real> > gp =
338  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(g)).getVector());
339  // SOLVE STATE EQUATION
340  std::vector<Real> u(this->nu_,0.0);
341  this->solve_state_equation(u,*zp);
342  // SOLVE ADJOINT EQUATION
343  std::vector<Real> p(this->nu_,0.0);
344  this->solve_adjoint_equation(p,u,*zp);
345  // Apply Transpose of Linearized Control Operator
347  // Build Gradient
348  for ( int i = 0; i < this->nz_; i++ ) {
349  (*gp)[i] += this->hz_*this->alpha_*(*zp)[i];
350  }
351  }
352 
353  void hessVec( ROL::Vector<Real> &hv, const ROL::Vector<Real> &v, const ROL::Vector<Real> &z, Real &tol ) {
354  if ( this->useCorrection_ ) {
355  this->hessVec_inertia(hv,v,z,tol);
356  }
357  else {
358  this->hessVec_true(hv,v,z,tol);
359  }
360  }
361 
362  void activateInertia(void) {
363  this->useCorrection_ = true;
364  }
365 
366  void deactivateInertia(void) {
367  this->useCorrection_ = false;
368  }
369 
370  void hessVec_true( ROL::Vector<Real> &hv, const ROL::Vector<Real> &v, const ROL::Vector<Real> &z, Real &tol ) {
371  Teuchos::RCP<const std::vector<Real> > vp =
372  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(v))).getVector();
373  Teuchos::RCP<const std::vector<Real> > zp =
374  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(z))).getVector();
375  Teuchos::RCP<std::vector<Real> > hvp =
376  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(hv)).getVector());
377  // SOLVE STATE EQUATION
378  std::vector<Real> u(this->nu_,0.0);
379  this->solve_state_equation(u,*zp);
380  // SOLVE ADJOINT EQUATION
381  std::vector<Real> p(this->nu_,0.0);
382  this->solve_adjoint_equation(p,u,*zp);
383  // SOLVE STATE SENSITIVITY EQUATION
384  std::vector<Real> w(this->nu_,0.0);
385  this->solve_state_sensitivity_equation(w,*vp,u,*zp);
386  // SOLVE ADJOINT SENSITIVITY EQUATION
387  std::vector<Real> q(this->nu_,0.0);
388  this->solve_adjoint_sensitivity_equation(q,w,*vp,p,u,*zp);
389  // Apply Transpose of Linearized Control Operator
391  // Apply Transpose of Linearized Control Operator
392  std::vector<Real> tmp(this->nz_,0.0);
393  this->apply_transposed_linearized_control_operator(tmp,*zp,w,p,false);
394  for (int i=0; i < this->nz_; i++) {
395  (*hvp)[i] -= tmp[i];
396  }
397  // Regularization hessVec
398  for (int i=0; i < this->nz_; i++) {
399  (*hvp)[i] += this->hz_*this->alpha_*(*vp)[i];
400  }
401  }
402 
403  void hessVec_inertia( ROL::Vector<Real> &hv, const ROL::Vector<Real> &v, const ROL::Vector<Real> &z, Real &tol ) {
404  Teuchos::RCP<std::vector<Real> > vp =
405  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(v))).getVector());
406  Teuchos::RCP<std::vector<Real> > hvp =
407  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(hv)).getVector());
408 
409  Teuchos::SerialDenseVector<int, Real> hv_teuchos(Teuchos::View, &((*hvp)[0]), this->nz_);
410  Teuchos::SerialDenseVector<int, Real> v_teuchos(Teuchos::View, &(( *vp)[0]), this->nz_);
411  hv_teuchos.multiply(Teuchos::NO_TRANS, Teuchos::NO_TRANS, 1.0, this->H_, v_teuchos, 0.0);
412  }
413 
414 };
415 
416 template<class Real>
418 private:
419  int dim_;
420  std::vector<Real> x_lo_;
421  std::vector<Real> x_up_;
422  Real min_diff_;
423 public:
424  BoundConstraint_PoissonInversion(std::vector<Real> &lo, std::vector<Real> &up) {
425  dim_ = lo.size();
426  x_lo_.clear();
427  x_lo_.assign(lo.begin(),lo.end());
428  x_up_.clear();
429  x_up_.assign(up.begin(),up.end());
430  for ( unsigned i = 0; i < (unsigned)dim_; i++ ) {
431  if ( i == 0 ) {
432  min_diff_ = x_up_[i]-x_lo_[i];
433  }
434  else {
435  min_diff_ = std::min(min_diff_,x_up_[i]-x_lo_[i]);
436  }
437  }
438  min_diff_ *= 0.5;
439  }
440  bool isFeasible( const ROL::Vector<Real> &x ) {
441  Teuchos::RCP<const std::vector<Real> > ex =
442  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(x))).getVector();
443  bool val = true;
444  int cnt = 1;
445  for ( int i = 0; i < this->dim_; i++ ) {
446  if ( (*ex)[i] >= this->x_lo_[i] && (*ex)[i] <= this->x_up_[i] ) { cnt *= 1; }
447  else { cnt *= 0; }
448  }
449  if ( cnt == 0 ) { val = false; }
450  return val;
451  }
453  Teuchos::RCP<std::vector<Real> > ex =
454  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(x)).getVector());
455  for ( int i = 0; i < this->dim_; i++ ) {
456  (*ex)[i] = std::max(this->x_lo_[i],std::min(this->x_up_[i],(*ex)[i]));
457  }
458  }
459  void pruneLowerActive(ROL::Vector<Real> &v, const ROL::Vector<Real> &x, Real eps) {
460  Teuchos::RCP<const std::vector<Real> > ex =
461  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(x))).getVector();
462  Teuchos::RCP<std::vector<Real> > ev =
463  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(v)).getVector());
464  Real epsn = std::min(eps,this->min_diff_);
465  for ( int i = 0; i < this->dim_; i++ ) {
466  if ( ((*ex)[i] <= this->x_lo_[i]+epsn) ) {
467  (*ev)[i] = 0.0;
468  }
469  }
470  }
471  void pruneUpperActive(ROL::Vector<Real> &v, const ROL::Vector<Real> &x, Real eps) {
472  Teuchos::RCP<const std::vector<Real> > ex =
473  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(x))).getVector();
474  Teuchos::RCP<std::vector<Real> > ev =
475  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(v)).getVector());
476  Real epsn = std::min(eps,this->min_diff_);
477  for ( int i = 0; i < this->dim_; i++ ) {
478  if ( ((*ex)[i] >= this->x_up_[i]-epsn) ) {
479  (*ev)[i] = 0.0;
480  }
481  }
482  }
484  Teuchos::RCP<const std::vector<Real> > ex =
485  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(x))).getVector();
486  Teuchos::RCP<const std::vector<Real> > eg =
487  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(g))).getVector();
488  Teuchos::RCP<std::vector<Real> > ev =
489  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<ROL::StdVector<Real> >(v)).getVector());
490  Real epsn = std::min(eps,this->min_diff_);
491  for ( int i = 0; i < this->dim_; i++ ) {
492  if ( ((*ex)[i] <= this->x_lo_[i]+epsn && (*eg)[i] > 0.0) ){
493  (*ev)[i] = 0.0;
494  }
495  }
496  }
498  Teuchos::RCP<const std::vector<Real> > ex =
499  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(x))).getVector();
500  Teuchos::RCP<const std::vector<Real> > eg =
501  (Teuchos::dyn_cast<ROL::StdVector<Real> >(const_cast<ROL::Vector<Real> &>(g))).getVector();
502  Teuchos::RCP<std::vector<Real> > ev =
503  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_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_up_[i]-epsn && (*eg)[i] < 0.0) ) {
507  (*ev)[i] = 0.0;
508  }
509  }
510  }
512  Teuchos::RCP<std::vector<Real> > us = Teuchos::rcp( new std::vector<Real>(this->dim_,0.0) );
513  us->assign(this->x_up_.begin(),this->x_up_.end());
514  Teuchos::RCP<ROL::Vector<Real> > up = Teuchos::rcp( new ROL::StdVector<Real>(us) );
515  u.set(*up);
516  }
518  Teuchos::RCP<std::vector<Real> > ls = Teuchos::rcp( new std::vector<Real>(this->dim_,0.0) );
519  ls->assign(this->x_lo_.begin(),this->x_lo_.end());
520  Teuchos::RCP<ROL::Vector<Real> > lp = Teuchos::rcp( new ROL::StdVector<Real>(ls) );
521  l.set(*lp);
522  }
523 };
524 
525 
526 
527 typedef double RealT;
528 
529 int main(int argc, char *argv[]) {
530 
531  Teuchos::GlobalMPISession mpiSession(&argc, &argv);
532 
533  // This little trick lets us print to std::cout only if a (dummy) command-line argument is provided.
534  int iprint = argc - 1;
535  Teuchos::RCP<std::ostream> outStream;
536  Teuchos::oblackholestream bhs; // outputs nothing
537  if (iprint > 0)
538  outStream = Teuchos::rcp(&std::cout, false);
539  else
540  outStream = Teuchos::rcp(&bhs, false);
541 
542  int errorFlag = 0;
543 
544  // *** Example body.
545 
546  try {
547 
548  int dim = 128; // Set problem dimension.
549  RealT alpha = 1.e-6;
550  Objective_PoissonInversion<RealT> obj(dim, alpha);
551 
552  // Iteration vector.
553  Teuchos::RCP<std::vector<RealT> > x_rcp = Teuchos::rcp( new std::vector<RealT> (dim, 0.0) );
554  Teuchos::RCP<std::vector<RealT> > y_rcp = Teuchos::rcp( new std::vector<RealT> (dim, 0.0) );
555  // Set initial guess.
556  for (int i=0; i<dim; i++) {
557  (*x_rcp)[i] = (RealT)rand()/(RealT)RAND_MAX + 1.e2;
558  (*y_rcp)[i] = (RealT)rand()/(RealT)RAND_MAX + 1.e2;
559  }
560  ROL::StdVector<RealT> x(x_rcp);
561  ROL::StdVector<RealT> y(y_rcp);
562  obj.checkGradient(x,y,true);
563  obj.checkHessVec(x,y,true);
564 
565  std::vector<RealT> lo(dim,1.0);
566  std::vector<RealT> up(dim,10.0);
568 
569  Teuchos::ParameterList parlist;
570  // Basic algorithm.
571  parlist.set("Trust-Region Subproblem Solver Type", "Truncated CG");
572  parlist.set("Initial Trust-Region Radius", 100.0);
573  // Secant parameters.
574  parlist.set("Secant Type", "Limited-Memory BFGS");
575  parlist.set("Maximum Secant Storage", 100);
576  // Krylov parameters.
577  parlist.set("Absolute Krylov Tolerance", 1.e-8);
578  parlist.set("Relative Krylov Tolerance", 1.e-4);
579  parlist.set("Maximum Number of Krylov Iterations", dim);
580  // PDAS parameters.
581  parlist.set("PDAS Relative Step Tolerance", 1.e-8);
582  parlist.set("PDAS Relative Gradient Tolerance", 1.e-6);
583  parlist.set("PDAS Maximum Number of Iterations", 10);
584  parlist.set("PDAS Dual Scaling", alpha);
585  // Define step.
586  parlist.set("Use Secant Hessian-Times-A-Vector", true);
588 
589  // Define status test.
590  RealT gtol = 1e-12; // norm of gradient tolerance
591  RealT stol = 1e-14; // norm of step tolerance
592  int maxit = 20000; // maximum number of iterations
593  ROL::StatusTest<RealT> status(gtol, stol, maxit);
594 
595  // Define algorithm.
596  ROL::DefaultAlgorithm<RealT> algo(step,status,false);
597 
598  x.zero();
599  obj.deactivateInertia();
600  algo.run(x,obj,icon,true);
601 
602  // Output control to file.
603  std::ofstream file;
604  file.open("control_PDAS.txt");
605  for ( unsigned i = 0; i < (unsigned)dim; i++ ) {
606  file << (*x_rcp)[i] << "\n";
607  }
608  file.close();
609 
610  // Projected Newtion.
611  // Define step.
612  parlist.set("Use Secant Hessian-Times-A-Vector", false);
613  ROL::TrustRegionStep<RealT> step_tr(parlist);
614  // Define algorithm.
615  ROL::DefaultAlgorithm<RealT> algo_tr(step_tr,status,false);
616  // Run Algorithm
617  y.zero();
618  obj.deactivateInertia();
619  algo_tr.run(y,obj,icon,true);
620 
621  std::ofstream file_tr;
622  file_tr.open("control_TR.txt");
623  for ( unsigned i = 0; i < (unsigned)dim; i++ ) {
624  file_tr << (*y_rcp)[i] << "\n";
625  }
626  file_tr.close();
627 
628  std::vector<RealT> u;
629  obj.solve_state_equation(u,*y_rcp);
630  std::ofstream file_u;
631  file_u.open("state.txt");
632  for ( unsigned i = 0; i < (unsigned)(dim-1); i++ ) {
633  file_u << u[i] << "\n";
634  }
635  file_u.close();
636 
637  Teuchos::RCP<ROL::Vector<RealT> > diff = x.clone();
638  diff->set(x);
639  diff->axpy(-1.0,y);
640  RealT error = diff->norm()/std::sqrt((RealT)dim-1.0);
641  std::cout << "\nError between PDAS solution and TR solution is " << error << "\n";
642  errorFlag = ((error > 1.e2*std::sqrt(ROL::ROL_EPSILON)) ? 1 : 0);
643 
644  }
645  catch (std::logic_error err) {
646  *outStream << err.what() << "\n";
647  errorFlag = -1000;
648  }; // end try
649 
650  if (errorFlag != 0)
651  std::cout << "End Result: TEST FAILED\n";
652  else
653  std::cout << "End Result: TEST PASSED\n";
654 
655  return 0;
656 
657 }
658 
Implements the computation of optimization steps with the Newton primal-dual active set method...
Provides the interface to evaluate objective functions.
BoundConstraint_PoissonInversion(std::vector< Real > &lo, std::vector< Real > &up)
void solve_state_equation(std::vector< Real > &u, const std::vector< Real > &z)
void apply_transposed_linearized_control_operator(std::vector< Real > &Bd, const std::vector< Real > &z, const std::vector< Real > &d, const std::vector< Real > &u, bool addBC=true)
Contains definitions of custom data types in ROL.
bool isFeasible(const ROL::Vector< Real > &x)
Check if the vector, v, is feasible.
virtual Teuchos::RCP< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
virtual void zero()
Set to zero vector.
Definition: ROL_Vector.hpp:155
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:72
void solve_adjoint_sensitivity_equation(std::vector< Real > &q, const std::vector< Real > &w, const std::vector< Real > &v, const std::vector< Real > &p, const std::vector< Real > &u, const std::vector< Real > &z)
virtual std::vector< std::vector< Real > > checkGradient(const Vector< Real > &x, const Vector< Real > &d, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
Finite-difference gradient check.
Teuchos::RCP< Vector< Real > > clone() const
Clone to make a new (uninitialized) vector.
void apply_mass(std::vector< Real > &Mz, const std::vector< Real > &z)
Provides the std::vector implementation of the ROL::Vector interface.
virtual std::vector< std::string > run(Vector< Real > &x, Objective< Real > &obj, bool print=false, std::ostream &outStream=std::cout)
Run algorithm on unconstrained problems (Type-U). This is the primary Type-U interface.
Real value(const ROL::Vector< Real > &z, Real &tol)
Compute value.
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 update(const ROL::Vector< Real > &z, bool flag, int iter)
Update objective function.
int main(int argc, char *argv[])
Provides an interface to check status of optimization algorithms.
void hessVec_true(ROL::Vector< Real > &hv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &z, Real &tol)
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 setVectorToLowerBound(ROL::Vector< Real > &l)
Set the input vector to the lower bound.
Provides the interface to apply upper and lower bound constraints.
void solve_adjoint_equation(std::vector< Real > &p, const std::vector< Real > &u, const std::vector< Real > &z)
void setVectorToUpperBound(ROL::Vector< Real > &u)
Set the input vector to the upper bound.
void solve_state_sensitivity_equation(std::vector< Real > &w, const std::vector< Real > &v, const std::vector< Real > &u, const std::vector< Real > &z)
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 project(ROL::Vector< Real > &x)
Project optimization variables onto the bounds.
void hessVec(ROL::Vector< Real > &hv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &z, Real &tol)
Apply Hessian approximation to vector.
virtual std::vector< std::vector< Real > > checkHessVec(const Vector< Real > &x, const Vector< Real > &v, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
Finite-difference Hessian-applied-to-vector check.
double RealT
virtual void set(const Vector &x)
Set where .
Definition: ROL_Vector.hpp:194
virtual Teuchos::RCP< Vector > basis(const int i) const
Return i-th basis vector.
Definition: ROL_Vector.hpp:170
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 apply_linearized_control_operator(std::vector< Real > &Bd, const std::vector< Real > &z, const std::vector< Real > &d, const std::vector< Real > &u, bool addBC=true)
void hessVec_inertia(ROL::Vector< Real > &hv, const ROL::Vector< Real > &v, const ROL::Vector< Real > &z, Real &tol)
Objective_PoissonInversion(int nz=32, Real alpha=1.e-4)
void gradient(ROL::Vector< Real > &g, const ROL::Vector< Real > &z, Real &tol)
Compute gradient.
Teuchos::SerialDenseMatrix< int, Real > H_
Provides the interface to compute optimization steps with trust regions.
static const double ROL_EPSILON
Platform-dependent machine epsilon.
Definition: ROL_Types.hpp:115