ROL
ROL_PoissonControl.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_POISSONCONTROL_HPP
54 #define ROL_POISSONCONTROL_HPP
55 
56 #include "ROL_StdVector.hpp"
57 #include "ROL_Objective.hpp"
58 
59 namespace ROL {
60 namespace ZOO {
61 
64  template<class Real>
65  class Objective_PoissonControl : public Objective<Real> {
66  private:
67  Real alpha_;
68 
69  public:
70 
71  Objective_PoissonControl(Real alpha = 1.e-4) : alpha_(alpha) {}
72 
73  void apply_mass(Vector<Real> &Mz, const Vector<Real> &z ) {
74  Teuchos::RCP<const std::vector<Real> > zp =
75  (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(z))).getVector();
76  Teuchos::RCP<std::vector<Real> > Mzp =
77  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(Mz)).getVector());
78 
79  int n = zp->size();
80  Real h = 1.0/((Real)n+1.0);
81  for (int i=0; i<n; i++) {
82  if ( i == 0 ) {
83  (*Mzp)[i] = h/6.0*(4.0*(*zp)[i] + (*zp)[i+1]);
84  }
85  else if ( i == n-1 ) {
86  (*Mzp)[i] = h/6.0*((*zp)[i-1] + 4.0*(*zp)[i]);
87  }
88  else {
89  (*Mzp)[i] = h/6.0*((*zp)[i-1] + 4.0*(*zp)[i] + (*zp)[i+1]);
90  }
91  }
92  }
93 
94  void solve_poisson(Vector<Real> & u, const Vector<Real> & z) {
95  Teuchos::RCP<std::vector<Real> > up =
96  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(u)).getVector());
97 
98  int n = up->size();
99  Real h = 1.0/((Real)n+1.0);
100  StdVector<Real> b( Teuchos::rcp( new std::vector<Real>(n,0.0) ) );
101  Teuchos::RCP<std::vector<Real> > bp =
102  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(b)).getVector());
103  this->apply_mass(b,z);
104 
105  Real d = 2.0/h;
106  Real o = -1.0/h;
107  Real m = 0.0;
108  std::vector<Real> c(n,o);
109  c[0] = c[0]/d;
110  (*up)[0] = (*bp)[0]/d;
111  for ( int i = 1; i < n; i++ ) {
112  m = 1.0/(d - o*c[i-1]);
113  c[i] = c[i]*m;
114  (*up)[i] = ( (*bp)[i] - o*(*up)[i-1] )*m;
115  }
116  for ( int i = n-1; i > 0; i-- ) {
117  (*up)[i-1] = (*up)[i-1] - c[i-1]*(*up)[i];
118  }
119  }
120 
121  Real evaluate_target(Real x) {
122  Real val = 1.0/3.0*std::pow(x,4.0) - 2.0/3.0*std::pow(x,3.0) + 1.0/3.0*x + 8.0*this->alpha_;
123  return val;
124  }
125 
126  Real value( const Vector<Real> &z, Real &tol ) {
127  Teuchos::RCP<const std::vector<Real> > zp =
128  (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(z))).getVector();
129  int n = zp->size();
130  Real h = 1.0/((Real)n+1.0);
131  // SOLVE STATE EQUATION
132  StdVector<Real> u( Teuchos::rcp( new std::vector<Real>(n,0.0) ) );
133  this->solve_poisson(u,z);
134  Teuchos::RCP<std::vector<Real> > up =
135  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(u)).getVector());
136 
137  Real val = 0.0;
138  Real res = 0.0;
139  Real res1 = 0.0;
140  Real res2 = 0.0;
141  Real res3 = 0.0;
142  for (int i=0; i<n; i++) {
143  res = this->alpha_*(*zp)[i];
144  if ( i == 0 ) {
145  res *= h/6.0*(4.0*(*zp)[i] + (*zp)[i+1]);
146  res1 = (*up)[i]-evaluate_target((Real)(i+1)*h);
147  res2 = (*up)[i+1]-evaluate_target((Real)(i+2)*h);
148  res += h/6.0*(4.0*res1 + res2)*res1;
149  }
150  else if ( i == n-1 ) {
151  res *= h/6.0*((*zp)[i-1] + 4.0*(*zp)[i]);
152  res1 = (*up)[i-1]-evaluate_target((Real)(i)*h);
153  res2 = (*up)[i]-evaluate_target((Real)(i+1)*h);
154  res += h/6.0*(res1 + 4.0*res2)*res2;
155  }
156  else {
157  res *= h/6.0*((*zp)[i-1] + 4.0*(*zp)[i] + (*zp)[i+1]);
158  res1 = (*up)[i-1]-evaluate_target((Real)(i)*h);
159  res2 = (*up)[i]-evaluate_target((Real)(i+1)*h);
160  res3 = (*up)[i+1]-evaluate_target((Real)(i+2)*h);
161  res += h/6.0*(res1 + 4.0*res2 + res3)*res2;
162  }
163  val += 0.5*res;
164  }
165  return val;
166  }
167 
168  void gradient( Vector<Real> &g, const Vector<Real> &z, Real &tol ) {
169  Teuchos::RCP<const std::vector<Real> > zp =
170  (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(z))).getVector();
171  Teuchos::RCP<std::vector<Real> > gp =
172  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(g)).getVector());
173  int n = zp->size();
174  Real h = 1.0/((Real)n+1.0);
175 
176  // SOLVE STATE EQUATION
177  StdVector<Real> u( Teuchos::rcp( new std::vector<Real>(n,0.0) ) );
178  this->solve_poisson(u,z);
179  Teuchos::RCP<std::vector<Real> > up =
180  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(u)).getVector());
181 
182  // SOLVE ADJOINT EQUATION
183  StdVector<Real> res( Teuchos::rcp( new std::vector<Real>(n,0.0) ) );
184  Teuchos::RCP<std::vector<Real> > rp =
185  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(res)).getVector());
186  for (int i=0; i<n; i++) {
187  (*rp)[i] = -((*up)[i]-evaluate_target((Real)(i+1)*h));
188  }
189  StdVector<Real> p( Teuchos::rcp( new std::vector<Real>(n,0.0) ) );
190  this->solve_poisson(p,res);
191  Teuchos::RCP<std::vector<Real> > pp =
192  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(p)).getVector());
193 
194  Real res1 = 0.0;
195  Real res2 = 0.0;
196  Real res3 = 0.0;
197  for (int i=0; i<n; i++) {
198  if ( i == 0 ) {
199  res1 = this->alpha_*(*zp)[i] - (*pp)[i];
200  res2 = this->alpha_*(*zp)[i+1] - (*pp)[i+1];
201  (*gp)[i] = h/6.0*(4.0*res1 + res2);
202  }
203  else if ( i == n-1 ) {
204  res1 = this->alpha_*(*zp)[i-1] - (*pp)[i-1];
205  res2 = this->alpha_*(*zp)[i] - (*pp)[i];
206  (*gp)[i] = h/6.0*(res1 + 4.0*res2);
207  }
208  else {
209  res1 = this->alpha_*(*zp)[i-1] - (*pp)[i-1];
210  res2 = this->alpha_*(*zp)[i] - (*pp)[i];
211  res3 = this->alpha_*(*zp)[i+1] - (*pp)[i+1];
212  (*gp)[i] = h/6.0*(res1 + 4.0*res2 + res3);
213  }
214  }
215  }
216 #if USE_HESSVEC
217  void hessVec( Vector<Real> &hv, const Vector<Real> &v, const Vector<Real> &z, Real &tol ) {
218  Teuchos::RCP<const std::vector<Real> > zp =
219  (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(z))).getVector();
220  Teuchos::RCP<const std::vector<Real> > vp =
221  (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(v))).getVector();
222  Teuchos::RCP<std::vector<Real> > hvp =
223  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(hv)).getVector());
224 
225  int n = zp->size();
226  Real h = 1.0/((Real)n+1.0);
227 
228  // SOLVE STATE EQUATION
229  StdVector<Real> u( Teuchos::rcp( new std::vector<Real>(n,0.0) ) );
230  this->solve_poisson(u,v);
231  Teuchos::RCP<std::vector<Real> > up =
232  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(u)).getVector());
233 
234  // SOLVE ADJOINT EQUATION
235  StdVector<Real> p( Teuchos::rcp( new std::vector<Real>(n,0.0) ) );
236  this->solve_poisson(p,u);
237  Teuchos::RCP<std::vector<Real> > pp =
238  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(p)).getVector());
239 
240  Real res1 = 0.0;
241  Real res2 = 0.0;
242  Real res3 = 0.0;
243  for (int i=0; i<n; i++) {
244  if ( i == 0 ) {
245  res1 = this->alpha_*(*vp)[i] + (*pp)[i];
246  res2 = this->alpha_*(*vp)[i+1] + (*pp)[i+1];
247  (*hvp)[i] = h/6.0*(4.0*res1 + res2);
248  }
249  else if ( i == n-1 ) {
250  res1 = this->alpha_*(*vp)[i-1] + (*pp)[i-1];
251  res2 = this->alpha_*(*vp)[i] + (*pp)[i];
252  (*hvp)[i] = h/6.0*(res1 + 4.0*res2);
253  }
254  else {
255  res1 = this->alpha_*(*vp)[i-1] + (*pp)[i-1];
256  res2 = this->alpha_*(*vp)[i] + (*pp)[i];
257  res3 = this->alpha_*(*vp)[i+1] + (*pp)[i+1];
258  (*hvp)[i] = h/6.0*(res1 + 4.0*res2 + res3);
259  }
260  }
261  }
262 #endif
263  };
264 
265  template<class Real>
266  void getPoissonControl( Teuchos::RCP<Objective<Real> > &obj, Vector<Real> &x0, Vector<Real> &x ) {
267  // Cast Initial Guess and Solution Vectors
268  Teuchos::RCP<std::vector<Real> > x0p =
269  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(x0)).getVector());
270  Teuchos::RCP<std::vector<Real> > xp =
271  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(x)).getVector());
272  int n = xp->size();
273  // Resize Vectors
274  n = 512;
275  x0p->resize(n);
276  xp->resize(n);
277  // Instantiate Objective Function
278  obj = Teuchos::rcp( new Objective_PoissonControl<Real> );
279  // Get Initial Guess
280  for (int i=0; i<n; i++) {
281  (*x0p)[i] = 0.0;
282  }
283  // Get Solution
284  Real h = 1.0/((Real)n+1.0);
285  Real pt = 0.0;
286  for( int i=0; i<n; i++ ) {
287  pt = (Real)(i+1)*h;
288  (*xp)[i] = 4.0*pt*(1.0-pt);
289  }
290  }
291 
292 } // End ZOO Namespace
293 } // End ROL Namespace
294 
295 #endif
Provides the interface to evaluate objective functions.
Real value(const Vector< Real > &z, Real &tol)
Compute value.
void getPoissonControl(Teuchos::RCP< Objective< Real > > &obj, Vector< Real > &x0, Vector< Real > &x)
virtual void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply Hessian approximation to vector.
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:72
void solve_poisson(Vector< Real > &u, const Vector< Real > &z)
Provides the std::vector implementation of the ROL::Vector interface.
void gradient(Vector< Real > &g, const Vector< Real > &z, Real &tol)
Compute gradient.
Poisson distributed control.
void apply_mass(Vector< Real > &Mz, const Vector< Real > &z)