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
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
61template<class Scalar> struct CellCornerData {
62 std::array<Scalar, 8> X;
63 std::array<Scalar, 8> Y;
64 std::array<Scalar, 8> Z;
65 CellCornerData() = default;
66
67 CellCornerData(const std::array<Scalar, 8>& x,
68 const std::array<Scalar, 8>& y,
69 const std::array<Scalar, 8>& z)
70 : X(x), Y(y), Z(z)
71 {}
72};
73
74template<class Scalar> class EquilReg;
75namespace Miscibility { template<class Scalar> class RsFunction; }
76
77namespace Details {
78template <class Scalar, class RHS>
79class RK4IVP
80{
81public:
82 RK4IVP(const RHS& f,
83 const std::array<Scalar,2>& span,
84 const Scalar y0,
85 const int N);
86
87 Scalar operator()(const Scalar x) const;
88
89private:
90 int N_;
91 std::array<Scalar,2> span_;
92 std::vector<Scalar> y_;
93 std::vector<Scalar> f_;
94
95 Scalar stepsize() const;
96};
97
98namespace PhasePressODE {
99template <class FluidSystem>
100class Water
101{
102 using Scalar = typename FluidSystem::Scalar;
103 using TabulatedFunction = Tabulated1DFunction<Scalar>;
104
105public:
106 Water(const TabulatedFunction& tempVdTable,
107 const TabulatedFunction& saltVdTable,
108 const int pvtRegionIdx,
109 const Scalar normGrav);
110
111 Scalar operator()(const Scalar depth,
112 const Scalar press) const;
113
114private:
115 const TabulatedFunction& tempVdTable_;
116 const TabulatedFunction& saltVdTable_;
117 const int pvtRegionIdx_;
118 const Scalar g_;
119
120 Scalar density(const Scalar depth,
121 const Scalar press) const;
122};
123
124template <class FluidSystem, class RS>
125class Oil
126{
127 using Scalar = typename FluidSystem::Scalar;
128 using TabulatedFunction = Tabulated1DFunction<Scalar>;
129
130public:
131 Oil(const TabulatedFunction& tempVdTable,
132 const RS& rs,
133 const int pvtRegionIdx,
134 const Scalar normGrav);
135
136 Scalar operator()(const Scalar depth,
137 const Scalar press) const;
138
139private:
140 const TabulatedFunction& tempVdTable_;
141 const RS& rs_;
142 const int pvtRegionIdx_;
143 const Scalar g_;
144
145 Scalar density(const Scalar depth,
146 const Scalar press) const;
147};
148
149template <class FluidSystem, class RV, class RVW>
150class Gas
151{
152 using Scalar = typename FluidSystem::Scalar;
153 using TabulatedFunction = Tabulated1DFunction<Scalar>;
154
155public:
156 Gas(const TabulatedFunction& tempVdTable,
157 const RV& rv,
158 const RVW& rvw,
159 const int pvtRegionIdx,
160 const Scalar normGrav);
161
162 Scalar operator()(const Scalar depth,
163 const Scalar press) const;
164
165private:
166 const TabulatedFunction& tempVdTable_;
167 const RV& rv_;
168 const RVW& rvw_;
169 const int pvtRegionIdx_;
170 const Scalar g_;
171
172 Scalar density(const Scalar depth,
173 const Scalar press) const;
174};
175
176} // namespace PhasePressODE
177
178template <class FluidSystem, class Region>
180{
181public:
182 using Scalar = typename FluidSystem::Scalar;
183 using VSpan = std::array<Scalar, 2>;
184
193 explicit PressureTable(const Scalar gravity,
194 const int samplePoints = 2000);
195
199 PressureTable(const PressureTable& rhs);
200
207
214
223
224 void equilibrate(const Region& reg,
225 const VSpan& span);
226
228 bool oilActive() const;
229
231 bool gasActive() const;
232
234 bool waterActive() const;
235
243 Scalar oil(const Scalar depth) const;
244
252 Scalar gas(const Scalar depth) const;
253
261 Scalar water(const Scalar depth) const;
262
263private:
264 template <class ODE>
265 class PressureFunction
266 {
267 public:
268 struct InitCond {
271 };
272
273 explicit PressureFunction(const ODE& ode,
274 const InitCond& ic,
275 const int nsample,
276 const VSpan& span);
277
278 PressureFunction(const PressureFunction& rhs);
279
280 PressureFunction(PressureFunction&& rhs) = default;
281
282 PressureFunction& operator=(const PressureFunction& rhs);
283
284 PressureFunction& operator=(PressureFunction&& rhs);
285
286 Scalar value(const Scalar depth) const;
287
288 private:
289 enum Direction : std::size_t { Up, Down, NumDir };
290
291 using Distribution = Details::RK4IVP<Scalar,ODE>;
292 using DistrPtr = std::unique_ptr<Distribution>;
293
294 InitCond initial_;
295 std::array<DistrPtr, Direction::NumDir> value_;
296 };
297
298 using OilPressODE = PhasePressODE::Oil<
299 FluidSystem, typename Region::CalcDissolution
300 >;
301
302 using GasPressODE = PhasePressODE::Gas<
303 FluidSystem, typename Region::CalcEvaporation, typename Region::CalcWaterEvaporation
304 >;
305
306 using WatPressODE = PhasePressODE::Water<FluidSystem>;
307
308 using OPress = PressureFunction<OilPressODE>;
309 using GPress = PressureFunction<GasPressODE>;
310 using WPress = PressureFunction<WatPressODE>;
311
312 using Strategy = void (PressureTable::*)
313 (const Region&, const VSpan&);
314
315 Scalar gravity_;
316 int nsample_;
317
318 std::unique_ptr<OPress> oil_{};
319 std::unique_ptr<GPress> gas_{};
320 std::unique_ptr<WPress> wat_{};
321
322 template <typename PressFunc>
323 void checkPtr(const PressFunc* phasePress,
324 const std::string& phaseName) const;
325
326 Strategy selectEquilibrationStrategy(const Region& reg) const;
327
328 void copyInPointers(const PressureTable& rhs);
329
330 void equil_WOG(const Region& reg, const VSpan& span);
331 void equil_GOW(const Region& reg, const VSpan& span);
332 void equil_OWG(const Region& reg, const VSpan& span);
333
334 void makeOilPressure(const typename OPress::InitCond& ic,
335 const Region& reg,
336 const VSpan& span);
337
338 void makeGasPressure(const typename GPress::InitCond& ic,
339 const Region& reg,
340 const VSpan& span);
341
342 void makeWatPressure(const typename WPress::InitCond& ic,
343 const Region& reg,
344 const VSpan& span);
345};
346
347// ===========================================================================
348
350template<class Scalar>
352 Scalar oil{0.0};
353 Scalar gas{0.0};
354 Scalar water{0.0};
355
356 PhaseQuantityValue& axpy(const PhaseQuantityValue& rhs, const Scalar a)
357 {
358 this->oil += a * rhs.oil;
359 this->gas += a * rhs.gas;
360 this->water += a * rhs.water;
361
362 return *this;
363 }
364
366 {
367 this->oil /= x;
368 this->gas /= x;
369 this->water /= x;
370
371 return *this;
372 }
373
374 void reset()
375 {
376 this->oil = this->gas = this->water = 0.0;
377 }
378};
379
397template <class MaterialLawManager, class FluidSystem, class Region, typename CellID>
399{
400public:
401 using Scalar = typename FluidSystem::Scalar;
405 struct Position {
406 CellID cell;
408 };
409
412
420 explicit PhaseSaturations(MaterialLawManager& matLawMgr,
421 const std::vector<Scalar>& swatInit);
422
427
430
433
448 const Region& reg,
449 const PTable& ptable);
450
456 {
457 return this->press_;
458 }
459
460private:
463 struct EvaluationPoint {
464 const Position* position{nullptr};
465 const Region* region {nullptr};
466 const PTable* ptable {nullptr};
467 };
468
472 using FluidState = ::Opm::
473 SimpleModularFluidState<Scalar, /*numPhases=*/3, /*numComponents=*/3,
474 FluidSystem,
475 /*storePressure=*/false,
476 /*storeTemperature=*/false,
477 /*storeComposition=*/false,
478 /*storeFugacity=*/false,
479 /*storeSaturation=*/true,
480 /*storeDensity=*/false,
481 /*storeViscosity=*/false,
482 /*storeEnthalpy=*/false>;
483
485 using MaterialLaw = typename MaterialLawManager::MaterialLaw;
486
488 using PhaseIdx = std::remove_cv_t<
489 std::remove_reference_t<decltype(FluidSystem::oilPhaseIdx)>
490 >;
491
493 MaterialLawManager& matLawMgr_;
494
496 const std::vector<Scalar>& swatInit_;
497
499 PhaseQuantityValue<Scalar> sat_;
500
502 PhaseQuantityValue<Scalar> press_;
503
505 EvaluationPoint evalPt_;
506
508 FluidState fluidState_;
509
511 std::array<Scalar, FluidSystem::numPhases> matLawCapPress_;
512
522 void setEvaluationPoint(const Position& x,
523 const Region& reg,
524 const PTable& ptable);
525
529 void initializePhaseQuantities();
530
534 void deriveOilSat();
535
540 void deriveGasSat();
541
548 void deriveWaterSat();
549
552 void fixUnphysicalTransition();
553
556 void accountForScaledSaturations();
557
558 // --------------------------------------------------------------------
559 // Note: Function 'applySwatInit' is non-const because the overload set
560 // needs to mutate the 'matLawMgr_'.
561 // --------------------------------------------------------------------
562
572 std::pair<Scalar, bool> applySwatInit(const Scalar pcow);
573
586 std::pair<Scalar, bool> applySwatInit(const Scalar pc, const Scalar sw);
587
590 void computeMaterialLawCapPress();
591
594 Scalar materialLawCapPressGasOil() const;
595
598 Scalar materialLawCapPressOilWater() const;
599
602 Scalar materialLawCapPressGasWater() const;
603
611 bool isConstCapPress(const PhaseIdx phaseIdx) const;
612
618 bool isOverlappingTransition() const;
619
639 Scalar fromDepthTable(const Scalar contactdepth,
640 const PhaseIdx phasePos,
641 const bool isincr) const;
642
660 Scalar invertCapPress(const Scalar pc,
661 const PhaseIdx phasePos,
662 const bool isincr) const;
663
665 PhaseIdx oilPos() const
666 {
667 return FluidSystem::oilPhaseIdx;
668 }
669
671 PhaseIdx gasPos() const
672 {
673 return FluidSystem::gasPhaseIdx;
674 }
675
677 PhaseIdx waterPos() const
678 {
679 return FluidSystem::waterPhaseIdx;
680 }
681};
682
683// ===========================================================================
684
685template <typename CellRange, class Scalar>
686void verticalExtent(const CellRange& cells,
687 const std::vector<std::pair<Scalar, Scalar>>& cellZMinMax,
688 const Parallel::Communication& comm,
689 std::array<Scalar,2>& span);
690
691template <class Scalar, class Element>
692std::pair<Scalar,Scalar> cellZMinMax(const Element& element);
693
694} // namespace Details
695
696namespace DeckDependent {
697
698template<class FluidSystem,
699 class Grid,
700 class GridView,
701 class ElementMapper,
702 class CartesianIndexMapper>
704{
705 using Element = typename GridView::template Codim<0>::Entity;
706 using Scalar = typename FluidSystem::Scalar;
707public:
708 template<class MaterialLawManager>
709 InitialStateComputer(MaterialLawManager& materialLawManager,
710 const EclipseState& eclipseState,
711 const Grid& grid,
712 const GridView& gridView,
713 const CartesianIndexMapper& cartMapper,
714 const Scalar grav,
715 const int num_pressure_points = 2000,
716 const bool applySwatInit = true);
717
718 using Vec = std::vector<Scalar>;
719 using PVec = std::vector<Vec>; // One per phase.
720
721 const Vec& temperature() const { return temperature_; }
722 const Vec& saltConcentration() const { return saltConcentration_; }
723 const Vec& saltSaturation() const { return saltSaturation_; }
724 const PVec& press() const { return pp_; }
725 const PVec& saturation() const { return sat_; }
726 const Vec& rs() const { return rs_; }
727 const Vec& rv() const { return rv_; }
728 const Vec& rvw() const { return rvw_; }
729
730private:
731 template <class RMap>
732 void updateInitialTemperature_(const EclipseState& eclState, const RMap& reg);
733
734 template <class RMap>
735 void updateInitialSaltConcentration_(const EclipseState& eclState, const RMap& reg);
736
737 template <class RMap>
738 void updateInitialSaltSaturation_(const EclipseState& eclState, const RMap& reg);
739
740 void updateCellProps_(const GridView& gridView,
741 const NumericalAquifers& aquifer);
742
743 void applyNumericalAquifers_(const GridView& gridView,
744 const NumericalAquifers& aquifer,
745 const bool co2store_or_h2store);
746
747 template<class RMap>
748 void setRegionPvtIdx(const EclipseState& eclState, const RMap& reg);
749
750 template <class RMap, class MaterialLawManager, class Comm>
751 void calcPressSatRsRv(const RMap& reg,
752 const std::vector<EquilRecord>& rec,
753 MaterialLawManager& materialLawManager,
754 const GridView& gridView,
755 const Comm& comm,
756 const Scalar grav);
757
758 template <class CellRange, class EquilibrationMethod>
759 void cellLoop(const CellRange& cells,
760 EquilibrationMethod&& eqmethod);
761
762 template <class CellRange, class PressTable, class PhaseSat>
763 void equilibrateCellCentres(const CellRange& cells,
764 const EquilReg<Scalar>& eqreg,
765 const PressTable& ptable,
766 PhaseSat& psat);
767
768 template <class CellRange, class PressTable, class PhaseSat>
769 void equilibrateHorizontal(const CellRange& cells,
770 const EquilReg<Scalar>& eqreg,
771 const int acc,
772 const PressTable& ptable,
773 PhaseSat& psat);
774
775 template<class CellRange, class PressTable, class PhaseSat>
776 void equilibrateTiltedFaultBlock(const CellRange& cells,
777 const EquilReg<Scalar>& eqreg,
778 const GridView& gridView, const int numLevels,
779 const PressTable& ptable, PhaseSat& psat);
780
781 template<class CellRange, class PressTable, class PhaseSat>
782 void equilibrateTiltedFaultBlockSimple(const CellRange& cells,
783 const EquilReg<Scalar>& eqreg,
784 const GridView& gridView, const int numLevels,
785 const PressTable& ptable, PhaseSat& psat);
786
787 std::vector< std::shared_ptr<Miscibility::RsFunction<Scalar>> > rsFunc_;
788 std::vector< std::shared_ptr<Miscibility::RsFunction<Scalar>> > rvFunc_;
789 std::vector< std::shared_ptr<Miscibility::RsFunction<Scalar>> > rvwFunc_;
790 using TabulatedFunction = Tabulated1DFunction<Scalar>;
791 std::vector<TabulatedFunction> tempVdTable_;
792 std::vector<TabulatedFunction> saltVdTable_;
793 std::vector<TabulatedFunction> saltpVdTable_;
794 std::vector<int> regionPvtIdx_;
795 Vec temperature_;
796 Vec saltConcentration_;
797 Vec saltSaturation_;
798 PVec pp_;
799 PVec sat_;
800 Vec rs_;
801 Vec rv_;
802 Vec rvw_;
803 const CartesianIndexMapper& cartesianIndexMapper_;
804 Vec swatInit_;
805 Vec cellCenterDepth_;
806 std::vector<std::pair<Scalar,Scalar>> cellCenterXY_;
807 std::vector<std::pair<Scalar,Scalar>> cellZSpan_;
808 std::vector<std::pair<Scalar,Scalar>> cellZMinMax_;
809 std::vector<CellCornerData<Scalar>> cellCorners_;
810 int num_pressure_points_;
811};
812
813} // namespace DeckDependent
814} // namespace EQUIL
815} // namespace Opm
816
817#endif // OPM_INIT_STATE_EQUIL_HPP
Dune::OwnerOverlapCopyCommunication< int, int > Comm
Definition: FlexibleSolver_impl.hpp:325
Definition: InitStateEquil.hpp:704
std::vector< Vec > PVec
Definition: InitStateEquil.hpp:719
const PVec & press() const
Definition: InitStateEquil.hpp:724
const Vec & rvw() const
Definition: InitStateEquil.hpp:728
InitialStateComputer(MaterialLawManager &materialLawManager, const EclipseState &eclipseState, const Grid &grid, const GridView &gridView, const CartesianIndexMapper &cartMapper, const Scalar grav, const int num_pressure_points=2000, const bool applySwatInit=true)
Definition: InitStateEquil_impl.hpp:1472
const Vec & rv() const
Definition: InitStateEquil.hpp:727
const Vec & saltSaturation() const
Definition: InitStateEquil.hpp:723
const Vec & saltConcentration() const
Definition: InitStateEquil.hpp:722
const Vec & temperature() const
Definition: InitStateEquil.hpp:721
std::vector< Scalar > Vec
Definition: InitStateEquil.hpp:718
const Vec & rs() const
Definition: InitStateEquil.hpp:726
const PVec & saturation() const
Definition: InitStateEquil.hpp:725
Definition: InitStateEquil.hpp:151
Gas(const TabulatedFunction &tempVdTable, const RV &rv, const RVW &rvw, const int pvtRegionIdx, const Scalar normGrav)
Definition: InitStateEquil_impl.hpp:472
Scalar operator()(const Scalar depth, const Scalar press) const
Definition: InitStateEquil_impl.hpp:488
Definition: InitStateEquil.hpp:126
Oil(const TabulatedFunction &tempVdTable, const RS &rs, const int pvtRegionIdx, const Scalar normGrav)
Definition: InitStateEquil_impl.hpp:424
Scalar operator()(const Scalar depth, const Scalar press) const
Definition: InitStateEquil_impl.hpp:438
Definition: InitStateEquil.hpp:101
Scalar operator()(const Scalar depth, const Scalar press) const
Definition: InitStateEquil_impl.hpp:398
Water(const TabulatedFunction &tempVdTable, const TabulatedFunction &saltVdTable, const int pvtRegionIdx, const Scalar normGrav)
Definition: InitStateEquil_impl.hpp:384
Definition: InitStateEquil.hpp:399
const PhaseQuantityValue< Scalar > & deriveSaturations(const Position &x, const Region &reg, const PTable &ptable)
Definition: InitStateEquil_impl.hpp:730
PhaseSaturations(MaterialLawManager &matLawMgr, const std::vector< Scalar > &swatInit)
Definition: InitStateEquil_impl.hpp:706
typename FluidSystem::Scalar Scalar
Definition: InitStateEquil.hpp:401
PressureTable< FluidSystem, Region > PTable
Convenience type alias.
Definition: InitStateEquil.hpp:411
PhaseSaturations & operator=(const PhaseSaturations &)=delete
Disabled assignment operator.
PhaseSaturations & operator=(PhaseSaturations &&)=delete
Disabled move-assignment operator.
const PhaseQuantityValue< Scalar > & correctedPhasePressures() const
Definition: InitStateEquil.hpp:455
Definition: InitStateEquil.hpp:180
PressureTable & operator=(const PressureTable &rhs)
Definition: InitStateEquil_impl.hpp:1151
Scalar water(const Scalar depth) const
Definition: InitStateEquil_impl.hpp:1231
Scalar gas(const Scalar depth) const
Definition: InitStateEquil_impl.hpp:1220
bool waterActive() const
Predicate for whether or not water is an active phase.
Definition: InitStateEquil_impl.hpp:1202
bool gasActive() const
Predicate for whether or not gas is an active phase.
Definition: InitStateEquil_impl.hpp:1195
Scalar oil(const Scalar depth) const
Definition: InitStateEquil_impl.hpp:1210
std::array< Scalar, 2 > VSpan
Definition: InitStateEquil.hpp:183
bool oilActive() const
Predicate for whether or not oil is an active phase.
Definition: InitStateEquil_impl.hpp:1188
typename FluidSystem::Scalar Scalar
Definition: InitStateEquil.hpp:182
void equilibrate(const Region &reg, const VSpan &span)
Definition: InitStateEquil_impl.hpp:1177
PressureTable(const Scalar gravity, const int samplePoints=2000)
Definition: InitStateEquil_impl.hpp:1121
Definition: InitStateEquil.hpp:80
Scalar operator()(const Scalar x) const
Definition: InitStateEquil_impl.hpp:350
RK4IVP(const RHS &f, const std::array< Scalar, 2 > &span, const Scalar y0, const int N)
Definition: InitStateEquil_impl.hpp:315
Definition: EquilibrationHelpers.hpp:627
std::pair< Scalar, Scalar > cellZMinMax(const Element &element)
Definition: InitStateEquil_impl.hpp:182
void verticalExtent(const CellRange &cells, const std::vector< std::pair< Scalar, Scalar > > &cellZMinMax, const Parallel::Communication &comm, std::array< Scalar, 2 > &span)
Definition: InitStateEquil_impl.hpp:65
Dune::Communication< MPIComm > Communication
Definition: ParallelCommunication.hpp:30
Definition: blackoilbioeffectsmodules.hh:43
The Opm property system, traits with inheritance.
Definition: InitStateEquil.hpp:61
std::array< Scalar, 8 > X
Definition: InitStateEquil.hpp:62
std::array< Scalar, 8 > Y
Definition: InitStateEquil.hpp:63
CellCornerData(const std::array< Scalar, 8 > &x, const std::array< Scalar, 8 > &y, const std::array< Scalar, 8 > &z)
Definition: InitStateEquil.hpp:67
std::array< Scalar, 8 > Z
Definition: InitStateEquil.hpp:64
Simple set of per-phase (named by primary component) quantities.
Definition: InitStateEquil.hpp:351
void reset()
Definition: InitStateEquil.hpp:374
Scalar gas
Definition: InitStateEquil.hpp:353
Scalar water
Definition: InitStateEquil.hpp:354
PhaseQuantityValue & operator/=(const Scalar x)
Definition: InitStateEquil.hpp:365
Scalar oil
Definition: InitStateEquil.hpp:352
PhaseQuantityValue & axpy(const PhaseQuantityValue &rhs, const Scalar a)
Definition: InitStateEquil.hpp:356
Definition: InitStateEquil.hpp:405
CellID cell
Definition: InitStateEquil.hpp:406
Scalar depth
Definition: InitStateEquil.hpp:407
Scalar depth
Definition: InitStateEquil.hpp:269
Scalar pressure
Definition: InitStateEquil.hpp:270