opm-common
Loading...
Searching...
No Matches
GasPvtMultiplexer.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 2 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*/
27#ifndef OPM_GAS_PVT_MULTIPLEXER_HPP
28#define OPM_GAS_PVT_MULTIPLEXER_HPP
29
37
38#include <functional>
39namespace Opm {
40
41#if HAVE_ECL_INPUT
42class EclipseState;
43class Schedule;
44#endif
45
46#define OPM_GAS_PVT_MULTIPLEXER_CALL(codeToCall, ...) \
47 switch (gasPvtApproach_) { \
48 case GasPvtApproach::DryGas: { \
49 auto& pvtImpl = getRealPvt<GasPvtApproach::DryGas>(); \
50 codeToCall; \
51 __VA_ARGS__; \
52 } \
53 case GasPvtApproach::DryHumidGas: { \
54 auto& pvtImpl = getRealPvt<GasPvtApproach::DryHumidGas>(); \
55 codeToCall; \
56 __VA_ARGS__; \
57 } \
58 case GasPvtApproach::WetHumidGas: { \
59 auto& pvtImpl = getRealPvt<GasPvtApproach::WetHumidGas>(); \
60 codeToCall; \
61 __VA_ARGS__; \
62 } \
63 case GasPvtApproach::WetGas: { \
64 auto& pvtImpl = getRealPvt<GasPvtApproach::WetGas>(); \
65 codeToCall; \
66 __VA_ARGS__; \
67 } \
68 case GasPvtApproach::ThermalGas: { \
69 auto& pvtImpl = getRealPvt<GasPvtApproach::ThermalGas>(); \
70 codeToCall; \
71 __VA_ARGS__; \
72 } \
73 case GasPvtApproach::Co2Gas: { \
74 auto& pvtImpl = getRealPvt<GasPvtApproach::Co2Gas>(); \
75 codeToCall; \
76 __VA_ARGS__; \
77 } \
78 case GasPvtApproach::H2Gas: { \
79 auto& pvtImpl = getRealPvt<GasPvtApproach::H2Gas>(); \
80 codeToCall; \
81 __VA_ARGS__; \
82 } \
83 default: \
84 case GasPvtApproach::NoGas: \
85 throw std::logic_error("Not implemented: Gas PVT of this deck!"); \
86 }
87
88enum class GasPvtApproach {
89 NoGas,
90 DryGas,
91 DryHumidGas,
92 WetHumidGas,
93 WetGas,
94 ThermalGas,
95 Co2Gas,
96 H2Gas
97};
98
109template <class Scalar, bool enableThermal = true>
110class GasPvtMultiplexer
111{
112public:
113 GasPvtMultiplexer()
114 : gasPvtApproach_(GasPvtApproach::NoGas)
115 , realGasPvt_(nullptr, [](void*){})
116 {
117 }
118
119 GasPvtMultiplexer(GasPvtApproach approach, void* realGasPvt)
120 : gasPvtApproach_(approach)
121 , realGasPvt_(realGasPvt, [this](void* ptr){ deleter(ptr); })
122 { }
123
124 GasPvtMultiplexer(const GasPvtMultiplexer<Scalar,enableThermal>& data)
125 {
126 *this = data;
127 }
128
129 ~GasPvtMultiplexer() = default;
130
131 bool mixingEnergy() const
132 {
133 return gasPvtApproach_ == GasPvtApproach::ThermalGas;
134 }
135
136 bool isActive() const {
137 return gasPvtApproach_ != GasPvtApproach::NoGas;
138 }
139
140#if HAVE_ECL_INPUT
146 void initFromState(const EclipseState& eclState, const Schedule& schedule);
147#endif // HAVE_ECL_INPUT
148
149 void setApproach(GasPvtApproach gasPvtAppr);
150
151 void initEnd();
152
156 unsigned numRegions() const;
157
158 void setVapPars(const Scalar par1, const Scalar par2);
159
163 Scalar gasReferenceDensity(unsigned regionIdx) const;
164
168 template <class Evaluation>
169 Evaluation internalEnergy(unsigned regionIdx,
170 const Evaluation& temperature,
171 const Evaluation& pressure,
172 const Evaluation& Rv,
173 const Evaluation& Rvw) const
174 { OPM_GAS_PVT_MULTIPLEXER_CALL(return pvtImpl.internalEnergy(regionIdx, temperature, pressure, Rv, Rvw)); }
175
176 Scalar hVap(unsigned regionIdx) const;
177
181 template <class Evaluation = Scalar>
182 Evaluation viscosity(unsigned regionIdx,
183 const Evaluation& temperature,
184 const Evaluation& pressure,
185 const Evaluation& Rv,
186 const Evaluation& Rvw ) const
187 { OPM_GAS_PVT_MULTIPLEXER_CALL(return pvtImpl.viscosity(regionIdx, temperature, pressure, Rv, Rvw)); }
188
192 template <class Evaluation = Scalar>
193 Evaluation saturatedViscosity(unsigned regionIdx,
194 const Evaluation& temperature,
195 const Evaluation& pressure) const
196 { OPM_GAS_PVT_MULTIPLEXER_CALL(return pvtImpl.saturatedViscosity(regionIdx, temperature, pressure)); }
197
201 template <class Evaluation = Scalar>
202 Evaluation inverseFormationVolumeFactor(unsigned regionIdx,
203 const Evaluation& temperature,
204 const Evaluation& pressure,
205 const Evaluation& Rv,
206 const Evaluation& Rvw) const
207 { OPM_GAS_PVT_MULTIPLEXER_CALL(return pvtImpl.inverseFormationVolumeFactor(regionIdx, temperature, pressure, Rv, Rvw)); }
208
212 template <class FluidState, class LhsEval = typename FluidState::Scalar>
213 std::pair<LhsEval, LhsEval>
214 inverseFormationVolumeFactorAndViscosity(const FluidState& fluidState, unsigned regionIdx)
215 { OPM_GAS_PVT_MULTIPLEXER_CALL(return pvtImpl.inverseFormationVolumeFactorAndViscosity(fluidState, regionIdx)); }
216
220 template <class Evaluation = Scalar>
221 Evaluation saturatedInverseFormationVolumeFactor(unsigned regionIdx,
222 const Evaluation& temperature,
223 const Evaluation& pressure) const
224 { OPM_GAS_PVT_MULTIPLEXER_CALL(return pvtImpl.saturatedInverseFormationVolumeFactor(regionIdx, temperature, pressure)); }
225
229 template <class Evaluation = Scalar>
230 Evaluation saturatedOilVaporizationFactor(unsigned regionIdx,
231 const Evaluation& temperature,
232 const Evaluation& pressure) const
233 { OPM_GAS_PVT_MULTIPLEXER_CALL(return pvtImpl.saturatedOilVaporizationFactor(regionIdx, temperature, pressure)); }
234
238 template <class Evaluation = Scalar>
239 Evaluation saturatedOilVaporizationFactor(unsigned regionIdx,
240 const Evaluation& temperature,
241 const Evaluation& pressure,
242 const Evaluation& oilSaturation,
243 const Evaluation& maxOilSaturation) const
244 { OPM_GAS_PVT_MULTIPLEXER_CALL(return pvtImpl.saturatedOilVaporizationFactor(regionIdx, temperature, pressure, oilSaturation, maxOilSaturation)); }
245
249 template <class Evaluation = Scalar>
250 Evaluation saturatedWaterVaporizationFactor(unsigned regionIdx,
251 const Evaluation& temperature,
252 const Evaluation& pressure) const
253 { OPM_GAS_PVT_MULTIPLEXER_CALL(return pvtImpl.saturatedWaterVaporizationFactor(regionIdx, temperature, pressure)); }
254
258 template <class Evaluation = Scalar>
259 Evaluation saturatedWaterVaporizationFactor(unsigned regionIdx,
260 const Evaluation& temperature,
261 const Evaluation& pressure,
262 const Evaluation& saltConcentration) const
263 { OPM_GAS_PVT_MULTIPLEXER_CALL(return pvtImpl.saturatedWaterVaporizationFactor(regionIdx, temperature, pressure, saltConcentration)); }
264
273 template <class Evaluation = Scalar>
274 Evaluation saturationPressure(unsigned regionIdx,
275 const Evaluation& temperature,
276 const Evaluation& Rv) const
277 { OPM_GAS_PVT_MULTIPLEXER_CALL(return pvtImpl.saturationPressure(regionIdx, temperature, Rv)); }
278
282 template <class Evaluation>
283 Evaluation diffusionCoefficient(const Evaluation& temperature,
284 const Evaluation& pressure,
285 unsigned compIdx) const
286 {
287 OPM_GAS_PVT_MULTIPLEXER_CALL(return pvtImpl.diffusionCoefficient(temperature, pressure, compIdx));
288 }
289
295 GasPvtApproach approach() const
296 { return gasPvtApproach_; }
297
298 // get the parameter object for the dry gas case
299 template <GasPvtApproach approachV>
300 typename std::enable_if<approachV == GasPvtApproach::DryGas, DryGasPvt<Scalar> >::type& getRealPvt()
301 {
302 assert(approach() == approachV);
303 return *static_cast<DryGasPvt<Scalar>* >(realGasPvt_.get());
304 }
305
306 template <GasPvtApproach approachV>
307 typename std::enable_if<approachV == GasPvtApproach::DryGas, const DryGasPvt<Scalar> >::type& getRealPvt() const
308 {
309 assert(approach() == approachV);
310 return *static_cast<const DryGasPvt<Scalar>* >(realGasPvt_.get());
311 }
312
313 // get the parameter object for the dry humid gas case
314 template <GasPvtApproach approachV>
315 typename std::enable_if<approachV == GasPvtApproach::DryHumidGas, DryHumidGasPvt<Scalar> >::type& getRealPvt()
316 {
317 assert(approach() == approachV);
318 return *static_cast<DryHumidGasPvt<Scalar>* >(realGasPvt_.get());
319 }
320
321 template <GasPvtApproach approachV>
322 typename std::enable_if<approachV == GasPvtApproach::DryHumidGas, const DryHumidGasPvt<Scalar> >::type& getRealPvt() const
323 {
324 assert(approach() == approachV);
325 return *static_cast<const DryHumidGasPvt<Scalar>* >(realGasPvt_.get());
326 }
327
328 // get the parameter object for the wet humid gas case
329 template <GasPvtApproach approachV>
330 typename std::enable_if<approachV == GasPvtApproach::WetHumidGas, WetHumidGasPvt<Scalar> >::type& getRealPvt()
331 {
332 assert(approach() == approachV);
333 return *static_cast<WetHumidGasPvt<Scalar>* >(realGasPvt_.get());
334 }
335
336 template <GasPvtApproach approachV>
337 typename std::enable_if<approachV == GasPvtApproach::WetHumidGas, const WetHumidGasPvt<Scalar> >::type& getRealPvt() const
338 {
339 assert(approach() == approachV);
340 return *static_cast<const WetHumidGasPvt<Scalar>* >(realGasPvt_.get());
341 }
342
343 // get the parameter object for the wet gas case
344 template <GasPvtApproach approachV>
345 typename std::enable_if<approachV == GasPvtApproach::WetGas, WetGasPvt<Scalar> >::type& getRealPvt()
346 {
347 assert(approach() == approachV);
348 return *static_cast<WetGasPvt<Scalar>* >(realGasPvt_.get());
349 }
350
351 template <GasPvtApproach approachV>
352 typename std::enable_if<approachV == GasPvtApproach::WetGas, const WetGasPvt<Scalar> >::type& getRealPvt() const
353 {
354 assert(approach() == approachV);
355 return *static_cast<const WetGasPvt<Scalar>* >(realGasPvt_.get());
356 }
357
358 // get the parameter object for the thermal gas case
359 template <GasPvtApproach approachV>
360 typename std::enable_if<approachV == GasPvtApproach::ThermalGas, GasPvtThermal<Scalar> >::type& getRealPvt()
361 {
362 assert(approach() == approachV);
363 return *static_cast<GasPvtThermal<Scalar>* >(realGasPvt_.get());
364 }
365 template <GasPvtApproach approachV>
366 typename std::enable_if<approachV == GasPvtApproach::ThermalGas, const GasPvtThermal<Scalar> >::type& getRealPvt() const
367 {
368 assert(approach() == approachV);
369 return *static_cast<const GasPvtThermal<Scalar>* >(realGasPvt_.get());
370 }
371
372 template <GasPvtApproach approachV>
373 typename std::enable_if<approachV == GasPvtApproach::Co2Gas, Co2GasPvt<Scalar> >::type& getRealPvt()
374 {
375 assert(approach() == approachV);
376 return *static_cast<Co2GasPvt<Scalar>* >(realGasPvt_.get());
377 }
378
379 template <GasPvtApproach approachV>
380 typename std::enable_if<approachV == GasPvtApproach::Co2Gas, const Co2GasPvt<Scalar> >::type& getRealPvt() const
381 {
382 assert(approach() == approachV);
383 return *static_cast<const Co2GasPvt<Scalar>* >(realGasPvt_.get());
384 }
385
386 template <GasPvtApproach approachV>
387 typename std::enable_if<approachV == GasPvtApproach::H2Gas, H2GasPvt<Scalar> >::type& getRealPvt()
388 {
389 assert(approach() == approachV);
390 return *static_cast<H2GasPvt<Scalar>* >(realGasPvt_.get());
391 }
392
393 template <GasPvtApproach approachV>
394 typename std::enable_if<approachV == GasPvtApproach::H2Gas, const H2GasPvt<Scalar> >::type& getRealPvt() const
395 {
396 assert(approach() == approachV);
397 return *static_cast<const H2GasPvt<Scalar>* >(realGasPvt_.get());
398 }
399
400 const void* realGasPvt() const { return realGasPvt_.get(); }
401
402 GasPvtMultiplexer<Scalar,enableThermal>&
403 operator=(const GasPvtMultiplexer<Scalar,enableThermal>& data);
404
405private:
406 using UniqueVoidPtrWithDeleter = std::unique_ptr<void, std::function<void(void*)>>;
407
408 template <class ConcreteGasPvt> UniqueVoidPtrWithDeleter makeGasPvt();
409
410 template <class ConcretePvt> UniqueVoidPtrWithDeleter copyPvt(const UniqueVoidPtrWithDeleter&);
411
412 GasPvtApproach gasPvtApproach_{GasPvtApproach::NoGas};
413 UniqueVoidPtrWithDeleter realGasPvt_;
414
415 void deleter(void* ptr);
416};
417
418} // namespace Opm
419
420#endif
This class represents the Pressure-Volume-Temperature relations of the gas phase for CO2.
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 with vaporized water...
This class implements temperature dependence of the PVT properties of gas.
This class represents the Pressure-Volume-Temperature relations of the gas phase for H2.
This class represents the Pressure-Volume-Temperature relations of the gas phas with vaporized oil.
This class represents the Pressure-Volume-Temperature relations of the gas phase with vaporized oil a...
This class represents the Pressure-Volume-Temperature relations of the gas phase without vaporized oi...
Definition DryGasPvt.hpp:50
Definition EclipseState.hpp:62
std::pair< LhsEval, LhsEval > inverseFormationVolumeFactorAndViscosity(const FluidState &fluidState, unsigned regionIdx)
Returns the formation volume factor [-] and viscosity [Pa s] of the fluid phase.
Definition GasPvtMultiplexer.hpp:214
Scalar gasReferenceDensity(unsigned regionIdx) const
Return the reference density which are considered by this PVT-object.
Definition GasPvtMultiplexer.cpp:57
Evaluation saturatedInverseFormationVolumeFactor(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure) const
Returns the formation volume factor [-] of oil saturated gas given a set of parameters.
Definition GasPvtMultiplexer.hpp:221
unsigned numRegions() const
Return the number of PVT regions which are considered by this PVT-object.
Definition GasPvtMultiplexer.cpp:41
Evaluation saturatedWaterVaporizationFactor(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure) const
Returns the water vaporization factor [m^3/m^3] of water saturated gas.
Definition GasPvtMultiplexer.hpp:250
Evaluation inverseFormationVolumeFactor(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure, const Evaluation &Rv, const Evaluation &Rvw) const
Returns the formation volume factor [-] of the fluid phase.
Definition GasPvtMultiplexer.hpp:202
GasPvtApproach approach() const
Returns the concrete approach for calculating the PVT relations.
Definition GasPvtMultiplexer.hpp:295
Evaluation internalEnergy(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure, const Evaluation &Rv, const Evaluation &Rvw) const
Returns the specific enthalpy [J/kg] of gas given a set of parameters.
Definition GasPvtMultiplexer.hpp:169
Evaluation viscosity(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure, const Evaluation &Rv, const Evaluation &Rvw) const
Returns the dynamic viscosity [Pa s] of the fluid phase given a set of parameters.
Definition GasPvtMultiplexer.hpp:182
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 GasPvtMultiplexer.hpp:283
Evaluation saturatedWaterVaporizationFactor(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure, const Evaluation &saltConcentration) const
Returns the water vaporization factor [m^3/m^3] of water saturated gas.
Definition GasPvtMultiplexer.hpp:259
Evaluation saturatedOilVaporizationFactor(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure) const
Returns the oil vaporization factor [m^3/m^3] of oil saturated gas.
Definition GasPvtMultiplexer.hpp:230
Evaluation saturationPressure(unsigned regionIdx, const Evaluation &temperature, const Evaluation &Rv) const
Returns the saturation pressure of the gas phase [Pa] depending on its mass fraction of the oil compo...
Definition GasPvtMultiplexer.hpp:274
Evaluation saturatedOilVaporizationFactor(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure, const Evaluation &oilSaturation, const Evaluation &maxOilSaturation) const
Returns the oil vaporization factor [m^3/m^3] of oil saturated gas.
Definition GasPvtMultiplexer.hpp:239
Evaluation saturatedViscosity(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure) const
Returns the dynamic viscosity [Pa s] of oil saturated gas given a set of parameters.
Definition GasPvtMultiplexer.hpp:193
Definition Schedule.hpp:101
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30