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
37
38namespace Opm {
39
40class EclipseState;
41class Schedule;
42
43#define OPM_OIL_PVT_MULTIPLEXER_CALL(codeToCall, ...) \
44 switch (approach_) { \
45 case OilPvtApproach::ConstantCompressibilityOil: { \
46 auto& pvtImpl = getRealPvt<OilPvtApproach::ConstantCompressibilityOil>(); \
47 codeToCall; \
48 __VA_ARGS__; \
49 } \
50 case OilPvtApproach::DeadOil: { \
51 auto& pvtImpl = getRealPvt<OilPvtApproach::DeadOil>(); \
52 codeToCall; \
53 __VA_ARGS__; \
54 } \
55 case OilPvtApproach::LiveOil: { \
56 auto& pvtImpl = getRealPvt<OilPvtApproach::LiveOil>(); \
57 codeToCall; \
58 __VA_ARGS__; \
59 } \
60 case OilPvtApproach::ThermalOil: { \
61 auto& pvtImpl = getRealPvt<OilPvtApproach::ThermalOil>(); \
62 codeToCall; \
63 __VA_ARGS__; \
64 } \
65 case OilPvtApproach::BrineCo2: { \
66 auto& pvtImpl = getRealPvt<OilPvtApproach::BrineCo2>(); \
67 codeToCall; \
68 __VA_ARGS__; \
69 } \
70 case OilPvtApproach::BrineH2: { \
71 auto& pvtImpl = getRealPvt<OilPvtApproach::BrineH2>(); \
72 codeToCall; \
73 __VA_ARGS__; \
74 } \
75 case OilPvtApproach::ConstantRsDeadOil: { \
76 auto& pvtImpl = getRealPvt<OilPvtApproach::ConstantRsDeadOil>(); \
77 codeToCall; \
78 __VA_ARGS__; \
79 } \
80 default: \
81 case OilPvtApproach::NoOil: \
82 throw std::logic_error("Not implemented: Oil PVT of this deck!"); \
83 }
84
85enum class OilPvtApproach {
86 NoOil,
87 LiveOil,
88 DeadOil,
89 ConstantCompressibilityOil,
90 ThermalOil,
91 BrineCo2,
92 BrineH2,
93 ConstantRsDeadOil
94};
95
108template <class Scalar, bool enableThermal = true>
109class OilPvtMultiplexer
110{
111public:
112 OilPvtMultiplexer()
113 : approach_(OilPvtApproach::NoOil)
114 , realOilPvt_(nullptr)
115 {
116 }
117
118 OilPvtMultiplexer(OilPvtApproach approach, void* realOilPvt)
119 : approach_(approach)
120 , realOilPvt_(realOilPvt)
121 { }
122
123 OilPvtMultiplexer(const OilPvtMultiplexer<Scalar,enableThermal>& data)
124 {
125 *this = data;
126 }
127
128 ~OilPvtMultiplexer();
129
130 bool mixingEnergy() const
131 {
132 return approach_ == OilPvtApproach::ThermalOil;
133 }
134
140 void initFromState(const EclipseState& eclState, const Schedule& schedule);
141
142 void initEnd();
143
144 bool isActive() const
145 {
146 return approach_ != OilPvtApproach::NoOil;
147 }
148
152 unsigned numRegions() const;
153
154 void setVapPars(const Scalar par1, const Scalar par2);
155
159 Scalar oilReferenceDensity(unsigned regionIdx) const;
160
164 template <class Evaluation>
165 Evaluation internalEnergy(unsigned regionIdx,
166 const Evaluation& temperature,
167 const Evaluation& pressure,
168 const Evaluation& Rs) const
169 { OPM_OIL_PVT_MULTIPLEXER_CALL(return pvtImpl.internalEnergy(regionIdx, temperature, pressure, Rs)); }
170
171 Scalar hVap(unsigned regionIdx) const
172 { OPM_OIL_PVT_MULTIPLEXER_CALL(return pvtImpl.hVap(regionIdx)); }
176 template <class Evaluation>
177 Evaluation viscosity(unsigned regionIdx,
178 const Evaluation& temperature,
179 const Evaluation& pressure,
180 const Evaluation& Rs) const
181 { OPM_OIL_PVT_MULTIPLEXER_CALL(return pvtImpl.viscosity(regionIdx, temperature, pressure, Rs)); }
182
186 template <class Evaluation>
187 Evaluation saturatedViscosity(unsigned regionIdx,
188 const Evaluation& temperature,
189 const Evaluation& pressure) const
190 { OPM_OIL_PVT_MULTIPLEXER_CALL(return pvtImpl.saturatedViscosity(regionIdx, temperature, pressure)); }
191
195 template <class Evaluation>
196 Evaluation inverseFormationVolumeFactor(unsigned regionIdx,
197 const Evaluation& temperature,
198 const Evaluation& pressure,
199 const Evaluation& Rs) const
200 { OPM_OIL_PVT_MULTIPLEXER_CALL(return pvtImpl.inverseFormationVolumeFactor(regionIdx, temperature, pressure, Rs)); }
201
205 template <class FluidState, class LhsEval = typename FluidState::ValueType>
206 std::pair<LhsEval, LhsEval>
207 inverseFormationVolumeFactorAndViscosity(const FluidState& fluidState, unsigned regionIdx)
208 { OPM_OIL_PVT_MULTIPLEXER_CALL(return pvtImpl.inverseFormationVolumeFactorAndViscosity(fluidState, regionIdx)); }
209
213 template <class Evaluation>
214 Evaluation saturatedInverseFormationVolumeFactor(unsigned regionIdx,
215 const Evaluation& temperature,
216 const Evaluation& pressure) const
217 { OPM_OIL_PVT_MULTIPLEXER_CALL(return pvtImpl.saturatedInverseFormationVolumeFactor(regionIdx, temperature, pressure)); }
218
222 template <class Evaluation>
223 Evaluation saturatedGasDissolutionFactor(unsigned regionIdx,
224 const Evaluation& temperature,
225 const Evaluation& pressure) const
226 { OPM_OIL_PVT_MULTIPLEXER_CALL(return pvtImpl.saturatedGasDissolutionFactor(regionIdx, temperature, pressure)); }
227
231 template <class Evaluation>
232 Evaluation saturatedGasDissolutionFactor(unsigned regionIdx,
233 const Evaluation& temperature,
234 const Evaluation& pressure,
235 const Evaluation& oilSaturation,
236 const Evaluation& maxOilSaturation) const
237 { OPM_OIL_PVT_MULTIPLEXER_CALL(return pvtImpl.saturatedGasDissolutionFactor(regionIdx, temperature, pressure, oilSaturation, maxOilSaturation)); }
238
246 template <class Evaluation>
247 Evaluation saturationPressure(unsigned regionIdx,
248 const Evaluation& temperature,
249 const Evaluation& Rs) const
250 { OPM_OIL_PVT_MULTIPLEXER_CALL(return pvtImpl.saturationPressure(regionIdx, temperature, Rs)); }
251
255 template <class Evaluation>
256 Evaluation diffusionCoefficient(const Evaluation& temperature,
257 const Evaluation& pressure,
258 unsigned compIdx) const
259 {
260 OPM_OIL_PVT_MULTIPLEXER_CALL(return pvtImpl.diffusionCoefficient(temperature, pressure, compIdx));
261 }
262
263 void setApproach(OilPvtApproach appr);
264
270 OilPvtApproach approach() const
271 { return approach_; }
272
273 // get the concrete parameter object for the oil phase
274 template <OilPvtApproach approachV>
275 typename std::enable_if<approachV == OilPvtApproach::LiveOil, LiveOilPvt<Scalar> >::type& getRealPvt()
276 {
277 assert(approach() == approachV);
278 return *static_cast<LiveOilPvt<Scalar>* >(realOilPvt_);
279 }
280
281 template <OilPvtApproach approachV>
282 typename std::enable_if<approachV == OilPvtApproach::LiveOil, const LiveOilPvt<Scalar> >::type& getRealPvt() const
283 {
284 assert(approach() == approachV);
285 return *static_cast<LiveOilPvt<Scalar>* >(realOilPvt_);
286 }
287
288 template <OilPvtApproach approachV>
289 typename std::enable_if<approachV == OilPvtApproach::DeadOil, DeadOilPvt<Scalar> >::type& getRealPvt()
290 {
291 assert(approach() == approachV);
292 return *static_cast<DeadOilPvt<Scalar>* >(realOilPvt_);
293 }
294
295 template <OilPvtApproach approachV>
296 typename std::enable_if<approachV == OilPvtApproach::DeadOil, const DeadOilPvt<Scalar> >::type& getRealPvt() const
297 {
298 assert(approach() == approachV);
299 return *static_cast<DeadOilPvt<Scalar>* >(realOilPvt_);
300 }
301
302 template <OilPvtApproach approachV>
303 typename std::enable_if<approachV == OilPvtApproach::ConstantCompressibilityOil, ConstantCompressibilityOilPvt<Scalar> >::type& getRealPvt()
304 {
305 assert(approach() == approachV);
306 return *static_cast<ConstantCompressibilityOilPvt<Scalar>* >(realOilPvt_);
307 }
308
309 template <OilPvtApproach approachV>
310 typename std::enable_if<approachV == OilPvtApproach::ConstantCompressibilityOil, const ConstantCompressibilityOilPvt<Scalar> >::type& getRealPvt() const
311 {
312 assert(approach() == approachV);
313 return *static_cast<ConstantCompressibilityOilPvt<Scalar>* >(realOilPvt_);
314 }
315
316 template <OilPvtApproach approachV>
317 typename std::enable_if<approachV == OilPvtApproach::ThermalOil, OilPvtThermal<Scalar> >::type& getRealPvt()
318 {
319 assert(approach() == approachV);
320 return *static_cast<OilPvtThermal<Scalar>* >(realOilPvt_);
321 }
322
323 template <OilPvtApproach approachV>
324 typename std::enable_if<approachV == OilPvtApproach::ThermalOil, const OilPvtThermal<Scalar> >::type& getRealPvt() const
325 {
326 assert(approach() == approachV);
327 return *static_cast<const OilPvtThermal<Scalar>* >(realOilPvt_);
328 }
329
330 template <OilPvtApproach approachV>
331 typename std::enable_if<approachV == OilPvtApproach::BrineCo2, BrineCo2Pvt<Scalar> >::type& getRealPvt()
332 {
333 assert(approach() == approachV);
334 return *static_cast<BrineCo2Pvt<Scalar>* >(realOilPvt_);
335 }
336
337 template <OilPvtApproach approachV>
338 typename std::enable_if<approachV == OilPvtApproach::BrineCo2, const BrineCo2Pvt<Scalar> >::type& getRealPvt() const
339 {
340 assert(approach() == approachV);
341 return *static_cast<const BrineCo2Pvt<Scalar>* >(realOilPvt_);
342 }
343
344 const void* realOilPvt() const { return realOilPvt_; }
345
346 template <OilPvtApproach approachV>
347 typename std::enable_if<approachV == OilPvtApproach::BrineH2, BrineH2Pvt<Scalar> >::type& getRealPvt()
348 {
349 assert(approach() == approachV);
350 return *static_cast<BrineH2Pvt<Scalar>* >(realOilPvt_);
351 }
352
353 template <OilPvtApproach approachV>
354 typename std::enable_if<approachV == OilPvtApproach::BrineH2, const BrineH2Pvt<Scalar> >::type& getRealPvt() const
355 {
356 assert(approach() == approachV);
357 return *static_cast<const BrineH2Pvt<Scalar>* >(realOilPvt_);
358 }
359
360 template <OilPvtApproach approachV>
361 typename std::enable_if<approachV == OilPvtApproach::ConstantRsDeadOil, ConstantRsDeadOilPvt<Scalar> >::type& getRealPvt()
362 {
363 assert(approach() == approachV);
364 return *static_cast<ConstantRsDeadOilPvt<Scalar>* >(realOilPvt_);
365 }
366
367 template <OilPvtApproach approachV>
368 typename std::enable_if<approachV == OilPvtApproach::ConstantRsDeadOil, const ConstantRsDeadOilPvt<Scalar> >::type& getRealPvt() const
369 {
370 assert(approach() == approachV);
371 return *static_cast<const ConstantRsDeadOilPvt<Scalar>* >(realOilPvt_);
372 }
373
374 OilPvtMultiplexer<Scalar,enableThermal>&
375 operator=(const OilPvtMultiplexer<Scalar,enableThermal>& data);
376
377private:
378 OilPvtApproach approach_{OilPvtApproach::NoOil};
379 void* realOilPvt_{nullptr};
380};
381
382} // namespace Opm
383
384#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 dead oil with constant dissolved g...
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:66
This class represents the Pressure-Volume-Temperature relations of the oil phas with dissolved gas.
Definition LiveOilPvt.hpp:49
unsigned numRegions() const
Return the number of PVT regions which are considered by this PVT-object.
Definition OilPvtMultiplexer.cpp:116
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:165
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:196
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:187
OilPvtApproach approach() const
Returns the concrete approach for calculating the PVT relations.
Definition OilPvtMultiplexer.hpp:270
Evaluation saturatedInverseFormationVolumeFactor(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure) const
Returns the formation volume factor [-] of the fluid phase.
Definition OilPvtMultiplexer.hpp:214
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:223
void initFromState(const EclipseState &eclState, const Schedule &schedule)
Initialize the parameters for water using an ECL state.
Definition OilPvtMultiplexer.cpp:74
Scalar oilReferenceDensity(unsigned regionIdx) const
Return the reference density which are considered by this PVT-object.
Definition OilPvtMultiplexer.cpp:130
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:177
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:247
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:232
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:207
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:256
Definition Schedule.hpp:101
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30