3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program 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 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program 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 this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_SPE5_FLUID_SYSTEM_HH
25#define DUMUX_SPE5_FLUID_SYSTEM_HH
26
27#include "spe5parametercache.hh"
28
31#include <dumux/io/name.hh>
32
33namespace Dumux {
34namespace FluidSystems {
35
55template <class Scalar>
56class Spe5
57{
59
62
63public:
65
66 /****************************************
67 * Fluid phase parameters
68 ****************************************/
69
71 static const int numPhases = 3;
72
74 static const int gPhaseIdx = 0;
76 static const int wPhaseIdx = 1;
78 static const int oPhaseIdx = 2;
79
82
87 static std::string phaseName(int phaseIdx)
88 {
89 assert(0 <= phaseIdx && phaseIdx < numPhases);
90 switch (phaseIdx)
91 {
92 case gPhaseIdx: return IOName::gaseousPhase();
93 case wPhaseIdx: return IOName::aqueousPhase();
94 case oPhaseIdx: return IOName::naplPhase();
95 }
96 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
97 }
98
102 static constexpr bool isMiscible()
103 { return true; }
104
109 static constexpr bool isGas(int phaseIdx)
110 {
111 assert(0 <= phaseIdx && phaseIdx < numPhases);
112 return phaseIdx == gPhaseIdx;
113 }
114
121 static constexpr bool isCompressible(int phaseIdx)
122 {
123 assert(0 <= phaseIdx && phaseIdx < numPhases);
124 return true;
125 }
126
140 static bool isIdealMixture(int phaseIdx)
141 {
142 // always use the reference oil for the fugacity coefficients,
143 // so they cannot be dependent on composition and they the
144 // phases thus always an ideal mixture
145 return phaseIdx == wPhaseIdx;
146 }
147
154 static bool isIdealGas(int phaseIdx)
155 {
156 assert(0 <= phaseIdx && phaseIdx < numPhases);
157 return false; // gas is not ideal here!
158 }
159
160 /****************************************
161 * Component related parameters
162 ****************************************/
163
165 static const int numComponents = 7;
166
167 static const int H2OIdx = 0;
168 static const int C1Idx = 1;
169 static const int C3Idx = 2;
170 static const int C6Idx = 3;
171 static const int C10Idx = 4;
172 static const int C15Idx = 5;
173 static const int C20Idx = 6;
174
179 static std::string componentName(int compIdx)
180 {
181 static std::string name[] = {
182 H2O::name(),
183 std::string("C1"),
184 std::string("C3"),
185 std::string("C6"),
186 std::string("C10"),
187 std::string("C15"),
188 std::string("C20")
189 };
190
191 assert(0 <= compIdx && compIdx < numComponents);
192 return name[compIdx];
193 }
194
199 static Scalar molarMass(int compIdx)
200 {
201 static const Scalar M[] = {
203 16.04e-3, // C1
204 44.10e-3, // C3
205 86.18e-3, // C6
206 142.29e-3, // C10
207 206.00e-3, // C15
208 282.00e-3 // C20
209 };
210
211 assert(0 <= compIdx && compIdx < numComponents);
212 return M[compIdx];
213 }
214
219 static Scalar criticalTemperature(int compIdx)
220 {
221 static const Scalar Tcrit[] = {
223 343.0*5/9, // C1
224 665.7*5/9, // C3
225 913.4*5/9, // C6
226 1111.8*5/9, // C10
227 1270.0*5/9, // C15
228 1380.0*5/9 // C20
229 };
230
231 assert(0 <= compIdx && compIdx < numComponents);
232 return Tcrit[compIdx];
233 }
234
239 static Scalar criticalPressure(int compIdx)
240 {
241 static const Scalar pcrit[] = {
242 H2O::criticalPressure(), // H2O
243 667.8 * 6894.7573, // C1
244 616.3 * 6894.7573, // C3
245 436.9 * 6894.7573, // C6
246 304.0 * 6894.7573, // C10
247 200.0 * 6894.7573, // C15
248 162.0 * 6894.7573 // C20
249 };
250
251 assert(0 <= compIdx && compIdx < numComponents);
252 return pcrit[compIdx];
253 }
254
259 static Scalar criticalMolarVolume(int compIdx)
260 {
261 static const Scalar R = 8.314472;
262 static const Scalar vcrit[] = {
264 0.290*R*criticalTemperature(1)/criticalPressure(1), // C1
265 0.277*R*criticalTemperature(2)/criticalPressure(2), // C3
266 0.264*R*criticalTemperature(3)/criticalPressure(3), // C6
267 0.257*R*criticalTemperature(4)/criticalPressure(4), // C10
268 0.245*R*criticalTemperature(5)/criticalPressure(5), // C15
269 0.235*R*criticalTemperature(6)/criticalPressure(6) // C20
270 };
271
272 assert(0 <= compIdx && compIdx < numComponents);
273 return vcrit[compIdx];
274 }
275
280 static Scalar acentricFactor(int compIdx)
281 {
282 static const Scalar accFac[] = {
283 0.344, // H2O (from Reid, et al.)
284 0.0130, // C1
285 0.1524, // C3
286 0.3007, // C6
287 0.4885, // C10
288 0.6500, // C15
289 0.8500 // C20
290 };
291
292 assert(0 <= compIdx && compIdx < numComponents);
293 return accFac[compIdx];
294 }
295
303 static Scalar interactionCoefficient(int comp1Idx, int comp2Idx)
304 {
305 using std::min;
306 using std::max;
307 int i = min(comp1Idx, comp2Idx);
308 int j = max(comp1Idx, comp2Idx);
309 if (i == C1Idx && (j == C15Idx || j == C20Idx))
310 return 0.05;
311 else if (i == C3Idx && (j == C15Idx || j == C20Idx))
312 return 0.005;
313 return 0;
314 }
315
316 /****************************************
317 * Methods which compute stuff
318 ****************************************/
319
323 static void init()
324 {
325 PengRobinsonParamsMixture<Scalar, ThisType, gPhaseIdx, /*useSpe5=*/true> prParams;
326
327 // find envelopes of the 'a' and 'b' parameters for the range
328 // 273.15K <= T <= 373.15K and 10e3 Pa <= p <= 100e6 Pa. For
329 // this we take advantage of the fact that 'a' and 'b' for
330 // mixtures is just a convex combination of the attractive and
331 // repulsive parameters of the pure components
332 Scalar minT = 273.15, maxT = 373.15;
333 Scalar minP = 1e4, maxP = 100e6;
334
335 Scalar minA = 1e100, maxA = -1e100;
336 Scalar minB = 1e100, maxB = -1e100;
337
338 prParams.updatePure(minT, minP);
339 using std::min;
340 using std::max;
341 for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
342 minA = min(prParams.pureParams(compIdx).a(), minA);
343 maxA = max(prParams.pureParams(compIdx).a(), maxA);
344 minB = min(prParams.pureParams(compIdx).b(), minB);
345 maxB = max(prParams.pureParams(compIdx).b(), maxB);
346 }
347
348 prParams.updatePure(maxT, minP);
349 for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
350 minA = min(prParams.pureParams(compIdx).a(), minA);
351 maxA = max(prParams.pureParams(compIdx).a(), maxA);
352 minB = min(prParams.pureParams(compIdx).b(), minB);
353 maxB = max(prParams.pureParams(compIdx).b(), maxB);
354 }
355
356 prParams.updatePure(minT, maxP);
357 for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
358 minA = min(prParams.pureParams(compIdx).a(), minA);
359 maxA = max(prParams.pureParams(compIdx).a(), maxA);
360 minB = min(prParams.pureParams(compIdx).b(), minB);
361 maxB = max(prParams.pureParams(compIdx).b(), maxB);
362 }
363
364 prParams.updatePure(maxT, maxP);
365 for (int compIdx = 0; compIdx < numComponents; ++compIdx) {
366 minA = min(prParams.pureParams(compIdx).a(), minA);
367 maxA = max(prParams.pureParams(compIdx).a(), maxA);
368 minB = min(prParams.pureParams(compIdx).b(), minB);
369 maxB = max(prParams.pureParams(compIdx).b(), maxB);
370 }
371
372 PengRobinson::init(/*aMin=*/minA, /*aMax=*/maxA, /*na=*/100,
373 /*bMin=*/minB, /*bMax=*/maxB, /*nb=*/200);
374 }
375
384 template <class FluidState>
385 static Scalar density(const FluidState &fluidState,
386 const ParameterCache &paramCache,
387 int phaseIdx)
388 {
389 assert(0 <= phaseIdx && phaseIdx < numPhases);
390
391 return fluidState.averageMolarMass(phaseIdx)/paramCache.molarVolume(phaseIdx);
392 }
393
403 template <class FluidState>
404 static Scalar molarDensity(const FluidState &fluidState, const ParameterCache &paramCache, int phaseIdx)
405 {
406 return 1.0/paramCache.molarVolume(phaseIdx);
407 }
408
415 template <class FluidState>
416 static Scalar viscosity(const FluidState &fs,
417 const ParameterCache &paramCache,
418 int phaseIdx)
419 {
420 assert(0 <= phaseIdx && phaseIdx <= numPhases);
421
422 if (phaseIdx == gPhaseIdx) {
423 // given by SPE-5 in table on page 64. we use a constant
424 // viscosity, though...
425 return 0.0170e-2 * 0.1;
426 }
427 else if (phaseIdx == wPhaseIdx)
428 // given by SPE-5: 0.7 centi-Poise = 0.0007 Pa s
429 return 0.7e-2 * 0.1;
430 else {
431 assert(phaseIdx == oPhaseIdx);
432 // given by SPE-5 in table on page 64. we use a constant
433 // viscosity, though...
434 return 0.208e-2 * 0.1;
435 }
436 }
437
453 template <class FluidState>
454 static Scalar fugacityCoefficient(const FluidState &fs,
455 const ParameterCache &paramCache,
456 int phaseIdx,
457 int compIdx)
458 {
459 assert(0 <= phaseIdx && phaseIdx <= numPhases);
460 assert(0 <= compIdx && compIdx <= numComponents);
461
462 if (phaseIdx == oPhaseIdx || phaseIdx == gPhaseIdx)
464 paramCache,
465 phaseIdx,
466 compIdx);
467 else {
468 assert(phaseIdx == wPhaseIdx);
469 return henryCoeffWater_(compIdx, fs.temperature(wPhaseIdx))
470 / fs.pressure(wPhaseIdx);
471 }
472 }
473
474
498 template <class FluidState>
499 static Scalar diffusionCoefficient(const FluidState &fs,
500 const ParameterCache &paramCache,
501 int phaseIdx,
502 int compIdx)
503 { DUNE_THROW(Dune::NotImplemented, "Diffusion coefficients"); }
504
516 template <class FluidState>
517 static Scalar binaryDiffusionCoefficient(const FluidState &fluidState,
518 const ParameterCache &paramCache,
519 int phaseIdx,
520 int compIIdx,
521 int compJIdx)
522 { DUNE_THROW(Dune::NotImplemented, "Binary diffusion coefficients"); }
523
531 template <class FluidState>
532 static Scalar enthalpy(const FluidState &fs,
533 const ParameterCache &paramCache,
534 int phaseIdx)
535 { DUNE_THROW(Dune::NotImplemented, "Enthalpies"); }
536
544 template <class FluidState>
545 static Scalar thermalConductivity(const FluidState &fluidState,
546 const ParameterCache &paramCache,
547 int phaseIdx)
548 { DUNE_THROW(Dune::NotImplemented, "Thermal conductivities"); }
549
557 template <class FluidState>
558 static Scalar heatCapacity(const FluidState &fluidState,
559 const ParameterCache &paramCache,
560 int phaseIdx)
561 { DUNE_THROW(Dune::NotImplemented, "Heat capacities"); }
562
563
564private:
565 static Scalar henryCoeffWater_(int compIdx, Scalar temperature)
566 {
567 // use henry's law for the solutes and the vapor pressure for
568 // the solvent.
569 switch (compIdx) {
571
572 // the values of the Henry constant for the solutes have
573 // are faked so far...
574 case C1Idx: return 5.57601e+09;
575 case C3Idx: return 1e10;
576 case C6Idx: return 1e10;
577 case C10Idx: return 1e10;
578 case C15Idx: return 1e10;
579 case C20Idx: return 1e10;
580 default: DUNE_THROW(Dune::InvalidStateException, "Unknown component index " << compIdx);
581 }
582 }
583};
584
585} // end namespace FluidSystems
586} // end namespace Dumux
587
588#endif
Provides 3rd order polynomial splines.
A collection of input/output field names for common physical quantities.
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...
Definition: adapt.hh:29
std::string temperature() noexcept
I/O name of temperature for equilibrium models.
Definition: name.hh:51
std::string gaseousPhase() noexcept
I/O name of gaseous phase.
Definition: name.hh:123
std::string naplPhase() noexcept
I/O name of napl phase.
Definition: name.hh:131
std::string aqueousPhase() noexcept
I/O name of aqueous phase.
Definition: name.hh:127
Material properties of pure water .
Definition: h2o.hh:61
static Scalar vaporPressure(Scalar T)
The vapor pressure in of pure water at a given temperature.
Definition: h2o.hh:131
static constexpr Scalar criticalTemperature()
Returns the critical temperature of water.
Definition: h2o.hh:93
static std::string name()
A human readable name for the water.
Definition: h2o.hh:75
static constexpr Scalar molarMass()
The molar mass in of water.
Definition: h2o.hh:81
static constexpr Scalar criticalPressure()
Returns the critical pressure of water.
Definition: h2o.hh:99
static constexpr Scalar criticalMolarVolume()
Returns the molar volume of water at the critical point.
Definition: h2o.hh:105
Implements the Peng-Robinson equation of state for liquids and gases.
Definition: pengrobinson.hh:57
static void init(Scalar aMin, Scalar aMax, int na, Scalar bMin, Scalar bMax, int nb)
Definition: pengrobinson.hh:62
Implements the Peng-Robinson equation of state for a mixture.
Definition: pengrobinsonmixture.hh:41
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:90
Scalar a() const
Returns the attractive parameter 'a' of the Peng-Robinson fluid.
Definition: pengrobinsonparams.hh:52
Scalar b() const
Returns the repulsive parameter 'b' of the Peng-Robinson fluid.
Definition: pengrobinsonparams.hh:60
The mixing rule for the oil and the gas phases of the SPE5 problem.
Definition: pengrobinsonparamsmixture.hh:63
const PureParams & pureParams(int compIdx) const
Return the Peng-Robinson parameters of a pure substance,.
Definition: pengrobinsonparamsmixture.hh:204
void updatePure(const FluidState &fluidState)
Update Peng-Robinson parameters for the pure components.
Definition: pengrobinsonparamsmixture.hh:76
The fluid system for the SPE-5 benchmark problem.
Definition: spe5.hh:57
static Scalar molarDensity(const FluidState &fluidState, const ParameterCache &paramCache, int phaseIdx)
The molar density of a fluid phase in .
Definition: spe5.hh:404
static bool isIdealMixture(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition: spe5.hh:140
static Scalar criticalTemperature(int compIdx)
Critical temperature of a component .
Definition: spe5.hh:219
static const int C15Idx
Definition: spe5.hh:172
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:517
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:558
static const int C3Idx
Definition: spe5.hh:169
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:454
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:545
static const int H2OIdx
Definition: spe5.hh:167
static std::string componentName(int compIdx)
Return the human readable name of a component.
Definition: spe5.hh:179
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:499
static const int C6Idx
Definition: spe5.hh:170
static const int C20Idx
Definition: spe5.hh:173
static Scalar molarMass(int compIdx)
Return the molar mass of a component in .
Definition: spe5.hh:199
static constexpr bool isMiscible()
Returns whether the fluids are miscible.
Definition: spe5.hh:102
static bool isIdealGas(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition: spe5.hh:154
static const int C1Idx
Definition: spe5.hh:168
static const int numComponents
Number of components in the fluid system.
Definition: spe5.hh:165
static const int oPhaseIdx
Index of the oil phase.
Definition: spe5.hh:78
static Scalar density(const FluidState &fluidState, const ParameterCache &paramCache, int phaseIdx)
Calculate the density of a fluid phase.
Definition: spe5.hh:385
static const int C10Idx
Definition: spe5.hh:171
static Scalar acentricFactor(int compIdx)
The acentric factor of a component .
Definition: spe5.hh:280
static const int gPhaseIdx
Index of the gas phase.
Definition: spe5.hh:74
static Scalar viscosity(const FluidState &fs, const ParameterCache &paramCache, int phaseIdx)
Calculate the dynamic viscosity of a fluid phase .
Definition: spe5.hh:416
static void init()
Initialize the fluid system's static parameters.
Definition: spe5.hh:323
static Scalar interactionCoefficient(int comp1Idx, int comp2Idx)
Returns the interaction coefficient for two components.
Definition: spe5.hh:303
static const int numPhases
Number of phases in the fluid system.
Definition: spe5.hh:71
static std::string phaseName(int phaseIdx)
Return the human readable name of a fluid phase.
Definition: spe5.hh:87
static const int wPhaseIdx
Index of the water phase.
Definition: spe5.hh:76
static constexpr bool isGas(int phaseIdx)
Return whether a phase is gaseous.
Definition: spe5.hh:109
static Scalar criticalPressure(int compIdx)
Critical pressure of a component .
Definition: spe5.hh:239
static Scalar criticalMolarVolume(int compIdx)
Molar volume of a component at the critical point .
Definition: spe5.hh:259
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:532
static constexpr bool isCompressible(int phaseIdx)
Return whether a phase is compressible.
Definition: spe5.hh:121
Specifies the parameters required by the SPE5 problem which are despondent on the thermodynamic state...
Definition: spe5parametercache.hh:45
Scalar molarVolume(int phaseIdx) const
Returns the molar volume of a phase .
Definition: spe5parametercache.hh:200