29#ifndef OPM_INIT_STATE_EQUIL_HPP
30#define OPM_INIT_STATE_EQUIL_HPP
34#include <opm/material/common/Tabulated1DFunction.hpp>
35#include <opm/material/fluidstates/SimpleModularFluidState.hpp>
37#include <opm/simulators/utils/ParallelCommunication.hpp>
50class NumericalAquifers;
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;
67 CellCornerData(
const std::array<Scalar, 8>& x,
68 const std::array<Scalar, 8>& y,
69 const std::array<Scalar, 8>& z)
74template<
class Scalar>
class EquilReg;
75namespace Miscibility {
template<
class Scalar>
class RsFunction; }
78template <
class Scalar,
class RHS>
83 const std::array<Scalar,2>& span,
87 Scalar operator()(
const Scalar x)
const;
91 std::array<Scalar,2> span_;
92 std::vector<Scalar> y_;
93 std::vector<Scalar> f_;
95 Scalar stepsize()
const;
98namespace PhasePressODE {
99template <
class Flu
idSystem>
102 using Scalar =
typename FluidSystem::Scalar;
103 using TabulatedFunction = Tabulated1DFunction<Scalar>;
106 Water(
const TabulatedFunction& tempVdTable,
107 const TabulatedFunction& saltVdTable,
108 const int pvtRegionIdx,
109 const Scalar normGrav);
111 Scalar operator()(
const Scalar depth,
112 const Scalar press)
const;
115 const TabulatedFunction& tempVdTable_;
116 const TabulatedFunction& saltVdTable_;
117 const int pvtRegionIdx_;
120 Scalar density(
const Scalar depth,
121 const Scalar press)
const;
124template <
class Flu
idSystem,
class RS>
127 using Scalar =
typename FluidSystem::Scalar;
128 using TabulatedFunction = Tabulated1DFunction<Scalar>;
131 Oil(
const TabulatedFunction& tempVdTable,
133 const int pvtRegionIdx,
134 const Scalar normGrav);
136 Scalar operator()(
const Scalar depth,
137 const Scalar press)
const;
140 const TabulatedFunction& tempVdTable_;
142 const int pvtRegionIdx_;
145 Scalar density(
const Scalar depth,
146 const Scalar press)
const;
149template <
class Flu
idSystem,
class RV,
class RVW>
152 using Scalar =
typename FluidSystem::Scalar;
153 using TabulatedFunction = Tabulated1DFunction<Scalar>;
156 Gas(
const TabulatedFunction& tempVdTable,
159 const int pvtRegionIdx,
160 const Scalar normGrav);
162 Scalar operator()(
const Scalar depth,
163 const Scalar press)
const;
166 const TabulatedFunction& tempVdTable_;
169 const int pvtRegionIdx_;
172 Scalar density(
const Scalar depth,
173 const Scalar press)
const;
178template <
class Flu
idSystem,
class Region>
182 using Scalar =
typename FluidSystem::Scalar;
183 using VSpan = std::array<Scalar, 2>;
194 const int samplePoints = 2000);
224 void equilibrate(
const Region& reg,
243 Scalar
oil(
const Scalar depth)
const;
252 Scalar
gas(
const Scalar depth)
const;
261 Scalar
water(
const Scalar depth)
const;
265 class PressureFunction
273 explicit PressureFunction(
const ODE& ode,
278 PressureFunction(
const PressureFunction& rhs);
280 PressureFunction(PressureFunction&& rhs) =
default;
282 PressureFunction& operator=(
const PressureFunction& rhs);
284 PressureFunction& operator=(PressureFunction&& rhs);
286 Scalar value(
const Scalar depth)
const;
289 enum Direction : std::size_t { Up, Down, NumDir };
292 using DistrPtr = std::unique_ptr<Distribution>;
295 std::array<DistrPtr, Direction::NumDir> value_;
299 FluidSystem,
typename Region::CalcDissolution
303 FluidSystem,
typename Region::CalcEvaporation,
typename Region::CalcWaterEvaporation
308 using OPress = PressureFunction<OilPressODE>;
309 using GPress = PressureFunction<GasPressODE>;
310 using WPress = PressureFunction<WatPressODE>;
313 (
const Region&,
const VSpan&);
318 std::unique_ptr<OPress> oil_{};
319 std::unique_ptr<GPress> gas_{};
320 std::unique_ptr<WPress> wat_{};
322 template <
typename PressFunc>
323 void checkPtr(
const PressFunc* phasePress,
324 const std::string& phaseName)
const;
326 Strategy selectEquilibrationStrategy(
const Region& reg)
const;
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);
334 void makeOilPressure(
const typename OPress::InitCond& ic,
338 void makeGasPressure(
const typename GPress::InitCond& ic,
342 void makeWatPressure(
const typename WPress::InitCond& ic,
350template<
class Scalar>
358 this->oil += a * rhs.oil;
359 this->gas += a * rhs.gas;
360 this->water += a * rhs.water;
376 this->oil = this->gas = this->water = 0.0;
397template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
401 using Scalar =
typename FluidSystem::Scalar;
421 const std::vector<Scalar>& swatInit);
463 struct EvaluationPoint {
464 const Position* position{
nullptr};
465 const Region* region {
nullptr};
466 const PTable* ptable {
nullptr};
472 using FluidState = ::Opm::
473 SimpleModularFluidState<Scalar, 3, 3,
485 using MaterialLaw =
typename MaterialLawManager::MaterialLaw;
488 using PhaseIdx = std::remove_cv_t<
489 std::remove_reference_t<
decltype(FluidSystem::oilPhaseIdx)>
493 MaterialLawManager& matLawMgr_;
496 const std::vector<Scalar>& swatInit_;
499 PhaseQuantityValue<Scalar> sat_;
502 PhaseQuantityValue<Scalar> press_;
505 EvaluationPoint evalPt_;
508 FluidState fluidState_;
511 std::array<Scalar, FluidSystem::numPhases> matLawCapPress_;
522 void setEvaluationPoint(
const Position& x,
529 void initializePhaseQuantities();
548 void deriveWaterSat();
552 void fixUnphysicalTransition();
556 void accountForScaledSaturations();
572 std::pair<Scalar, bool> applySwatInit(
const Scalar pcow);
586 std::pair<Scalar, bool> applySwatInit(
const Scalar pc,
const Scalar sw);
590 void computeMaterialLawCapPress();
594 Scalar materialLawCapPressGasOil()
const;
598 Scalar materialLawCapPressOilWater()
const;
602 Scalar materialLawCapPressGasWater()
const;
611 bool isConstCapPress(
const PhaseIdx phaseIdx)
const;
618 bool isOverlappingTransition()
const;
639 Scalar fromDepthTable(
const Scalar contactdepth,
640 const PhaseIdx phasePos,
641 const bool isincr)
const;
660 Scalar invertCapPress(
const Scalar pc,
661 const PhaseIdx phasePos,
662 const bool isincr)
const;
665 PhaseIdx oilPos()
const
667 return FluidSystem::oilPhaseIdx;
671 PhaseIdx gasPos()
const
673 return FluidSystem::gasPhaseIdx;
677 PhaseIdx waterPos()
const
679 return FluidSystem::waterPhaseIdx;
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);
691template <
class Scalar,
class Element>
692std::pair<Scalar,Scalar> cellZMinMax(
const Element& element);
696namespace DeckDependent {
698template<
class FluidSystem,
702 class CartesianIndexMapper>
703class InitialStateComputer
705 using Element =
typename GridView::template Codim<0>::Entity;
706 using Scalar =
typename FluidSystem::Scalar;
708 template<
class MaterialLawManager>
709 InitialStateComputer(MaterialLawManager& materialLawManager,
710 const EclipseState& eclipseState,
712 const GridView& gridView,
713 const CartesianIndexMapper& cartMapper,
715 const int num_pressure_points = 2000,
716 const bool applySwatInit =
true);
718 using Vec = std::vector<Scalar>;
719 using PVec = std::vector<Vec>;
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_; }
731 template <
class RMap>
732 void updateInitialTemperature_(
const EclipseState& eclState,
const RMap& reg);
734 template <
class RMap>
735 void updateInitialSaltConcentration_(
const EclipseState& eclState,
const RMap& reg);
737 template <
class RMap>
738 void updateInitialSaltSaturation_(
const EclipseState& eclState,
const RMap& reg);
740 void updateCellProps_(
const GridView& gridView,
741 const NumericalAquifers& aquifer);
743 void applyNumericalAquifers_(
const GridView& gridView,
744 const NumericalAquifers& aquifer,
745 const bool co2store_or_h2store);
748 void setRegionPvtIdx(
const EclipseState& eclState,
const RMap& reg);
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,
758 template <
class CellRange,
class EquilibrationMethod>
759 void cellLoop(
const CellRange& cells,
760 EquilibrationMethod&& eqmethod);
762 template <
class CellRange,
class PressTable,
class PhaseSat>
763 void equilibrateCellCentres(
const CellRange& cells,
765 const PressTable& ptable,
768 template <
class CellRange,
class PressTable,
class PhaseSat>
769 void equilibrateHorizontal(
const CellRange& cells,
772 const PressTable& ptable,
775 template<
class CellRange,
class PressTable,
class PhaseSat>
776 void equilibrateTiltedFaultBlock(
const CellRange& cells,
778 const GridView& gridView,
const int numLevels,
779 const PressTable& ptable, PhaseSat& psat);
781 template<
class CellRange,
class PressTable,
class PhaseSat>
782 void equilibrateTiltedFaultBlockSimple(
const CellRange& cells,
784 const GridView& gridView,
const int numLevels,
785 const PressTable& ptable, PhaseSat& psat);
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_;
796 Vec saltConcentration_;
803 const CartesianIndexMapper& cartesianIndexMapper_;
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_;
Definition InitStateEquil.hpp:151
Definition InitStateEquil.hpp:126
Definition InitStateEquil.hpp:101
const PhaseQuantityValue< Scalar > & deriveSaturations(const Position &x, const Region ®, const PTable &ptable)
Calculate phase saturations at particular point of the simulation model geometry.
Definition InitStateEquil_impl.hpp:734
PhaseSaturations(MaterialLawManager &matLawMgr, const std::vector< Scalar > &swatInit)
Constructor.
Definition InitStateEquil_impl.hpp:710
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
Retrieve saturation-corrected phase pressures.
Definition InitStateEquil.hpp:455
Definition InitStateEquil.hpp:180
PressureTable & operator=(const PressureTable &rhs)
Assignment operator.
Definition InitStateEquil_impl.hpp:1155
Scalar water(const Scalar depth) const
Evaluate water phase pressure at specified depth.
Definition InitStateEquil_impl.hpp:1235
Scalar gas(const Scalar depth) const
Evaluate gas phase pressure at specified depth.
Definition InitStateEquil_impl.hpp:1224
bool waterActive() const
Predicate for whether or not water is an active phase.
Definition InitStateEquil_impl.hpp:1206
bool gasActive() const
Predicate for whether or not gas is an active phase.
Definition InitStateEquil_impl.hpp:1199
Scalar oil(const Scalar depth) const
Evaluate oil phase pressure at specified depth.
Definition InitStateEquil_impl.hpp:1214
bool oilActive() const
Predicate for whether or not oil is an active phase.
Definition InitStateEquil_impl.hpp:1192
PressureTable(const Scalar gravity, const int samplePoints=2000)
Constructor.
Definition InitStateEquil_impl.hpp:1125
Definition InitStateEquil.hpp:80
Aggregate information base of an equilibration region.
Definition EquilibrationHelpers.hpp:676
Types and routines relating to phase mixing in equilibration calculations.
Definition EquilibrationHelpers.hpp:94
Types and routines that collectively implement a basic ECLIPSE-style equilibration-based initialisati...
Definition EquilibrationHelpers.cpp:28
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:45
The Opm property system, traits with inheritance.
Simple set of per-phase (named by primary component) quantities.
Definition InitStateEquil.hpp:351
Evaluation point within a model geometry.
Definition InitStateEquil.hpp:405
Definition InitStateEquil.hpp:268