ROL
ROL_PoissonInversion.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 #ifndef USE_HESSVEC
50 #define USE_HESSVEC 1
51 #endif
52 
53 #ifndef ROL_POISSONINVERSION_HPP
54 #define ROL_POISSONINVERSION_HPP
55 
56 #include "ROL_StdVector.hpp"
57 #include "ROL_Objective.hpp"
58 #include "ROL_HelperFunctions.hpp"
59 
60 #include "Teuchos_LAPACK.hpp"
61 
62 namespace ROL {
63 namespace ZOO {
64 
67  template<class Real>
68  class Objective_PoissonInversion : public Objective<Real> {
69  private:
70  int nu_;
71  int nz_;
72 
73  Real hu_;
74  Real hz_;
75 
76  Real alpha_;
77 
78  Real eps_;
79  int reg_type_;
80 
81  public:
82 
83  /* CONSTRUCTOR */
84  Objective_PoissonInversion(int nz = 32, Real alpha = 1.e-4) : nz_(nz), alpha_(alpha) {
85  nu_ = nz_-1;
86  hu_ = 1.0/((Real)nu_+1.0);
87  hz_ = hu_;
88  eps_ = 1.e-4;
89  reg_type_ = 2;
90  }
91 
92  /* REGULARIZATION DEFINITIONS */
93  Real reg_value(const Vector<Real> &z) {
94  Teuchos::RCP<const std::vector<Real> > zp = ROL::StdVector_Helper::constDownCast(z);
95 // Teuchos::RCP<const std::vector<Real> > zp =
96 // (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(z))).getVector();
97  Real val = 0.0;
98  for (int i = 0; i < this->nz_; i++) {
99  if ( this->reg_type_ == 2 ) {
100  val += this->alpha_/2.0 * this->hz_ * (*zp)[i]*(*zp)[i];
101  }
102  else if ( this->reg_type_ == 1 ) {
103  val += this->alpha_ * this->hz_ * std::sqrt((*zp)[i]*(*zp)[i] + this->eps_);
104  }
105  else if ( this->reg_type_ == 0 ) {
106  if ( i < this->nz_-1 ) {
107  val += this->alpha_ * std::sqrt(std::pow((*zp)[i]-(*zp)[i+1],2.0)+this->eps_);
108  }
109  }
110  }
111  return val;
112  }
113 
115  if ( this->reg_type_ == 2 ) {
116  g.set(z);
117  g.scale(this->alpha_*this->hz_);
118  }
119  else if ( this->reg_type_ == 1 ) {
120  Teuchos::RCP<const std::vector<Real> > zp = ROL::StdVector_Helper::constDownCast(z);
121  Teuchos::RCP<std::vector<Real> > gp = ROL::StdVector_Helper::downCast(g);
122 // Teuchos::RCP<const std::vector<Real> > zp =
123 // (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(z))).getVector();
124 // Teuchos::RCP<std::vector<Real> > gp =
125 // Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(g)).getVector());
126  for (int i = 0; i < this->nz_; i++) {
127  (*gp)[i] = this->alpha_ * this->hz_ * (*zp)[i]/std::sqrt(std::pow((*zp)[i],2.0)+this->eps_);
128  }
129  }
130  else if ( this->reg_type_ == 0 ) {
131  Teuchos::RCP<const std::vector<Real> > zp = ROL::StdVector_Helper::constDownCast(z);
132  Teuchos::RCP<std::vector<Real> > gp = ROL::StdVector_Helper::downCast(g);
133 // Teuchos::RCP<const std::vector<Real> > zp =
134 // (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(z))).getVector();
135 // Teuchos::RCP<std::vector<Real> > gp =
136 // Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(g)).getVector());
137  Real diff = 0.0;
138  for (int i = 0; i < this->nz_; i++) {
139  if ( i == 0 ) {
140  diff = (*zp)[i]-(*zp)[i+1];
141  (*gp)[i] = this->alpha_ * diff/std::sqrt(std::pow(diff,2.0)+this->eps_);
142  }
143  else if ( i == this->nz_-1 ) {
144  diff = (*zp)[i-1]-(*zp)[i];
145  (*gp)[i] = -this->alpha_ * diff/std::sqrt(std::pow(diff,2.0)+this->eps_);
146  }
147  else {
148  diff = (*zp)[i]-(*zp)[i+1];
149  (*gp)[i] = this->alpha_ * diff/std::sqrt(std::pow(diff,2.0)+this->eps_);
150  diff = (*zp)[i-1]-(*zp)[i];
151  (*gp)[i] -= this->alpha_ * diff/std::sqrt(std::pow(diff,2.0)+this->eps_);
152  }
153  }
154  }
155  }
156 
157  void reg_hessVec(Vector<Real> &hv, const Vector<Real> &v, const Vector<Real> &z) {
158  if ( this->reg_type_ == 2 ) {
159  hv.set(v);
160  hv.scale(this->alpha_*this->hz_);
161  }
162  else if ( this->reg_type_ == 1 ) {
163  Teuchos::RCP<const std::vector<Real> > zp = ROL::StdVector_Helper::constDownCast(z);
164  Teuchos::RCP<const std::vector<Real> > vp = ROL::StdVector_Helper::constDownCast(v);
165  Teuchos::RCP<std::vector<Real> > hvp = ROL::StdVector_Helper::downCast(hv);
166 // Teuchos::RCP<const std::vector<Real> > zp =
167 // (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(z))).getVector();
168 // Teuchos::RCP<const std::vector<Real> > vp =
169 // (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(v))).getVector();
170 // Teuchos::RCP<std::vector<Real> > hvp =
171 // Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(hv)).getVector());
172  for (int i = 0; i < this->nz_; i++) {
173  (*hvp)[i] = this->alpha_*this->hz_*(*vp)[i]*this->eps_/std::pow(std::pow((*zp)[i],2.0)+this->eps_,3.0/2.0);
174  }
175  }
176  else if ( this->reg_type_ == 0 ) {
177  Teuchos::RCP<const std::vector<Real> > zp = ROL::StdVector_Helper::constDownCast(z);
178  Teuchos::RCP<const std::vector<Real> > vp = ROL::StdVector_Helper::constDownCast(v);
179  Teuchos::RCP<std::vector<Real> > hvp = ROL::StdVector_Helper::downCast(hv);
180 // Teuchos::RCP<const std::vector<Real> > zp =
181 // (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(z))).getVector();
182 // Teuchos::RCP<const std::vector<Real> > vp =
183 // (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(v))).getVector();
184 // Teuchos::RCP<std::vector<Real> > hvp =
185 // Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(hv)).getVector());
186  Real diff1 = 0.0;
187  Real diff2 = 0.0;
188  for (int i = 0; i < this->nz_; i++) {
189  if ( i == 0 ) {
190  diff1 = (*zp)[i]-(*zp)[i+1];
191  diff1 = this->eps_/std::pow(std::pow(diff1,2.0)+this->eps_,3.0/2.0);
192  (*hvp)[i] = this->alpha_* ((*vp)[i]*diff1 - (*vp)[i+1]*diff1);
193  }
194  else if ( i == this->nz_-1 ) {
195  diff2 = (*zp)[i-1]-(*zp)[i];
196  diff2 = this->eps_/std::pow(std::pow(diff2,2.0)+this->eps_,3.0/2.0);
197  (*hvp)[i] = this->alpha_* (-(*vp)[i-1]*diff2 + (*vp)[i]*diff2);
198  }
199  else {
200  diff1 = (*zp)[i]-(*zp)[i+1];
201  diff1 = this->eps_/std::pow(std::pow(diff1,2.0)+this->eps_,3.0/2.0);
202  diff2 = (*zp)[i-1]-(*zp)[i];
203  diff2 = this->eps_/std::pow(std::pow(diff2,2.0)+this->eps_,3.0/2.0);
204  (*hvp)[i] = this->alpha_* (-(*vp)[i-1]*diff2 + (*vp)[i]*(diff1 + diff2) - (*vp)[i+1]*diff1);
205  }
206  }
207  }
208  }
209 
210  /* FINITE ELEMENT DEFINTIONS */
211  void apply_mass(Vector<Real> &Mf, const Vector<Real> &f ) {
212  Teuchos::RCP<const std::vector<Real> > fp = ROL::StdVector_Helper::constDownCast(f);
213  Teuchos::RCP<std::vector<Real> > Mfp = ROL::StdVector_Helper::downCast(Mf);
214 // Teuchos::RCP<const std::vector<Real> > fp =
215 // (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(f))).getVector();
216 // Teuchos::RCP<std::vector<Real> > Mfp =
217 // Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(Mf)).getVector());
218 
219  for (int i = 0; i < this->nu_; i++) {
220  if ( i == 0 ) {
221  (*Mfp)[i] = this->hu_/6.0*(4.0*(*fp)[i] + (*fp)[i+1]);
222  }
223  else if ( i == this->nu_-1 ) {
224  (*Mfp)[i] = this->hu_/6.0*((*fp)[i-1] + 4.0*(*fp)[i]);
225  }
226  else {
227  (*Mfp)[i] = this->hu_/6.0*((*fp)[i-1] + 4.0*(*fp)[i] + (*fp)[i+1]);
228  }
229  }
230  }
231 
233  Teuchos::RCP<const std::vector<Real> > zp = ROL::StdVector_Helper::constDownCast(z);
234  Teuchos::RCP<std::vector<Real> > up = ROL::StdVector_Helper::downCast(u);
235  Teuchos::RCP<std::vector<Real> > bp = ROL::StdVector_Helper::downCast(b);
236 // Teuchos::RCP<const std::vector<Real> > zp =
237 // (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(z))).getVector();
238 // Teuchos::RCP<std::vector<Real> > up =
239 // Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(u)).getVector());
240 // Teuchos::RCP<std::vector<Real> > bp =
241 // Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(b)).getVector());
242 
243  // Get Diagonal and Off-Diagonal Entries of PDE Jacobian
244  std::vector<Real> d(this->nu_,1.0);
245  std::vector<Real> o(this->nu_-1,1.0);
246  for ( int i = 0; i < this->nu_; i++ ) {
247  d[i] = (std::exp((*zp)[i]) + std::exp((*zp)[i+1]))/this->hu_;
248  if ( i < this->nu_-1 ) {
249  o[i] *= -std::exp((*zp)[i+1])/this->hu_;
250  }
251  }
252 
253  // Solve Tridiagonal System Using LAPACK's SPD Tridiagonal Solver
254  Teuchos::LAPACK<int,Real> lp;
255  int info;
256  int ldb = this->nu_;
257  int nhrs = 1;
258  lp.PTTRF(this->nu_,&d[0],&o[0],&info);
259  lp.PTTRS(this->nu_,nhrs,&d[0],&o[0],&(*bp)[0],ldb,&info);
260  u.set(b);
261  }
262 
263  Real evaluate_target(Real x) {
264  return x*(1.0-x);
265  }
266 
268  const Vector<Real> &d, const Vector<Real> &u ) {
269  Teuchos::RCP<const std::vector<Real> > zp = ROL::StdVector_Helper::constDownCast(z);
270  Teuchos::RCP<const std::vector<Real> > up = ROL::StdVector_Helper::constDownCast(u);
271  Teuchos::RCP<const std::vector<Real> > dp = ROL::StdVector_Helper::constDownCast(d);
272  Teuchos::RCP<std::vector<Real> > Bdp = ROL::StdVector_Helper::downCast(Bd);
273 // Teuchos::RCP<const std::vector<Real> > zp =
274 // (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(z))).getVector();
275 // Teuchos::RCP<const std::vector<Real> > up =
276 // (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(u))).getVector();
277 // Teuchos::RCP<const std::vector<Real> > dp =
278 // (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(d))).getVector();
279 // Teuchos::RCP<std::vector<Real> > Bdp =
280 // Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(Bd)).getVector());
281 
282  for (int i = 0; i < this->nu_; i++) {
283  if ( i == 0 ) {
284  (*Bdp)[i] = 1.0/this->hu_*( std::exp((*zp)[i])*(*up)[i]*(*dp)[i]
285  + std::exp((*zp)[i+1])*((*up)[i]-(*up)[i+1])*(*dp)[i+1] );
286  }
287  else if ( i == this->nu_-1 ) {
288  (*Bdp)[i] = 1.0/this->hu_*( std::exp((*zp)[i])*((*up)[i]-(*up)[i-1])*(*dp)[i]
289  + std::exp((*zp)[i+1])*(*up)[i]*(*dp)[i+1] );
290  }
291  else {
292  (*Bdp)[i] = 1.0/this->hu_*( std::exp((*zp)[i])*((*up)[i]-(*up)[i-1])*(*dp)[i]
293  + std::exp((*zp)[i+1])*((*up)[i]-(*up)[i+1])*(*dp)[i+1] );
294  }
295  }
296  }
297 
299  const Vector<Real> &d, const Vector<Real> &u ) {
300  Teuchos::RCP<const std::vector<Real> > zp = ROL::StdVector_Helper::constDownCast(z);
301  Teuchos::RCP<const std::vector<Real> > up = ROL::StdVector_Helper::constDownCast(u);
302  Teuchos::RCP<const std::vector<Real> > dp = ROL::StdVector_Helper::constDownCast(d);
303  Teuchos::RCP<std::vector<Real> > Bdp = ROL::StdVector_Helper::downCast(Bd);
304 // Teuchos::RCP<const std::vector<Real> > zp =
305 // (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(z))).getVector();
306 // Teuchos::RCP<const std::vector<Real> > up =
307 // (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(u))).getVector();
308 // Teuchos::RCP<const std::vector<Real> > dp =
309 // (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(d))).getVector();
310 // Teuchos::RCP<std::vector<Real> > Bdp =
311 // Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(Bd)).getVector());
312 
313  for (int i = 0; i < this->nz_; i++) {
314  if ( i == 0 ) {
315  (*Bdp)[i] = std::exp((*zp)[i])/this->hu_*(*up)[i]*(*dp)[i];
316  }
317  else if ( i == this->nz_-1 ) {
318  (*Bdp)[i] = std::exp((*zp)[i])/this->hu_*(*up)[i-1]*(*dp)[i-1];
319  }
320  else {
321  (*Bdp)[i] = std::exp((*zp)[i])/this->hu_*( ((*up)[i]-(*up)[i-1])*((*dp)[i]-(*dp)[i-1]) );
322  }
323  }
324  }
325 
327  const Vector<Real> &d, const Vector<Real> &u ) {
328  Teuchos::RCP<const std::vector<Real> > zp = ROL::StdVector_Helper::constDownCast(z);
329  Teuchos::RCP<const std::vector<Real> > vp = ROL::StdVector_Helper::constDownCast(v);
330  Teuchos::RCP<const std::vector<Real> > up = ROL::StdVector_Helper::constDownCast(u);
331  Teuchos::RCP<const std::vector<Real> > dp = ROL::StdVector_Helper::constDownCast(d);
332  Teuchos::RCP<std::vector<Real> > Bdp = ROL::StdVector_Helper::downCast(Bd);
333 // Teuchos::RCP<const std::vector<Real> > zp =
334 // (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(z))).getVector();
335 // Teuchos::RCP<const std::vector<Real> > vp =
336 // (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(v))).getVector();
337 // Teuchos::RCP<const std::vector<Real> > up =
338 // (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(u))).getVector();
339 // Teuchos::RCP<const std::vector<Real> > dp =
340 // (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(d))).getVector();
341 // Teuchos::RCP<std::vector<Real> > Bdp =
342 // Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(Bd)).getVector());
343 
344  for (int i = 0; i < this->nz_; i++) {
345  if ( i == 0 ) {
346  (*Bdp)[i] = (*vp)[i]*std::exp((*zp)[i])/this->hu_*(*up)[i]*(*dp)[i];
347  }
348  else if ( i == this->nz_-1 ) {
349  (*Bdp)[i] = (*vp)[i]*std::exp((*zp)[i])/this->hu_*(*up)[i-1]*(*dp)[i-1];
350  }
351  else {
352  (*Bdp)[i] = (*vp)[i]*std::exp((*zp)[i])/this->hu_*( ((*up)[i]-(*up)[i-1])*((*dp)[i]-(*dp)[i-1]) );
353  }
354  }
355  }
356 
357  /* STATE AND ADJOINT EQUATION DEFINTIONS */
359  Real k1 = 1.0;
360  Real k2 = 2.0;
361  // Right Hand Side
362  Teuchos::RCP<std::vector<Real> > bp = Teuchos::rcp( new std::vector<Real> (this->nu_, 0.0) );
363  for ( int i = 0; i < this->nu_; i++ ) {
364  if ( (Real)(i+1)*this->hu_ < 0.5 ) {
365  (*bp)[i] = 2.0*k1*this->hu_;
366  }
367  else if ( std::abs((Real)(i+1)*this->hu_ - 0.5) < ROL_EPSILON ) {
368  (*bp)[i] = (k1+k2)*this->hu_;
369  }
370  else if ( (Real)(i+1)*this->hu_ > 0.5 ) {
371  (*bp)[i] = 2.0*k2*this->hu_;
372  }
373  }
374 
375  StdVector<Real> b(bp);
376  // Solve Equation
377  this->solve_poisson(u,z,b);
378  }
379 
381  Teuchos::RCP<const std::vector<Real> > up = ROL::StdVector_Helper::constDownCast(u);
382  StdVector<Real> res( Teuchos::rcp( new std::vector<Real>(this->nu_,0.0) ) );
383  Teuchos::RCP<std::vector<Real> > rp = ROL::StdVector_Helper::downCast(res);
384 // Teuchos::RCP<const std::vector<Real> > up =
385 // (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(u))).getVector();
386 // StdVector<Real> res( Teuchos::rcp( new std::vector<Real>(this->nu_,0.0) ) );
387 // Teuchos::RCP<std::vector<Real> > rp =
388 // Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(res)).getVector());
389  for (int i = 0; i < this->nu_; i++) {
390  (*rp)[i] = -((*up)[i]-this->evaluate_target((Real)(i+1)*this->hu_));
391  }
392  StdVector<Real> Mres( Teuchos::rcp( new std::vector<Real>(this->nu_,0.0) ) );
393  this->apply_mass(Mres,res);
394  this->solve_poisson(p,z,Mres);
395  }
396 
398  const Vector<Real> &u, const Vector<Real> &z) {
399  StdVector<Real> b( Teuchos::rcp( new std::vector<Real>(this->nu_,0.0) ) );
400  this->apply_linearized_control_operator(b,z,v,u);
401  this->solve_poisson(w,z,b);
402  }
403 
405  const Vector<Real> &p, const Vector<Real> &u, const Vector<Real> &z) {
406  StdVector<Real> res( Teuchos::rcp( new std::vector<Real>(this->nu_,0.0) ) );
407  this->apply_mass(res,w);
408  StdVector<Real> res1( Teuchos::rcp( new std::vector<Real>(this->nu_,0.0) ) );
409  this->apply_linearized_control_operator(res1,z,v,p);
410  res.axpy(-1.0,res1);
411  this->solve_poisson(q,z,res);
412  }
413 
414  /* OBJECTIVE FUNCTION DEFINITIONS */
415  Real value( const Vector<Real> &z, Real &tol ) {
416  // SOLVE STATE EQUATION
417  StdVector<Real> u( Teuchos::rcp( new std::vector<Real>(this->nu_,0.0) ) );
418  Teuchos::RCP<std::vector<Real> > up = ROL::StdVector_Helper::downCast(u);
419 // Teuchos::RCP<std::vector<Real> > up =
420 // Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(u)).getVector());
421  this->solve_state_equation(u,z);
422 
423  // COMPUTE MISFIT
424  StdVector<Real> res( Teuchos::rcp( new std::vector<Real>(this->nu_,0.0) ) );
425  Teuchos::RCP<std::vector<Real> > rp = ROL::StdVector_Helper::downCast(res);
426 // Teuchos::RCP<std::vector<Real> > rp =
427 // Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(res)).getVector());
428  for (int i = 0; i < this->nu_; i++) {
429  (*rp)[i] = ((*up)[i]-this->evaluate_target((Real)(i+1)*this->hu_));
430  }
431  StdVector<Real> Mres( Teuchos::rcp( new std::vector<Real>(this->nu_,0.0) ) );
432  this->apply_mass(Mres,res);
433  return 0.5*Mres.dot(res) + this->reg_value(z);
434  }
435 
436  void gradient( Vector<Real> &g, const Vector<Real> &z, Real &tol ) {
437  // SOLVE STATE EQUATION
438  StdVector<Real> u( Teuchos::rcp( new std::vector<Real>(this->nu_,0.0) ) );
439  this->solve_state_equation(u,z);
440 
441  // SOLVE ADJOINT EQUATION
442  StdVector<Real> p( Teuchos::rcp( new std::vector<Real>(this->nu_,0.0) ) );
443  this->solve_adjoint_equation(p,u,z);
444 
445  // Apply Transpose of Linearized Control Operator
447 
448  // Regularization gradient
449  StdVector<Real> g_reg( Teuchos::rcp( new std::vector<Real>(this->nz_,0.0) ) );
450  this->reg_gradient(g_reg,z);
451 
452  // Build Gradient
453  g.plus(g_reg);
454  }
455 #if USE_HESSVEC
456  void hessVec( Vector<Real> &hv, const Vector<Real> &v, const Vector<Real> &z, Real &tol ) {
457  // SOLVE STATE EQUATION
458  StdVector<Real> u( Teuchos::rcp( new std::vector<Real>(this->nu_,0.0) ) );
459  this->solve_state_equation(u,z);
460 
461  // SOLVE ADJOINT EQUATION
462  StdVector<Real> p( Teuchos::rcp( new std::vector<Real>(this->nu_,0.0) ) );
463  this->solve_adjoint_equation(p,u,z);
464 
465  // SOLVE STATE SENSITIVITY EQUATION
466  StdVector<Real> w( Teuchos::rcp( new std::vector<Real>(this->nu_,0.0) ) );
467  this->solve_state_sensitivity_equation(w,v,u,z);
468 
469  // SOLVE ADJOINT SENSITIVITY EQUATION
470  StdVector<Real> q( Teuchos::rcp( new std::vector<Real>(this->nu_,0.0) ) );
471  this->solve_adjoint_sensitivity_equation(q,w,v,p,u,z);
472 
473  // Apply Transpose of Linearized Control Operator
475 
476  // Apply Transpose of Linearized Control Operator
477  StdVector<Real> tmp( Teuchos::rcp( new std::vector<Real>(this->nz_,0.0) ) );
479  hv.axpy(-1.0,tmp);
480 
481  // Apply Transpose of 2nd Derivative of Control Operator
482  tmp.zero();
484  hv.plus(tmp);
485 
486  // Regularization hessVec
487  StdVector<Real> hv_reg( Teuchos::rcp( new std::vector<Real>(this->nz_,0.0) ) );
488  this->reg_hessVec(hv_reg,v,z);
489 
490  // Build hessVec
491  hv.plus(hv_reg);
492  }
493 #endif
494 
495  void invHessVec( Vector<Real> &hv, const Vector<Real> &v, const Vector<Real> &x, Real &tol ) {
496 
497  // Cast hv and v vectors to std::vector.
498  Teuchos::RCP<std::vector<Real> > hvp = ROL::StdVector_Helper::downCast(hv);
499  std::vector<Real> vp(*ROL::StdVector_Helper::constDownCast(v));
500 // Teuchos::RCP<std::vector<Real> > hvp =
501 // Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(hv)).getVector());
502 // Teuchos::RCP<std::vector<Real> > vp =
503 // Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(v))).getVector());
504 
505  int dim = vp.size();
506 
507  // Compute dense Hessian.
508  Teuchos::SerialDenseMatrix<int, Real> H(dim, dim);
509  Objective_PoissonInversion<Real> & obj = *this;
510  H = computeDenseHessian<Real>(obj, x);
511 
512  // Compute eigenvalues, sort real part.
513  std::vector<std::vector<Real> > eigenvals = computeEigenvalues<Real>(H);
514  std::sort((eigenvals[0]).begin(), (eigenvals[0]).end());
515 
516  // Perform 'inertia' correction.
517  Real inertia = (eigenvals[0])[0];
518  Real correction = 0.0;
519  if ( inertia <= 0.0 ) {
520  correction = 2.0*std::abs(inertia);
521  if ( inertia == 0.0 ) {
522  int cnt = 0;
523  while ( eigenvals[0][cnt] == 0.0 ) {
524  cnt++;
525  }
526  correction = 0.5*eigenvals[0][cnt];
527  if ( cnt == dim-1 ) {
528  correction = 1.0;
529  }
530  }
531  for (int i=0; i<dim; i++) {
532  H(i,i) += correction;
533  }
534  }
535 
536  // Compute dense inverse Hessian.
537  Teuchos::SerialDenseMatrix<int, Real> invH = computeInverse<Real>(H);
538 
539  // Apply dense inverse Hessian.
540  Teuchos::SerialDenseVector<int, Real> hv_teuchos(Teuchos::View, &((*hvp)[0]), dim);
541  const Teuchos::SerialDenseVector<int, Real> v_teuchos(Teuchos::View, &(vp[0]), dim);
542  hv_teuchos.multiply(Teuchos::NO_TRANS, Teuchos::NO_TRANS, 1.0, invH, v_teuchos, 0.0);
543  }
544 
545  };
546 
547  template<class Real>
548  void getPoissonInversion( Teuchos::RCP<Objective<Real> > &obj, Vector<Real> &x0, Vector<Real> &x ) {
549  // Cast Initial Guess and Solution Vectors
550  Teuchos::RCP<std::vector<Real> > x0p = ROL::StdVector_Helper::downCast(x0);
551  Teuchos::RCP<std::vector<Real> > xp = ROL::StdVector_Helper::downCast(x);
552 // Teuchos::RCP<std::vector<Real> > x0p =
553 // Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(x0)).getVector());
554 // Teuchos::RCP<std::vector<Real> > xp =
555 // Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(x)).getVector());
556  int n = xp->size();
557  // Resize Vectors
558  n = 128;
559  x0p->resize(n);
560  xp->resize(n);
561  // Instantiate Objective Function
562  obj = Teuchos::rcp( new Objective_PoissonInversion<Real>(n,1.e-6) );
563  // Get Initial Guess
564  for (int i=0; i<n; i++) {
565  (*x0p)[i] = 1.5;
566  }
567  // Get Solution
568  Real h = 1.0/((Real)n+1);
569  Real pt = 0.0;
570  Real k1 = 1.0;
571  Real k2 = 2.0;
572  for( int i=0; i<n; i++ ) {
573  pt = (Real)(i+1)*h;
574  if ( pt >= 0.0 && pt < 0.5 ) {
575  (*xp)[i] = std::log(k1);
576  }
577  else if ( pt >= 0.5 && pt < 1.0 ) {
578  (*xp)[i] = std::log(k2);
579  }
580  }
581  }
582 
583 } // End ZOO Namespace
584 } // End ROL Namespace
585 
586 #endif
Provides the interface to evaluate objective functions.
void apply_transposed_linearized_control_operator(Vector< Real > &Bd, const Vector< Real > &z, const Vector< Real > &d, const Vector< Real > &u)
virtual void scale(const Real alpha)=0
Compute where .
void getPoissonInversion(Teuchos::RCP< Objective< Real > > &obj, Vector< Real > &x0, Vector< Real > &x)
void solve_state_sensitivity_equation(Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z)
virtual void plus(const Vector &x)=0
Compute , where .
virtual void axpy(const Real alpha, const Vector &x)
Compute where .
Definition: ROL_Vector.hpp:141
void gradient(Vector< Real > &g, const Vector< Real > &z, Real &tol)
Compute gradient.
Real dot(const Vector< Real > &x) const
Compute where .
virtual void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply Hessian approximation to vector.
Contains definitions for helper functions in ROL.
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:72
void solve_poisson(Vector< Real > &u, const Vector< Real > &z, Vector< Real > &b)
void solve_adjoint_sensitivity_equation(Vector< Real > &q, const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &p, const Vector< Real > &u, const Vector< Real > &z)
void apply_transposed_linearized_control_operator_2(Vector< Real > &Bd, const Vector< Real > &z, const Vector< Real > &v, const Vector< Real > &d, const Vector< Real > &u)
Provides the std::vector implementation of the ROL::Vector interface.
void invHessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply inverse Hessian approximation to vector.
void solve_state_equation(Vector< Real > &u, const Vector< Real > &z)
void reg_hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &z)
Teuchos::RCP< std::vector< Real > > downCast(Vector< Real > &x)
void apply_mass(Vector< Real > &Mf, const Vector< Real > &f)
void solve_adjoint_equation(Vector< Real > &p, const Vector< Real > &u, const Vector< Real > &z)
Objective_PoissonInversion(int nz=32, Real alpha=1.e-4)
void apply_linearized_control_operator(Vector< Real > &Bd, const Vector< Real > &z, const Vector< Real > &d, const Vector< Real > &u)
Teuchos::RCP< const std::vector< Real > > constDownCast(const Vector< Real > &x)
Real value(const Vector< Real > &z, Real &tol)
Compute value.
virtual void set(const Vector &x)
Set where .
Definition: ROL_Vector.hpp:194
Real reg_value(const Vector< Real > &z)
void reg_gradient(Vector< Real > &g, const Vector< Real > &z)
static const double ROL_EPSILON
Platform-dependent machine epsilon.
Definition: ROL_Types.hpp:115