opm-common
Loading...
Searching...
No Matches
OilPvtMultiplexer.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_OIL_PVT_MULTIPLEXER_HPP
28#define OPM_OIL_PVT_MULTIPLEXER_HPP
29
36
37namespace Opm {
38
39#if HAVE_ECL_INPUT
40class EclipseState;
41class Schedule;
42#endif
43
44#define OPM_OIL_PVT_MULTIPLEXER_CALL(codeToCall, ...) \
45 switch (approach_) { \
46 case OilPvtApproach::ConstantCompressibilityOil: { \
47 auto& pvtImpl = getRealPvt<OilPvtApproach::ConstantCompressibilityOil>(); \
48 codeToCall; \
49 __VA_ARGS__; \
50 } \
51 case OilPvtApproach::DeadOil: { \
52 auto& pvtImpl = getRealPvt<OilPvtApproach::DeadOil>(); \
53 codeToCall; \
54 __VA_ARGS__; \
55 } \
56 case OilPvtApproach::LiveOil: { \
57 auto& pvtImpl = getRealPvt<OilPvtApproach::LiveOil>(); \
58 codeToCall; \
59 __VA_ARGS__; \
60 } \
61 case OilPvtApproach::ThermalOil: { \
62 auto& pvtImpl = getRealPvt<OilPvtApproach::ThermalOil>(); \
63 codeToCall; \
64 __VA_ARGS__; \
65 } \
66 case OilPvtApproach::BrineCo2: { \
67 auto& pvtImpl = getRealPvt<OilPvtApproach::BrineCo2>(); \
68 codeToCall; \
69 __VA_ARGS__; \
70 } \
71 case OilPvtApproach::BrineH2: { \
72 auto& pvtImpl = getRealPvt<OilPvtApproach::BrineH2>(); \
73 codeToCall; \
74 __VA_ARGS__; \
75 } \
76 default: \
77 case OilPvtApproach::NoOil: \
78 throw std::logic_error("Not implemented: Oil PVT of this deck!"); \
79 }
80
81enum class OilPvtApproach {
82 NoOil,
83 LiveOil,
84 DeadOil,
85 ConstantCompressibilityOil,
86 ThermalOil,
87 BrineCo2,
88 BrineH2
89};
90
103template <class Scalar, bool enableThermal = true>
104class OilPvtMultiplexer
105{
106public:
107 OilPvtMultiplexer()
108 : approach_(OilPvtApproach::NoOil)
109 , realOilPvt_(nullptr)
110 {
111 }
112
113 OilPvtMultiplexer(OilPvtApproach approach, void* realOilPvt)
114 : approach_(approach)
115 , realOilPvt_(realOilPvt)
116 { }
117
118 OilPvtMultiplexer(const OilPvtMultiplexer<Scalar,enableThermal>& data)
119 {
120 *this = data;
121 }
122
123 ~OilPvtMultiplexer();
124
125 bool mixingEnergy() const
126 {
127 return approach_ == OilPvtApproach::ThermalOil;
128 }
129#if HAVE_ECL_INPUT
135 void initFromState(const EclipseState& eclState, const Schedule& schedule);
136#endif // HAVE_ECL_INPUT
137
138
139 void initEnd();
140
141 bool isActive() const
142 {
143 return approach_ != OilPvtApproach::NoOil;
144 }
145
149 unsigned numRegions() const;
150
151 void setVapPars(const Scalar par1, const Scalar par2);
152
156 Scalar oilReferenceDensity(unsigned regionIdx) const;
157
161 template <class Evaluation>
162 Evaluation internalEnergy(unsigned regionIdx,
163 const Evaluation& temperature,
164 const Evaluation& pressure,
165 const Evaluation& Rs) const
166 { OPM_OIL_PVT_MULTIPLEXER_CALL(return pvtImpl.internalEnergy(regionIdx, temperature, pressure, Rs)); }
167
168 Scalar hVap(unsigned regionIdx) const
169 { OPM_OIL_PVT_MULTIPLEXER_CALL(return pvtImpl.hVap(regionIdx)); }
173 template <class Evaluation>
174 Evaluation viscosity(unsigned regionIdx,
175 const Evaluation& temperature,
176 const Evaluation& pressure,
177 const Evaluation& Rs) const
178 { OPM_OIL_PVT_MULTIPLEXER_CALL(return pvtImpl.viscosity(regionIdx, temperature, pressure, Rs)); }
179
183 template <class Evaluation>
184 Evaluation saturatedViscosity(unsigned regionIdx,
185 const Evaluation& temperature,
186 const Evaluation& pressure) const
187 { OPM_OIL_PVT_MULTIPLEXER_CALL(return pvtImpl.saturatedViscosity(regionIdx, temperature, pressure)); }
188
192 template <class Evaluation>
193 Evaluation inverseFormationVolumeFactor(unsigned regionIdx,
194 const Evaluation& temperature,
195 const Evaluation& pressure,
196 const Evaluation& Rs) const
197 { OPM_OIL_PVT_MULTIPLEXER_CALL(return pvtImpl.inverseFormationVolumeFactor(regionIdx, temperature, pressure, Rs)); }
198
202 template <class FluidState, class LhsEval = typename FluidState::Scalar>
203 std::pair<LhsEval, LhsEval>
204 inverseFormationVolumeFactorAndViscosity(const FluidState& fluidState, unsigned regionIdx)
205 { OPM_OIL_PVT_MULTIPLEXER_CALL(return pvtImpl.inverseFormationVolumeFactorAndViscosity(fluidState, regionIdx)); }
206
210 template <class Evaluation>
211 Evaluation saturatedInverseFormationVolumeFactor(unsigned regionIdx,
212 const Evaluation& temperature,
213 const Evaluation& pressure) const
214 { OPM_OIL_PVT_MULTIPLEXER_CALL(return pvtImpl.saturatedInverseFormationVolumeFactor(regionIdx, temperature, pressure)); }
215
219 template <class Evaluation>
220 Evaluation saturatedGasDissolutionFactor(unsigned regionIdx,
221 const Evaluation& temperature,
222 const Evaluation& pressure) const
223 { OPM_OIL_PVT_MULTIPLEXER_CALL(return pvtImpl.saturatedGasDissolutionFactor(regionIdx, temperature, pressure)); }
224
228 template <class Evaluation>
229 Evaluation saturatedGasDissolutionFactor(unsigned regionIdx,
230 const Evaluation& temperature,
231 const Evaluation& pressure,
232 const Evaluation& oilSaturation,
233 const Evaluation& maxOilSaturation) const
234 { OPM_OIL_PVT_MULTIPLEXER_CALL(return pvtImpl.saturatedGasDissolutionFactor(regionIdx, temperature, pressure, oilSaturation, maxOilSaturation)); }
235
243 template <class Evaluation>
244 Evaluation saturationPressure(unsigned regionIdx,
245 const Evaluation& temperature,
246 const Evaluation& Rs) const
247 { OPM_OIL_PVT_MULTIPLEXER_CALL(return pvtImpl.saturationPressure(regionIdx, temperature, Rs)); }
248
252 template <class Evaluation>
253 Evaluation diffusionCoefficient(const Evaluation& temperature,
254 const Evaluation& pressure,
255 unsigned compIdx) const
256 {
257 OPM_OIL_PVT_MULTIPLEXER_CALL(return pvtImpl.diffusionCoefficient(temperature, pressure, compIdx));
258 }
259
260 void setApproach(OilPvtApproach appr);
261
267 OilPvtApproach approach() const
268 { return approach_; }
269
270 // get the concrete parameter object for the oil phase
271 template <OilPvtApproach approachV>
272 typename std::enable_if<approachV == OilPvtApproach::LiveOil, LiveOilPvt<Scalar> >::type& getRealPvt()
273 {
274 assert(approach() == approachV);
275 return *static_cast<LiveOilPvt<Scalar>* >(realOilPvt_);
276 }
277
278 template <OilPvtApproach approachV>
279 typename std::enable_if<approachV == OilPvtApproach::LiveOil, const LiveOilPvt<Scalar> >::type& getRealPvt() const
280 {
281 assert(approach() == approachV);
282 return *static_cast<LiveOilPvt<Scalar>* >(realOilPvt_);
283 }
284
285 template <OilPvtApproach approachV>
286 typename std::enable_if<approachV == OilPvtApproach::DeadOil, DeadOilPvt<Scalar> >::type& getRealPvt()
287 {
288 assert(approach() == approachV);
289 return *static_cast<DeadOilPvt<Scalar>* >(realOilPvt_);
290 }
291
292 template <OilPvtApproach approachV>
293 typename std::enable_if<approachV == OilPvtApproach::DeadOil, const DeadOilPvt<Scalar> >::type& getRealPvt() const
294 {
295 assert(approach() == approachV);
296 return *static_cast<DeadOilPvt<Scalar>* >(realOilPvt_);
297 }
298
299 template <OilPvtApproach approachV>
300 typename std::enable_if<approachV == OilPvtApproach::ConstantCompressibilityOil, ConstantCompressibilityOilPvt<Scalar> >::type& getRealPvt()
301 {
302 assert(approach() == approachV);
303 return *static_cast<ConstantCompressibilityOilPvt<Scalar>* >(realOilPvt_);
304 }
305
306 template <OilPvtApproach approachV>
307 typename std::enable_if<approachV == OilPvtApproach::ConstantCompressibilityOil, const ConstantCompressibilityOilPvt<Scalar> >::type& getRealPvt() const
308 {
309 assert(approach() == approachV);
310 return *static_cast<ConstantCompressibilityOilPvt<Scalar>* >(realOilPvt_);
311 }
312
313 template <OilPvtApproach approachV>
314 typename std::enable_if<approachV == OilPvtApproach::ThermalOil, OilPvtThermal<Scalar> >::type& getRealPvt()
315 {
316 assert(approach() == approachV);
317 return *static_cast<OilPvtThermal<Scalar>* >(realOilPvt_);
318 }
319
320 template <OilPvtApproach approachV>
321 typename std::enable_if<approachV == OilPvtApproach::ThermalOil, const OilPvtThermal<Scalar> >::type& getRealPvt() const
322 {
323 assert(approach() == approachV);
324 return *static_cast<const OilPvtThermal<Scalar>* >(realOilPvt_);
325 }
326
327 template <OilPvtApproach approachV>
328 typename std::enable_if<approachV == OilPvtApproach::BrineCo2, BrineCo2Pvt<Scalar> >::type& getRealPvt()
329 {
330 assert(approach() == approachV);
331 return *static_cast<BrineCo2Pvt<Scalar>* >(realOilPvt_);
332 }
333
334 template <OilPvtApproach approachV>
335 typename std::enable_if<approachV == OilPvtApproach::BrineCo2, const BrineCo2Pvt<Scalar> >::type& getRealPvt() const
336 {
337 assert(approach() == approachV);
338 return *static_cast<const BrineCo2Pvt<Scalar>* >(realOilPvt_);
339 }
340
341 const void* realOilPvt() const { return realOilPvt_; }
342
343 template <OilPvtApproach approachV>
344 typename std::enable_if<approachV == OilPvtApproach::BrineH2, BrineH2Pvt<Scalar> >::type& getRealPvt()
345 {
346 assert(approach() == approachV);
347 return *static_cast<BrineH2Pvt<Scalar>* >(realOilPvt_);
348 }
349
350 template <OilPvtApproach approachV>
351 typename std::enable_if<approachV == OilPvtApproach::BrineH2, const BrineH2Pvt<Scalar> >::type& getRealPvt() const
352 {
353 assert(approach() == approachV);
354 return *static_cast<const BrineH2Pvt<Scalar>* >(realOilPvt_);
355 }
356
357 OilPvtMultiplexer<Scalar,enableThermal>&
358 operator=(const OilPvtMultiplexer<Scalar,enableThermal>& data);
359
360private:
361 OilPvtApproach approach_{OilPvtApproach::NoOil};
362 void* realOilPvt_{nullptr};
363};
364
365} // namespace Opm
366
367#endif
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 oil phase without dissolved ga...
This class represents the Pressure-Volume-Temperature relations of the oil phase without dissolved ga...
This class represents the Pressure-Volume-Temperature relations of the oil phas with dissolved gas.
This class implements temperature dependence of the PVT properties of oil.
Definition EclipseState.hpp:62
This class represents the Pressure-Volume-Temperature relations of the oil phas with dissolved gas.
Definition LiveOilPvt.hpp:51
unsigned numRegions() const
Return the number of PVT regions which are considered by this PVT-object.
Definition OilPvtMultiplexer.cpp:102
Evaluation internalEnergy(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure, const Evaluation &Rs) const
Returns the specific enthalpy [J/kg] oil given a set of parameters.
Definition OilPvtMultiplexer.hpp:162
Evaluation inverseFormationVolumeFactor(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure, const Evaluation &Rs) const
Returns the formation volume factor [-] of the fluid phase.
Definition OilPvtMultiplexer.hpp:193
Evaluation saturatedViscosity(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure) const
Returns the dynamic viscosity [Pa s] of the fluid phase given a set of parameters.
Definition OilPvtMultiplexer.hpp:184
OilPvtApproach approach() const
Returns the concrete approach for calculating the PVT relations.
Definition OilPvtMultiplexer.hpp:267
Evaluation saturatedInverseFormationVolumeFactor(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure) const
Returns the formation volume factor [-] of the fluid phase.
Definition OilPvtMultiplexer.hpp:211
Evaluation saturatedGasDissolutionFactor(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure) const
Returns the gas dissolution factor [m^3/m^3] of saturated oil.
Definition OilPvtMultiplexer.hpp:220
Scalar oilReferenceDensity(unsigned regionIdx) const
Return the reference density which are considered by this PVT-object.
Definition OilPvtMultiplexer.cpp:116
Evaluation viscosity(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure, const Evaluation &Rs) const
Returns the dynamic viscosity [Pa s] of the fluid phase given a set of parameters.
Definition OilPvtMultiplexer.hpp:174
Evaluation saturationPressure(unsigned regionIdx, const Evaluation &temperature, const Evaluation &Rs) const
Returns the saturation pressure [Pa] of oil given the mass fraction of the gas component in the oil p...
Definition OilPvtMultiplexer.hpp:244
Evaluation saturatedGasDissolutionFactor(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure, const Evaluation &oilSaturation, const Evaluation &maxOilSaturation) const
Returns the gas dissolution factor [m^3/m^3] of saturated oil.
Definition OilPvtMultiplexer.hpp:229
std::pair< LhsEval, LhsEval > inverseFormationVolumeFactorAndViscosity(const FluidState &fluidState, unsigned regionIdx)
Returns the formation volume factor [-] and viscosity [Pa s] of the fluid phase.
Definition OilPvtMultiplexer.hpp:204
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 OilPvtMultiplexer.hpp:253
Definition Schedule.hpp:101
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30