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> class EquilReg;
62namespace Miscibility { template<class Scalar> class RsFunction; }
63
64namespace Details {
65template <class Scalar, class RHS>
66class RK4IVP
67{
68public:
69 RK4IVP(const RHS& f,
70 const std::array<Scalar,2>& span,
71 const Scalar y0,
72 const int N);
73
74 Scalar operator()(const Scalar x) const;
75
76private:
77 int N_;
78 std::array<Scalar,2> span_;
79 std::vector<Scalar> y_;
80 std::vector<Scalar> f_;
81
82 Scalar stepsize() const;
83};
84
85namespace PhasePressODE {
86template <class FluidSystem>
87class Water
88{
89 using Scalar = typename FluidSystem::Scalar;
90 using TabulatedFunction = Tabulated1DFunction<Scalar>;
91
92public:
93 Water(const TabulatedFunction& tempVdTable,
94 const TabulatedFunction& saltVdTable,
95 const int pvtRegionIdx,
96 const Scalar normGrav);
97
98 Scalar operator()(const Scalar depth,
99 const Scalar press) const;
100
101private:
102 const TabulatedFunction& tempVdTable_;
103 const TabulatedFunction& saltVdTable_;
104 const int pvtRegionIdx_;
105 const Scalar g_;
106
107 Scalar density(const Scalar depth,
108 const Scalar press) const;
109};
110
111template <class FluidSystem, class RS>
112class Oil
113{
114 using Scalar = typename FluidSystem::Scalar;
115 using TabulatedFunction = Tabulated1DFunction<Scalar>;
116
117public:
118 Oil(const TabulatedFunction& tempVdTable,
119 const RS& rs,
120 const int pvtRegionIdx,
121 const Scalar normGrav);
122
123 Scalar operator()(const Scalar depth,
124 const Scalar press) const;
125
126private:
127 const TabulatedFunction& tempVdTable_;
128 const RS& rs_;
129 const int pvtRegionIdx_;
130 const Scalar g_;
131
132 Scalar density(const Scalar depth,
133 const Scalar press) const;
134};
135
136template <class FluidSystem, class RV, class RVW>
137class Gas
138{
139 using Scalar = typename FluidSystem::Scalar;
140 using TabulatedFunction = Tabulated1DFunction<Scalar>;
141
142public:
143 Gas(const TabulatedFunction& tempVdTable,
144 const RV& rv,
145 const RVW& rvw,
146 const int pvtRegionIdx,
147 const Scalar normGrav);
148
149 Scalar operator()(const Scalar depth,
150 const Scalar press) const;
151
152private:
153 const TabulatedFunction& tempVdTable_;
154 const RV& rv_;
155 const RVW& rvw_;
156 const int pvtRegionIdx_;
157 const Scalar g_;
158
159 Scalar density(const Scalar depth,
160 const Scalar press) const;
161};
162
163} // namespace PhasePressODE
164
165template <class FluidSystem, class Region>
167{
168public:
169 using Scalar = typename FluidSystem::Scalar;
170 using VSpan = std::array<Scalar, 2>;
171
180 explicit PressureTable(const Scalar gravity,
181 const int samplePoints = 2000);
182
186 PressureTable(const PressureTable& rhs);
187
194
201
210
211 void equilibrate(const Region& reg,
212 const VSpan& span);
213
215 bool oilActive() const;
216
218 bool gasActive() const;
219
221 bool waterActive() const;
222
230 Scalar oil(const Scalar depth) const;
231
239 Scalar gas(const Scalar depth) const;
240
248 Scalar water(const Scalar depth) const;
249
250private:
251 template <class ODE>
252 class PressureFunction
253 {
254 public:
255 struct InitCond {
258 };
259
260 explicit PressureFunction(const ODE& ode,
261 const InitCond& ic,
262 const int nsample,
263 const VSpan& span);
264
265 PressureFunction(const PressureFunction& rhs);
266
267 PressureFunction(PressureFunction&& rhs) = default;
268
269 PressureFunction& operator=(const PressureFunction& rhs);
270
271 PressureFunction& operator=(PressureFunction&& rhs);
272
273 Scalar value(const Scalar depth) const;
274
275 private:
276 enum Direction : std::size_t { Up, Down, NumDir };
277
278 using Distribution = Details::RK4IVP<Scalar,ODE>;
279 using DistrPtr = std::unique_ptr<Distribution>;
280
281 InitCond initial_;
282 std::array<DistrPtr, Direction::NumDir> value_;
283 };
284
285 using OilPressODE = PhasePressODE::Oil<
286 FluidSystem, typename Region::CalcDissolution
287 >;
288
289 using GasPressODE = PhasePressODE::Gas<
290 FluidSystem, typename Region::CalcEvaporation, typename Region::CalcWaterEvaporation
291 >;
292
293 using WatPressODE = PhasePressODE::Water<FluidSystem>;
294
295 using OPress = PressureFunction<OilPressODE>;
296 using GPress = PressureFunction<GasPressODE>;
297 using WPress = PressureFunction<WatPressODE>;
298
299 using Strategy = void (PressureTable::*)
300 (const Region&, const VSpan&);
301
302 Scalar gravity_;
303 int nsample_;
304
305 std::unique_ptr<OPress> oil_{};
306 std::unique_ptr<GPress> gas_{};
307 std::unique_ptr<WPress> wat_{};
308
309 template <typename PressFunc>
310 void checkPtr(const PressFunc* phasePress,
311 const std::string& phaseName) const;
312
313 Strategy selectEquilibrationStrategy(const Region& reg) const;
314
315 void copyInPointers(const PressureTable& rhs);
316
317 void equil_WOG(const Region& reg, const VSpan& span);
318 void equil_GOW(const Region& reg, const VSpan& span);
319 void equil_OWG(const Region& reg, const VSpan& span);
320
321 void makeOilPressure(const typename OPress::InitCond& ic,
322 const Region& reg,
323 const VSpan& span);
324
325 void makeGasPressure(const typename GPress::InitCond& ic,
326 const Region& reg,
327 const VSpan& span);
328
329 void makeWatPressure(const typename WPress::InitCond& ic,
330 const Region& reg,
331 const VSpan& span);
332};
333
334// ===========================================================================
335
337template<class Scalar>
339 Scalar oil{0.0};
340 Scalar gas{0.0};
341 Scalar water{0.0};
342
343 PhaseQuantityValue& axpy(const PhaseQuantityValue& rhs, const Scalar a)
344 {
345 this->oil += a * rhs.oil;
346 this->gas += a * rhs.gas;
347 this->water += a * rhs.water;
348
349 return *this;
350 }
351
353 {
354 this->oil /= x;
355 this->gas /= x;
356 this->water /= x;
357
358 return *this;
359 }
360
361 void reset()
362 {
363 this->oil = this->gas = this->water = 0.0;
364 }
365};
366
384template <class MaterialLawManager, class FluidSystem, class Region, typename CellID>
386{
387public:
388 using Scalar = typename FluidSystem::Scalar;
392 struct Position {
393 CellID cell;
395 };
396
399
407 explicit PhaseSaturations(MaterialLawManager& matLawMgr,
408 const std::vector<Scalar>& swatInit);
409
414
417
420
435 const Region& reg,
436 const PTable& ptable);
437
443 {
444 return this->press_;
445 }
446
447private:
450 struct EvaluationPoint {
451 const Position* position{nullptr};
452 const Region* region {nullptr};
453 const PTable* ptable {nullptr};
454 };
455
459 using FluidState = ::Opm::
460 SimpleModularFluidState<Scalar, /*numPhases=*/3, /*numComponents=*/3,
461 FluidSystem,
462 /*storePressure=*/false,
463 /*storeTemperature=*/false,
464 /*storeComposition=*/false,
465 /*storeFugacity=*/false,
466 /*storeSaturation=*/true,
467 /*storeDensity=*/false,
468 /*storeViscosity=*/false,
469 /*storeEnthalpy=*/false>;
470
472 using MaterialLaw = typename MaterialLawManager::MaterialLaw;
473
475 using PhaseIdx = std::remove_cv_t<
476 std::remove_reference_t<decltype(FluidSystem::oilPhaseIdx)>
477 >;
478
480 MaterialLawManager& matLawMgr_;
481
483 const std::vector<Scalar>& swatInit_;
484
486 PhaseQuantityValue<Scalar> sat_;
487
489 PhaseQuantityValue<Scalar> press_;
490
492 EvaluationPoint evalPt_;
493
495 FluidState fluidState_;
496
498 std::array<Scalar, FluidSystem::numPhases> matLawCapPress_;
499
509 void setEvaluationPoint(const Position& x,
510 const Region& reg,
511 const PTable& ptable);
512
516 void initializePhaseQuantities();
517
521 void deriveOilSat();
522
527 void deriveGasSat();
528
535 void deriveWaterSat();
536
539 void fixUnphysicalTransition();
540
543 void accountForScaledSaturations();
544
545 // --------------------------------------------------------------------
546 // Note: Function 'applySwatInit' is non-const because the overload set
547 // needs to mutate the 'matLawMgr_'.
548 // --------------------------------------------------------------------
549
559 std::pair<Scalar, bool> applySwatInit(const Scalar pcow);
560
573 std::pair<Scalar, bool> applySwatInit(const Scalar pc, const Scalar sw);
574
577 void computeMaterialLawCapPress();
578
581 Scalar materialLawCapPressGasOil() const;
582
585 Scalar materialLawCapPressOilWater() const;
586
589 Scalar materialLawCapPressGasWater() const;
590
598 bool isConstCapPress(const PhaseIdx phaseIdx) const;
599
605 bool isOverlappingTransition() const;
606
626 Scalar fromDepthTable(const Scalar contactdepth,
627 const PhaseIdx phasePos,
628 const bool isincr) const;
629
647 Scalar invertCapPress(const Scalar pc,
648 const PhaseIdx phasePos,
649 const bool isincr) const;
650
652 PhaseIdx oilPos() const
653 {
654 return FluidSystem::oilPhaseIdx;
655 }
656
658 PhaseIdx gasPos() const
659 {
660 return FluidSystem::gasPhaseIdx;
661 }
662
664 PhaseIdx waterPos() const
665 {
666 return FluidSystem::waterPhaseIdx;
667 }
668};
669
670// ===========================================================================
671
672template <typename CellRange, class Scalar>
673void verticalExtent(const CellRange& cells,
674 const std::vector<std::pair<Scalar, Scalar>>& cellZMinMax,
675 const Parallel::Communication& comm,
676 std::array<Scalar,2>& span);
677
678template <class Scalar, class Element>
679std::pair<Scalar,Scalar> cellZMinMax(const Element& element);
680
681} // namespace Details
682
683namespace DeckDependent {
684
685template<class FluidSystem,
686 class Grid,
687 class GridView,
688 class ElementMapper,
689 class CartesianIndexMapper>
691{
692 using Element = typename GridView::template Codim<0>::Entity;
693 using Scalar = typename FluidSystem::Scalar;
694public:
695 template<class MaterialLawManager>
696 InitialStateComputer(MaterialLawManager& materialLawManager,
697 const EclipseState& eclipseState,
698 const Grid& grid,
699 const GridView& gridView,
700 const CartesianIndexMapper& cartMapper,
701 const Scalar grav,
702 const int num_pressure_points = 2000,
703 const bool applySwatInit = true);
704
705 using Vec = std::vector<Scalar>;
706 using PVec = std::vector<Vec>; // One per phase.
707
708 const Vec& temperature() const { return temperature_; }
709 const Vec& saltConcentration() const { return saltConcentration_; }
710 const Vec& saltSaturation() const { return saltSaturation_; }
711 const PVec& press() const { return pp_; }
712 const PVec& saturation() const { return sat_; }
713 const Vec& rs() const { return rs_; }
714 const Vec& rv() const { return rv_; }
715 const Vec& rvw() const { return rvw_; }
716
717private:
718 template <class RMap>
719 void updateInitialTemperature_(const EclipseState& eclState, const RMap& reg);
720
721 template <class RMap>
722 void updateInitialSaltConcentration_(const EclipseState& eclState, const RMap& reg);
723
724 template <class RMap>
725 void updateInitialSaltSaturation_(const EclipseState& eclState, const RMap& reg);
726
727 void updateCellProps_(const GridView& gridView,
728 const NumericalAquifers& aquifer);
729
730 void applyNumericalAquifers_(const GridView& gridView,
731 const NumericalAquifers& aquifer,
732 const bool co2store_or_h2store);
733
734 template<class RMap>
735 void setRegionPvtIdx(const EclipseState& eclState, const RMap& reg);
736
737 template <class RMap, class MaterialLawManager, class Comm>
738 void calcPressSatRsRv(const RMap& reg,
739 const std::vector<EquilRecord>& rec,
740 MaterialLawManager& materialLawManager,
741 const Comm& comm,
742 const Scalar grav);
743
744 template <class CellRange, class EquilibrationMethod>
745 void cellLoop(const CellRange& cells,
746 EquilibrationMethod&& eqmethod);
747
748 template <class CellRange, class PressTable, class PhaseSat>
749 void equilibrateCellCentres(const CellRange& cells,
750 const EquilReg<Scalar>& eqreg,
751 const PressTable& ptable,
752 PhaseSat& psat);
753
754 template <class CellRange, class PressTable, class PhaseSat>
755 void equilibrateHorizontal(const CellRange& cells,
756 const EquilReg<Scalar>& eqreg,
757 const int acc,
758 const PressTable& ptable,
759 PhaseSat& psat);
760
761 std::vector< std::shared_ptr<Miscibility::RsFunction<Scalar>> > rsFunc_;
762 std::vector< std::shared_ptr<Miscibility::RsFunction<Scalar>> > rvFunc_;
763 std::vector< std::shared_ptr<Miscibility::RsFunction<Scalar>> > rvwFunc_;
764 using TabulatedFunction = Tabulated1DFunction<Scalar>;
765 std::vector<TabulatedFunction> tempVdTable_;
766 std::vector<TabulatedFunction> saltVdTable_;
767 std::vector<TabulatedFunction> saltpVdTable_;
768 std::vector<int> regionPvtIdx_;
769 Vec temperature_;
770 Vec saltConcentration_;
771 Vec saltSaturation_;
772 PVec pp_;
773 PVec sat_;
774 Vec rs_;
775 Vec rv_;
776 Vec rvw_;
777 const CartesianIndexMapper& cartesianIndexMapper_;
778 Vec swatInit_;
779 Vec cellCenterDepth_;
780 std::vector<std::pair<Scalar,Scalar>> cellZSpan_;
781 std::vector<std::pair<Scalar,Scalar>> cellZMinMax_;
782 int num_pressure_points_;
783};
784
785} // namespace DeckDependent
786} // namespace EQUIL
787} // namespace Opm
788
789#endif // OPM_INIT_STATE_EQUIL_HPP
Dune::OwnerOverlapCopyCommunication< int, int > Comm
Definition: FlexibleSolver_impl.hpp:285
Definition: InitStateEquil.hpp:691
std::vector< Vec > PVec
Definition: InitStateEquil.hpp:706
const PVec & press() const
Definition: InitStateEquil.hpp:711
const Vec & rvw() const
Definition: InitStateEquil.hpp:715
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:1338
const Vec & rv() const
Definition: InitStateEquil.hpp:714
const Vec & saltSaturation() const
Definition: InitStateEquil.hpp:710
const Vec & saltConcentration() const
Definition: InitStateEquil.hpp:709
const Vec & temperature() const
Definition: InitStateEquil.hpp:708
std::vector< Scalar > Vec
Definition: InitStateEquil.hpp:705
const Vec & rs() const
Definition: InitStateEquil.hpp:713
const PVec & saturation() const
Definition: InitStateEquil.hpp:712
Definition: InitStateEquil.hpp:138
Gas(const TabulatedFunction &tempVdTable, const RV &rv, const RVW &rvw, const int pvtRegionIdx, const Scalar normGrav)
Definition: InitStateEquil_impl.hpp:338
Scalar operator()(const Scalar depth, const Scalar press) const
Definition: InitStateEquil_impl.hpp:354
Definition: InitStateEquil.hpp:113
Oil(const TabulatedFunction &tempVdTable, const RS &rs, const int pvtRegionIdx, const Scalar normGrav)
Definition: InitStateEquil_impl.hpp:290
Scalar operator()(const Scalar depth, const Scalar press) const
Definition: InitStateEquil_impl.hpp:304
Definition: InitStateEquil.hpp:88
Scalar operator()(const Scalar depth, const Scalar press) const
Definition: InitStateEquil_impl.hpp:264
Water(const TabulatedFunction &tempVdTable, const TabulatedFunction &saltVdTable, const int pvtRegionIdx, const Scalar normGrav)
Definition: InitStateEquil_impl.hpp:250
Definition: InitStateEquil.hpp:386
const PhaseQuantityValue< Scalar > & deriveSaturations(const Position &x, const Region &reg, const PTable &ptable)
Definition: InitStateEquil_impl.hpp:596
PhaseSaturations(MaterialLawManager &matLawMgr, const std::vector< Scalar > &swatInit)
Definition: InitStateEquil_impl.hpp:572
typename FluidSystem::Scalar Scalar
Definition: InitStateEquil.hpp:388
PressureTable< FluidSystem, Region > PTable
Convenience type alias.
Definition: InitStateEquil.hpp:398
PhaseSaturations & operator=(const PhaseSaturations &)=delete
Disabled assignment operator.
PhaseSaturations & operator=(PhaseSaturations &&)=delete
Disabled move-assignment operator.
const PhaseQuantityValue< Scalar > & correctedPhasePressures() const
Definition: InitStateEquil.hpp:442
Definition: InitStateEquil.hpp:167
PressureTable & operator=(const PressureTable &rhs)
Definition: InitStateEquil_impl.hpp:1017
Scalar water(const Scalar depth) const
Definition: InitStateEquil_impl.hpp:1097
Scalar gas(const Scalar depth) const
Definition: InitStateEquil_impl.hpp:1086
bool waterActive() const
Predicate for whether or not water is an active phase.
Definition: InitStateEquil_impl.hpp:1068
bool gasActive() const
Predicate for whether or not gas is an active phase.
Definition: InitStateEquil_impl.hpp:1061
Scalar oil(const Scalar depth) const
Definition: InitStateEquil_impl.hpp:1076
std::array< Scalar, 2 > VSpan
Definition: InitStateEquil.hpp:170
bool oilActive() const
Predicate for whether or not oil is an active phase.
Definition: InitStateEquil_impl.hpp:1054
typename FluidSystem::Scalar Scalar
Definition: InitStateEquil.hpp:169
void equilibrate(const Region &reg, const VSpan &span)
Definition: InitStateEquil_impl.hpp:1043
PressureTable(const Scalar gravity, const int samplePoints=2000)
Definition: InitStateEquil_impl.hpp:987
Definition: InitStateEquil.hpp:67
Scalar operator()(const Scalar x) const
Definition: InitStateEquil_impl.hpp:216
RK4IVP(const RHS &f, const std::array< Scalar, 2 > &span, const Scalar y0, const int N)
Definition: InitStateEquil_impl.hpp:181
Definition: EquilibrationHelpers.hpp:617
std::pair< Scalar, Scalar > cellZMinMax(const Element &element)
Definition: InitStateEquil_impl.hpp:162
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:64
Dune::Communication< MPIComm > Communication
Definition: ParallelCommunication.hpp:30
Definition: blackoilboundaryratevector.hh:37
The Opm property system, traits with inheritance.
Simple set of per-phase (named by primary component) quantities.
Definition: InitStateEquil.hpp:338
void reset()
Definition: InitStateEquil.hpp:361
Scalar gas
Definition: InitStateEquil.hpp:340
Scalar water
Definition: InitStateEquil.hpp:341
PhaseQuantityValue & operator/=(const Scalar x)
Definition: InitStateEquil.hpp:352
Scalar oil
Definition: InitStateEquil.hpp:339
PhaseQuantityValue & axpy(const PhaseQuantityValue &rhs, const Scalar a)
Definition: InitStateEquil.hpp:343
Definition: InitStateEquil.hpp:392
CellID cell
Definition: InitStateEquil.hpp:393
Scalar depth
Definition: InitStateEquil.hpp:394
Scalar depth
Definition: InitStateEquil.hpp:256
Scalar pressure
Definition: InitStateEquil.hpp:257