27#ifndef OPM_WATER_PVT_MULTIPLEXER_HPP
28#define OPM_WATER_PVT_MULTIPLEXER_HPP
36#define OPM_WATER_PVT_MULTIPLEXER_CALL(codeToCall, ...) \
37 switch (approach_) { \
38 case WaterPvtApproach::ConstantCompressibilityWater: { \
39 auto& pvtImpl = getRealPvt<WaterPvtApproach::ConstantCompressibilityWater>(); \
43 case WaterPvtApproach::ConstantCompressibilityBrine: { \
44 auto& pvtImpl = getRealPvt<WaterPvtApproach::ConstantCompressibilityBrine>(); \
48 case WaterPvtApproach::ThermalWater: { \
49 auto& pvtImpl = getRealPvt<WaterPvtApproach::ThermalWater>(); \
53 case WaterPvtApproach::BrineCo2: { \
54 auto& pvtImpl = getRealPvt<WaterPvtApproach::BrineCo2>(); \
58 case WaterPvtApproach::BrineH2: { \
59 auto& pvtImpl = getRealPvt<WaterPvtApproach::BrineH2>(); \
64 case WaterPvtApproach::NoWater: \
65 throw std::logic_error("Not implemented: Water PVT of this deck!"); \
70enum class WaterPvtApproach {
72 ConstantCompressibilityBrine,
73 ConstantCompressibilityWater,
88template <
class Scalar,
bool enableThermal = true,
bool enableBrine = true>
89class WaterPvtMultiplexer
93 : approach_(WaterPvtApproach::NoWater)
94 , realWaterPvt_(
nullptr)
98 WaterPvtMultiplexer(WaterPvtApproach
approach,
void* realWaterPvt)
100 , realWaterPvt_(realWaterPvt)
103 WaterPvtMultiplexer(
const WaterPvtMultiplexer<Scalar,enableThermal,enableBrine>& data)
108 ~WaterPvtMultiplexer();
110 bool mixingEnergy()
const
112 return approach_ == WaterPvtApproach::ThermalWater;
131 void setVapPars(
const Scalar par1,
const Scalar par2);
141 template <
class Evaluation>
143 const Evaluation& temperature,
144 const Evaluation& pressure,
145 const Evaluation& Rsw,
146 const Evaluation& saltconcentration)
const
147 { OPM_WATER_PVT_MULTIPLEXER_CALL(
return pvtImpl.internalEnergy(regionIdx, temperature, pressure, Rsw, saltconcentration)); }
149 Scalar hVap(
unsigned regionIdx)
const;
154 template <
class Evaluation>
156 const Evaluation& temperature,
157 const Evaluation& pressure,
158 const Evaluation& Rsw,
159 const Evaluation& saltconcentration)
const
161 OPM_WATER_PVT_MULTIPLEXER_CALL(
return pvtImpl.viscosity(regionIdx, temperature, pressure, Rsw, saltconcentration));
164 bool isActive()
const
166 return approach_ != WaterPvtApproach::NoWater;
172 template <
class Evaluation>
174 const Evaluation& temperature,
175 const Evaluation& pressure,
176 const Evaluation& saltconcentration)
const
178 OPM_WATER_PVT_MULTIPLEXER_CALL(
return pvtImpl.saturatedViscosity(regionIdx, temperature, pressure, saltconcentration));
184 template <
class Evaluation>
186 const Evaluation& temperature,
187 const Evaluation& pressure,
188 const Evaluation& Rsw,
189 const Evaluation& saltconcentration)
const
191 OPM_WATER_PVT_MULTIPLEXER_CALL(
return pvtImpl.inverseFormationVolumeFactor(regionIdx, temperature, pressure, Rsw, saltconcentration));
197 template <
class Flu
idState,
class LhsEval =
typename Flu
idState::Scalar>
198 std::pair<LhsEval, LhsEval>
200 { OPM_WATER_PVT_MULTIPLEXER_CALL(
return pvtImpl.inverseFormationVolumeFactorAndViscosity(fluidState, regionIdx)); }
205 template <
class Evaluation>
207 const Evaluation& temperature,
208 const Evaluation& pressure,
209 const Evaluation& saltconcentration)
const
211 OPM_WATER_PVT_MULTIPLEXER_CALL(
return pvtImpl.saturatedInverseFormationVolumeFactor(regionIdx, temperature, pressure, saltconcentration));
217 template <
class Evaluation>
219 const Evaluation& temperature,
220 const Evaluation& pressure,
221 const Evaluation& saltconcentration)
const
223 OPM_WATER_PVT_MULTIPLEXER_CALL(
return pvtImpl.saturatedGasDissolutionFactor(regionIdx, temperature, pressure, saltconcentration));
233 template <
class Evaluation>
235 const Evaluation& temperature,
236 const Evaluation& Rs,
237 const Evaluation& saltconcentration)
const
238 { OPM_WATER_PVT_MULTIPLEXER_CALL(
return pvtImpl.saturationPressure(regionIdx, temperature, Rs, saltconcentration)); }
243 template <
class Evaluation>
245 const Evaluation& pressure,
246 unsigned compIdx)
const
248 OPM_WATER_PVT_MULTIPLEXER_CALL(
return pvtImpl.diffusionCoefficient(temperature, pressure, compIdx));
251 void setApproach(WaterPvtApproach appr);
259 {
return approach_; }
262 template <WaterPvtApproach approachV>
263 typename std::enable_if<approachV == WaterPvtApproach::ConstantCompressibilityWater, ConstantCompressibilityWaterPvt<Scalar> >::type& getRealPvt()
269 template <WaterPvtApproach approachV>
270 typename std::enable_if<approachV == WaterPvtApproach::ConstantCompressibilityWater, const ConstantCompressibilityWaterPvt<Scalar> >::type& getRealPvt()
const
273 return *
static_cast<ConstantCompressibilityWaterPvt<Scalar>*
>(realWaterPvt_);
276 template <WaterPvtApproach approachV>
277 typename std::enable_if<approachV == WaterPvtApproach::ConstantCompressibilityBrine, ConstantCompressibilityBrinePvt<Scalar> >::type& getRealPvt()
280 return *
static_cast<ConstantCompressibilityBrinePvt<Scalar>*
>(realWaterPvt_);
283 template <WaterPvtApproach approachV>
284 typename std::enable_if<approachV == WaterPvtApproach::ConstantCompressibilityBrine, const ConstantCompressibilityBrinePvt<Scalar> >::type& getRealPvt()
const
287 return *
static_cast<ConstantCompressibilityBrinePvt<Scalar>*
>(realWaterPvt_);
290 template <WaterPvtApproach approachV>
291 typename std::enable_if<approachV == WaterPvtApproach::ThermalWater, WaterPvtThermal<Scalar, enableBrine> >::type& getRealPvt()
294 return *
static_cast<WaterPvtThermal<Scalar, enableBrine>*
>(realWaterPvt_);
297 template <WaterPvtApproach approachV>
298 typename std::enable_if<approachV == WaterPvtApproach::ThermalWater, const WaterPvtThermal<Scalar, enableBrine> >::type& getRealPvt()
const
301 return *
static_cast<WaterPvtThermal<Scalar, enableBrine>*
>(realWaterPvt_);
304 template <WaterPvtApproach approachV>
305 typename std::enable_if<approachV == WaterPvtApproach::BrineCo2, BrineCo2Pvt<Scalar> >::type& getRealPvt()
308 return *
static_cast<BrineCo2Pvt<Scalar>*
>(realWaterPvt_);
311 template <WaterPvtApproach approachV>
312 typename std::enable_if<approachV == WaterPvtApproach::BrineCo2, const BrineCo2Pvt<Scalar> >::type& getRealPvt()
const
315 return *
static_cast<const BrineCo2Pvt<Scalar>*
>(realWaterPvt_);
318 template <WaterPvtApproach approachV>
319 typename std::enable_if<approachV == WaterPvtApproach::BrineH2, BrineH2Pvt<Scalar> >::type& getRealPvt()
322 return *
static_cast<BrineH2Pvt<Scalar>*
>(realWaterPvt_);
325 template <WaterPvtApproach approachV>
326 typename std::enable_if<approachV == WaterPvtApproach::BrineH2, const BrineH2Pvt<Scalar> >::type& getRealPvt()
const
329 return *
static_cast<const BrineH2Pvt<Scalar>*
>(realWaterPvt_);
332 const void* realWaterPvt()
const {
return realWaterPvt_; }
334 WaterPvtMultiplexer<Scalar,enableThermal,enableBrine>&
335 operator=(
const WaterPvtMultiplexer<Scalar,enableThermal,enableBrine>& data);
338 WaterPvtApproach approach_{WaterPvtApproach::NoWater};
339 void* realWaterPvt_{
nullptr};
This class represents the Pressure-Volume-Temperature relations of the liquid phase for a CO2-Brine s...
This class represents the Pressure-Volume-Temperature relations of the liquid phase for a H2-Brine sy...
This class represents the Pressure-Volume-Temperature relations of the gas phase without vaporized oi...
This class represents the Pressure-Volume-Temperature relations of the gas phase without vaporized oi...
This class implements temperature dependence of the PVT properties of water.
This class represents the Pressure-Volume-Temperature relations of the gas phase without vaporized oi...
Definition ConstantCompressibilityWaterPvt.hpp:49
Definition EclipseState.hpp:62
Definition Schedule.hpp:101
Evaluation viscosity(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure, const Evaluation &Rsw, const Evaluation &saltconcentration) const
Returns the dynamic viscosity [Pa s] of the fluid phase given a set of parameters.
Definition WaterPvtMultiplexer.hpp:155
Evaluation saturatedInverseFormationVolumeFactor(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure, const Evaluation &saltconcentration) const
Returns the formation volume factor [-] of the fluid phase.
Definition WaterPvtMultiplexer.hpp:206
Evaluation saturatedViscosity(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure, const Evaluation &saltconcentration) const
Returns the dynamic viscosity [Pa s] of the fluid phase given a set of parameters.
Definition WaterPvtMultiplexer.hpp:173
unsigned numRegions() const
Return the number of PVT regions which are considered by this PVT-object.
Definition WaterPvtMultiplexer.cpp:97
std::pair< LhsEval, LhsEval > inverseFormationVolumeFactorAndViscosity(const FluidState &fluidState, unsigned regionIdx)
Returns the formation volume factor [-] and viscosity [Pa s] of the fluid phase.
Definition WaterPvtMultiplexer.hpp:199
Evaluation saturatedGasDissolutionFactor(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure, const Evaluation &saltconcentration) const
Returns the gas dissolution factor [m^3/m^3] of saturated water.
Definition WaterPvtMultiplexer.hpp:218
Evaluation internalEnergy(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure, const Evaluation &Rsw, const Evaluation &saltconcentration) const
Returns the specific enthalpy [J/kg] of gas given a set of parameters.
Definition WaterPvtMultiplexer.hpp:142
Evaluation saturationPressure(unsigned regionIdx, const Evaluation &temperature, const Evaluation &Rs, const Evaluation &saltconcentration) const
Returns the saturation pressure [Pa] of water given the mass fraction of the gas component in the wat...
Definition WaterPvtMultiplexer.hpp:234
WaterPvtApproach approach() const
Returns the concrete approach for calculating the PVT relations.
Definition WaterPvtMultiplexer.hpp:258
Scalar waterReferenceDensity(unsigned regionIdx) const
Return the reference density which are considered by this PVT-object.
Definition WaterPvtMultiplexer.cpp:111
Evaluation inverseFormationVolumeFactor(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure, const Evaluation &Rsw, const Evaluation &saltconcentration) const
Returns the formation volume factor [-] of the fluid phase.
Definition WaterPvtMultiplexer.hpp:185
Evaluation diffusionCoefficient(const Evaluation &temperature, const Evaluation &pressure, unsigned compIdx) const
Calculate the binary molecular diffusion coefficient for a component in a fluid phase [mol^2 * s / (k...
Definition WaterPvtMultiplexer.hpp:244
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30