InitStateEquil.hpp
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 3 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
29#ifndef OPM_INIT_STATE_EQUIL_HPP
30#define OPM_INIT_STATE_EQUIL_HPP
31
32#include <opm/models/utils/propertysystem.hh>
33
34#include <opm/material/common/Tabulated1DFunction.hpp>
35#include <opm/material/fluidstates/SimpleModularFluidState.hpp>
36
38
39#include <array>
40#include <cstddef>
41#include <memory>
42#include <utility>
43#include <vector>
44#include <string>
45
46namespace Opm {
47
48class EclipseState;
49class EquilRecord;
50class NumericalAquifers;
51
59namespace EQUIL {
60
61class EquilReg;
62namespace Miscibility { class RsFunction; }
63
64namespace Details {
65template <class RHS>
66class RK4IVP {
67public:
68 RK4IVP(const RHS& f,
69 const std::array<double,2>& span,
70 const double y0,
71 const int N);
72
73 double
74 operator()(const double x) const;
75
76private:
77 int N_;
78 std::array<double,2> span_;
79 std::vector<double> y_;
80 std::vector<double> f_;
81
82 double stepsize() const;
83};
84
85namespace PhasePressODE {
86template <class FluidSystem>
87class Water
88{
89using TabulatedFunction = Tabulated1DFunction<double>;
90public:
91 Water(const TabulatedFunction& tempVdTable,
92 const TabulatedFunction& saltVdTable,
93 const int pvtRegionIdx,
94 const double normGrav);
95
96 double operator()(const double depth,
97 const double press) const;
98
99private:
100 const TabulatedFunction& tempVdTable_;
101 const TabulatedFunction& saltVdTable_;
102 const int pvtRegionIdx_;
103 const double g_;
104
105 double density(const double depth,
106 const double press) const;
107};
108
109template <class FluidSystem, class RS>
110class Oil
111{
112using TabulatedFunction = Tabulated1DFunction<double>;
113public:
114 Oil(const TabulatedFunction& tempVdTable,
115 const RS& rs,
116 const int pvtRegionIdx,
117 const double normGrav);
118
119 double operator()(const double depth,
120 const double press) const;
121
122private:
123 const TabulatedFunction& tempVdTable_;
124 const RS& rs_;
125 const int pvtRegionIdx_;
126 const double g_;
127
128 double density(const double depth,
129 const double press) const;
130};
131
132template <class FluidSystem, class RV, class RVW>
133class Gas
134{
135using TabulatedFunction = Tabulated1DFunction<double>;
136public:
137 Gas(const TabulatedFunction& tempVdTable,
138 const RV& rv,
139 const RVW& rvw,
140 const int pvtRegionIdx,
141 const double normGrav);
142
143 double operator()(const double depth,
144 const double press) const;
145
146private:
147 const TabulatedFunction& tempVdTable_;
148 const RV& rv_;
149 const RVW& rvw_;
150 const int pvtRegionIdx_;
151 const double g_;
152
153 double density(const double depth,
154 const double press) const;
155};
156
157} // namespace PhasePressODE
158
159template <class FluidSystem, class Region>
161{
162public:
163 using VSpan = std::array<double, 2>;
164
173 explicit PressureTable(const double gravity,
174 const int samplePoints = 2000);
175
179 PressureTable(const PressureTable& rhs);
180
187
194
203
204 void equilibrate(const Region& reg,
205 const VSpan& span);
206
208 bool oilActive() const;
209
211 bool gasActive() const;
212
214 bool waterActive() const;
215
223 double oil(const double depth) const;
224
232 double gas(const double depth) const;
233
241 double water(const double depth) const;
242
243private:
244 template <class ODE>
245 class PressureFunction
246 {
247 public:
248 struct InitCond {
249 double depth;
250 double pressure;
251 };
252
253 explicit PressureFunction(const ODE& ode,
254 const InitCond& ic,
255 const int nsample,
256 const VSpan& span);
257
258 PressureFunction(const PressureFunction& rhs);
259
260 PressureFunction(PressureFunction&& rhs) = default;
261
262 PressureFunction& operator=(const PressureFunction& rhs);
263
264 PressureFunction& operator=(PressureFunction&& rhs);
265
266 double value(const double depth) const;
267
268 private:
269 enum Direction : std::size_t { Up, Down, NumDir };
270
271 using Distribution = Details::RK4IVP<ODE>;
272 using DistrPtr = std::unique_ptr<Distribution>;
273
274 InitCond initial_;
275 std::array<DistrPtr, Direction::NumDir> value_;
276 };
277
278 using OilPressODE = PhasePressODE::Oil<
279 FluidSystem, typename Region::CalcDissolution
280 >;
281
282 using GasPressODE = PhasePressODE::Gas<
283 FluidSystem, typename Region::CalcEvaporation, typename Region::CalcWaterEvaporation
284 >;
285
286 using WatPressODE = PhasePressODE::Water<FluidSystem>;
287
288 using OPress = PressureFunction<OilPressODE>;
289 using GPress = PressureFunction<GasPressODE>;
290 using WPress = PressureFunction<WatPressODE>;
291
292 using Strategy = void (PressureTable::*)
293 (const Region&, const VSpan&);
294
295 double gravity_;
296 int nsample_;
297
298 std::unique_ptr<OPress> oil_{};
299 std::unique_ptr<GPress> gas_{};
300 std::unique_ptr<WPress> wat_{};
301
302 template <typename PressFunc>
303 void checkPtr(const PressFunc* phasePress,
304 const std::string& phaseName) const;
305
306 Strategy selectEquilibrationStrategy(const Region& reg) const;
307
308 void copyInPointers(const PressureTable& rhs);
309
310 void equil_WOG(const Region& reg, const VSpan& span);
311 void equil_GOW(const Region& reg, const VSpan& span);
312 void equil_OWG(const Region& reg, const VSpan& span);
313
314 void makeOilPressure(const typename OPress::InitCond& ic,
315 const Region& reg,
316 const VSpan& span);
317
318 void makeGasPressure(const typename GPress::InitCond& ic,
319 const Region& reg,
320 const VSpan& span);
321
322 void makeWatPressure(const typename WPress::InitCond& ic,
323 const Region& reg,
324 const VSpan& span);
325};
326
327// ===========================================================================
328
331 double oil{0.0};
332 double gas{0.0};
333 double water{0.0};
334
335 PhaseQuantityValue& axpy(const PhaseQuantityValue& rhs, const double a)
336 {
337 this->oil += a * rhs.oil;
338 this->gas += a * rhs.gas;
339 this->water += a * rhs.water;
340
341 return *this;
342 }
343
345 {
346 this->oil /= x;
347 this->gas /= x;
348 this->water /= x;
349
350 return *this;
351 }
352
353 void reset()
354 {
355 this->oil = this->gas = this->water = 0.0;
356 }
357};
358
376template <class MaterialLawManager, class FluidSystem, class Region, typename CellID>
378{
379public:
383 struct Position {
384 CellID cell;
385 double depth;
386 };
387
390
398 explicit PhaseSaturations(MaterialLawManager& matLawMgr,
399 const std::vector<double>& swatInit);
400
405
408
411
424 const PhaseQuantityValue&
426 const Region& reg,
427 const PTable& ptable);
428
434 {
435 return this->press_;
436 }
437
438private:
441 struct EvaluationPoint {
442 const Position* position{nullptr};
443 const Region* region {nullptr};
444 const PTable* ptable {nullptr};
445 };
446
450 using FluidState = ::Opm::
451 SimpleModularFluidState<double, /*numPhases=*/3, /*numComponents=*/3,
452 FluidSystem,
453 /*storePressure=*/false,
454 /*storeTemperature=*/false,
455 /*storeComposition=*/false,
456 /*storeFugacity=*/false,
457 /*storeSaturation=*/true,
458 /*storeDensity=*/false,
459 /*storeViscosity=*/false,
460 /*storeEnthalpy=*/false>;
461
463 using MaterialLaw = typename MaterialLawManager::MaterialLaw;
464
466 using PhaseIdx = std::remove_cv_t<
467 std::remove_reference_t<decltype(FluidSystem::oilPhaseIdx)>
468 >;
469
471 MaterialLawManager& matLawMgr_;
472
474 const std::vector<double>& swatInit_;
475
477 PhaseQuantityValue sat_;
478
480 PhaseQuantityValue press_;
481
483 EvaluationPoint evalPt_;
484
486 FluidState fluidState_;
487
489 std::array<double, FluidSystem::numPhases> matLawCapPress_;
490
500 void setEvaluationPoint(const Position& x,
501 const Region& reg,
502 const PTable& ptable);
503
507 void initializePhaseQuantities();
508
512 void deriveOilSat();
513
518 void deriveGasSat();
519
526 void deriveWaterSat();
527
530 void fixUnphysicalTransition();
531
534 void accountForScaledSaturations();
535
536 // --------------------------------------------------------------------
537 // Note: Function 'applySwatInit' is non-const because the overload set
538 // needs to mutate the 'matLawMgr_'.
539 // --------------------------------------------------------------------
540
550 std::pair<double, bool> applySwatInit(const double pcow);
551
564 std::pair<double, bool> applySwatInit(const double pc, const double sw);
565
568 void computeMaterialLawCapPress();
569
572 double materialLawCapPressGasOil() const;
573
576 double materialLawCapPressOilWater() const;
577
580 double materialLawCapPressGasWater() const;
581
589 bool isConstCapPress(const PhaseIdx phaseIdx) const;
590
596 bool isOverlappingTransition() const;
597
617 double fromDepthTable(const double contactdepth,
618 const PhaseIdx phasePos,
619 const bool isincr) const;
620
638 double invertCapPress(const double pc,
639 const PhaseIdx phasePos,
640 const bool isincr) const;
641
643 PhaseIdx oilPos() const
644 {
645 return FluidSystem::oilPhaseIdx;
646 }
647
649 PhaseIdx gasPos() const
650 {
651 return FluidSystem::gasPhaseIdx;
652 }
653
655 PhaseIdx waterPos() const
656 {
657 return FluidSystem::waterPhaseIdx;
658 }
659};
660
661// ===========================================================================
662
663template <typename CellRange>
664void verticalExtent(const CellRange& cells,
665 const std::vector<std::pair<double, double>>& cellZMinMax,
666 const Parallel::Communication& comm,
667 std::array<double,2>& span);
668
669template <class Element>
670std::pair<double,double> cellZMinMax(const Element& element);
671
672} // namespace Details
673
674namespace DeckDependent {
675
676template<class FluidSystem,
677 class Grid,
678 class GridView,
679 class ElementMapper,
680 class CartesianIndexMapper>
682{
683 using Element = typename GridView::template Codim<0>::Entity;
684public:
685 template<class MaterialLawManager>
686 InitialStateComputer(MaterialLawManager& materialLawManager,
687 const EclipseState& eclipseState,
688 const Grid& grid,
689 const GridView& gridView,
690 const CartesianIndexMapper& cartMapper,
691 const double grav,
692 const int num_pressure_points = 2000,
693 const bool applySwatInit = true);
694
695 using Vec = std::vector<double>;
696 using PVec = std::vector<Vec>; // One per phase.
697
698 const Vec& temperature() const { return temperature_; }
699 const Vec& saltConcentration() const { return saltConcentration_; }
700 const Vec& saltSaturation() const { return saltSaturation_; }
701 const PVec& press() const { return pp_; }
702 const PVec& saturation() const { return sat_; }
703 const Vec& rs() const { return rs_; }
704 const Vec& rv() const { return rv_; }
705 const Vec& rvw() const { return rvw_; }
706
707private:
708 template <class RMap>
709 void updateInitialTemperature_(const EclipseState& eclState, const RMap& reg);
710
711 template <class RMap>
712 void updateInitialSaltConcentration_(const EclipseState& eclState, const RMap& reg);
713
714 template <class RMap>
715 void updateInitialSaltSaturation_(const EclipseState& eclState, const RMap& reg);
716
717 void updateCellProps_(const GridView& gridView,
718 const NumericalAquifers& aquifer);
719
720 void applyNumericalAquifers_(const GridView& gridView,
721 const NumericalAquifers& aquifer,
722 const bool co2store_or_h2store);
723
724 template<class RMap>
725 void setRegionPvtIdx(const EclipseState& eclState, const RMap& reg);
726
727 template <class RMap, class MaterialLawManager, class Comm>
728 void calcPressSatRsRv(const RMap& reg,
729 const std::vector<EquilRecord>& rec,
730 MaterialLawManager& materialLawManager,
731 const Comm& comm,
732 const double grav);
733
734 template <class CellRange, class EquilibrationMethod>
735 void cellLoop(const CellRange& cells,
736 EquilibrationMethod&& eqmethod);
737
738 template <class CellRange, class PressTable, class PhaseSat>
739 void equilibrateCellCentres(const CellRange& cells,
740 const EquilReg& eqreg,
741 const PressTable& ptable,
742 PhaseSat& psat);
743
744 template <class CellRange, class PressTable, class PhaseSat>
745 void equilibrateHorizontal(const CellRange& cells,
746 const EquilReg& eqreg,
747 const int acc,
748 const PressTable& ptable,
749 PhaseSat& psat);
750
751 std::vector< std::shared_ptr<Miscibility::RsFunction> > rsFunc_;
752 std::vector< std::shared_ptr<Miscibility::RsFunction> > rvFunc_;
753 std::vector< std::shared_ptr<Miscibility::RsFunction> > rvwFunc_;
754 using TabulatedFunction = Tabulated1DFunction<double>;
755 std::vector<TabulatedFunction> tempVdTable_;
756 std::vector<TabulatedFunction> saltVdTable_;
757 std::vector<TabulatedFunction> saltpVdTable_;
758 std::vector<int> regionPvtIdx_;
759 Vec temperature_;
760 Vec saltConcentration_;
761 Vec saltSaturation_;
762 PVec pp_;
763 PVec sat_;
764 Vec rs_;
765 Vec rv_;
766 Vec rvw_;
767 const CartesianIndexMapper& cartesianIndexMapper_;
768 Vec swatInit_;
769 Vec cellCenterDepth_;
770 std::vector<std::pair<double,double>> cellZSpan_;
771 std::vector<std::pair<double,double>> cellZMinMax_;
772 int num_pressure_points_;
773};
774
775} // namespace DeckDependent
776} // namespace EQUIL
777} // namespace Opm
778
779#endif // OPM_INIT_STATE_EQUIL_HPP
Dune::OwnerOverlapCopyCommunication< int, int > Comm
Definition: FlexibleSolver_impl.hpp:274
Definition: InitStateEquil.hpp:682
std::vector< Vec > PVec
Definition: InitStateEquil.hpp:696
const PVec & press() const
Definition: InitStateEquil.hpp:701
const Vec & rvw() const
Definition: InitStateEquil.hpp:705
const Vec & rv() const
Definition: InitStateEquil.hpp:704
const Vec & saltSaturation() const
Definition: InitStateEquil.hpp:700
const Vec & saltConcentration() const
Definition: InitStateEquil.hpp:699
const Vec & temperature() const
Definition: InitStateEquil.hpp:698
const Vec & rs() const
Definition: InitStateEquil.hpp:703
std::vector< double > Vec
Definition: InitStateEquil.hpp:695
const PVec & saturation() const
Definition: InitStateEquil.hpp:702
InitialStateComputer(MaterialLawManager &materialLawManager, const EclipseState &eclipseState, const Grid &grid, const GridView &gridView, const CartesianIndexMapper &cartMapper, const double grav, const int num_pressure_points=2000, const bool applySwatInit=true)
Definition: InitStateEquil_impl.hpp:1308
Definition: InitStateEquil.hpp:134
Gas(const TabulatedFunction &tempVdTable, const RV &rv, const RVW &rvw, const int pvtRegionIdx, const double normGrav)
Definition: InitStateEquil_impl.hpp:329
double operator()(const double depth, const double press) const
Definition: InitStateEquil_impl.hpp:344
Definition: InitStateEquil.hpp:111
double operator()(const double depth, const double press) const
Definition: InitStateEquil_impl.hpp:296
Oil(const TabulatedFunction &tempVdTable, const RS &rs, const int pvtRegionIdx, const double normGrav)
Definition: InitStateEquil_impl.hpp:283
Definition: InitStateEquil.hpp:88
Water(const TabulatedFunction &tempVdTable, const TabulatedFunction &saltVdTable, const int pvtRegionIdx, const double normGrav)
Definition: InitStateEquil_impl.hpp:249
double operator()(const double depth, const double press) const
Definition: InitStateEquil_impl.hpp:262
Definition: InitStateEquil.hpp:378
const PhaseQuantityValue & correctedPhasePressures() const
Definition: InitStateEquil.hpp:433
PressureTable< FluidSystem, Region > PTable
Convenience type alias.
Definition: InitStateEquil.hpp:389
PhaseSaturations & operator=(const PhaseSaturations &)=delete
Disabled assignment operator.
const PhaseQuantityValue & deriveSaturations(const Position &x, const Region &reg, const PTable &ptable)
Definition: InitStateEquil_impl.hpp:574
PhaseSaturations & operator=(PhaseSaturations &&)=delete
Disabled move-assignment operator.
PhaseSaturations(MaterialLawManager &matLawMgr, const std::vector< double > &swatInit)
Definition: InitStateEquil_impl.hpp:550
Definition: InitStateEquil.hpp:161
PressureTable & operator=(const PressureTable &rhs)
Definition: InitStateEquil_impl.hpp:990
double gas(const double depth) const
Definition: InitStateEquil_impl.hpp:1057
bool waterActive() const
Predicate for whether or not water is an active phase.
Definition: InitStateEquil_impl.hpp:1041
bool gasActive() const
Predicate for whether or not gas is an active phase.
Definition: InitStateEquil_impl.hpp:1034
bool oilActive() const
Predicate for whether or not oil is an active phase.
Definition: InitStateEquil_impl.hpp:1027
PressureTable(const double gravity, const int samplePoints=2000)
Definition: InitStateEquil_impl.hpp:960
void equilibrate(const Region &reg, const VSpan &span)
Definition: InitStateEquil_impl.hpp:1016
double water(const double depth) const
Definition: InitStateEquil_impl.hpp:1067
std::array< double, 2 > VSpan
Definition: InitStateEquil.hpp:163
double oil(const double depth) const
Definition: InitStateEquil_impl.hpp:1048
Definition: InitStateEquil.hpp:66
double operator()(const double x) const
Definition: InitStateEquil_impl.hpp:215
RK4IVP(const RHS &f, const std::array< double, 2 > &span, const double y0, const int N)
Definition: InitStateEquil_impl.hpp:180
Definition: EquilibrationHelpers.hpp:600
std::pair< double, double > cellZMinMax(const Element &element)
Definition: InitStateEquil_impl.hpp:161
void verticalExtent(const CellRange &cells, const std::vector< std::pair< double, double > > &cellZMinMax, const Parallel::Communication &comm, std::array< double, 2 > &span)
Definition: InitStateEquil_impl.hpp:64
Dune::Communication< MPIComm > Communication
Definition: ParallelCommunication.hpp:30
Definition: BlackoilPhases.hpp:27
Simple set of per-phase (named by primary component) quantities.
Definition: InitStateEquil.hpp:330
double water
Definition: InitStateEquil.hpp:333
PhaseQuantityValue & operator/=(const double x)
Definition: InitStateEquil.hpp:344
PhaseQuantityValue & axpy(const PhaseQuantityValue &rhs, const double a)
Definition: InitStateEquil.hpp:335
double gas
Definition: InitStateEquil.hpp:332
double oil
Definition: InitStateEquil.hpp:331
void reset()
Definition: InitStateEquil.hpp:353
Definition: InitStateEquil.hpp:383
CellID cell
Definition: InitStateEquil.hpp:384
double depth
Definition: InitStateEquil.hpp:385
double pressure
Definition: InitStateEquil.hpp:250
double depth
Definition: InitStateEquil.hpp:249