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
37#include <array>
38#include <cstddef>
39#include <memory>
40#include <utility>
41#include <vector>
42#include <string>
43
44namespace Opm {
45
46class EclipseState;
47class EquilRecord;
48class NumericalAquifers;
49
57namespace EQUIL {
58
59class EquilReg;
60namespace Miscibility { class RsFunction; }
61
62namespace Details {
63template <class RHS>
64class RK4IVP {
65public:
66 RK4IVP(const RHS& f,
67 const std::array<double,2>& span,
68 const double y0,
69 const int N);
70
71 double
72 operator()(const double x) const;
73
74private:
75 int N_;
76 std::array<double,2> span_;
77 std::vector<double> y_;
78 std::vector<double> f_;
79
80 double stepsize() const;
81};
82
83namespace PhasePressODE {
84template <class FluidSystem>
85class Water
86{
87using TabulatedFunction = Tabulated1DFunction<double>;
88public:
89 Water(const TabulatedFunction& tempVdTable,
90 const TabulatedFunction& saltVdTable,
91 const int pvtRegionIdx,
92 const double normGrav);
93
94 double operator()(const double depth,
95 const double press) const;
96
97private:
98 const TabulatedFunction& tempVdTable_;
99 const TabulatedFunction& saltVdTable_;
100 const int pvtRegionIdx_;
101 const double g_;
102
103 double density(const double depth,
104 const double press) const;
105};
106
107template <class FluidSystem, class RS>
108class Oil
109{
110using TabulatedFunction = Tabulated1DFunction<double>;
111public:
112 Oil(const TabulatedFunction& tempVdTable,
113 const RS& rs,
114 const int pvtRegionIdx,
115 const double normGrav);
116
117 double operator()(const double depth,
118 const double press) const;
119
120private:
121 const TabulatedFunction& tempVdTable_;
122 const RS& rs_;
123 const int pvtRegionIdx_;
124 const double g_;
125
126 double density(const double depth,
127 const double press) const;
128};
129
130template <class FluidSystem, class RV, class RVW>
131class Gas
132{
133using TabulatedFunction = Tabulated1DFunction<double>;
134public:
135 Gas(const TabulatedFunction& tempVdTable,
136 const RV& rv,
137 const RVW& rvw,
138 const int pvtRegionIdx,
139 const double normGrav);
140
141 double operator()(const double depth,
142 const double press) const;
143
144private:
145 const TabulatedFunction& tempVdTable_;
146 const RV& rv_;
147 const RVW& rvw_;
148 const int pvtRegionIdx_;
149 const double g_;
150
151 double density(const double depth,
152 const double press) const;
153};
154
155} // namespace PhasePressODE
156
157template <class FluidSystem, class Region>
159{
160public:
161 using VSpan = std::array<double, 2>;
162
171 explicit PressureTable(const double gravity,
172 const int samplePoints = 2000);
173
177 PressureTable(const PressureTable& rhs);
178
185
192
201
202 void equilibrate(const Region& reg,
203 const VSpan& span);
204
206 bool oilActive() const;
207
209 bool gasActive() const;
210
212 bool waterActive() const;
213
221 double oil(const double depth) const;
222
230 double gas(const double depth) const;
231
239 double water(const double depth) const;
240
241private:
242 template <class ODE>
243 class PressureFunction
244 {
245 public:
246 struct InitCond {
247 double depth;
248 double pressure;
249 };
250
251 explicit PressureFunction(const ODE& ode,
252 const InitCond& ic,
253 const int nsample,
254 const VSpan& span);
255
256 PressureFunction(const PressureFunction& rhs);
257
258 PressureFunction(PressureFunction&& rhs) = default;
259
260 PressureFunction& operator=(const PressureFunction& rhs);
261
262 PressureFunction& operator=(PressureFunction&& rhs);
263
264 double value(const double depth) const;
265
266 private:
267 enum Direction : std::size_t { Up, Down, NumDir };
268
269 using Distribution = Details::RK4IVP<ODE>;
270 using DistrPtr = std::unique_ptr<Distribution>;
271
272 InitCond initial_;
273 std::array<DistrPtr, Direction::NumDir> value_;
274 };
275
276 using OilPressODE = PhasePressODE::Oil<
277 FluidSystem, typename Region::CalcDissolution
278 >;
279
280 using GasPressODE = PhasePressODE::Gas<
281 FluidSystem, typename Region::CalcEvaporation, typename Region::CalcWaterEvaporation
282 >;
283
284 using WatPressODE = PhasePressODE::Water<FluidSystem>;
285
286 using OPress = PressureFunction<OilPressODE>;
287 using GPress = PressureFunction<GasPressODE>;
288 using WPress = PressureFunction<WatPressODE>;
289
290 using Strategy = void (PressureTable::*)
291 (const Region&, const VSpan&);
292
293 double gravity_;
294 int nsample_;
295
296 std::unique_ptr<OPress> oil_{};
297 std::unique_ptr<GPress> gas_{};
298 std::unique_ptr<WPress> wat_{};
299
300 template <typename PressFunc>
301 void checkPtr(const PressFunc* phasePress,
302 const std::string& phaseName) const;
303
304 Strategy selectEquilibrationStrategy(const Region& reg) const;
305
306 void copyInPointers(const PressureTable& rhs);
307
308 void equil_WOG(const Region& reg, const VSpan& span);
309 void equil_GOW(const Region& reg, const VSpan& span);
310 void equil_OWG(const Region& reg, const VSpan& span);
311
312 void makeOilPressure(const typename OPress::InitCond& ic,
313 const Region& reg,
314 const VSpan& span);
315
316 void makeGasPressure(const typename GPress::InitCond& ic,
317 const Region& reg,
318 const VSpan& span);
319
320 void makeWatPressure(const typename WPress::InitCond& ic,
321 const Region& reg,
322 const VSpan& span);
323};
324
325// ===========================================================================
326
329 double oil{0.0};
330 double gas{0.0};
331 double water{0.0};
332
333 PhaseQuantityValue& axpy(const PhaseQuantityValue& rhs, const double a)
334 {
335 this->oil += a * rhs.oil;
336 this->gas += a * rhs.gas;
337 this->water += a * rhs.water;
338
339 return *this;
340 }
341
343 {
344 this->oil /= x;
345 this->gas /= x;
346 this->water /= x;
347
348 return *this;
349 }
350
351 void reset()
352 {
353 this->oil = this->gas = this->water = 0.0;
354 }
355};
356
374template <class MaterialLawManager, class FluidSystem, class Region, typename CellID>
376{
377public:
381 struct Position {
382 CellID cell;
383 double depth;
384 };
385
388
396 explicit PhaseSaturations(MaterialLawManager& matLawMgr,
397 const std::vector<double>& swatInit);
398
403
406
409
422 const PhaseQuantityValue&
424 const Region& reg,
425 const PTable& ptable);
426
432 {
433 return this->press_;
434 }
435
436private:
439 struct EvaluationPoint {
440 const Position* position{nullptr};
441 const Region* region {nullptr};
442 const PTable* ptable {nullptr};
443 };
444
448 using FluidState = ::Opm::
449 SimpleModularFluidState<double, /*numPhases=*/3, /*numComponents=*/3,
450 FluidSystem,
451 /*storePressure=*/false,
452 /*storeTemperature=*/false,
453 /*storeComposition=*/false,
454 /*storeFugacity=*/false,
455 /*storeSaturation=*/true,
456 /*storeDensity=*/false,
457 /*storeViscosity=*/false,
458 /*storeEnthalpy=*/false>;
459
461 using MaterialLaw = typename MaterialLawManager::MaterialLaw;
462
464 using PhaseIdx = std::remove_cv_t<
465 std::remove_reference_t<decltype(FluidSystem::oilPhaseIdx)>
466 >;
467
469 MaterialLawManager& matLawMgr_;
470
472 const std::vector<double>& swatInit_;
473
475 PhaseQuantityValue sat_;
476
478 PhaseQuantityValue press_;
479
481 EvaluationPoint evalPt_;
482
484 FluidState fluidState_;
485
487 std::array<double, FluidSystem::numPhases> matLawCapPress_;
488
498 void setEvaluationPoint(const Position& x,
499 const Region& reg,
500 const PTable& ptable);
501
505 void initializePhaseQuantities();
506
510 void deriveOilSat();
511
516 void deriveGasSat();
517
524 void deriveWaterSat();
525
528 void fixUnphysicalTransition();
529
532 void accountForScaledSaturations();
533
534 // --------------------------------------------------------------------
535 // Note: Function 'applySwatInit' is non-const because the overload set
536 // needs to mutate the 'matLawMgr_'.
537 // --------------------------------------------------------------------
538
548 std::pair<double, bool> applySwatInit(const double pcow);
549
562 std::pair<double, bool> applySwatInit(const double pc, const double sw);
563
566 void computeMaterialLawCapPress();
567
570 double materialLawCapPressGasOil() const;
571
574 double materialLawCapPressOilWater() const;
575
578 double materialLawCapPressGasWater() const;
579
587 bool isConstCapPress(const PhaseIdx phaseIdx) const;
588
594 bool isOverlappingTransition() const;
595
615 double fromDepthTable(const double contactdepth,
616 const PhaseIdx phasePos,
617 const bool isincr) const;
618
636 double invertCapPress(const double pc,
637 const PhaseIdx phasePos,
638 const bool isincr) const;
639
641 PhaseIdx oilPos() const
642 {
643 return FluidSystem::oilPhaseIdx;
644 }
645
647 PhaseIdx gasPos() const
648 {
649 return FluidSystem::gasPhaseIdx;
650 }
651
653 PhaseIdx waterPos() const
654 {
655 return FluidSystem::waterPhaseIdx;
656 }
657};
658
659// ===========================================================================
660
661template <typename CellRange, typename Comm>
662void verticalExtent(const CellRange& cells,
663 const std::vector<std::pair<double, double>>& cellZMinMax,
664 const Comm& comm,
665 std::array<double,2>& span);
666
667template <class Element>
668std::pair<double,double> cellZMinMax(const Element& element);
669
670} // namespace Details
671
672namespace DeckDependent {
673
674template<class FluidSystem,
675 class Grid,
676 class GridView,
677 class ElementMapper,
678 class CartesianIndexMapper>
680{
681 using Element = typename GridView::template Codim<0>::Entity;
682public:
683 template<class MaterialLawManager>
684 InitialStateComputer(MaterialLawManager& materialLawManager,
685 const EclipseState& eclipseState,
686 const Grid& grid,
687 const GridView& gridView,
688 const CartesianIndexMapper& cartMapper,
689 const double grav,
690 const int num_pressure_points = 2000,
691 const bool applySwatInit = true);
692
693 using Vec = std::vector<double>;
694 using PVec = std::vector<Vec>; // One per phase.
695
696 const Vec& temperature() const { return temperature_; }
697 const Vec& saltConcentration() const { return saltConcentration_; }
698 const Vec& saltSaturation() const { return saltSaturation_; }
699 const PVec& press() const { return pp_; }
700 const PVec& saturation() const { return sat_; }
701 const Vec& rs() const { return rs_; }
702 const Vec& rv() const { return rv_; }
703 const Vec& rvw() const { return rvw_; }
704
705private:
706 template <class RMap>
707 void updateInitialTemperature_(const EclipseState& eclState, const RMap& reg);
708
709 template <class RMap>
710 void updateInitialSaltConcentration_(const EclipseState& eclState, const RMap& reg);
711
712 template <class RMap>
713 void updateInitialSaltSaturation_(const EclipseState& eclState, const RMap& reg);
714
715 void updateCellProps_(const GridView& gridView,
716 const NumericalAquifers& aquifer);
717
718 void applyNumericalAquifers_(const GridView& gridView,
719 const NumericalAquifers& aquifer,
720 const bool co2store_or_h2store);
721
722 template<class RMap>
723 void setRegionPvtIdx(const EclipseState& eclState, const RMap& reg);
724
725 template <class RMap, class MaterialLawManager, class Comm>
726 void calcPressSatRsRv(const RMap& reg,
727 const std::vector<EquilRecord>& rec,
728 MaterialLawManager& materialLawManager,
729 const Comm& comm,
730 const double grav);
731
732 template <class CellRange, class EquilibrationMethod>
733 void cellLoop(const CellRange& cells,
734 EquilibrationMethod&& eqmethod);
735
736 template <class CellRange, class PressTable, class PhaseSat>
737 void equilibrateCellCentres(const CellRange& cells,
738 const EquilReg& eqreg,
739 const PressTable& ptable,
740 PhaseSat& psat);
741
742 template <class CellRange, class PressTable, class PhaseSat>
743 void equilibrateHorizontal(const CellRange& cells,
744 const EquilReg& eqreg,
745 const int acc,
746 const PressTable& ptable,
747 PhaseSat& psat);
748
749 std::vector< std::shared_ptr<Miscibility::RsFunction> > rsFunc_;
750 std::vector< std::shared_ptr<Miscibility::RsFunction> > rvFunc_;
751 std::vector< std::shared_ptr<Miscibility::RsFunction> > rvwFunc_;
752 using TabulatedFunction = Tabulated1DFunction<double>;
753 std::vector<TabulatedFunction> tempVdTable_;
754 std::vector<TabulatedFunction> saltVdTable_;
755 std::vector<TabulatedFunction> saltpVdTable_;
756 std::vector<int> regionPvtIdx_;
757 Vec temperature_;
758 Vec saltConcentration_;
759 Vec saltSaturation_;
760 PVec pp_;
761 PVec sat_;
762 Vec rs_;
763 Vec rv_;
764 Vec rvw_;
765 const CartesianIndexMapper& cartesianIndexMapper_;
766 Vec swatInit_;
767 Vec cellCenterDepth_;
768 std::vector<std::pair<double,double>> cellZSpan_;
769 std::vector<std::pair<double,double>> cellZMinMax_;
770 int num_pressure_points_;
771};
772
773} // namespace DeckDependent
774} // namespace EQUIL
775} // namespace Opm
776
777#endif // OPM_INIT_STATE_EQUIL_HPP
Dune::OwnerOverlapCopyCommunication< int, int > Comm
Definition: FlexibleSolver_impl.hpp:270
Definition: InitStateEquil.hpp:680
std::vector< Vec > PVec
Definition: InitStateEquil.hpp:694
const PVec & press() const
Definition: InitStateEquil.hpp:699
const Vec & rvw() const
Definition: InitStateEquil.hpp:703
const Vec & rv() const
Definition: InitStateEquil.hpp:702
const Vec & saltSaturation() const
Definition: InitStateEquil.hpp:698
const Vec & saltConcentration() const
Definition: InitStateEquil.hpp:697
const Vec & temperature() const
Definition: InitStateEquil.hpp:696
const Vec & rs() const
Definition: InitStateEquil.hpp:701
std::vector< double > Vec
Definition: InitStateEquil.hpp:693
const PVec & saturation() const
Definition: InitStateEquil.hpp:700
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:132
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:109
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:86
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:376
const PhaseQuantityValue & correctedPhasePressures() const
Definition: InitStateEquil.hpp:431
PressureTable< FluidSystem, Region > PTable
Convenience type alias.
Definition: InitStateEquil.hpp:387
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:159
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:161
double oil(const double depth) const
Definition: InitStateEquil_impl.hpp:1048
Definition: InitStateEquil.hpp:64
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 Comm &comm, std::array< double, 2 > &span)
Definition: InitStateEquil_impl.hpp:64
Definition: BlackoilPhases.hpp:27
Simple set of per-phase (named by primary component) quantities.
Definition: InitStateEquil.hpp:328
double water
Definition: InitStateEquil.hpp:331
PhaseQuantityValue & operator/=(const double x)
Definition: InitStateEquil.hpp:342
PhaseQuantityValue & axpy(const PhaseQuantityValue &rhs, const double a)
Definition: InitStateEquil.hpp:333
double gas
Definition: InitStateEquil.hpp:330
double oil
Definition: InitStateEquil.hpp:329
void reset()
Definition: InitStateEquil.hpp:351
Definition: InitStateEquil.hpp:381
CellID cell
Definition: InitStateEquil.hpp:382
double depth
Definition: InitStateEquil.hpp:383
double pressure
Definition: InitStateEquil.hpp:248
double depth
Definition: InitStateEquil.hpp:247