initStateEquil.hpp
Go to the documentation of this file.
1 /*
2  Copyright 2014 SINTEF ICT, Applied Mathematics.
3  Copyright 2015 Dr. Blatt - HPC-Simulation-Software & Services
4  Copyright 2015 NTNU
5 
6  This file is part of the Open Porous Media project (OPM).
7 
8  OPM is free software: you can redistribute it and/or modify
9  it under the terms of the GNU General Public License as published by
10  the Free Software Foundation, either version 3 of the License, or
11  (at your option) any later version.
12 
13  OPM is distributed in the hope that it will be useful,
14  but WITHOUT ANY WARRANTY; without even the implied warranty of
15  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16  GNU General Public License for more details.
17 
18  You should have received a copy of the GNU General Public License
19  along with OPM. If not, see <http://www.gnu.org/licenses/>.
20 */
21 
22 #ifndef OPM_INITSTATEEQUIL_HEADER_INCLUDED
23 #define OPM_INITSTATEEQUIL_HEADER_INCLUDED
24 
32 #include <opm/parser/eclipse/Utility/EquilWrapper.hpp>
33 #include <opm/parser/eclipse/EclipseState/EclipseState.hpp>
34 
35 #include <array>
36 #include <cassert>
37 #include <utility>
38 #include <vector>
39 
45 struct UnstructuredGrid;
46 
47 namespace Opm
48 {
49 
65  template<class Grid>
66  void initStateEquil(const Grid& grid,
67  const BlackoilPropertiesInterface& props,
68  const Opm::DeckConstPtr deck,
69  const Opm::EclipseStateConstPtr eclipseState,
70  const double gravity,
71  BlackoilState& state);
72 
73 
81  namespace Equil {
82 
121  template <class Grid, class Region, class CellRange>
122  std::vector< std::vector<double> >
123  phasePressures(const Grid& G,
124  const Region& reg,
125  const CellRange& cells,
126  const double grav = unit::gravity);
127 
128 
129 
154  template <class Grid, class Region, class CellRange>
155  std::vector< std::vector<double> >
156  phaseSaturations(const Grid& grid,
157  const Region& reg,
158  const CellRange& cells,
159  BlackoilPropertiesInterface& props,
160  const std::vector<double> swat_init,
161  std::vector< std::vector<double> >& phase_pressures);
162 
163 
164 
183  template <class Grid, class CellRangeType>
184  std::vector<double> computeRs(const Grid& grid,
185  const CellRangeType& cells,
186  const std::vector<double> oil_pressure,
187  const std::vector<double>& temperature,
188  const Miscibility::RsFunction& rs_func,
189  const std::vector<double> gas_saturation);
190 
191  namespace DeckDependent {
192  inline
193  std::vector<EquilRecord>
194  getEquil(const Opm::DeckConstPtr deck)
195  {
196  if (deck->hasKeyword("EQUIL")) {
197 
198  Opm::EquilWrapper eql(deck->getKeyword("EQUIL"));
199 
200  const int nrec = eql.numRegions();
201 
202  std::vector<EquilRecord> ret;
203  ret.reserve(nrec);
204  for (int r = 0; r < nrec; ++r) {
205 
206  EquilRecord record =
207  {
208  { eql.datumDepth(r) ,
209  eql.datumDepthPressure(r) }
210  ,
211  { eql.waterOilContactDepth(r) ,
212  eql.waterOilContactCapillaryPressure(r) }
213  ,
214  { eql.gasOilContactDepth(r) ,
215  eql.gasOilContactCapillaryPressure(r) }
216  ,
217  eql.liveOilInitProceedure(r)
218  ,
219  eql.wetGasInitProceedure(r)
220  ,
221  eql.initializationTargetAccuracy(r)
222  };
223  if (record.N != 0) {
224  OPM_THROW(std::domain_error,
225  "kw EQUIL, item 9: Only N=0 supported.");
226  }
227  ret.push_back(record);
228  }
229 
230  return ret;
231  }
232  else {
233  OPM_THROW(std::domain_error,
234  "Deck does not provide equilibration data.");
235  }
236  }
237 
238  template<class Grid>
239  inline
240  std::vector<int>
241  equilnum(const Opm::DeckConstPtr deck,
242  const Opm::EclipseStateConstPtr eclipseState,
243  const Grid& G )
244  {
245  std::vector<int> eqlnum;
246  if (deck->hasKeyword("EQLNUM")) {
247  const int nc = UgGridHelpers::numCells(G);
248  eqlnum.resize(nc);
249  const std::vector<int>& e =
250  eclipseState->getIntGridProperty("EQLNUM")->getData();
251  const int* gc = UgGridHelpers::globalCell(G);
252  for (int cell = 0; cell < nc; ++cell) {
253  const int deck_pos = (gc == NULL) ? cell : gc[cell];
254  eqlnum[cell] = e[deck_pos] - 1;
255  }
256  }
257  else {
258  // No explicit equilibration region.
259  // All cells in region zero.
260  eqlnum.assign(UgGridHelpers::numCells(G), 0);
261  }
262 
263  return eqlnum;
264  }
265 
266 
267  class InitialStateComputer {
268  public:
269  template<class Grid>
270  InitialStateComputer(BlackoilPropertiesInterface& props,
271  const Opm::DeckConstPtr deck,
272  const Opm::EclipseStateConstPtr eclipseState,
273  const Grid& G ,
274  const double grav = unit::gravity)
275  : pp_(props.numPhases(),
276  std::vector<double>(UgGridHelpers::numCells(G))),
277  sat_(props.numPhases(),
278  std::vector<double>(UgGridHelpers::numCells(G))),
279  rs_(UgGridHelpers::numCells(G)),
280  rv_(UgGridHelpers::numCells(G))
281  {
282  // Get the equilibration records.
283  const std::vector<EquilRecord> rec = getEquil(deck);
284  std::shared_ptr<const TableManager> tables = eclipseState->getTableManager();
285  // Create (inverse) region mapping.
286  const RegionMapping<> eqlmap(equilnum(deck, eclipseState, G));
287 
288  // Create Rs functions.
289  rs_func_.reserve(rec.size());
290  if (deck->hasKeyword("DISGAS")) {
291  const TableContainer& rsvdTables = tables->getRsvdTables();
292  for (size_t i = 0; i < rec.size(); ++i) {
293  const int cell = *(eqlmap.cells(i).begin());
294  if (rec[i].live_oil_table_index > 0) {
295  if (rsvdTables.size() > 0 && size_t(rec[i].live_oil_table_index) <= rsvdTables.size()) {
296  const RsvdTable& rsvdTable = rsvdTables.getTable<RsvdTable>(i);
297  rs_func_.push_back(std::make_shared<Miscibility::RsVD>(props,
298  cell,
299  rsvdTable.getDepthColumn(),
300  rsvdTable.getRsColumn()));
301  } else {
302  OPM_THROW(std::runtime_error, "Cannot initialise: RSVD table " << (rec[i].live_oil_table_index) << " not available.");
303  }
304  } else {
305  if (rec[i].goc.depth != rec[i].main.depth) {
306  OPM_THROW(std::runtime_error,
307  "Cannot initialise: when no explicit RSVD table is given, \n"
308  "datum depth must be at the gas-oil-contact. "
309  "In EQUIL region " << (i + 1) << " (counting from 1), this does not hold.");
310  }
311  const double p_contact = rec[i].main.press;
312  const double T_contact = 273.15 + 20; // standard temperature for now
313  rs_func_.push_back(std::make_shared<Miscibility::RsSatAtContact>(props, cell, p_contact, T_contact));
314  }
315  }
316  } else {
317  for (size_t i = 0; i < rec.size(); ++i) {
318  rs_func_.push_back(std::make_shared<Miscibility::NoMixing>());
319  }
320  }
321 
322  rv_func_.reserve(rec.size());
323  if (deck->hasKeyword("VAPOIL")) {
324  const TableContainer& rvvdTables = tables->getRvvdTables();
325  for (size_t i = 0; i < rec.size(); ++i) {
326  const int cell = *(eqlmap.cells(i).begin());
327  if (rec[i].wet_gas_table_index > 0) {
328  if (rvvdTables.size() > 0 && size_t(rec[i].wet_gas_table_index) <= rvvdTables.size()) {
329  const RvvdTable& rvvdTable = rvvdTables.getTable<RvvdTable>(i);
330  rv_func_.push_back(std::make_shared<Miscibility::RvVD>(props,
331  cell,
332  rvvdTable.getDepthColumn(),
333  rvvdTable.getRvColumn()));
334  } else {
335  OPM_THROW(std::runtime_error, "Cannot initialise: RVVD table " << (rec[i].wet_gas_table_index) << " not available.");
336  }
337  } else {
338  if (rec[i].goc.depth != rec[i].main.depth) {
339  OPM_THROW(std::runtime_error,
340  "Cannot initialise: when no explicit RVVD table is given, \n"
341  "datum depth must be at the gas-oil-contact. "
342  "In EQUIL region " << (i + 1) << " (counting from 1), this does not hold.");
343  }
344  const double p_contact = rec[i].main.press + rec[i].goc.press;
345  const double T_contact = 273.15 + 20; // standard temperature for now
346  rv_func_.push_back(std::make_shared<Miscibility::RvSatAtContact>(props, cell, p_contact, T_contact));
347  }
348  }
349  } else {
350  for (size_t i = 0; i < rec.size(); ++i) {
351  rv_func_.push_back(std::make_shared<Miscibility::NoMixing>());
352  }
353  }
354 
355 
356  // Check for presence of kw SWATINIT
357  if (deck->hasKeyword("SWATINIT")) {
358  const std::vector<double>& swat_init = eclipseState->getDoubleGridProperty("SWATINIT")->getData();
359  const int nc = UgGridHelpers::numCells(G);
360  swat_init_.resize(nc);
361  const int* gc = UgGridHelpers::globalCell(G);
362  for (int c = 0; c < nc; ++c) {
363  const int deck_pos = (gc == NULL) ? c : gc[c];
364  swat_init_[c] = swat_init[deck_pos];
365  }
366  }
367 
368  // Compute pressures, saturations, rs and rv factors.
369  calcPressSatRsRv(eqlmap, rec, props, G, grav);
370 
371  // Modify oil pressure in no-oil regions so that the pressures of present phases can
372  // be recovered from the oil pressure and capillary relations.
373  }
374 
375  typedef std::vector<double> Vec;
376  typedef std::vector<Vec> PVec; // One per phase.
377 
378  const PVec& press() const { return pp_; }
379  const PVec& saturation() const { return sat_; }
380  const Vec& rs() const { return rs_; }
381  const Vec& rv() const { return rv_; }
382 
383  private:
384  typedef DensityCalculator<BlackoilPropertiesInterface> RhoCalc;
385  typedef EquilReg<RhoCalc> EqReg;
386 
387  std::vector< std::shared_ptr<Miscibility::RsFunction> > rs_func_;
388  std::vector< std::shared_ptr<Miscibility::RsFunction> > rv_func_;
389 
390  PVec pp_;
391  PVec sat_;
392  Vec rs_;
393  Vec rv_;
394  Vec swat_init_;
395 
396  template <class RMap, class Grid>
397  void
398  calcPressSatRsRv(const RMap& reg ,
399  const std::vector< EquilRecord >& rec ,
401  const Grid& G ,
402  const double grav)
403  {
404  for (const auto& r : reg.activeRegions()) {
405  const auto& cells = reg.cells(r);
406  const int repcell = *cells.begin();
407 
408  const RhoCalc calc(props, repcell);
409  const EqReg eqreg(rec[r], calc,
410  rs_func_[r], rv_func_[r],
411  props.phaseUsage());
412 
413  PVec pressures = phasePressures(G, eqreg, cells, grav);
414  const std::vector<double>& temp = temperature(G, eqreg, cells);
415 
416  const PVec sat = phaseSaturations(G, eqreg, cells, props, swat_init_, pressures);
417 
418  const int np = props.numPhases();
419  for (int p = 0; p < np; ++p) {
420  copyFromRegion(pressures[p], cells, pp_[p]);
421  copyFromRegion(sat[p], cells, sat_[p]);
422  }
425  const int oilpos = props.phaseUsage().phase_pos[BlackoilPhases::Liquid];
426  const int gaspos = props.phaseUsage().phase_pos[BlackoilPhases::Vapour];
427  const Vec rs_vals = computeRs(G, cells, pressures[oilpos], temp, *(rs_func_[r]), sat[gaspos]);
428  const Vec rv_vals = computeRs(G, cells, pressures[gaspos], temp, *(rv_func_[r]), sat[oilpos]);
429  copyFromRegion(rs_vals, cells, rs_);
430  copyFromRegion(rv_vals, cells, rv_);
431  }
432  }
433  }
434 
435  template <class CellRangeType>
436  void copyFromRegion(const Vec& source,
437  const CellRangeType& cells,
438  Vec& destination)
439  {
440  auto s = source.begin();
441  auto c = cells.begin();
442  const auto e = cells.end();
443  for (; c != e; ++c, ++s) {
444  destination[*c] = *s;
445  }
446  }
447 
448  };
449  } // namespace DeckDependent
450  } // namespace Equil
451 } // namespace Opm
452 
454 
455 #endif // OPM_INITSTATEEQUIL_HEADER_INCLUDED
Definition: grid.h:98
void initStateEquil(const Grid &grid, const BlackoilPropertiesInterface &props, const Opm::DeckConstPtr deck, const Opm::EclipseStateConstPtr eclipseState, const double gravity, BlackoilState &state)
Definition: AnisotropicEikonal.hpp:43
Definition: BlackoilPhases.hpp:32
int phase_used[MaxNumPhases]
Definition: BlackoilPhases.hpp:39
virtual int numPhases() const =0
int numCells(const UnstructuredGrid &grid)
Get the number of cells of a grid.
Definition: BlackoilPropertiesInterface.hpp:37
int phase_pos[MaxNumPhases]
Definition: BlackoilPhases.hpp:40
STL namespace.
const int * globalCell(const UnstructuredGrid &grid)
Get the local to global index mapping.
std::vector< std::vector< double > > phasePressures(const Grid &G, const Region &reg, const CellRange &cells, const double grav)
Definition: initStateEquil_impl.hpp:594
virtual PhaseUsage phaseUsage() const =0
const double gravity
Definition: Units.hpp:120
std::vector< double > computeRs(const Grid &grid, const CellRangeType &cells, const std::vector< double > oil_pressure, const std::vector< double > &temperature, const Miscibility::RsFunction &rs_func, const std::vector< double > gas_saturation)
Definition: initStateEquil_impl.hpp:821
Definition: BlackoilPhases.hpp:32
const UnstructuredGrid & grid
Definition: ColumnExtract.hpp:31
std::vector< double > temperature(const Grid &, const Region &, const CellRange &cells)
Definition: initStateEquil_impl.hpp:667
std::vector< std::vector< double > > phaseSaturations(const Grid &G, const Region &reg, const CellRange &cells, BlackoilPropertiesInterface &props, const std::vector< double > swat_init, std::vector< std::vector< double > > &phase_pressures)
Definition: initStateEquil_impl.hpp:677