WaterPvtMultiplexer.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_WATER_PVT_MULTIPLEXER_HPP
28#define OPM_WATER_PVT_MULTIPLEXER_HPP
29
32#include "WaterPvtThermal.hpp"
33
34#if HAVE_ECL_INPUT
35#include <opm/input/eclipse/EclipseState/EclipseState.hpp>
36#include <opm/input/eclipse/EclipseState/Runspec.hpp>
37#endif
38
39#define OPM_WATER_PVT_MULTIPLEXER_CALL(codeToCall) \
40 switch (approach_) { \
41 case WaterPvtApproach::ConstantCompressibilityWaterPvt: { \
42 auto& pvtImpl = getRealPvt<WaterPvtApproach::ConstantCompressibilityWaterPvt>(); \
43 codeToCall; \
44 break; \
45 } \
46 case WaterPvtApproach::ConstantCompressibilityBrinePvt: { \
47 auto& pvtImpl = getRealPvt<WaterPvtApproach::ConstantCompressibilityBrinePvt>(); \
48 codeToCall; \
49 break; \
50 } \
51 case WaterPvtApproach::ThermalWaterPvt: { \
52 auto& pvtImpl = getRealPvt<WaterPvtApproach::ThermalWaterPvt>(); \
53 codeToCall; \
54 break; \
55 } \
56 case WaterPvtApproach::NoWaterPvt: \
57 throw std::logic_error("Not implemented: Water PVT of this deck!"); \
58 }
59
60namespace Opm {
61
62enum class WaterPvtApproach {
67};
68
73template <class Scalar, bool enableThermal = true, bool enableBrine = true>
75{
76public:
78 {
80 realWaterPvt_ = nullptr;
81 }
82
84 : approach_(approach)
85 , realWaterPvt_(realWaterPvt)
86 { }
87
89 {
90 *this = data;
91 }
92
94 {
95 switch (approach_) {
97 delete &getRealPvt<WaterPvtApproach::ConstantCompressibilityWaterPvt>();
98 break;
99 }
101 delete &getRealPvt<WaterPvtApproach::ConstantCompressibilityBrinePvt>();
102 break;
103 }
105 delete &getRealPvt<WaterPvtApproach::ThermalWaterPvt>();
106 break;
107 }
109 break;
110 }
111 }
112
113#if HAVE_ECL_INPUT
119 void initFromState(const EclipseState& eclState, const Schedule& schedule)
120 {
121 if (!eclState.runspec().phases().active(Phase::WATER))
122 return;
123
124 if (enableThermal && eclState.getSimulationConfig().isThermal())
126 else if (!eclState.getTableManager().getPvtwTable().empty())
128 else if (enableBrine && !eclState.getTableManager().getPvtwSaltTables().empty())
130
131 OPM_WATER_PVT_MULTIPLEXER_CALL(pvtImpl.initFromState(eclState, schedule));
132 }
133#endif // HAVE_ECL_INPUT
134
135 void initEnd()
136 { OPM_WATER_PVT_MULTIPLEXER_CALL(pvtImpl.initEnd()); }
137
141 unsigned numRegions() const
142 { OPM_WATER_PVT_MULTIPLEXER_CALL(return pvtImpl.numRegions()); return 1; }
143
147 const Scalar waterReferenceDensity(unsigned regionIdx)
148 { OPM_WATER_PVT_MULTIPLEXER_CALL(return pvtImpl.waterReferenceDensity(regionIdx)); return 1000.; }
149
153 template <class Evaluation>
154 Evaluation internalEnergy(unsigned regionIdx,
155 const Evaluation& temperature,
156 const Evaluation& pressure,
157 const Evaluation& saltconcentration) const
158 { OPM_WATER_PVT_MULTIPLEXER_CALL(return pvtImpl.internalEnergy(regionIdx, temperature, pressure, saltconcentration)); return 0; }
159
163 template <class Evaluation>
164 Evaluation viscosity(unsigned regionIdx,
165 const Evaluation& temperature,
166 const Evaluation& pressure,
167 const Evaluation& saltconcentration) const
168 {
169 OPM_WATER_PVT_MULTIPLEXER_CALL(return pvtImpl.viscosity(regionIdx, temperature, pressure, saltconcentration));
170 return 0;
171 }
172
176 template <class Evaluation>
177 Evaluation inverseFormationVolumeFactor(unsigned regionIdx,
178 const Evaluation& temperature,
179 const Evaluation& pressure,
180 const Evaluation& saltconcentration) const
181 { OPM_WATER_PVT_MULTIPLEXER_CALL(return pvtImpl.inverseFormationVolumeFactor(regionIdx, temperature, pressure, saltconcentration));
182 return 0;
183 }
184
186 {
187 switch (appr) {
189 realWaterPvt_ = new ConstantCompressibilityWaterPvt<Scalar>;
190 break;
191
193 realWaterPvt_ = new ConstantCompressibilityBrinePvt<Scalar>;
194 break;
195
197 realWaterPvt_ = new WaterPvtThermal<Scalar, enableBrine>;
198 break;
199
201 throw std::logic_error("Not implemented: Water PVT of this deck!");
202 }
203
204 approach_ = appr;
205 }
206
213 { return approach_; }
214
215 // get the concrete parameter object for the water phase
216 template <WaterPvtApproach approachV>
217 typename std::enable_if<approachV == WaterPvtApproach::ConstantCompressibilityWaterPvt, ConstantCompressibilityWaterPvt<Scalar> >::type& getRealPvt()
218 {
219 assert(approach() == approachV);
220 return *static_cast<ConstantCompressibilityWaterPvt<Scalar>* >(realWaterPvt_);
221 }
222
223 template <WaterPvtApproach approachV>
224 typename std::enable_if<approachV == WaterPvtApproach::ConstantCompressibilityWaterPvt, const ConstantCompressibilityWaterPvt<Scalar> >::type& getRealPvt() const
225 {
226 assert(approach() == approachV);
227 return *static_cast<ConstantCompressibilityWaterPvt<Scalar>* >(realWaterPvt_);
228 }
229
230 template <WaterPvtApproach approachV>
231 typename std::enable_if<approachV == WaterPvtApproach::ConstantCompressibilityBrinePvt, ConstantCompressibilityBrinePvt<Scalar> >::type& getRealPvt()
232 {
233 assert(approach() == approachV);
234 return *static_cast<ConstantCompressibilityBrinePvt<Scalar>* >(realWaterPvt_);
235 }
236
237 template <WaterPvtApproach approachV>
238 typename std::enable_if<approachV == WaterPvtApproach::ConstantCompressibilityBrinePvt, const ConstantCompressibilityBrinePvt<Scalar> >::type& getRealPvt() const
239 {
240 assert(approach() == approachV);
241 return *static_cast<ConstantCompressibilityBrinePvt<Scalar>* >(realWaterPvt_);
242 }
243
244 template <WaterPvtApproach approachV>
245 typename std::enable_if<approachV == WaterPvtApproach::ThermalWaterPvt, WaterPvtThermal<Scalar, enableBrine> >::type& getRealPvt()
246 {
247 assert(approach() == approachV);
248 return *static_cast<WaterPvtThermal<Scalar, enableBrine>* >(realWaterPvt_);
249 }
250
251 template <WaterPvtApproach approachV>
252 typename std::enable_if<approachV == WaterPvtApproach::ThermalWaterPvt, const WaterPvtThermal<Scalar, enableBrine> >::type& getRealPvt() const
253 {
254 assert(approach() == approachV);
255 return *static_cast<WaterPvtThermal<Scalar, enableBrine>* >(realWaterPvt_);
256 }
257
258 const void* realWaterPvt() const { return realWaterPvt_; }
259
261 {
262 if (this->approach() != data.approach())
263 return false;
264
265 switch (approach_) {
267 return *static_cast<const ConstantCompressibilityWaterPvt<Scalar>*>(realWaterPvt_) ==
268 *static_cast<const ConstantCompressibilityWaterPvt<Scalar>*>(data.realWaterPvt_);
270 return *static_cast<const ConstantCompressibilityBrinePvt<Scalar>*>(realWaterPvt_) ==
271 *static_cast<const ConstantCompressibilityBrinePvt<Scalar>*>(data.realWaterPvt_);
273 return *static_cast<const WaterPvtThermal<Scalar, enableBrine>*>(realWaterPvt_) ==
274 *static_cast<const WaterPvtThermal<Scalar, enableBrine>*>(data.realWaterPvt_);
275 default:
276 return true;
277 }
278 }
279
281 {
282 approach_ = data.approach_;
283 switch (approach_) {
285 realWaterPvt_ = new ConstantCompressibilityWaterPvt<Scalar>(*static_cast<const ConstantCompressibilityWaterPvt<Scalar>*>(data.realWaterPvt_));
286 break;
288 realWaterPvt_ = new ConstantCompressibilityBrinePvt<Scalar>(*static_cast<const ConstantCompressibilityBrinePvt<Scalar>*>(data.realWaterPvt_));
289 break;
291 realWaterPvt_ = new WaterPvtThermal<Scalar, enableBrine>(*static_cast<const WaterPvtThermal<Scalar, enableBrine>*>(data.realWaterPvt_));
292 break;
293 default:
294 break;
295 }
296
297 return *this;
298 }
299
300private:
301 WaterPvtApproach approach_;
302 void* realWaterPvt_;
303};
304
305#undef OPM_WATER_PVT_MULTIPLEXER_CALL
306
307} // namespace Opm
308
309#endif
#define OPM_WATER_PVT_MULTIPLEXER_CALL(codeToCall)
Definition: WaterPvtMultiplexer.hpp:39
This class represents the Pressure-Volume-Temperature relations of the gas phase without vaporized oi...
Definition: ConstantCompressibilityBrinePvt.hpp:49
This class represents the Pressure-Volume-Temperature relations of the gas phase without vaporized oi...
Definition: ConstantCompressibilityWaterPvt.hpp:43
This class represents the Pressure-Volume-Temperature relations of the water phase in the black-oil m...
Definition: WaterPvtMultiplexer.hpp:75
std::enable_if< approachV==WaterPvtApproach::ConstantCompressibilityWaterPvt, constConstantCompressibilityWaterPvt< Scalar > >::type & getRealPvt() const
Definition: WaterPvtMultiplexer.hpp:224
std::enable_if< approachV==WaterPvtApproach::ConstantCompressibilityWaterPvt, ConstantCompressibilityWaterPvt< Scalar > >::type & getRealPvt()
Definition: WaterPvtMultiplexer.hpp:217
WaterPvtMultiplexer(const WaterPvtMultiplexer< Scalar, enableThermal, enableBrine > &data)
Definition: WaterPvtMultiplexer.hpp:88
unsigned numRegions() const
Return the number of PVT regions which are considered by this PVT-object.
Definition: WaterPvtMultiplexer.hpp:141
const Scalar waterReferenceDensity(unsigned regionIdx)
Return the reference density which are considered by this PVT-object.
Definition: WaterPvtMultiplexer.hpp:147
Evaluation inverseFormationVolumeFactor(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure, const Evaluation &saltconcentration) const
Returns the formation volume factor [-] of the fluid phase.
Definition: WaterPvtMultiplexer.hpp:177
WaterPvtMultiplexer(WaterPvtApproach approach, void *realWaterPvt)
Definition: WaterPvtMultiplexer.hpp:83
std::enable_if< approachV==WaterPvtApproach::ThermalWaterPvt, constWaterPvtThermal< Scalar, enableBrine > >::type & getRealPvt() const
Definition: WaterPvtMultiplexer.hpp:252
void initEnd()
Definition: WaterPvtMultiplexer.hpp:135
std::enable_if< approachV==WaterPvtApproach::ConstantCompressibilityBrinePvt, constConstantCompressibilityBrinePvt< Scalar > >::type & getRealPvt() const
Definition: WaterPvtMultiplexer.hpp:238
WaterPvtMultiplexer()
Definition: WaterPvtMultiplexer.hpp:77
bool operator==(const WaterPvtMultiplexer< Scalar, enableThermal, enableBrine > &data) const
Definition: WaterPvtMultiplexer.hpp:260
std::enable_if< approachV==WaterPvtApproach::ThermalWaterPvt, WaterPvtThermal< Scalar, enableBrine > >::type & getRealPvt()
Definition: WaterPvtMultiplexer.hpp:245
~WaterPvtMultiplexer()
Definition: WaterPvtMultiplexer.hpp:93
Evaluation viscosity(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:164
void setApproach(WaterPvtApproach appr)
Definition: WaterPvtMultiplexer.hpp:185
WaterPvtApproach approach() const
Returns the concrete approach for calculating the PVT relations.
Definition: WaterPvtMultiplexer.hpp:212
Evaluation internalEnergy(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure, const Evaluation &saltconcentration) const
Returns the specific enthalpy [J/kg] of gas given a set of parameters.
Definition: WaterPvtMultiplexer.hpp:154
std::enable_if< approachV==WaterPvtApproach::ConstantCompressibilityBrinePvt, ConstantCompressibilityBrinePvt< Scalar > >::type & getRealPvt()
Definition: WaterPvtMultiplexer.hpp:231
WaterPvtMultiplexer< Scalar, enableThermal, enableBrine > & operator=(const WaterPvtMultiplexer< Scalar, enableThermal, enableBrine > &data)
Definition: WaterPvtMultiplexer.hpp:280
const void * realWaterPvt() const
Definition: WaterPvtMultiplexer.hpp:258
This class implements temperature dependence of the PVT properties of water.
Definition: WaterPvtThermal.hpp:50
Definition: Air_Mesitylene.hpp:34
WaterPvtApproach
Definition: WaterPvtMultiplexer.hpp:62