3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
combustionfluidsystem.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_PURE_WATER_FLUID_SYSTEM_HH
25#define DUMUX_PURE_WATER_FLUID_SYSTEM_HH
26
27#include <cassert>
28
30
35
38
39namespace Dumux {
40namespace FluidSystems {
41
49template <class Scalar>
51 : public Base<Scalar, CombustionFluidsystem<Scalar> >
52{
55
56 // convenience using declarations
60
61public:
62 /****************************************
63 * Fluid phase related static parameters
64 ****************************************/
65
67 static constexpr int numPhases = 2;
68
69 static constexpr int wPhaseIdx = 0; // index of the wetting phase
70 static constexpr int nPhaseIdx = 1; // index of the non-wetting phase
71 static constexpr int phase0Idx = 0; // index of the wetting phase
72 static constexpr int phase1Idx = 1; // index of the non-wetting phase
73
74 // export component indices to indicate the main component
75 // of the corresponding phase at atmospheric pressure 1 bar
76 // and room temperature 20°C:
77 static constexpr int wCompIdx = wPhaseIdx;
78 static constexpr int nCompIdx = nPhaseIdx;
79 static constexpr int comp0Idx = 0; // index of the wetting phase
80 static constexpr int comp1Idx = 1; // index of the non-wetting phase
81
87 static std::string phaseName(int phaseIdx)
88 {
89 assert(0 <= phaseIdx && phaseIdx < numPhases);
90 switch (phaseIdx)
91 {
92 case wPhaseIdx: return "liq";
93 case nPhaseIdx: return "gas";
94 }
95 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
96 }
97
103 static constexpr bool isGas(int phaseIdx)
104 {
105 assert(0 <= phaseIdx && phaseIdx < numPhases);
106 return phaseIdx == nPhaseIdx;
107 }
108
123 static bool isIdealMixture(int phaseIdx)
124 {
125 assert(0 <= phaseIdx && phaseIdx < numPhases);
126 // we assume Henry's and Raoult's laws for the water phase and
127 // and no interaction between gas molecules of different
128 // components, so all phases are ideal mixtures!
129 return true;
130 }
131
141 static constexpr bool isCompressible(int phaseIdx)
142 {
143 assert(0 <= phaseIdx && phaseIdx < numPhases);
144 // gases are always compressible
145 if (phaseIdx == nPhaseIdx)
146 return true;
147 // the water component decides for the liquid phase...
149 }
150
157 static bool isIdealGas(int phaseIdx)
158 {
159 assert(0 <= phaseIdx && phaseIdx < numPhases);
160
161 if (phaseIdx == nPhaseIdx)
162 // let the components decide
163 return H2O::gasIsIdeal() && N2::gasIsIdeal();
164 return false; // not a gas
165 }
166
167 /****************************************
168 * Component related static parameters
169 ****************************************/
170
172 static constexpr int numComponents = 2;
173
174 static constexpr int H2OIdx = wCompIdx;
175 static constexpr int N2Idx = nCompIdx;
176
178 using H2O = SimpleH2O;
179
181 using N2 = SimpleN2;
182
188 static std::string componentName(int compIdx)
189 {
190 static std::string name[] = {
191 H2O::name(),
192 N2::name()
193 };
194
195 assert(0 <= compIdx && compIdx < numComponents);
196 return name[compIdx];
197 }
198
204 static Scalar molarMass(int compIdx)
205 {
206 static const Scalar M[] = {
209 };
210
211 assert(0 <= compIdx && compIdx < numComponents);
212 return M[compIdx];
213 }
214
220 static Scalar criticalTemperature(int compIdx)
221 {
222 static const Scalar Tcrit[] = {
225 };
226
227 assert(0 <= compIdx && compIdx < numComponents);
228 return Tcrit[compIdx];
229 }
230
236 static Scalar criticalPressure(int compIdx)
237 {
238 static const Scalar pcrit[] = {
241 };
242
243 assert(0 <= compIdx && compIdx < numComponents);
244 return pcrit[compIdx];
245 }
246
252 static Scalar criticalMolarVolume(int compIdx)
253 {
254 DUNE_THROW(Dune::NotImplemented, "criticalMolarVolume()");
255 }
256
262 static Scalar acentricFactor(int compIdx)
263 {
264 static const Scalar accFac[] = {
265 H2O::acentricFactor(), // H2O (from Reid, et al.)
266 N2::acentricFactor()
267 };
268
269 assert(0 <= compIdx && compIdx < numComponents);
270 return accFac[compIdx];
271 }
272
273 /****************************************
274 * thermodynamic relations
275 ****************************************/
276
283 static void init()
284 {
285 init(/*tempMin=*/273.15,
286 /*tempMax=*/623.15,
287 /*numTemp=*/100,
288 /*pMin=*/0.0,
289 /*pMax=*/20e6,
290 /*numP=*/200);
291 }
292
304 static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp,
305 Scalar pressMin, Scalar pressMax, unsigned nPress)
306 {
307 std::cout << "Using very simple pure water fluid system\n";
308 }
309
310 using Base::density;
317 template <class FluidState>
318 static Scalar density(const FluidState &fluidState,
319 int phaseIdx)
320 {
321 assert(0 <= phaseIdx && phaseIdx < numPhases);
322
323 // liquid phase
324 if (phaseIdx == wPhaseIdx)
325 {
326 return 1044.0;
327 }
328 else if (phaseIdx == nPhaseIdx)// gas phase
329 {
330 return 1.679;
331 }
332 else DUNE_THROW(Dune::NotImplemented, "Wrong phase index");
333 }
334
335 using Base::molarDensity;
342 template <class FluidState>
343 static Scalar molarDensity(const FluidState &fluidState, int phaseIdx)
344 {
345 if (phaseIdx == wPhaseIdx)
346 {
347 return density(fluidState, phaseIdx)/fluidState.averageMolarMass(phaseIdx);
348 }
349 else if (phaseIdx == nPhaseIdx)
350 {
351 return density(fluidState, phaseIdx)/fluidState.averageMolarMass(phaseIdx);
352 }
353 else DUNE_THROW(Dune::NotImplemented, "Wrong phase index");
354 }
355
356 using Base::viscosity;
363 template <class FluidState>
364 static Scalar viscosity(const FluidState &fluidState,
365 int phaseIdx)
366 {
367 assert(0 <= phaseIdx && phaseIdx < numPhases);
368
369 // liquid phase
370 if (phaseIdx == wPhaseIdx)
371 {
372 return 2.694e-7 * density(fluidState, phaseIdx);
373 }
374 else if (phaseIdx == nPhaseIdx) // gas phase
375 {
376 return 7.16e-6 * density(fluidState, phaseIdx);
377 }
378 else DUNE_THROW(Dune::NotImplemented, "Wrong phase index");
379 }
380
387 template <class FluidState>
388 static Scalar vaporTemperature(const FluidState &fluidState,
389 const unsigned int phaseIdx)
390 {
391 assert(0 <= phaseIdx && phaseIdx < numPhases);
392 Scalar pressure = fluidState.pressure(nPhaseIdx);
393
395 }
396
425 template <class FluidState>
426 static Scalar fugacityCoefficient(const FluidState &fluidState,
427 int phaseIdx,
428 int compIdx)
429 {
430 assert(0 <= phaseIdx && phaseIdx < numPhases);
431 assert(0 <= compIdx && compIdx < numComponents);
432
433 Scalar T = fluidState.temperature(phaseIdx);
434 Scalar p = fluidState.pressure(phaseIdx);
435
436 // liquid phase
437 if (phaseIdx == wPhaseIdx)
438 {
439 if (compIdx == H2OIdx)
440 return H2O::vaporPressure(T)/p;
441 return BinaryCoeff::H2O_N2::henry(T)/p;
442 }
443
444 // for the gas phase, assume an ideal gas when it comes to
445 // fugacity (-> fugacity == partial pressure)
446 return 1.0;
447 }
448
458 template <class FluidState>
459 static Scalar diffusionCoefficient(const FluidState &fluidState,
460 int phaseIdx,
461 int compIdx)
462 {
463 DUNE_THROW(Dune::NotImplemented, "Diffusion coefficients");
464 }
465
477 template <class FluidState>
478 static Scalar binaryDiffusionCoefficient(const FluidState &fluidState,
479 int phaseIdx,
480 int compIIdx,
481 int compJIdx)
482
483 {
484 DUNE_THROW(Dune::NotImplemented, "Binary Diffusion coefficients");
485 }
486
487 using Base::enthalpy;
494 template <class FluidState>
495 static Scalar enthalpy(const FluidState &fluidState,
496 int phaseIdx)
497 {
498 assert(0 <= phaseIdx && phaseIdx < numPhases);
499 Scalar temperature = fluidState.temperature(phaseIdx);
500
501 const Scalar cp = heatCapacity(fluidState, phaseIdx);
502
503 // liquid phase
504 if (phaseIdx == wPhaseIdx)
505 {
506 return cp * (temperature - 373.15);
507 }
508 else if (phaseIdx == nPhaseIdx) // gas phase
509 {
510 return cp * (temperature - 373.15) + 2.257e6;
511 }
512 else DUNE_THROW(Dune::NotImplemented, "Wrong phase index");
513 }
514
524 template <class FluidState>
525 static Scalar thermalConductivity(const FluidState &fluidState,
526 const int phaseIdx)
527 {
528 assert(0 <= phaseIdx && phaseIdx < numPhases);
529 // liquid phase
530 if (phaseIdx == wPhaseIdx)
531 {
532 return 0.68;
533 }
534 else if (phaseIdx == nPhaseIdx) // gas phase
535 {
536 return 0.0248;
537 }
538 else DUNE_THROW(Dune::NotImplemented, "Wrong phase index");
539 }
540
541 using Base::heatCapacity;
549 template <class FluidState>
550 static Scalar heatCapacity(const FluidState &fluidState,
551 int phaseIdx)
552 {
553 assert(0 <= phaseIdx && phaseIdx < numPhases);
554 // liquid phase
555 if (phaseIdx == wPhaseIdx)
556 {
557 return 4.217e3;
558 }
559 else if (phaseIdx == nPhaseIdx) // gas phase
560 {
561 return 2.029e3;
562 }
563 else DUNE_THROW(Dune::NotImplemented, "Wrong phase index");
564 }
565};
566
567} // end namespace FluidSystems
568} // end namespace Dumux
569
570#endif
Properties of pure molecular nitrogen .
A much simpler (and thus potentially less buggy) version of pure water.
Relations valid for an ideal gas.
Binary coefficients for water and nitrogen.
Some exceptions thrown in DuMux
Some templates to wrap the valgrind macros.
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
std::string temperature() noexcept
I/O name of temperature for equilibrium models.
Definition: name.hh:51
std::string pressure(int phaseIdx) noexcept
I/O name of pressure for multiphase systems.
Definition: name.hh:34
static Scalar henry(Scalar temperature)
Henry coefficient for molecular nitrogen in liquid water.
Definition: h2o_n2.hh:48
static Scalar vaporTemperature(Scalar pressure)
Returns the saturation temperature in of pure water at a given pressure.
Definition: region4.hh:94
Properties of pure molecular nitrogen .
Definition: n2.hh:47
static Scalar criticalTemperature()
Returns the critical temperature of molecular nitrogen.
Definition: n2.hh:66
static Scalar criticalPressure()
Returns the critical pressure of molecular nitrogen.
Definition: n2.hh:72
static std::string name()
A human readable name for nitrogen.
Definition: n2.hh:54
static constexpr Scalar molarMass()
The molar mass in of molecular nitrogen.
Definition: n2.hh:60
static constexpr bool gasIsIdeal()
Returns true if the gas phase is assumed to be ideal.
Definition: n2.hh:155
A much simpler (and thus potentially less buggy) version of pure water.
Definition: simpleh2o.hh:51
static Scalar criticalTemperature()
Returns the critical temperature of water.
Definition: simpleh2o.hh:72
static std::string name()
A human readable name for the water.
Definition: simpleh2o.hh:60
static Scalar criticalPressure()
Returns the critical pressure of water.
Definition: simpleh2o.hh:78
static constexpr Scalar molarMass()
The molar mass in of water.
Definition: simpleh2o.hh:66
static constexpr bool liquidIsCompressible()
Returns true if the liquid phase is assumed to be compressible.
Definition: simpleh2o.hh:198
static Scalar vaporPressure(Scalar T)
The vapor pressure in of pure water at a given temperature.
Definition: simpleh2o.hh:105
static constexpr bool gasIsIdeal()
Returns true if the gas phase is assumed to be ideal.
Definition: simpleh2o.hh:238
Fluid system base class.
Definition: fluidsystems/base.hh:45
Scalar Scalar
export the scalar type
Definition: fluidsystems/base.hh:48
static Scalar density(const FluidState &fluidState, int phaseIdx)
Calculate the density of a fluid phase.
Definition: fluidsystems/base.hh:134
static Scalar thermalConductivity(const FluidState &fluidState, int phaseIdx)
Thermal conductivity of a fluid phase .
Definition: fluidsystems/base.hh:390
static Scalar fugacityCoefficient(const FluidState &fluidState, int phaseIdx, int compIdx)
Calculate the fugacity coefficient of an individual component in a fluid phase.
Definition: fluidsystems/base.hh:197
static Scalar diffusionCoefficient(const FluidState &fluidState, int phaseIdx, int compIdx)
Calculate the binary molecular diffusion coefficient for a component in a fluid phase .
Definition: fluidsystems/base.hh:278
static Scalar binaryDiffusionCoefficient(const FluidState &fluidState, int phaseIdx, int compIIdx, int compJIdx)
Given a phase's composition, temperature and pressure, return the binary diffusion coefficient for c...
Definition: fluidsystems/base.hh:326
static Scalar enthalpy(const FluidState &fluidState, int phaseIdx)
Given a phase's composition, temperature, pressure and density, calculate its specific enthalpy .
Definition: fluidsystems/base.hh:363
static Scalar molarDensity(const FluidState &fluidState, int phaseIdx)
Calculate the molar density of a fluid phase.
Definition: fluidsystems/base.hh:160
static Scalar viscosity(const FluidState &fluidState, int phaseIdx)
Calculate the dynamic viscosity of a fluid phase .
Definition: fluidsystems/base.hh:236
static Scalar heatCapacity(const FluidState &fluidState, int phaseIdx)
Specific isobaric heat capacity of a fluid phase .
Definition: fluidsystems/base.hh:424
Relations valid for an ideal gas.
Definition: idealgas.hh:37
A two-phase fluid system with water as sole component.
Definition: combustionfluidsystem.hh:52
static std::string phaseName(int phaseIdx)
Returns the human readable name of a fluid phase.
Definition: combustionfluidsystem.hh:87
static Scalar criticalPressure(int compIdx)
Critical pressure of a component .
Definition: combustionfluidsystem.hh:236
static Scalar heatCapacity(const FluidState &fluidState, int phaseIdx)
Specific isobaric heat capacity of a fluid phase. .
Definition: combustionfluidsystem.hh:550
static void init()
Initializes the fluid system's static parameters generically.
Definition: combustionfluidsystem.hh:283
static constexpr int numPhases
Number of phases in the fluid system.
Definition: combustionfluidsystem.hh:67
static constexpr int nPhaseIdx
Definition: combustionfluidsystem.hh:70
static constexpr int comp0Idx
Definition: combustionfluidsystem.hh:79
static constexpr int nCompIdx
Definition: combustionfluidsystem.hh:78
static bool isIdealMixture(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition: combustionfluidsystem.hh:123
static constexpr int H2OIdx
Definition: combustionfluidsystem.hh:174
static constexpr int comp1Idx
Definition: combustionfluidsystem.hh:80
static constexpr int phase0Idx
Definition: combustionfluidsystem.hh:71
static constexpr bool isCompressible(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition: combustionfluidsystem.hh:141
static Scalar acentricFactor(int compIdx)
The acentric factor of a component .
Definition: combustionfluidsystem.hh:262
static Scalar criticalMolarVolume(int compIdx)
Molar volume of a component at the critical point .
Definition: combustionfluidsystem.hh:252
static Scalar fugacityCoefficient(const FluidState &fluidState, int phaseIdx, int compIdx)
Calculates the fugacity coefficient of an individual component in a fluid phase.
Definition: combustionfluidsystem.hh:426
static Scalar diffusionCoefficient(const FluidState &fluidState, int phaseIdx, int compIdx)
Calculates the molecular diffusion coefficient for a component in a fluid phase .
Definition: combustionfluidsystem.hh:459
static Scalar binaryDiffusionCoefficient(const FluidState &fluidState, int phaseIdx, int compIIdx, int compJIdx)
Given a phase's composition, temperature and pressure, returns the binary diffusion coefficient for ...
Definition: combustionfluidsystem.hh:478
static Scalar criticalTemperature(int compIdx)
Critical temperature of a component .
Definition: combustionfluidsystem.hh:220
static Scalar enthalpy(const FluidState &fluidState, int phaseIdx)
Calculates specific enthalpy .
Definition: combustionfluidsystem.hh:495
static constexpr int phase1Idx
Definition: combustionfluidsystem.hh:72
static constexpr int N2Idx
Definition: combustionfluidsystem.hh:175
static Scalar thermalConductivity(const FluidState &fluidState, const int phaseIdx)
Thermal conductivity of a fluid phase .
Definition: combustionfluidsystem.hh:525
static constexpr int wPhaseIdx
Definition: combustionfluidsystem.hh:69
static Scalar molarMass(int compIdx)
Returns the molar mass of a component in .
Definition: combustionfluidsystem.hh:204
static constexpr bool isGas(int phaseIdx)
Returns whether a phase is gaseous.
Definition: combustionfluidsystem.hh:103
static Scalar density(const FluidState &fluidState, int phaseIdx)
Calculates the density of a fluid phase.
Definition: combustionfluidsystem.hh:318
static constexpr int wCompIdx
Definition: combustionfluidsystem.hh:77
static std::string componentName(int compIdx)
Returns the human readable name of a component.
Definition: combustionfluidsystem.hh:188
static constexpr int numComponents
Number of components in the fluid system.
Definition: combustionfluidsystem.hh:172
static Scalar vaporTemperature(const FluidState &fluidState, const unsigned int phaseIdx)
Calculates the temperature of vapor at a given pressure on the vapor pressure curve.
Definition: combustionfluidsystem.hh:388
static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp, Scalar pressMin, Scalar pressMax, unsigned nPress)
Initializes the fluid system's static parameters using problem specific temperature and pressure rang...
Definition: combustionfluidsystem.hh:304
static Scalar molarDensity(const FluidState &fluidState, int phaseIdx)
The molar density of a fluid phase in .
Definition: combustionfluidsystem.hh:343
static bool isIdealGas(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition: combustionfluidsystem.hh:157
static Scalar viscosity(const FluidState &fluidState, int phaseIdx)
Calculates the dynamic viscosity of a fluid phase .
Definition: combustionfluidsystem.hh:364
Fluid system base class.