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
31#include "DeadOilPvt.hpp"
32#include "LiveOilPvt.hpp"
33#include "OilPvtThermal.hpp"
34#include "BrineCo2Pvt.hpp"
35
36#if HAVE_ECL_INPUT
37#include <opm/input/eclipse/EclipseState/EclipseState.hpp>
38#include <opm/input/eclipse/EclipseState/Runspec.hpp>
39#endif
40
41namespace Opm {
42#define OPM_OIL_PVT_MULTIPLEXER_CALL(codeToCall) \
43 switch (approach_) { \
44 case OilPvtApproach::ConstantCompressibilityOilPvt: { \
45 auto& pvtImpl = getRealPvt<OilPvtApproach::ConstantCompressibilityOilPvt>(); \
46 codeToCall; \
47 break; \
48 } \
49 case OilPvtApproach::DeadOilPvt: { \
50 auto& pvtImpl = getRealPvt<OilPvtApproach::DeadOilPvt>(); \
51 codeToCall; \
52 break; \
53 } \
54 case OilPvtApproach::LiveOilPvt: { \
55 auto& pvtImpl = getRealPvt<OilPvtApproach::LiveOilPvt>(); \
56 codeToCall; \
57 break; \
58 } \
59 case OilPvtApproach::ThermalOilPvt: { \
60 auto& pvtImpl = getRealPvt<OilPvtApproach::ThermalOilPvt>(); \
61 codeToCall; \
62 break; \
63 } \
64 case OilPvtApproach::BrineCo2Pvt: { \
65 auto& pvtImpl = getRealPvt<OilPvtApproach::BrineCo2Pvt>(); \
66 codeToCall; \
67 break; \
68 } \
69 case OilPvtApproach::NoOilPvt: \
70 throw std::logic_error("Not implemented: Oil PVT of this deck!"); \
71 } \
72
73enum class OilPvtApproach {
80};
81
94template <class Scalar, bool enableThermal = true>
96{
97public:
99 {
100 approach_ = OilPvtApproach::NoOilPvt;
101 realOilPvt_ = nullptr;
102 }
103
105 : approach_(approach)
106 , realOilPvt_(realOilPvt)
107 { }
108
110 {
111 *this = data;
112 }
113
115 {
116 switch (approach_) {
118 delete &getRealPvt<OilPvtApproach::LiveOilPvt>();
119 break;
120 }
122 delete &getRealPvt<OilPvtApproach::DeadOilPvt>();
123 break;
124 }
126 delete &getRealPvt<OilPvtApproach::ConstantCompressibilityOilPvt>();
127 break;
128 }
130 delete &getRealPvt<OilPvtApproach::ThermalOilPvt>();
131 break;
132 }
134 delete &getRealPvt<OilPvtApproach::BrineCo2Pvt>();
135 break;
136 }
137
139 break;
140 }
141 }
142
143#if HAVE_ECL_INPUT
149 void initFromState(const EclipseState& eclState, const Schedule& schedule)
150 {
151 if (!eclState.runspec().phases().active(Phase::OIL))
152 return;
153 // TODO move the BrineCo2 approach to the waterPvtMultiplexer
154 // when a proper gas-water simulator is supported
155 if (eclState.runspec().co2Storage())
157 else if (enableThermal && eclState.getSimulationConfig().isThermal())
159 else if (!eclState.getTableManager().getPvcdoTable().empty())
161 else if (eclState.getTableManager().hasTables("PVDO"))
163 else if (!eclState.getTableManager().getPvtoTables().empty())
165
166 OPM_OIL_PVT_MULTIPLEXER_CALL(pvtImpl.initFromState(eclState, schedule));
167 }
168#endif // HAVE_ECL_INPUT
169
170
171 void initEnd()
172 { OPM_OIL_PVT_MULTIPLEXER_CALL(pvtImpl.initEnd()); }
173
177 unsigned numRegions() const
178 { OPM_OIL_PVT_MULTIPLEXER_CALL(return pvtImpl.numRegions()); return 1; }
179
183 const Scalar oilReferenceDensity(unsigned regionIdx) const
184 { OPM_OIL_PVT_MULTIPLEXER_CALL(return pvtImpl.oilReferenceDensity(regionIdx)); return 700.; }
185
189 template <class Evaluation>
190 Evaluation internalEnergy(unsigned regionIdx,
191 const Evaluation& temperature,
192 const Evaluation& pressure,
193 const Evaluation& Rs) const
194 { OPM_OIL_PVT_MULTIPLEXER_CALL(return pvtImpl.internalEnergy(regionIdx, temperature, pressure, Rs)); return 0; }
195
199 template <class Evaluation>
200 Evaluation viscosity(unsigned regionIdx,
201 const Evaluation& temperature,
202 const Evaluation& pressure,
203 const Evaluation& Rs) const
204 { OPM_OIL_PVT_MULTIPLEXER_CALL(return pvtImpl.viscosity(regionIdx, temperature, pressure, Rs)); return 0; }
205
209 template <class Evaluation>
210 Evaluation saturatedViscosity(unsigned regionIdx,
211 const Evaluation& temperature,
212 const Evaluation& pressure) const
213 { OPM_OIL_PVT_MULTIPLEXER_CALL(return pvtImpl.saturatedViscosity(regionIdx, temperature, pressure)); return 0; }
214
218 template <class Evaluation>
219 Evaluation inverseFormationVolumeFactor(unsigned regionIdx,
220 const Evaluation& temperature,
221 const Evaluation& pressure,
222 const Evaluation& Rs) const
223 { OPM_OIL_PVT_MULTIPLEXER_CALL(return pvtImpl.inverseFormationVolumeFactor(regionIdx, temperature, pressure, Rs)); return 0; }
224
228 template <class Evaluation>
229 Evaluation saturatedInverseFormationVolumeFactor(unsigned regionIdx,
230 const Evaluation& temperature,
231 const Evaluation& pressure) const
232 { OPM_OIL_PVT_MULTIPLEXER_CALL(return pvtImpl.saturatedInverseFormationVolumeFactor(regionIdx, temperature, pressure)); return 0; }
233
237 template <class Evaluation>
238 Evaluation saturatedGasDissolutionFactor(unsigned regionIdx,
239 const Evaluation& temperature,
240 const Evaluation& pressure) const
241 { OPM_OIL_PVT_MULTIPLEXER_CALL(return pvtImpl.saturatedGasDissolutionFactor(regionIdx, temperature, pressure)); return 0; }
242
246 template <class Evaluation>
247 Evaluation saturatedGasDissolutionFactor(unsigned regionIdx,
248 const Evaluation& temperature,
249 const Evaluation& pressure,
250 const Evaluation& oilSaturation,
251 const Evaluation& maxOilSaturation) const
252 { OPM_OIL_PVT_MULTIPLEXER_CALL(return pvtImpl.saturatedGasDissolutionFactor(regionIdx, temperature, pressure, oilSaturation, maxOilSaturation)); return 0; }
253
261 template <class Evaluation>
262 Evaluation saturationPressure(unsigned regionIdx,
263 const Evaluation& temperature,
264 const Evaluation& Rs) const
265 { OPM_OIL_PVT_MULTIPLEXER_CALL(return pvtImpl.saturationPressure(regionIdx, temperature, Rs)); return 0; }
266
270 template <class Evaluation>
271 Evaluation diffusionCoefficient(const Evaluation& temperature,
272 const Evaluation& pressure,
273 unsigned compIdx) const
274 {
275 OPM_OIL_PVT_MULTIPLEXER_CALL(return pvtImpl.diffusionCoefficient(temperature, pressure, compIdx)); return 0;
276 }
277
279 {
280 switch (appr) {
282 realOilPvt_ = new LiveOilPvt<Scalar>;
283 break;
284
286 realOilPvt_ = new DeadOilPvt<Scalar>;
287 break;
288
291 break;
292
294 realOilPvt_ = new OilPvtThermal<Scalar>;
295 break;
296
298 realOilPvt_ = new BrineCo2Pvt<Scalar>;
299 break;
300
302 throw std::logic_error("Not implemented: Oil PVT of this deck!");
303 }
304
305 approach_ = appr;
306 }
307
314 { return approach_; }
315
316 // get the concrete parameter object for the oil phase
317 template <OilPvtApproach approachV>
318 typename std::enable_if<approachV == OilPvtApproach::LiveOilPvt, LiveOilPvt<Scalar> >::type& getRealPvt()
319 {
320 assert(approach() == approachV);
321 return *static_cast<LiveOilPvt<Scalar>* >(realOilPvt_);
322 }
323
324 template <OilPvtApproach approachV>
325 typename std::enable_if<approachV == OilPvtApproach::LiveOilPvt, const LiveOilPvt<Scalar> >::type& getRealPvt() const
326 {
327 assert(approach() == approachV);
328 return *static_cast<LiveOilPvt<Scalar>* >(realOilPvt_);
329 }
330
331 template <OilPvtApproach approachV>
332 typename std::enable_if<approachV == OilPvtApproach::DeadOilPvt, DeadOilPvt<Scalar> >::type& getRealPvt()
333 {
334 assert(approach() == approachV);
335 return *static_cast<DeadOilPvt<Scalar>* >(realOilPvt_);
336 }
337
338 template <OilPvtApproach approachV>
339 typename std::enable_if<approachV == OilPvtApproach::DeadOilPvt, const DeadOilPvt<Scalar> >::type& getRealPvt() const
340 {
341 assert(approach() == approachV);
342 return *static_cast<DeadOilPvt<Scalar>* >(realOilPvt_);
343 }
344
345 template <OilPvtApproach approachV>
346 typename std::enable_if<approachV == OilPvtApproach::ConstantCompressibilityOilPvt, ConstantCompressibilityOilPvt<Scalar> >::type& getRealPvt()
347 {
348 assert(approach() == approachV);
349 return *static_cast<ConstantCompressibilityOilPvt<Scalar>* >(realOilPvt_);
350 }
351
352 template <OilPvtApproach approachV>
353 typename std::enable_if<approachV == OilPvtApproach::ConstantCompressibilityOilPvt, const ConstantCompressibilityOilPvt<Scalar> >::type& getRealPvt() const
354 {
355 assert(approach() == approachV);
356 return *static_cast<ConstantCompressibilityOilPvt<Scalar>* >(realOilPvt_);
357 }
358
359 template <OilPvtApproach approachV>
360 typename std::enable_if<approachV == OilPvtApproach::ThermalOilPvt, OilPvtThermal<Scalar> >::type& getRealPvt()
361 {
362 assert(approach() == approachV);
363 return *static_cast<OilPvtThermal<Scalar>* >(realOilPvt_);
364 }
365
366 template <OilPvtApproach approachV>
367 typename std::enable_if<approachV == OilPvtApproach::ThermalOilPvt, const OilPvtThermal<Scalar> >::type& getRealPvt() const
368 {
369 assert(approach() == approachV);
370 return *static_cast<const OilPvtThermal<Scalar>* >(realOilPvt_);
371 }
372
373 template <OilPvtApproach approachV>
374 typename std::enable_if<approachV == OilPvtApproach::BrineCo2Pvt, BrineCo2Pvt<Scalar> >::type& getRealPvt()
375 {
376 assert(approach() == approachV);
377 return *static_cast<BrineCo2Pvt<Scalar>* >(realOilPvt_);
378 }
379
380 template <OilPvtApproach approachV>
381 typename std::enable_if<approachV == OilPvtApproach::BrineCo2Pvt, const BrineCo2Pvt<Scalar> >::type& getRealPvt() const
382 {
383 assert(approach() == approachV);
384 return *static_cast<const BrineCo2Pvt<Scalar>* >(realOilPvt_);
385 }
386
387 const void* realOilPvt() const { return realOilPvt_; }
388
390 {
391 if (this->approach() != data.approach())
392 return false;
393
394 switch (approach_) {
396 return *static_cast<const ConstantCompressibilityOilPvt<Scalar>*>(realOilPvt_) ==
397 *static_cast<const ConstantCompressibilityOilPvt<Scalar>*>(data.realOilPvt_);
399 return *static_cast<const DeadOilPvt<Scalar>*>(realOilPvt_) ==
400 *static_cast<const DeadOilPvt<Scalar>*>(data.realOilPvt_);
402 return *static_cast<const LiveOilPvt<Scalar>*>(realOilPvt_) ==
403 *static_cast<const LiveOilPvt<Scalar>*>(data.realOilPvt_);
405 return *static_cast<const OilPvtThermal<Scalar>*>(realOilPvt_) ==
406 *static_cast<const OilPvtThermal<Scalar>*>(data.realOilPvt_);
408 return *static_cast<const BrineCo2Pvt<Scalar>*>(realOilPvt_) ==
409 *static_cast<const BrineCo2Pvt<Scalar>*>(data.realOilPvt_);
410 default:
411 return true;
412 }
413 }
414
416 {
417 approach_ = data.approach_;
418 switch (approach_) {
420 realOilPvt_ = new ConstantCompressibilityOilPvt<Scalar>(*static_cast<const ConstantCompressibilityOilPvt<Scalar>*>(data.realOilPvt_));
421 break;
423 realOilPvt_ = new DeadOilPvt<Scalar>(*static_cast<const DeadOilPvt<Scalar>*>(data.realOilPvt_));
424 break;
426 realOilPvt_ = new LiveOilPvt<Scalar>(*static_cast<const LiveOilPvt<Scalar>*>(data.realOilPvt_));
427 break;
429 realOilPvt_ = new OilPvtThermal<Scalar>(*static_cast<const OilPvtThermal<Scalar>*>(data.realOilPvt_));
430 break;
432 realOilPvt_ = new BrineCo2Pvt<Scalar>(*static_cast<const BrineCo2Pvt<Scalar>*>(data.realOilPvt_));
433 break;
434 default:
435 break;
436 }
437
438 return *this;
439 }
440
441private:
442 OilPvtApproach approach_;
443 void* realOilPvt_;
444};
445
446#undef OPM_OIL_PVT_MULTIPLEXER_CALL
447
448} // namespace Opm
449
450#endif
#define OPM_OIL_PVT_MULTIPLEXER_CALL(codeToCall)
Definition: OilPvtMultiplexer.hpp:42
This class represents the Pressure-Volume-Temperature relations of the liquid phase for a CO2-Brine s...
Definition: BrineCo2Pvt.hpp:57
This class represents the Pressure-Volume-Temperature relations of the oil phase without dissolved ga...
Definition: ConstantCompressibilityOilPvt.hpp:42
This class represents the Pressure-Volume-Temperature relations of the oil phase without dissolved ga...
Definition: DeadOilPvt.hpp:45
This class represents the Pressure-Volume-Temperature relations of the oil phas with dissolved gas.
Definition: LiveOilPvt.hpp:52
This class represents the Pressure-Volume-Temperature relations of the oil phase in the black-oil mod...
Definition: OilPvtMultiplexer.hpp:96
std::enable_if< approachV==OilPvtApproach::DeadOilPvt, constDeadOilPvt< Scalar > >::type & getRealPvt() const
Definition: OilPvtMultiplexer.hpp:339
unsigned numRegions() const
Return the number of PVT regions which are considered by this PVT-object.
Definition: OilPvtMultiplexer.hpp:177
const void * realOilPvt() const
Definition: OilPvtMultiplexer.hpp:387
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:190
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:219
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:210
std::enable_if< approachV==OilPvtApproach::LiveOilPvt, LiveOilPvt< Scalar > >::type & getRealPvt()
Definition: OilPvtMultiplexer.hpp:318
OilPvtMultiplexer(OilPvtApproach approach, void *realOilPvt)
Definition: OilPvtMultiplexer.hpp:104
OilPvtApproach approach() const
Returns the concrete approach for calculating the PVT relations.
Definition: OilPvtMultiplexer.hpp:313
std::enable_if< approachV==OilPvtApproach::LiveOilPvt, constLiveOilPvt< Scalar > >::type & getRealPvt() const
Definition: OilPvtMultiplexer.hpp:325
OilPvtMultiplexer()
Definition: OilPvtMultiplexer.hpp:98
Evaluation saturatedInverseFormationVolumeFactor(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure) const
Returns the formation volume factor [-] of the fluid phase.
Definition: OilPvtMultiplexer.hpp:229
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:238
OilPvtMultiplexer< Scalar, enableThermal > & operator=(const OilPvtMultiplexer< Scalar, enableThermal > &data)
Definition: OilPvtMultiplexer.hpp:415
std::enable_if< approachV==OilPvtApproach::BrineCo2Pvt, constBrineCo2Pvt< Scalar > >::type & getRealPvt() const
Definition: OilPvtMultiplexer.hpp:381
std::enable_if< approachV==OilPvtApproach::DeadOilPvt, DeadOilPvt< Scalar > >::type & getRealPvt()
Definition: OilPvtMultiplexer.hpp:332
OilPvtMultiplexer(const OilPvtMultiplexer< Scalar, enableThermal > &data)
Definition: OilPvtMultiplexer.hpp:109
void setApproach(OilPvtApproach appr)
Definition: OilPvtMultiplexer.hpp:278
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:200
void initEnd()
Definition: OilPvtMultiplexer.hpp:171
std::enable_if< approachV==OilPvtApproach::ThermalOilPvt, OilPvtThermal< Scalar > >::type & getRealPvt()
Definition: OilPvtMultiplexer.hpp:360
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:262
bool operator==(const OilPvtMultiplexer< Scalar, enableThermal > &data) const
Definition: OilPvtMultiplexer.hpp:389
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:247
std::enable_if< approachV==OilPvtApproach::BrineCo2Pvt, BrineCo2Pvt< Scalar > >::type & getRealPvt()
Definition: OilPvtMultiplexer.hpp:374
std::enable_if< approachV==OilPvtApproach::ConstantCompressibilityOilPvt, ConstantCompressibilityOilPvt< Scalar > >::type & getRealPvt()
Definition: OilPvtMultiplexer.hpp:346
~OilPvtMultiplexer()
Definition: OilPvtMultiplexer.hpp:114
const Scalar oilReferenceDensity(unsigned regionIdx) const
Return the reference density which are considered by this PVT-object.
Definition: OilPvtMultiplexer.hpp:183
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:271
std::enable_if< approachV==OilPvtApproach::ConstantCompressibilityOilPvt, constConstantCompressibilityOilPvt< Scalar > >::type & getRealPvt() const
Definition: OilPvtMultiplexer.hpp:353
std::enable_if< approachV==OilPvtApproach::ThermalOilPvt, constOilPvtThermal< Scalar > >::type & getRealPvt() const
Definition: OilPvtMultiplexer.hpp:367
This class implements temperature dependence of the PVT properties of oil.
Definition: OilPvtThermal.hpp:50
Definition: Air_Mesitylene.hpp:34
OilPvtApproach
Definition: OilPvtMultiplexer.hpp:73