RedistributeDataHandles.hpp
Go to the documentation of this file.
1 /*
2  Copyright 2015 Dr. Blatt - HPC-Simulation-Software & Services.
3  Coypright 2015 NTNU
4  Copyright 2015 Statoil AS
5  Copyright 2015 IRIS AS
6 
7  This file is part of the Open Porous Media project (OPM).
8 
9  OPM is free software: you can redistribute it and/or modify
10  it under the terms of the GNU General Public License as published by
11  the Free Software Foundation, either version 3 of the License, or
12  (at your option) any later version.
13 
14  OPM is distributed in the hope that it will be useful,
15  but WITHOUT ANY WARRANTY; without even the implied warranty of
16  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17  GNU General Public License for more details.
18 
19  You should have received a copy of the GNU General Public License
20  along with OPM. If not, see <http://www.gnu.org/licenses/>.
21 */
22 #ifndef OPM_REDISTRIBUTEDATAHANDLES_HEADER
23 #define OPM_REDISTRIBUTEDATAHANDLES_HEADER
24 
25 #include <opm/core/simulator/BlackoilState.hpp>
26 
29 
30 #include<boost/any.hpp>
31 
32 namespace Opm
33 {
34 
35 template <class Grid>
36 inline void distributeGridAndData( Grid& ,
37  EclipseStateConstPtr ,
38  BlackoilState& ,
41  boost::any& ,
42  const bool )
43 {
44 }
45 
46 #if HAVE_DUNE_CORNERPOINT
47 class GeologyDataHandle
49 {
50 public:
52  typedef double DataType;
58  GeologyDataHandle(const Dune::CpGrid& sendGrid,
59  const Dune::CpGrid& recvGrid,
60  const DerivedGeology& sendGeology,
61  DerivedGeology& recvGeology)
62  : sendGrid_(sendGrid), recvGrid_(recvGrid), sendGeology_(sendGeology),
63  recvGeology_(recvGeology)
64  {}
65 
66  bool fixedsize(int /*dim*/, int /*codim*/)
67  {
68  return false;
69  }
70  template<class T>
71  std::size_t size(const T& e)
72  {
73  if ( T::codimension == 0)
74  {
75  return 1 + sendGrid_.numCellFaces(e.index());
76  }
77  else
78  {
79  OPM_THROW(std::logic_error, "Data handle can only be used for elements");
80  }
81  }
82  template<class B, class T>
83  void gather(B& buffer, const T& e)
84  {
85  assert( T::codimension == 0);
86  buffer.write(sendGeology_.poreVolume()[e.index()]);
87  for ( int i=0; i< sendGrid_.numCellFaces(e.index()); ++i )
88  {
89  buffer.write(sendGeology_.transmissibility()[sendGrid_.cellFace(e.index(), i)]);
90  }
91  }
92  template<class B, class T>
93  void scatter(B& buffer, const T& e, std::size_t /* size */)
94  {
95  assert( T::codimension == 0);
96  double val;
97  buffer.read(val);
98  recvGeology_.poreVolume()[e.index()]=val;
99  for ( int i=0; i< recvGrid_.numCellFaces(e.index()); ++i )
100  {
101  buffer.read(val);
102  recvGeology_.transmissibility()[recvGrid_.cellFace(e.index(), i)]=val;
103  }
104  }
105  bool contains(int dim, int codim)
106  {
107  return dim==3 && codim==0;
108  }
109 private:
111  const Dune::CpGrid& sendGrid_;
113  const Dune::CpGrid& recvGrid_;
115  const DerivedGeology& sendGeology_;
117  DerivedGeology& recvGeology_;
118 };
119 
121 class BlackoilStateDataHandle
122 {
123 public:
125  typedef double DataType;
131  BlackoilStateDataHandle(const Dune::CpGrid& sendGrid,
132  const Dune::CpGrid& recvGrid,
133  const BlackoilState& sendState,
134  BlackoilState& recvState)
135  : sendGrid_(sendGrid), recvGrid_(recvGrid), sendState_(sendState), recvState_(recvState)
136  {}
137 
138  bool fixedsize(int /*dim*/, int /*codim*/)
139  {
140  return false;
141  }
142 
143  template<class T>
144  std::size_t size(const T& e)
145  {
146  if ( T::codimension == 0)
147  {
148  return 2 * sendState_.numPhases() +4+2*sendGrid_.numCellFaces(e.index());
149  }
150  else
151  {
152  OPM_THROW(std::logic_error, "Data handle can only be used for elements");
153  }
154  }
155 
156  template<class B, class T>
157  void gather(B& buffer, const T& e)
158  {
159  assert( T::codimension == 0);
160 
161  for ( int i=0; i<sendState_.numPhases(); ++i )
162  {
163  buffer.write(sendState_.surfacevol()[e.index()*sendState_.numPhases()+i]);
164  }
165  buffer.write(sendState_.gasoilratio()[e.index()]);
166  buffer.write(sendState_.rv()[e.index()]);
167  buffer.write(sendState_.pressure()[e.index()]);
168  buffer.write(sendState_.temperature()[e.index()]);
169  for ( int i=0; i<sendState_.numPhases(); ++i )
170  {
171  buffer.write(sendState_.saturation()[e.index()*sendState_.numPhases()+i]);
172  }
173  for ( int i=0; i<sendGrid_.numCellFaces(e.index()); ++i )
174  {
175  buffer.write(sendState_.facepressure()[sendGrid_.cellFace(e.index(), i)]);
176  }
177  for ( int i=0; i<sendGrid_.numCellFaces(e.index()); ++i )
178  {
179  buffer.write(sendState_.faceflux()[sendGrid_.cellFace(e.index(), i)]);
180  }
181  }
182  template<class B, class T>
183  void scatter(B& buffer, const T& e, std::size_t size_arg)
184  {
185  assert( T::codimension == 0);
186  assert( int(size_arg) == 2 * recvState_.numPhases() +4+2*recvGrid_.numCellFaces(e.index()));
187  static_cast<void>(size_arg);
188 
189  double val;
190  for ( int i=0; i<recvState_.numPhases(); ++i )
191  {
192  buffer.read(val);
193  recvState_.surfacevol()[e.index()*sendState_.numPhases()+i]=val;
194  }
195  buffer.read(val);
196  recvState_.gasoilratio()[e.index()]=val;
197  buffer.read(val);
198  recvState_.rv()[e.index()]=val;
199  buffer.read(val);
200  recvState_.pressure()[e.index()]=val;
201  buffer.read(val);
202  recvState_.temperature()[e.index()]=val;
203  for ( int i=0; i<recvState_.numPhases(); ++i )
204  {
205  buffer.read(val);
206  recvState_.saturation()[e.index()*sendState_.numPhases()+i]=val;
207  }
208  for ( int i=0; i<recvGrid_.numCellFaces(e.index()); ++i )
209  {
210  buffer.read(val);
211  recvState_.facepressure()[recvGrid_.cellFace(e.index(), i)]=val;
212  }
213  for ( int i=0; i<recvGrid_.numCellFaces(e.index()); ++i )
214  {
215  buffer.read(val);
216  recvState_.faceflux()[recvGrid_.cellFace(e.index(), i)]=val;
217  }
218  }
219  bool contains(int dim, int codim)
220  {
221  return dim==3 && codim==0;
222  }
223 private:
225  const Dune::CpGrid& sendGrid_;
227  const Dune::CpGrid& recvGrid_;
229  const BlackoilState& sendState_;
230  // \brief The state where we will store the received values.
231  BlackoilState& recvState_;
232 };
233 
234 class BlackoilPropsDataHandle
235 {
236 public:
238  typedef double DataType;
242  BlackoilPropsDataHandle(const BlackoilPropsAdFromDeck& sendProps,
243  BlackoilPropsAdFromDeck& recvProps)
244  : sendProps_(sendProps), recvProps_(recvProps),
245  size_(1)
246  {
247  // satOilMax might be non empty. In this case we will need to send it, too.
248  // It has to have the same size as the cellPvtRegionIdx_
249  if ( sendProps.satOilMax_.size()>0 )
250  {
251  recvProps_.satOilMax_.resize(recvProps_.cellPvtRegionIdx_.size(),
252  -std::numeric_limits<double>::max());
253  ++size_;
254  }
255  }
256 
257  bool fixedsize(int /*dim*/, int /*codim*/)
258  {
259  return true;
260  }
261 
262  template<class T>
263  std::size_t size(const T&)
264  {
265  if ( T::codimension == 0)
266  {
267  // We only send cellPvtRegionIdx_, and maybe satOilMax_
268  return size_;
269  }
270  else
271  {
272  OPM_THROW(std::logic_error, "Data handle can only be used for elements");
273  }
274  }
275  template<class B, class T>
276  void gather(B& buffer, const T& e)
277  {
278  assert( T::codimension == 0);
279 
280  buffer.write(sendProps_.cellPvtRegionIndex()[e.index()]);
281  if ( size_ > 1 ) {
282  buffer.write(sendProps_.satOilMax_[e.index()]);
283  }
284  }
285  template<class B, class T>
286  void scatter(B& buffer, const T& e, std::size_t size_arg)
287  {
288  assert( T::codimension == 0);
289  assert( size_arg==size_ ); (void) size_arg;
290  double val;
291  buffer.read(val);
292  recvProps_.cellPvtRegionIdx_[e.index()]=val;
293  if ( size_ > 1 ) {
294  buffer.read(val);
295  recvProps_.satOilMax_[e.index()]=val;
296  }
297  }
298  bool contains(int dim, int codim)
299  {
300  return dim==3 && codim==0;
301  }
302 private:
304  const BlackoilPropsAdFromDeck& sendProps_;
305  // \brief The properties where we will store the received values.
306  BlackoilPropsAdFromDeck& recvProps_;
307  std::size_t size_;
308 };
309 
310 inline
311 void distributeGridAndData( Dune::CpGrid& grid,
312  EclipseStateConstPtr eclipseState,
313  BlackoilState& state,
314  BlackoilPropsAdFromDeck& properties,
315  DerivedGeology& geology,
316  boost::any& parallelInformation,
317  const bool useLocalPerm)
318 {
319  BlackoilState distributed_state;
320 
321  Dune::CpGrid global_grid ( grid );
322  global_grid.switchToGlobalView();
323 
324  // distribute the grid and switch to the distributed view
325  grid.loadBalance(eclipseState);
326  grid.switchToDistributedView();
327 
328  BlackoilPropsAdFromDeck distributed_props(properties, grid.numCells());
329  distributed_state.init(grid.numCells(), grid.numFaces(), state.numPhases());
330  // init does not resize surfacevol. Do it manually.
331  distributed_state.surfacevol().resize(grid.numCells()*state.numPhases(),
332  std::numeric_limits<double>::max());
333  BlackoilStateDataHandle state_handle(global_grid, grid,
334  state, distributed_state);
335  BlackoilPropsDataHandle props_handle(properties,
336  distributed_props);
337  grid.scatterData(state_handle);
338  grid.scatterData(props_handle);
339  // Create a distributed Geology. Some values will be updated using communication
340  // below
341  DerivedGeology distributed_geology(grid,
342  distributed_props, eclipseState,
343  useLocalPerm, geology.gravity());
344  GeologyDataHandle geo_handle(global_grid, grid,
345  geology, distributed_geology);
346  grid.scatterData(geo_handle);
347 
348  // copy states
349  properties = distributed_props;
350  geology = distributed_geology;
351  state = distributed_state;
352 
353  extractParallelGridInformationToISTL(grid, parallelInformation);
354 }
355 #endif
356 
357 } // end namespace Opm
358 #endif
Definition: GeoProps.hpp:53
Definition: AdditionalObjectDeleter.hpp:22
Definition: BlackoilPropsAdFromDeck.hpp:58
void distributeGridAndData(Grid &, EclipseStateConstPtr, BlackoilState &, BlackoilPropsAdFromDeck &, DerivedGeology &, boost::any &, const bool)
Definition: RedistributeDataHandles.hpp:36