version 3.8
spe5.hh
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_SPE5_FLUID_SYSTEM_HH
13#define DUMUX_SPE5_FLUID_SYSTEM_HH
14
15#include "spe5parametercache.hh"
16
19#include <dumux/io/name.hh>
20
21namespace Dumux {
22namespace FluidSystems {
23
43template <class Scalar>
44class Spe5
45{
47
50
51public:
53
54 /****************************************
55 * Fluid phase parameters
56 ****************************************/
57
59 static const int numPhases = 3;
60
62 static const int gPhaseIdx = 0;
64 static const int wPhaseIdx = 1;
66 static const int oPhaseIdx = 2;
67
70
75 static std::string phaseName(int phaseIdx)
76 {
77 assert(0 <= phaseIdx && phaseIdx < numPhases);
78 switch (phaseIdx)
79 {
80 case gPhaseIdx: return IOName::gaseousPhase();
81 case wPhaseIdx: return IOName::aqueousPhase();
82 case oPhaseIdx: return IOName::naplPhase();
83 }
84 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
85 }
86
90 static constexpr bool isMiscible()
91 { return true; }
92
97 static constexpr bool isGas(int phaseIdx)
98 {
99 assert(0 <= phaseIdx && phaseIdx < numPhases);
100 return phaseIdx == gPhaseIdx;
101 }
102
109 static constexpr bool isCompressible(int phaseIdx)
110 {
111 assert(0 <= phaseIdx && phaseIdx < numPhases);
112 return true;
113 }
114
128 static bool isIdealMixture(int phaseIdx)
129 {
130 // always use the reference oil for the fugacity coefficients,
131 // so they cannot be dependent on composition and they the
132 // phases thus always an ideal mixture
133 return phaseIdx == wPhaseIdx;
134 }
135
142 static bool isIdealGas(int phaseIdx)
143 {
144 assert(0 <= phaseIdx && phaseIdx < numPhases);
145 return false; // gas is not ideal here!
146 }
147
148 /****************************************
149 * Component related parameters
150 ****************************************/
151
153 static const int numComponents = 7;
154
155 static const int H2OIdx = 0;
156 static const int C1Idx = 1;
157 static const int C3Idx = 2;
158 static const int C6Idx = 3;
159 static const int C10Idx = 4;
160 static const int C15Idx = 5;
161 static const int C20Idx = 6;
162
167 static std::string componentName(int compIdx)
168 {
169 static std::string name[] = {
170 H2O::name(),
171 std::string("C1"),
172 std::string("C3"),
173 std::string("C6"),
174 std::string("C10"),
175 std::string("C15"),
176 std::string("C20")
177 };
178
179 assert(0 <= compIdx && compIdx < numComponents);
180 return name[compIdx];
181 }
182
187 static Scalar molarMass(int compIdx)
188 {
189 static const Scalar M[] = {
191 16.04e-3, // C1
192 44.10e-3, // C3
193 86.18e-3, // C6
194 142.29e-3, // C10
195 206.00e-3, // C15
196 282.00e-3 // C20
197 };
198
199 assert(0 <= compIdx && compIdx < numComponents);
200 return M[compIdx];
201 }
202
207 static Scalar criticalTemperature(int compIdx)
208 {
209 static const Scalar Tcrit[] = {
211 343.0*5/9, // C1
212 665.7*5/9, // C3
213 913.4*5/9, // C6
214 1111.8*5/9, // C10
215 1270.0*5/9, // C15
216 1380.0*5/9 // C20
217 };
218
219 assert(0 <= compIdx && compIdx < numComponents);
220 return Tcrit[compIdx];
221 }
222
227 static Scalar criticalPressure(int compIdx)
228 {
229 static const Scalar pcrit[] = {
230 H2O::criticalPressure(), // H2O
231 667.8 * 6894.7573, // C1
232 616.3 * 6894.7573, // C3
233 436.9 * 6894.7573, // C6
234 304.0 * 6894.7573, // C10
235 200.0 * 6894.7573, // C15
236 162.0 * 6894.7573 // C20
237 };
238
239 assert(0 <= compIdx && compIdx < numComponents);
240 return pcrit[compIdx];
241 }
242
247 static Scalar criticalMolarVolume(int compIdx)
248 {
249 static const Scalar R = 8.314472;
250 static const Scalar vcrit[] = {
252 0.290*R*criticalTemperature(1)/criticalPressure(1), // C1
253 0.277*R*criticalTemperature(2)/criticalPressure(2), // C3
254 0.264*R*criticalTemperature(3)/criticalPressure(3), // C6
255 0.257*R*criticalTemperature(4)/criticalPressure(4), // C10
256 0.245*R*criticalTemperature(5)/criticalPressure(5), // C15
257 0.235*R*criticalTemperature(6)/criticalPressure(6) // C20
258 };
259
260 assert(0 <= compIdx && compIdx < numComponents);
261 return vcrit[compIdx];
262 }
263
268 static Scalar acentricFactor(int compIdx)
269 {
270 static const Scalar accFac[] = {
271 0.344, // H2O (from Reid, et al.)
272 0.0130, // C1
273 0.1524, // C3
274 0.3007, // C6
275 0.4885, // C10
276 0.6500, // C15
277 0.8500 // C20
278 };
279
280 assert(0 <= compIdx && compIdx < numComponents);
281 return accFac[compIdx];
282 }
283
291 static Scalar interactionCoefficient(int comp1Idx, int comp2Idx)
292 {
293 using std::min;
294 using std::max;
295 int i = min(comp1Idx, comp2Idx);
296 int j = max(comp1Idx, comp2Idx);
297 if (i == C1Idx && (j == C15Idx || j == C20Idx))
298 return 0.05;
299 else if (i == C3Idx && (j == C15Idx || j == C20Idx))
300 return 0.005;
301 return 0;
302 }
303
304 /****************************************
305 * Methods which compute stuff
306 ****************************************/
307
311 static void init()
312 {
313 PengRobinsonParamsMixture<Scalar, ThisType, gPhaseIdx, /*useSpe5=*/true> prParams;
314
315 // find envelopes of the 'a' and 'b' parameters for the range
316 // 273.15K <= T <= 373.15K and 10e3 Pa <= p <= 100e6 Pa. For
317 // this we take advantage of the fact that 'a' and 'b' for
318 // mixtures is just a convex combination of the attractive and
319 // repulsive parameters of the pure components
320 Scalar minT = 273.15, maxT = 373.15;
321 Scalar minP = 1e4, maxP = 100e6;
322
323 Scalar minA = 1e100, maxA = -1e100;
324 Scalar minB = 1e100, maxB = -1e100;
325
326 prParams.updatePure(minT, minP);
327 using std::min;
328 using std::max;
329 for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
330 minA = min(prParams.pureParams(compIdx).a(), minA);
331 maxA = max(prParams.pureParams(compIdx).a(), maxA);
332 minB = min(prParams.pureParams(compIdx).b(), minB);
333 maxB = max(prParams.pureParams(compIdx).b(), maxB);
334 }
335
336 prParams.updatePure(maxT, minP);
337 for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
338 minA = min(prParams.pureParams(compIdx).a(), minA);
339 maxA = max(prParams.pureParams(compIdx).a(), maxA);
340 minB = min(prParams.pureParams(compIdx).b(), minB);
341 maxB = max(prParams.pureParams(compIdx).b(), maxB);
342 }
343
344 prParams.updatePure(minT, maxP);
345 for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
346 minA = min(prParams.pureParams(compIdx).a(), minA);
347 maxA = max(prParams.pureParams(compIdx).a(), maxA);
348 minB = min(prParams.pureParams(compIdx).b(), minB);
349 maxB = max(prParams.pureParams(compIdx).b(), maxB);
350 }
351
352 prParams.updatePure(maxT, maxP);
353 for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
354 minA = min(prParams.pureParams(compIdx).a(), minA);
355 maxA = max(prParams.pureParams(compIdx).a(), maxA);
356 minB = min(prParams.pureParams(compIdx).b(), minB);
357 maxB = max(prParams.pureParams(compIdx).b(), maxB);
358 }
359
360 PengRobinson::init(/*aMin=*/minA, /*aMax=*/maxA, /*na=*/100,
361 /*bMin=*/minB, /*bMax=*/maxB, /*nb=*/200);
362 }
363
372 template <class FluidState>
373 static Scalar density(const FluidState &fluidState,
374 const ParameterCache &paramCache,
375 int phaseIdx)
376 {
377 assert(0 <= phaseIdx && phaseIdx < numPhases);
378
379 return fluidState.averageMolarMass(phaseIdx)/paramCache.molarVolume(phaseIdx);
380 }
381
391 template <class FluidState>
392 static Scalar molarDensity(const FluidState &fluidState, const ParameterCache &paramCache, int phaseIdx)
393 {
394 return 1.0/paramCache.molarVolume(phaseIdx);
395 }
396
403 template <class FluidState>
404 static Scalar viscosity(const FluidState &fs,
405 const ParameterCache &paramCache,
406 int phaseIdx)
407 {
408 assert(0 <= phaseIdx && phaseIdx <= numPhases);
409
410 if (phaseIdx == gPhaseIdx) {
411 // given by SPE-5 in table on page 64. we use a constant
412 // viscosity, though...
413 return 0.0170e-2 * 0.1;
414 }
415 else if (phaseIdx == wPhaseIdx)
416 // given by SPE-5: 0.7 centi-Poise = 0.0007 Pa s
417 return 0.7e-2 * 0.1;
418 else {
419 assert(phaseIdx == oPhaseIdx);
420 // given by SPE-5 in table on page 64. we use a constant
421 // viscosity, though...
422 return 0.208e-2 * 0.1;
423 }
424 }
425
441 template <class FluidState>
442 static Scalar fugacityCoefficient(const FluidState &fs,
443 const ParameterCache &paramCache,
444 int phaseIdx,
445 int compIdx)
446 {
447 assert(0 <= phaseIdx && phaseIdx <= numPhases);
448 assert(0 <= compIdx && compIdx <= numComponents);
449
450 if (phaseIdx == oPhaseIdx || phaseIdx == gPhaseIdx)
452 paramCache,
453 phaseIdx,
454 compIdx);
455 else {
456 assert(phaseIdx == wPhaseIdx);
457 return henryCoeffWater_(compIdx, fs.temperature(wPhaseIdx))
458 / fs.pressure(wPhaseIdx);
459 }
460 }
461
462
486 template <class FluidState>
487 static Scalar diffusionCoefficient(const FluidState &fs,
488 const ParameterCache &paramCache,
489 int phaseIdx,
490 int compIdx)
491 { DUNE_THROW(Dune::NotImplemented, "Diffusion coefficients"); }
492
504 template <class FluidState>
505 static Scalar binaryDiffusionCoefficient(const FluidState &fluidState,
506 const ParameterCache &paramCache,
507 int phaseIdx,
508 int compIIdx,
509 int compJIdx)
510 { DUNE_THROW(Dune::NotImplemented, "Binary diffusion coefficients"); }
511
519 template <class FluidState>
520 static Scalar enthalpy(const FluidState &fs,
521 const ParameterCache &paramCache,
522 int phaseIdx)
523 { DUNE_THROW(Dune::NotImplemented, "Enthalpies"); }
524
532 template <class FluidState>
533 static Scalar thermalConductivity(const FluidState &fluidState,
534 const ParameterCache &paramCache,
535 int phaseIdx)
536 { DUNE_THROW(Dune::NotImplemented, "Thermal conductivities"); }
537
545 template <class FluidState>
546 static Scalar heatCapacity(const FluidState &fluidState,
547 const ParameterCache &paramCache,
548 int phaseIdx)
549 { DUNE_THROW(Dune::NotImplemented, "Heat capacities"); }
550
551
552private:
553 static Scalar henryCoeffWater_(int compIdx, Scalar temperature)
554 {
555 // use henry's law for the solutes and the vapor pressure for
556 // the solvent.
557 switch (compIdx) {
559
560 // the values of the Henry constant for the solutes have
561 // are faked so far...
562 case C1Idx: return 5.57601e+09;
563 case C3Idx: return 1e10;
564 case C6Idx: return 1e10;
565 case C10Idx: return 1e10;
566 case C15Idx: return 1e10;
567 case C20Idx: return 1e10;
568 default: DUNE_THROW(Dune::InvalidStateException, "Unknown component index " << compIdx);
569 }
570 }
571};
572
573} // end namespace FluidSystems
574} // end namespace Dumux
575
576#endif
Material properties of pure water .
Definition: h2o.hh:49
static Scalar vaporPressure(Scalar T)
The vapor pressure in of pure water at a given temperature.
Definition: h2o.hh:119
static constexpr Scalar criticalTemperature()
Returns the critical temperature of water.
Definition: h2o.hh:81
static std::string name()
A human readable name for the water.
Definition: h2o.hh:63
static constexpr Scalar molarMass()
The molar mass in of water.
Definition: h2o.hh:69
static constexpr Scalar criticalPressure()
Returns the critical pressure of water.
Definition: h2o.hh:87
static constexpr Scalar criticalMolarVolume()
Returns the molar volume of water at the critical point.
Definition: h2o.hh:93
The fluid system for the SPE-5 benchmark problem.
Definition: spe5.hh:45
static Scalar molarDensity(const FluidState &fluidState, const ParameterCache &paramCache, int phaseIdx)
The molar density of a fluid phase in .
Definition: spe5.hh:392
static bool isIdealMixture(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition: spe5.hh:128
static Scalar criticalTemperature(int compIdx)
Critical temperature of a component .
Definition: spe5.hh:207
static const int C15Idx
Definition: spe5.hh:160
static Scalar binaryDiffusionCoefficient(const FluidState &fluidState, const ParameterCache &paramCache, int phaseIdx, int compIIdx, int compJIdx)
Given a phase's composition, temperature and pressure, return the binary diffusion coefficient for c...
Definition: spe5.hh:505
static Scalar heatCapacity(const FluidState &fluidState, const ParameterCache &paramCache, int phaseIdx)
Given a phase's composition, temperature and pressure, calculate its heat capacity .
Definition: spe5.hh:546
static const int C3Idx
Definition: spe5.hh:157
static Scalar fugacityCoefficient(const FluidState &fs, const ParameterCache &paramCache, int phaseIdx, int compIdx)
Calculate the fugacity coefficient of an individual component in a fluid phase.
Definition: spe5.hh:442
static Scalar thermalConductivity(const FluidState &fluidState, const ParameterCache &paramCache, int phaseIdx)
Given a phase's composition, temperature and pressure, calculate its thermal conductivity .
Definition: spe5.hh:533
static const int H2OIdx
Definition: spe5.hh:155
static std::string componentName(int compIdx)
Return the human readable name of a component.
Definition: spe5.hh:167
static Scalar diffusionCoefficient(const FluidState &fs, const ParameterCache &paramCache, int phaseIdx, int compIdx)
Calculate the binary molecular diffusion coefficient for a component in a fluid phase .
Definition: spe5.hh:487
static const int C6Idx
Definition: spe5.hh:158
static const int C20Idx
Definition: spe5.hh:161
static Scalar molarMass(int compIdx)
Return the molar mass of a component in .
Definition: spe5.hh:187
static constexpr bool isMiscible()
Returns whether the fluids are miscible.
Definition: spe5.hh:90
static bool isIdealGas(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition: spe5.hh:142
static const int C1Idx
Definition: spe5.hh:156
static const int numComponents
Number of components in the fluid system.
Definition: spe5.hh:153
static const int oPhaseIdx
Index of the oil phase.
Definition: spe5.hh:66
static Scalar density(const FluidState &fluidState, const ParameterCache &paramCache, int phaseIdx)
Calculate the density of a fluid phase.
Definition: spe5.hh:373
static const int C10Idx
Definition: spe5.hh:159
static Scalar acentricFactor(int compIdx)
The acentric factor of a component .
Definition: spe5.hh:268
static const int gPhaseIdx
Index of the gas phase.
Definition: spe5.hh:62
static Scalar viscosity(const FluidState &fs, const ParameterCache &paramCache, int phaseIdx)
Calculate the dynamic viscosity of a fluid phase .
Definition: spe5.hh:404
static void init()
Initialize the fluid system's static parameters.
Definition: spe5.hh:311
static Scalar interactionCoefficient(int comp1Idx, int comp2Idx)
Returns the interaction coefficient for two components.
Definition: spe5.hh:291
static const int numPhases
Number of phases in the fluid system.
Definition: spe5.hh:59
static std::string phaseName(int phaseIdx)
Return the human readable name of a fluid phase.
Definition: spe5.hh:75
static const int wPhaseIdx
Index of the water phase.
Definition: spe5.hh:64
static constexpr bool isGas(int phaseIdx)
Return whether a phase is gaseous.
Definition: spe5.hh:97
static Scalar criticalPressure(int compIdx)
Critical pressure of a component .
Definition: spe5.hh:227
static Scalar criticalMolarVolume(int compIdx)
Molar volume of a component at the critical point .
Definition: spe5.hh:247
static Scalar enthalpy(const FluidState &fs, const ParameterCache &paramCache, int phaseIdx)
Given a phase's composition, temperature and pressure, calculate its specific enthalpy .
Definition: spe5.hh:520
static constexpr bool isCompressible(int phaseIdx)
Return whether a phase is compressible.
Definition: spe5.hh:109
Implements the Peng-Robinson equation of state for liquids and gases.
Definition: pengrobinson.hh:48
static void init(Scalar aMin, Scalar aMax, int na, Scalar bMin, Scalar bMax, int nb)
Definition: pengrobinson.hh:53
Implements the Peng-Robinson equation of state for a mixture.
Definition: pengrobinsonmixture.hh:29
static Scalar computeFugacityCoefficient(const FluidState &fs, const Params &params, int phaseIdx, int compIdx)
Returns the fugacity coefficient of an individual component in the phase.
Definition: pengrobinsonmixture.hh:61
Scalar a() const
Returns the attractive parameter 'a' of the Peng-Robinson fluid.
Definition: pengrobinsonparams.hh:38
Scalar b() const
Returns the repulsive parameter 'b' of the Peng-Robinson fluid.
Definition: pengrobinsonparams.hh:46
The mixing rule for the oil and the gas phases of the SPE5 problem.
Definition: pengrobinsonparamsmixture.hh:51
const PureParams & pureParams(int compIdx) const
Return the Peng-Robinson parameters of a pure substance,.
Definition: pengrobinsonparamsmixture.hh:180
void updatePure(const FluidState &fluidState)
Update Peng-Robinson parameters for the pure components.
Definition: pengrobinsonparamsmixture.hh:64
Specifies the parameters required by the SPE5 problem which are despondent on the thermodynamic state...
Definition: spe5parametercache.hh:33
Scalar molarVolume(int phaseIdx) const
Returns the molar volume of a phase .
Definition: spe5parametercache.hh:186
A collection of input/output field names for common physical quantities.
std::string temperature() noexcept
I/O name of temperature for equilibrium models.
Definition: name.hh:39
std::string gaseousPhase() noexcept
I/O name of gaseous phase.
Definition: name.hh:111
std::string naplPhase() noexcept
I/O name of napl phase.
Definition: name.hh:119
std::string aqueousPhase() noexcept
I/O name of aqueous phase.
Definition: name.hh:115
Definition: adapt.hh:17
Implements the Peng-Robinson equation of state for a mixture.
Specifies the parameters required by the SPE5 problem which are despondent on the thermodynamic state...
Provides 3rd order polynomial splines.