29#ifndef OPM_EQUILIBRATION_HELPERS_HPP
30#define OPM_EQUILIBRATION_HELPERS_HPP
32#include <opm/material/common/Tabulated1DFunction.hpp>
34#include <opm/input/eclipse/EclipseState/InitConfig/Equil.hpp>
126 const Scalar sat = 0.0)
const = 0;
133template<
class Scalar>
148 const Scalar = 0.0)
const override
160template <
class Flu
idSystem>
164 using Scalar =
typename FluidSystem::Scalar;
172 RsVD(
const int pvtRegionIdx,
173 const std::vector<Scalar>& depth,
174 const std::vector<Scalar>& rs);
197 const Scalar satGas = 0.0)
const override;
200 using RsVsDepthFunc = Tabulated1DFunction<Scalar>;
202 const int pvtRegionIdx_;
203 RsVsDepthFunc rsVsDepth_;
205 Scalar satRs(
const Scalar press,
const Scalar temp)
const;
214template <
class Flu
idSystem>
218 using Scalar =
typename FluidSystem::Scalar;
226 PBVD(
const int pvtRegionIdx,
227 const std::vector<Scalar>& depth,
228 const std::vector<Scalar>& pbub);
248 const Scalar cellPress,
250 const Scalar satGas = 0.0)
const override;
253 using PbubVsDepthFunc = Tabulated1DFunction<Scalar>;
255 const int pvtRegionIdx_;
256 PbubVsDepthFunc pbubVsDepth_;
258 Scalar satRs(
const Scalar press,
const Scalar temp)
const;
267template <
class Flu
idSystem>
270 using Scalar =
typename FluidSystem::Scalar;
279 PDVD(
const int pvtRegionIdx,
280 const std::vector<Scalar>& depth,
281 const std::vector<Scalar>& pdew);
301 const Scalar cellPress,
303 const Scalar satOil = 0.0)
const override;
306 using PdewVsDepthFunc = Tabulated1DFunction<Scalar>;
308 const int pvtRegionIdx_;
309 PdewVsDepthFunc pdewVsDepth_;
311 Scalar satRv(
const Scalar press,
const Scalar temp)
const;
320template <
class Flu
idSystem>
323 using Scalar =
typename FluidSystem::Scalar;
332 RvVD(
const int pvtRegionIdx,
333 const std::vector<Scalar>& depth,
334 const std::vector<Scalar>& rv);
357 const Scalar satOil = 0.0)
const override;
360 using RvVsDepthFunc = Tabulated1DFunction<Scalar>;
362 const int pvtRegionIdx_;
363 RvVsDepthFunc rvVsDepth_;
365 Scalar satRv(
const Scalar press,
const Scalar temp)
const;
374template <
class Flu
idSystem>
377 using Scalar =
typename FluidSystem::Scalar;
386 RvwVD(
const int pvtRegionIdx,
387 const std::vector<Scalar>& depth,
388 const std::vector<Scalar>& rvw);
411 const Scalar satWat = 0.0)
const override;
414 using RvwVsDepthFunc = Tabulated1DFunction<Scalar>;
416 const int pvtRegionIdx_;
417 RvwVsDepthFunc rvwVsDepth_;
419 Scalar satRvw(
const Scalar press,
const Scalar temp)
const;
437template <
class Flu
idSystem>
440 using Scalar =
typename FluidSystem::Scalar;
450 const Scalar pContact,
451 const Scalar T_contact);
474 const Scalar satGas = 0.0)
const override;
477 const int pvtRegionIdx_;
478 Scalar rsSatContact_;
480 Scalar satRs(
const Scalar press,
const Scalar temp)
const;
498template <
class Flu
idSystem>
501 using Scalar =
typename FluidSystem::Scalar;
511 const Scalar pContact,
512 const Scalar T_contact);
535 const Scalar satOil = 0.0)
const override;
538 const int pvtRegionIdx_;
539 Scalar rvSatContact_;
541 Scalar satRv(
const Scalar press,
const Scalar temp)
const;
558template <
class Flu
idSystem>
561 using Scalar =
typename FluidSystem::Scalar;
571 const Scalar pContact,
572 const Scalar T_contact);
595 const Scalar satWat = 0.0)
const override;
598 const int pvtRegionIdx_;
599 Scalar rvwSatContact_;
601 Scalar satRvw(
const Scalar press,
const Scalar temp)
const;
610template <
class Flu
idSystem>
614 using Scalar =
typename FluidSystem::Scalar;
622 RsConst(
const Scalar rs_constant,
623 const Scalar pb_constant);
632 const Scalar satGas = 0.0)
const override;
644 bool isConstantRs()
const {
return true; }
646 Scalar satRs(
const Scalar press,
const Scalar temp)
const;
649 Scalar rs_constant_{};
650 Scalar pb_constant_{};
674template<
class Scalar>
677 using TabulatedFunction = Tabulated1DFunction<Scalar>;
695 const TabulatedFunction& tempVdTable,
696 const TabulatedFunction& saltVdTable,
718 Scalar
datum()
const;
776 const TabulatedFunction& saltVdTable()
const;
777 const TabulatedFunction& tempVdTable()
const;
786 std::shared_ptr<Miscibility::RsFunction<Scalar>> rs_;
787 std::shared_ptr<Miscibility::RsFunction<Scalar>> rv_;
788 std::shared_ptr<Miscibility::RsFunction<Scalar>> rvw_;
789 const TabulatedFunction& tempVdTable_;
790 const TabulatedFunction& saltVdTable_;
799template <
class Flu
idSystem,
class MaterialLawManager>
802 using Scalar =
typename FluidSystem::Scalar;
803 PcEq(
const MaterialLawManager& materialLawManager,
806 const Scalar targetPc);
808 Scalar operator()(Scalar s)
const;
811 const MaterialLawManager& materialLawManager_;
814 const Scalar targetPc_;
817template <
class Flu
idSystem,
class MaterialLawManager>
818typename FluidSystem::Scalar
819minSaturations(
const MaterialLawManager& materialLawManager,
820 const int phase,
const int cell);
822template <
class Flu
idSystem,
class MaterialLawManager>
823typename FluidSystem::Scalar
824maxSaturations(
const MaterialLawManager& materialLawManager,
825 const int phase,
const int cell);
829template <
class Flu
idSystem,
class MaterialLawManager>
830typename FluidSystem::Scalar
831satFromPc(
const MaterialLawManager& materialLawManager,
834 const typename FluidSystem::Scalar targetPc,
835 const bool increasing =
false);
840template <
class Flu
idSystem,
class MaterialLawManager>
843 using Scalar =
typename FluidSystem::Scalar;
844 PcEqSum(
const MaterialLawManager& materialLawManager,
848 const Scalar targetPc);
850 Scalar operator()(Scalar s)
const;
853 const MaterialLawManager& materialLawManager_;
857 const Scalar targetPc_;
863template <
class Flu
idSystem,
class MaterialLawManager>
864typename FluidSystem::Scalar
869 const typename FluidSystem::Scalar targetPc);
872template <
class Flu
idSystem,
class MaterialLawManager>
873typename FluidSystem::Scalar
874satFromDepth(
const MaterialLawManager& materialLawManager,
875 const typename FluidSystem::Scalar cellDepth,
876 const typename FluidSystem::Scalar contactDepth,
879 const bool increasing =
false);
882template <
class Flu
idSystem,
class MaterialLawManager>
883bool isConstPc(
const MaterialLawManager& materialLawManager,
Scalar datum() const
Datum depth in current region.
Definition EquilibrationHelpers_impl.hpp:426
int pvtIdx() const
Retrieve pvtIdx of the region.
Definition EquilibrationHelpers_impl.hpp:503
Miscibility::RsFunction< Scalar > CalcEvaporation
Type of vapourised oil-gas ratio calculator.
Definition EquilibrationHelpers.hpp:707
EquilReg(const EquilRecord &rec, std::shared_ptr< Miscibility::RsFunction< Scalar > > rs, std::shared_ptr< Miscibility::RsFunction< Scalar > > rv, std::shared_ptr< Miscibility::RsFunction< Scalar > > rvw, const TabulatedFunction &tempVdTable, const TabulatedFunction &saltVdTable, const int pvtIdx)
Constructor.
Definition EquilibrationHelpers_impl.hpp:408
Scalar pcgoGoc() const
Gas-oil capillary pressure at gas-oil contact.
Definition EquilibrationHelpers_impl.hpp:456
int equilibrationAccuracy() const
Accuracy/strategy for initial fluid-in-place calculation.
Definition EquilibrationHelpers_impl.hpp:462
const CalcEvaporation & evaporationCalculator() const
Retrieve vapourised oil-gas ratio calculator of current region.
Definition EquilibrationHelpers_impl.hpp:476
Miscibility::RsFunction< Scalar > CalcWaterEvaporation
Type of vapourised water-gas ratio calculator.
Definition EquilibrationHelpers.hpp:712
Scalar zgoc() const
Depth of gas-oil contact.
Definition EquilibrationHelpers_impl.hpp:450
Miscibility::RsFunction< Scalar > CalcDissolution
Type of dissolved gas-oil ratio calculator.
Definition EquilibrationHelpers.hpp:702
Scalar pressure() const
Pressure at datum depth in current region.
Definition EquilibrationHelpers_impl.hpp:432
Scalar pcowWoc() const
water-oil capillary pressure at water-oil contact.
Definition EquilibrationHelpers_impl.hpp:444
Scalar zwoc() const
Depth of water-oil contact.
Definition EquilibrationHelpers_impl.hpp:438
const CalcDissolution & dissolutionCalculator() const
Retrieve dissolved gas-oil ratio calculator of current region.
Definition EquilibrationHelpers_impl.hpp:469
const CalcWaterEvaporation & waterEvaporationCalculator() const
Retrieve vapourised water-gas ratio calculator of current region.
Definition EquilibrationHelpers_impl.hpp:483
Type that implements "no phase mixing" policy.
Definition EquilibrationHelpers.hpp:135
Scalar operator()(const Scalar, const Scalar, const Scalar, const Scalar=0.0) const override
Function call.
Definition EquilibrationHelpers.hpp:145
Scalar operator()(const Scalar depth, const Scalar cellPress, const Scalar temp, const Scalar satGas=0.0) const override
Function call.
Definition EquilibrationHelpers_impl.hpp:113
PBVD(const int pvtRegionIdx, const std::vector< Scalar > &depth, const std::vector< Scalar > &pbub)
Constructor.
Definition EquilibrationHelpers_impl.hpp:102
Scalar operator()(const Scalar depth, const Scalar cellPress, const Scalar temp, const Scalar satOil=0.0) const override
Function call.
Definition EquilibrationHelpers_impl.hpp:150
PDVD(const int pvtRegionIdx, const std::vector< Scalar > &depth, const std::vector< Scalar > &pdew)
Constructor.
Definition EquilibrationHelpers_impl.hpp:139
RsConst(const Scalar rs_constant, const Scalar pb_constant)
Constructor.
Definition EquilibrationHelpers_impl.hpp:359
Scalar operator()(const Scalar depth, const Scalar press, const Scalar temp, const Scalar satGas=0.0) const override
Function call - returns constant Rs if pressure above bubble point.
Definition EquilibrationHelpers_impl.hpp:385
Scalar constantRs() const
Get the constant Rs value.
Definition EquilibrationHelpers.hpp:642
Scalar bubblePoint() const
Get the constant bubble point pressure.
Definition EquilibrationHelpers.hpp:637
Base class for phase mixing functions.
Definition EquilibrationHelpers.hpp:101
virtual Scalar operator()(const Scalar depth, const Scalar press, const Scalar temp, const Scalar sat=0.0) const =0
Function call operator.
Scalar operator()(const Scalar depth, const Scalar press, const Scalar temp, const Scalar satGas=0.0) const override
Function call.
Definition EquilibrationHelpers_impl.hpp:75
RsVD(const int pvtRegionIdx, const std::vector< Scalar > &depth, const std::vector< Scalar > &rs)
Constructor.
Definition EquilibrationHelpers_impl.hpp:64
Scalar operator()(const Scalar depth, const Scalar press, const Scalar temp, const Scalar satOil=0.0) const override
Function call.
Definition EquilibrationHelpers_impl.hpp:187
RvVD(const int pvtRegionIdx, const std::vector< Scalar > &depth, const std::vector< Scalar > &rv)
Constructor.
Definition EquilibrationHelpers_impl.hpp:176
Scalar operator()(const Scalar depth, const Scalar press, const Scalar temp, const Scalar satWat=0.0) const override
Function call.
Definition EquilibrationHelpers_impl.hpp:230
RvwVD(const int pvtRegionIdx, const std::vector< Scalar > &depth, const std::vector< Scalar > &rvw)
Constructor.
Definition EquilibrationHelpers_impl.hpp:219
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
FluidSystem::Scalar satFromDepth(const MaterialLawManager &materialLawManager, const typename FluidSystem::Scalar cellDepth, const typename FluidSystem::Scalar contactDepth, const int phase, const int cell, const bool increasing=false)
Compute saturation from depth. Used for constant capillary pressure function.
Definition EquilibrationHelpers_impl.hpp:703
FluidSystem::Scalar satFromPc(const MaterialLawManager &materialLawManager, const int phase, const int cell, const typename FluidSystem::Scalar targetPc, const bool increasing=false)
Compute saturation of some phase corresponding to a given capillary pressure.
Definition EquilibrationHelpers_impl.hpp:639
bool isConstPc(const MaterialLawManager &materialLawManager, const int phase, const int cell)
Return true if capillary pressure function is constant.
Definition EquilibrationHelpers_impl.hpp:723
FluidSystem::Scalar satFromSumOfPcs(const MaterialLawManager &materialLawManager, const int phase1, const int phase2, const int cell, const typename FluidSystem::Scalar targetPc)
Compute saturation of some phase corresponding to a given capillary pressure, where the capillary pre...
Definition EquilibrationHelpers_impl.hpp:672
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:45