BlackoilPropsAdFromDeck.hpp
Go to the documentation of this file.
1 /*
2  Copyright 2013 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_BLACKOILPROPSADFROMDECK_HEADER_INCLUDED
23 #define OPM_BLACKOILPROPSADFROMDECK_HEADER_INCLUDED
24 
27 
28 #include <opm/core/props/BlackoilPhases.hpp>
29 #include <opm/core/props/satfunc/SaturationPropsFromDeck.hpp>
30 #include <opm/core/props/rock/RockFromDeck.hpp>
31 
32 #include <opm/parser/eclipse/Deck/Deck.hpp>
33 #include <opm/parser/eclipse/EclipseState/EclipseState.hpp>
34 
35 #include <memory>
36 #include <array>
37 #include <vector>
38 
39 #ifdef HAVE_DUNE_CORNERPOINT
40 #include <opm/common/utility/platform_dependent/disable_warnings.h>
41 #include <dune/grid/CpGrid.hpp>
42 #include <opm/common/utility/platform_dependent/reenable_warnings.h>
43 #endif
44 
45 namespace Opm
46 {
47  class PvtInterface;
48 
59  {
61  public:
62  typedef typename SaturationPropsFromDeck::MaterialLawManager MaterialLawManager;
63 
78  BlackoilPropsAdFromDeck(Opm::DeckConstPtr deck,
79  Opm::EclipseStateConstPtr eclState,
80  std::shared_ptr<MaterialLawManager> materialLawManager,
81  const UnstructuredGrid& grid,
82  const bool init_rock = true );
83 
84 #ifdef HAVE_DUNE_CORNERPOINT
85  BlackoilPropsAdFromDeck(Opm::DeckConstPtr deck,
100  Opm::EclipseStateConstPtr eclState,
101  std::shared_ptr<MaterialLawManager> materialLawManager,
102  const Dune::CpGrid& grid,
103  const bool init_rock = true );
104 #endif
105 
113  BlackoilPropsAdFromDeck(Opm::DeckConstPtr deck,
114  Opm::EclipseStateConstPtr eclState,
115  const UnstructuredGrid& grid,
116  const bool init_rock = true );
117 
118 #ifdef HAVE_DUNE_CORNERPOINT
119  BlackoilPropsAdFromDeck(Opm::DeckConstPtr deck,
127  Opm::EclipseStateConstPtr eclState,
128  const Dune::CpGrid& grid,
129  const bool init_rock = true );
130 #endif
131 
143  const int number_of_cells);
144 
145 
147  // Rock interface //
149 
151  int numDimensions() const;
152 
154  int numCells() const;
155 
158  virtual const int* cellPvtRegionIndex() const
159  { return &cellPvtRegionIdx_[0]; }
160 
162  const double* porosity() const;
163 
167  const double* permeability() const;
168 
169 
171  // Fluid interface //
173 
175  typedef ADB::V V;
176  typedef std::vector<int> Cells;
177 
179  int numPhases() const;
180 
182  PhaseUsage phaseUsage() const;
183 
184  // ------ Density ------
185 
188  const double* surfaceDensity(const int cellIdx = 0) const;
189 
190 
191  // ------ Viscosity ------
192 
198  ADB muWat(const ADB& pw,
199  const ADB& T,
200  const Cells& cells) const;
201 
209  ADB muOil(const ADB& po,
210  const ADB& T,
211  const ADB& rs,
212  const std::vector<PhasePresence>& cond,
213  const Cells& cells) const;
214 
222  ADB muGas(const ADB& pg,
223  const ADB& T,
224  const ADB& rv,
225  const std::vector<PhasePresence>& cond,
226  const Cells& cells) const;
227 
228  // ------ Formation volume factor (b) ------
229 
235  ADB bWat(const ADB& pw,
236  const ADB& T,
237  const Cells& cells) const;
238 
246  ADB bOil(const ADB& po,
247  const ADB& T,
248  const ADB& rs,
249  const std::vector<PhasePresence>& cond,
250  const Cells& cells) const;
251 
259  ADB bGas(const ADB& pg,
260  const ADB& T,
261  const ADB& rv,
262  const std::vector<PhasePresence>& cond,
263  const Cells& cells) const;
264 
265  // ------ Rs bubble point curve ------
266 
271  V rsSat(const V& po,
272  const Cells& cells) const;
273 
279  V rsSat(const V& po,
280  const V& so,
281  const Cells& cells) const;
282 
287  ADB rsSat(const ADB& po,
288  const Cells& cells) const;
289 
295  ADB rsSat(const ADB& po,
296  const ADB& so,
297  const Cells& cells) const;
298 
299  // ------ Rv condensation curve ------
300 
305  ADB rvSat(const ADB& po,
306  const Cells& cells) const;
307 
313  ADB rvSat(const ADB& po,
314  const ADB& so,
315  const Cells& cells) const;
316 
317  // ------ Relative permeability ------
318 
326  std::vector<ADB> relperm(const ADB& sw,
327  const ADB& so,
328  const ADB& sg,
329  const Cells& cells) const;
330 
339  std::vector<ADB> capPress(const ADB& sw,
340  const ADB& so,
341  const ADB& sg,
342  const Cells& cells) const;
343 
346  void updateSatHyst(const std::vector<double>& saturation,
347  const std::vector<int>& cells);
348 
350  void updateSatOilMax(const std::vector<double>& saturation);
351 
355  void setSwatInitScaling(const std::vector<double>& saturation,
356  const std::vector<double>& pc);
357 
358 
359  private:
361  void init(Opm::DeckConstPtr deck,
362  Opm::EclipseStateConstPtr eclState,
363  std::shared_ptr<MaterialLawManager> materialLawManager,
364  int number_of_cells,
365  const int* global_cell,
366  const int* cart_dims,
367  const bool init_rock);
368 
370  void applyVap(V& r,
371  const V& so,
372  const std::vector<int>& cells,
373  const double vap) const;
374 
375  void applyVap(ADB& r,
376  const ADB& so,
377  const std::vector<int>& cells,
378  const double vap) const;
379 
380  // Fills pvt_region_ with cellPvtRegionIdx_[cells].
381  void mapPvtRegions(const std::vector<int>& cells) const;
382 
383  RockFromDeck rock_;
384 
385  // This has to be a shared pointer as we must
386  // be able to make a copy of *this in the parallel case.
387  std::shared_ptr<MaterialLawManager> materialLawManager_;
388  std::shared_ptr<SaturationPropsFromDeck> satprops_;
389 
390  PhaseUsage phase_usage_;
391  // bool has_vapoil_;
392  // bool has_disgas_;
393 
394  // The PVT region which is to be used for each cell
395  std::vector<int> cellPvtRegionIdx_;
396 
397  // Used for storing the region-per-cell array computed in calls
398  // to pvt functions.
399  mutable std::vector<int> pvt_region_;
400 
401  // The PVT properties. One object per active fluid phase.
402  std::vector<std::shared_ptr<Opm::PvtInterface> > props_;
403 
404  // Densities, one std::array per PVT region.
405  std::vector<std::array<double, BlackoilPhases::MaxNumPhases> > densities_;
406 
407  // VAPPARS
408  double vap1_;
409  double vap2_;
410  std::vector<double> satOilMax_;
411  double vap_satmax_guard_; //Threshold value to promote stability
412 
413  };
414 } // namespace Opm
415 
416 #endif // OPM_BLACKOILPROPSADFROMDECK_HEADER_INCLUDED
std::vector< ADB > relperm(const ADB &sw, const ADB &so, const ADB &sg, const Cells &cells) const
ADB muGas(const ADB &pg, const ADB &T, const ADB &rv, const std::vector< PhasePresence > &cond, const Cells &cells) const
ADB::V V
Definition: BlackoilPropsAdFromDeck.hpp:175
Definition: BlackoilPropsAdInterface.hpp:38
Eigen::Array< double, Eigen::Dynamic, 1 > V
Underlying type for values.
Definition: AutoDiffBlock.hpp:98
const double * surfaceDensity(const int cellIdx=0) const
Definition: AdditionalObjectDeleter.hpp:22
void setSwatInitScaling(const std::vector< double > &saturation, const std::vector< double > &pc)
ADB rvSat(const ADB &po, const Cells &cells) const
ADB bOil(const ADB &po, const ADB &T, const ADB &rs, const std::vector< PhasePresence > &cond, const Cells &cells) const
Definition: BlackoilPropsAdFromDeck.hpp:58
ADB bWat(const ADB &pw, const ADB &T, const Cells &cells) const
AutoDiffBlock< double > ADB
Definition: BlackoilPropsAdFromDeck.hpp:174
ADB muOil(const ADB &po, const ADB &T, const ADB &rs, const std::vector< PhasePresence > &cond, const Cells &cells) const
V rsSat(const V &po, const Cells &cells) const
friend class BlackoilPropsDataHandle
Definition: BlackoilPropsAdFromDeck.hpp:60
void updateSatOilMax(const std::vector< double > &saturation)
Update for max oil saturation.
ADB bGas(const ADB &pg, const ADB &T, const ADB &rv, const std::vector< PhasePresence > &cond, const Cells &cells) const
ADB muWat(const ADB &pw, const ADB &T, const Cells &cells) const
BlackoilPropsAdFromDeck(Opm::DeckConstPtr deck, Opm::EclipseStateConstPtr eclState, std::shared_ptr< MaterialLawManager > materialLawManager, const UnstructuredGrid &grid, const bool init_rock=true)
PhaseUsage phaseUsage() const
virtual const int * cellPvtRegionIndex() const
Definition: BlackoilPropsAdFromDeck.hpp:158
void updateSatHyst(const std::vector< double > &saturation, const std::vector< int > &cells)
std::vector< int > Cells
Definition: BlackoilPropsAdFromDeck.hpp:176
const double * porosity() const
const double * permeability() const
std::vector< ADB > capPress(const ADB &sw, const ADB &so, const ADB &sg, const Cells &cells) const
SaturationPropsFromDeck::MaterialLawManager MaterialLawManager
Definition: BlackoilPropsAdFromDeck.hpp:62