opm-common
Loading...
Searching...
No Matches
GenericOilGasWaterFluidSystem.hpp
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 Copyright 2023 SINTEF Digital, Mathematics and Cybernetics.
5
6 This file is part of the Open Porous Media project (OPM).
7
8 OPM is free software: you can redistribute it and/or modify
9 it under the terms of the GNU General Public License as published by
10 the Free Software Foundation, either version 2 of the License, or
11 (at your option) any later version.
12
13 OPM is distributed in the hope that it will be useful,
14 but WITHOUT ANY WARRANTY; without even the implied warranty of
15 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16 GNU General Public License for more details.
17
18 You should have received a copy of the GNU General Public License
19 along with OPM. If not, see <http://www.gnu.org/licenses/>.
20
21 Consult the COPYING file in the top-level source directory of this
22 module for the precise wording of the license and the list of
23 copyright holders.
24*/
25
26#ifndef OPM_GENERIC_OIL_GAS_WATER_FLUIDSYSTEM_HPP
27#define OPM_GENERIC_OIL_GAS_WATER_FLUIDSYSTEM_HPP
28
29#include <opm/common/OpmLog/OpmLog.hpp>
30
31#include <opm/input/eclipse/EclipseState/EclipseState.hpp>
32#include <opm/input/eclipse/Schedule/Schedule.hpp>
33
34#include <opm/material/eos/CubicEOS.hpp>
37#include <opm/material/fluidsystems/PTFlashParameterCache.hpp> // TODO: this is something else need to check
39
40#include <cassert>
41#include <cstddef>
42#include <string>
43#include <string_view>
44
45#include <fmt/format.h>
46
47namespace Opm {
48
57 template<class Scalar, int NumComp, bool enableWater>
59 : public BaseFluidSystem<Scalar, GenericOilGasWaterFluidSystem<Scalar, NumComp, enableWater> > {
60 public:
61 // TODO: I do not think these should be constant in fluidsystem, will try to make it non-constant later
62 static constexpr bool waterEnabled = enableWater;
63 static constexpr int numPhases = enableWater ? 3 : 2;
64 static constexpr int numComponents = NumComp;
65 static constexpr int numMisciblePhases = 2;
66 // \Note: not totally sure when we should distinguish numMiscibleComponents and numComponents.
67 // Possibly when with a dummy phase like water?
68 static constexpr int numMiscibleComponents = NumComp;
69 // TODO: phase location should be more general
70 static constexpr int oilPhaseIdx = 0;
71 static constexpr int gasPhaseIdx = 1;
72 static constexpr int waterPhaseIdx = 2;
73
74 static constexpr int waterCompIdx = -1;
75 static constexpr int oilCompIdx = 0;
76 static constexpr int gasCompIdx = 1;
77 static constexpr int compositionSwitchIdx = -1; // equil initializer
78
79 template <class ValueType>
83 using WaterPvt = WaterPvtMultiplexer<Scalar>;
84
85 struct ComponentParam {
86 std::string name;
87 Scalar molar_mass; // unit: g/mol
88 Scalar critic_temp; // unit: K
89 Scalar critic_pres; // unit: parscal
90 Scalar critic_vol; // unit: m^3/kmol
91 Scalar acentric_factor; // unit: dimension less
92
93 ComponentParam(const std::string_view name_, const Scalar molar_mass_, const Scalar critic_temp_,
94 const Scalar critic_pres_, const Scalar critic_vol_, const Scalar acentric_factor_)
95 : name(name_),
96 molar_mass(molar_mass_),
97 critic_temp(critic_temp_),
98 critic_pres(critic_pres_),
99 critic_vol(critic_vol_),
100 acentric_factor(acentric_factor_)
101 {}
102 };
103
104 static bool phaseIsActive(unsigned phaseIdx)
105 {
106 if constexpr (enableWater) {
107 assert(phaseIdx < numPhases);
108 return true;
109 }
110 else {
111 if (phaseIdx == waterPhaseIdx)
112 return false;
113 assert(phaseIdx < numPhases);
114 return true;
115 }
116 }
117
118 template<typename Param>
119 static void addComponent(const Param& param)
120 {
121 // Check if the current size is less than the maximum allowed components.
122 if (component_param_.size() < numComponents) {
123 component_param_.push_back(param);
124 } else {
125 // Adding another component would exceed the limit.
126 const std::string msg = fmt::format("The fluid system has reached its maximum capacity of {} components,"
127 "the component '{}' will not be added.", NumComp, param.name);
128 OpmLog::note(msg);
129 // Optionally, throw an exception?
130 }
131 }
132
136 static void initFromState(const EclipseState& eclState, const Schedule& schedule)
137 {
138 // TODO: we are not considering the EOS region for now
139 const auto& comp_config = eclState.compositionalConfig();
140 // how should we utilize the numComps from the CompositionalConfig?
142 const std::size_t num_comps = comp_config.numComps();
143 // const std::size_t num_eos_region = comp_config.
144 assert(num_comps == NumComp);
145 const auto& names = comp_config.compName();
146 // const auto& eos_type = comp_config.eosType(0);
147 const auto& molar_weight = comp_config.molecularWeights(0);
148 const auto& acentric_factor = comp_config.acentricFactors(0);
149 const auto& critic_pressure = comp_config.criticalPressure(0);
150 const auto& critic_temp = comp_config.criticalTemperature(0);
151 const auto& critic_volume = comp_config.criticalVolume(0);
152 FluidSystem::init();
153 using CompParm = typename FluidSystem::ComponentParam;
154 for (std::size_t c = 0; c < num_comps; ++c) {
155 // we use m^3/kmol for the critic volume in the flash calculation, so we multiply 1.e3 for the critic volume
156 FluidSystem::addComponent(CompParm{names[c], molar_weight[c], critic_temp[c], critic_pressure[c],
157 critic_volume[c] * 1.e3, acentric_factor[c]});
158 }
159 FluidSystem::printComponentParams();
160 interaction_coefficients_ = comp_config.binaryInteractionCoefficient(0);
161
162 // Init. water pvt from deck
163 waterPvt_->initFromState(eclState, schedule);
164
165 }
166
167 static void init()
168 {
169 waterPvt_ = std::make_shared<WaterPvt>();
170 component_param_.reserve(numComponents);
171 }
172
176 static void setWaterPvt(std::shared_ptr<WaterPvt> pvtObj)
177 { waterPvt_ = std::move(pvtObj); }
178
184 static Scalar acentricFactor(unsigned compIdx)
185 {
186 assert(isConsistent());
187 assert(compIdx < numComponents);
188
189 return component_param_[compIdx].acentric_factor;
190 }
191
196 static Scalar criticalTemperature(unsigned compIdx)
197 {
198 assert(isConsistent());
199 assert(compIdx < numComponents);
200
201 return component_param_[compIdx].critic_temp;
202 }
203
208 static Scalar criticalPressure(unsigned compIdx)
209 {
210 assert(isConsistent());
211 assert(compIdx < numComponents);
212
213 return component_param_[compIdx].critic_pres;
214 }
215
220 static Scalar criticalVolume(unsigned compIdx)
221 {
222 assert(isConsistent());
223 assert(compIdx < numComponents);
224
225 return component_param_[compIdx].critic_vol;
226 }
227
229 static Scalar molarMass(unsigned compIdx)
230 {
231 assert(isConsistent());
232 assert(compIdx < numComponents);
233
234 return component_param_[compIdx].molar_mass;
235 }
236
241 static Scalar interactionCoefficient(unsigned comp1Idx, unsigned comp2Idx)
242 {
243 assert(isConsistent());
244 assert(comp1Idx < numComponents);
245 assert(comp2Idx < numComponents);
246 if (interaction_coefficients_.empty() || comp2Idx == comp1Idx) {
247 return 0.0;
248 }
249 // make sure row is the bigger value compared to column number
250 const auto [column, row] = std::minmax(comp1Idx, comp2Idx);
251 const unsigned index = (row * (row - 1) / 2 + column); // it is the current understanding
252 return interaction_coefficients_[index];
253 }
254
256 static std::string_view phaseName(unsigned phaseIdx)
257 {
258 static const std::string_view name[] = {"o", // oleic phase
259 "g", // gas phase
260 "w"}; // aqueous phase
261
262 assert(phaseIdx < numPhases);
263 return name[phaseIdx];
264 }
265
267 static std::string_view componentName(unsigned compIdx)
268 {
269 assert(isConsistent());
270 assert(compIdx < numComponents);
271
272 return component_param_[compIdx].name;
273 }
274
278 template <class FluidState, class LhsEval = typename FluidState::ValueType, class ParamCacheEval = LhsEval>
279 static LhsEval density(const FluidState& fluidState,
280 const ParameterCache<ParamCacheEval>& paramCache,
281 unsigned phaseIdx)
282 {
283 assert(isConsistent());
284 assert(phaseIdx < numPhases);
285
286 if (phaseIdx == oilPhaseIdx || phaseIdx == gasPhaseIdx) {
287 return decay<LhsEval>(fluidState.averageMolarMass(phaseIdx) / paramCache.molarVolume(phaseIdx));
288 }
289 else {
290 const LhsEval& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
291 const LhsEval& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
292 const Scalar& rho_w_ref = waterPvt_->waterReferenceDensity(0);
293 const LhsEval& bw = waterPvt_->inverseFormationVolumeFactor(0, T, p, LhsEval(0.0), LhsEval(0.0));
294 return rho_w_ref * bw;
295 }
296
297 return {};
298 }
299
301 template <class FluidState, class LhsEval = typename FluidState::ValueType, class ParamCacheEval = LhsEval>
302 static LhsEval viscosity(const FluidState& fluidState,
303 const ParameterCache<ParamCacheEval>& paramCache,
304 unsigned phaseIdx)
305 {
306 assert(isConsistent());
307 assert(phaseIdx < numPhases);
308
309 if (phaseIdx == oilPhaseIdx || phaseIdx == gasPhaseIdx) {
310 // Use LBC method to calculate viscosity
311 return decay<LhsEval>(ViscosityModel::LBC(fluidState, paramCache, phaseIdx));
312 }
313 else {
314 const LhsEval& p = decay<LhsEval>(fluidState.pressure(phaseIdx));
315 const LhsEval& T = decay<LhsEval>(fluidState.temperature(phaseIdx));
316 return waterPvt_->viscosity(0, T, p, LhsEval(0.0), LhsEval(0.0));
317 }
318 }
319
321 template <class FluidState, class LhsEval = typename FluidState::ValueType, class ParamCacheEval = LhsEval>
322 static LhsEval fugacityCoefficient(const FluidState& fluidState,
323 const ParameterCache<ParamCacheEval>& paramCache,
324 unsigned phaseIdx,
325 unsigned compIdx)
326 {
327 if (waterEnabled && phaseIdx == static_cast<unsigned int>(waterPhaseIdx))
328 return LhsEval(0.0);
329
330 assert(isConsistent());
331 assert(phaseIdx < numPhases);
332 assert(compIdx < numComponents);
333
334 return decay<LhsEval>(CubicEOS::computeFugacityCoefficient(fluidState, paramCache, phaseIdx, compIdx));
335 }
336
337 // TODO: the following interfaces are needed by function checkFluidSystem()
339 static bool isCompressible([[maybe_unused]] unsigned phaseIdx)
340 {
341 assert(phaseIdx < numPhases);
342
343 return true;
344 }
345
347 static bool isIdealMixture([[maybe_unused]] unsigned phaseIdx)
348 {
349 assert(phaseIdx < numPhases);
350
351 return false;
352 }
353
355 static bool isLiquid(unsigned phaseIdx)
356 {
357 assert(phaseIdx < numPhases);
358
359 return (phaseIdx == oilPhaseIdx);
360 }
361
363 static bool isIdealGas(unsigned phaseIdx)
364 {
365 assert(phaseIdx < numPhases);
366
367 return (phaseIdx == gasPhaseIdx);
368 }
369 // the following funcitons are needed to compile the GenericOutputBlackoilModule
370 // not implemented for this FluidSystem yet
371 template <class LhsEval>
372 static LhsEval convertXwGToxwG(const LhsEval&, unsigned)
373 {
374 assert(false && "convertXwGToxwG not implemented for GenericOilGasWaterFluidSystem!");
375 return 0.;
376 }
377 template <class LhsEval>
378 static LhsEval convertXoGToxoG(const LhsEval&, unsigned)
379 {
380 assert(false && "convertXoGToxoG not implemented for GenericOilGasWaterFluidSystem!");
381 return 0.;
382 }
383
384 template <class LhsEval>
385 static LhsEval convertxoGToXoG(const LhsEval&, unsigned)
386 {
387 assert(false && "convertxoGToXoG not implemented for GenericOilGasWaterFluidSystem!");
388 return 0.;
389 }
390
391 template <class LhsEval>
392 static LhsEval convertXgOToxgO(const LhsEval&, unsigned)
393 {
394 assert(false && "convertXgOToxgO not implemented for GenericOilGasWaterFluidSystem!");
395 return 0.;
396 }
397
398 template <class LhsEval>
399 static LhsEval convertRswToXwG(const LhsEval&, unsigned)
400 {
401 assert(false && "convertRswToXwG not implemented for GenericOilGasWaterFluidSystem!");
402 return 0.;
403 }
404
405 template <class LhsEval>
406 static LhsEval convertRvwToXgW(const LhsEval&, unsigned)
407 {
408 assert(false && "convertRvwToXgW not implemented for GenericOilGasWaterFluidSystem!");
409 return 0.;
410 }
411
412 template <class LhsEval>
413 static LhsEval convertXgWToxgW(const LhsEval&, unsigned)
414 {
415 assert(false && "convertXgWToxgW not implemented for GenericOilGasWaterFluidSystem!");
416 return 0.;
417 }
418
419 static bool enableDissolvedGas()
420 {
421 return false;
422 }
423
424 static bool enableDissolvedGasInWater()
425 {
426 return false;
427 }
428
429 static bool enableVaporizedWater()
430 {
431 return false;
432 }
433
434 static bool enableVaporizedOil()
435 {
436 return false;
437 }
438
439 private:
440 static bool isConsistent() {
441 return component_param_.size() == NumComp;
442 }
443
444 static std::vector<ComponentParam> component_param_;
445 static std::vector<Scalar> interaction_coefficients_;
446 static std::shared_ptr<WaterPvt> waterPvt_;
447
448 public:
449 static std::string printComponentParams() {
450 std::string result = "Components Information:\n";
451 for (const auto& param : component_param_) {
452 result += fmt::format("Name: {}\n", param.name);
453 result += fmt::format("Molar Mass: {} g/mol\n", param.molar_mass);
454 result += fmt::format("Critical Temperature: {} K\n", param.critic_temp);
455 result += fmt::format("Critical Pressure: {} Pascal\n", param.critic_pres);
456 result += fmt::format("Critical Volume: {} m^3/kmol\n", param.critic_vol);
457 result += fmt::format("Acentric Factor: {}\n", param.acentric_factor);
458 result += "---------------------------------\n";
459 }
460 return result;
461 }
462 };
463
464 template <class Scalar, int NumComp, bool enableWater>
465 std::vector<typename GenericOilGasWaterFluidSystem<Scalar, NumComp, enableWater>::ComponentParam>
466 GenericOilGasWaterFluidSystem<Scalar, NumComp, enableWater>::component_param_;
467
468 template <class Scalar, int NumComp, bool enableWater>
469 std::vector<Scalar>
470 GenericOilGasWaterFluidSystem<Scalar, NumComp, enableWater>::interaction_coefficients_;
471
472 template <class Scalar, int NumComp, bool enableWater>
473 std::shared_ptr<WaterPvtMultiplexer<Scalar> >
474 GenericOilGasWaterFluidSystem<Scalar, NumComp, enableWater>::waterPvt_;
475
476}
477#endif // OPM_GENERIC_OIL_GAS_WATER_FLUIDSYSTEM_HPP
The base class for all fluid systems.
Specifies the parameter cache used by the SPE-5 fluid system.
This class represents the Pressure-Volume-Temperature relations of the water phase in the black-oil m...
The base class for all fluid systems.
Definition BaseFluidSystem.hpp:43
Definition CubicEOS.hpp:34
Definition EclipseState.hpp:66
A two phase system that can contain NumComp components.
Definition GenericOilGasWaterFluidSystem.hpp:59
static bool isCompressible(unsigned phaseIdx)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition GenericOilGasWaterFluidSystem.hpp:339
static LhsEval fugacityCoefficient(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &paramCache, unsigned phaseIdx, unsigned compIdx)
Calculate the fugacity coefficient [Pa] of an individual component in a fluid phase.
Definition GenericOilGasWaterFluidSystem.hpp:322
static Scalar acentricFactor(unsigned compIdx)
The acentric factor of a component [].
Definition GenericOilGasWaterFluidSystem.hpp:184
static void initFromState(const EclipseState &eclState, const Schedule &schedule)
Initialize the fluid system using an ECL deck object.
Definition GenericOilGasWaterFluidSystem.hpp:136
static std::string_view phaseName(unsigned phaseIdx)
Return the human readable name of a fluid phase.
Definition GenericOilGasWaterFluidSystem.hpp:256
static bool isIdealGas(unsigned phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition GenericOilGasWaterFluidSystem.hpp:363
static Scalar criticalTemperature(unsigned compIdx)
Critical temperature of a component [K].
Definition GenericOilGasWaterFluidSystem.hpp:196
static LhsEval density(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &paramCache, unsigned phaseIdx)
Calculate the density [kg/m^3] of a fluid phase.
Definition GenericOilGasWaterFluidSystem.hpp:279
static Scalar interactionCoefficient(unsigned comp1Idx, unsigned comp2Idx)
Returns the interaction coefficient for two components.
Definition GenericOilGasWaterFluidSystem.hpp:241
static Scalar criticalPressure(unsigned compIdx)
Critical pressure of a component [Pa].
Definition GenericOilGasWaterFluidSystem.hpp:208
static bool isIdealMixture(unsigned phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition GenericOilGasWaterFluidSystem.hpp:347
static Scalar molarMass(unsigned compIdx)
Return the molar mass of a component in [kg/mol].
Definition GenericOilGasWaterFluidSystem.hpp:229
static std::string_view componentName(unsigned compIdx)
Return the human readable name of a component.
Definition GenericOilGasWaterFluidSystem.hpp:267
static bool isLiquid(unsigned phaseIdx)
Return whether a phase is liquid.
Definition GenericOilGasWaterFluidSystem.hpp:355
static Scalar criticalVolume(unsigned compIdx)
Critical volume of a component [m3].
Definition GenericOilGasWaterFluidSystem.hpp:220
static void setWaterPvt(std::shared_ptr< WaterPvt > pvtObj)
Set the pressure-volume-saturation (PVT) relations for the water phase.
Definition GenericOilGasWaterFluidSystem.hpp:176
static LhsEval viscosity(const FluidState &fluidState, const ParameterCache< ParamCacheEval > &paramCache, unsigned phaseIdx)
Calculate the dynamic viscosity of a fluid phase [Pa*s].
Definition GenericOilGasWaterFluidSystem.hpp:302
Specifies the parameter cache used by the SPE-5 fluid system.
Definition PTFlashParameterCache.hpp:51
Scalar molarVolume(unsigned phaseIdx) const
Returns the molar volume of a phase [m^3/mol].
Definition PTFlashParameterCache.hpp:283
Definition Schedule.hpp:101
Definition LBC.hpp:40
This class represents the Pressure-Volume-Temperature relations of the water phase in the black-oil m...
Definition WaterPvtMultiplexer.hpp:88
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30