3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
brineair.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_BRINE_AIR_FLUID_SYSTEM_HH
25#define DUMUX_BRINE_AIR_FLUID_SYSTEM_HH
26
27#include <array>
28#include <cassert>
29#include <iomanip>
30
40
42
43#include <dumux/io/name.hh>
44
45#include "brine.hh"
46
47namespace Dumux {
48namespace FluidSystems {
49
54template<bool fastButSimplifiedRelations = false>
56{
57 static constexpr bool useBrineDensityAsLiquidMixtureDensity() { return fastButSimplifiedRelations;}
58 static constexpr bool useIdealGasDensity() { return fastButSimplifiedRelations; }
59 static constexpr bool useKelvinVaporPressure() { return false; }
60};
61
70template <class Scalar,
72 class Policy = BrineAirDefaultPolicy<>>
74: public Base<Scalar, BrineAir<Scalar, H2Otype, Policy>>
75{
79
80public:
82 using H2O = H2Otype;
85
88
91
94
95 /****************************************
96 * Fluid phase related static parameters
97 ****************************************/
98 static constexpr int numPhases = 2; // one liquid and one gas phase
99 static constexpr int numComponents = 3; // H2O, Air, NaCl
100
101 static constexpr int liquidPhaseIdx = 0; // index of the liquid phase
102 static constexpr int gasPhaseIdx = 1; // index of the gas phase
103
104 static constexpr int phase0Idx = liquidPhaseIdx; // index of the first phase
105 static constexpr int phase1Idx = gasPhaseIdx; // index of the second phase
106
107 // export component indices to indicate the main component
108 // of the corresponding phase at atmospheric pressure 1 bar
109 // and room temperature 20°C:
110 static constexpr int H2OIdx = 0;
111 static constexpr int AirIdx = 1;
112 static constexpr int NaClIdx = 2;
113 static constexpr int comp0Idx = H2OIdx;
114 static constexpr int comp1Idx = AirIdx;
115 static constexpr int comp2Idx = NaClIdx;
116
117private:
118 struct BrineAdapterPolicy
119 {
120 using FluidSystem = Brine;
121
122 static constexpr int phaseIdx(int brinePhaseIdx) { return liquidPhaseIdx; }
123 static constexpr int compIdx(int brineCompIdx)
124 {
125 switch (brineCompIdx)
126 {
127 case Brine::H2OIdx: return H2OIdx;
128 case Brine::NaClIdx: return NaClIdx;
129 default: return 0; // this will never be reached, only needed to suppress compiler warning
130 }
131 }
132 };
133
134 template<class FluidState>
136
137public:
138
139 /****************************************
140 * phase related static parameters
141 ****************************************/
142
147 static std::string phaseName(int phaseIdx)
148 {
149 assert(0 <= phaseIdx && phaseIdx < numPhases);
150 switch (phaseIdx)
151 {
152 case liquidPhaseIdx: return IOName::liquidPhase();
153 case gasPhaseIdx: return IOName::gaseousPhase();
154 }
155 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
156 }
157
161 static constexpr bool isMiscible()
162 { return true; }
163
168 static constexpr bool isGas(int phaseIdx)
169 {
170 assert(0 <= phaseIdx && phaseIdx < numPhases);
171 return phaseIdx == gasPhaseIdx;
172 }
173
188 static bool isIdealMixture(int phaseIdx)
189 {
190 assert(0 <= phaseIdx && phaseIdx < numPhases);
191 // we assume Henry's and Raoult's laws for the water phase and
192 // and no interaction between gas molecules of different
193 // components, so all phases are ideal mixtures!
194 return true;
195 }
196
206 static constexpr bool isCompressible(int phaseIdx)
207 {
208 assert(0 <= phaseIdx && phaseIdx < numPhases);
209 // ideal gases are always compressible
210 if (phaseIdx == gasPhaseIdx)
211 return true;
212 // let brine decide for the liquid phase...
214 }
215
221 static bool isIdealGas(int phaseIdx)
222 {
223 assert(0 <= phaseIdx && phaseIdx < numPhases);
224 // let the fluids decide
225 if (phaseIdx == gasPhaseIdx)
226 return H2O::gasIsIdeal() && Air::gasIsIdeal();
227 return false; // not a gas
228 }
229
234 static constexpr int getMainComponent(int phaseIdx)
235 {
236 assert(0 <= phaseIdx && phaseIdx < numPhases);
237 if (phaseIdx == liquidPhaseIdx)
238 return H2OIdx;
239 else
240 return AirIdx;
241 }
242
243 /****************************************
244 * Component related static parameters
245 ****************************************/
250 static std::string componentName(int compIdx)
251 {
252 assert(0 <= compIdx && compIdx < numComponents);
253 switch (compIdx)
254 {
255 case H2OIdx: return H2O::name();
256 case AirIdx: return Air::name();
257 case NaClIdx: return NaCl::name();
258 }
259 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << compIdx);
260 }
261
266 static Scalar molarMass(int compIdx)
267 {
268 assert(0 <= compIdx && compIdx < numComponents);
269 switch (compIdx)
270 {
271 case H2OIdx: return H2O::molarMass();
272 case AirIdx: return Air::molarMass();
273 case NaClIdx: return NaCl::molarMass();
274 }
275 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << compIdx);
276 }
277
284 template <class FluidState>
285 static Scalar vaporPressure(const FluidState& fluidState, int compIdx)
286 {
287 // The vapor pressure of the water is affected by the
288 // salinity, thus, we forward to the interface of Brine here
289 if (compIdx == H2OIdx)
290 {
291 if (!Policy::useKelvinVaporPressure())
293 else
294 {
295 //additionally taking the influence of the capillary pressure on the vapor pressure, described by Kelvin's equation, into account.
296 const auto t = fluidState.temperature(liquidPhaseIdx);
297 const auto pc = (fluidState.wettingPhase() == (int) liquidPhaseIdx)
298 ? fluidState.pressure(gasPhaseIdx)-fluidState.pressure(liquidPhaseIdx)
299 : fluidState.pressure(liquidPhaseIdx)-fluidState.pressure(gasPhaseIdx);
300 //influence of the capillary pressure on the vapor pressure, the influence of the salt concentration is already taken into account in the brine vapour pressure
302 *exp( -pc / (Brine::molarDensity(fluidState, liquidPhaseIdx) * (Dumux::Constants<Scalar>::R*t)));
303 }
304 }
305 else if (compIdx == NaClIdx)
306 DUNE_THROW(Dune::NotImplemented, "NaCl::vaporPressure(t)");
307 else
308 DUNE_THROW(Dune::NotImplemented, "Invalid component index " << compIdx);
309 }
310
311 /****************************************
312 * thermodynamic relations
313 ****************************************/
320 static void init()
321 {
322 init(/*tempMin=*/273.15,
323 /*tempMax=*/800.0,
324 /*numTemptempSteps=*/200,
325 /*startPressure=*/-10,
326 /*endPressure=*/20e6,
327 /*pressureSteps=*/200);
328 }
329
341 static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp,
342 Scalar pressMin, Scalar pressMax, unsigned nPress)
343 {
344 std::cout << "The brine-air fluid system was configured with the following policy:\n";
345 std::cout << " - use brine density as liquid mixture density: " << std::boolalpha << Policy::useBrineDensityAsLiquidMixtureDensity() << "\n";
346 std::cout << " - use ideal gas density: " << std::boolalpha << Policy::useIdealGasDensity() << std::endl;
347
348 if (H2O::isTabulated)
349 H2O::init(tempMin, tempMax, nTemp, pressMin, pressMax, nPress);
350 }
351
352 using Base::density;
366 template <class FluidState>
367 static Scalar density(const FluidState &fluidState, int phaseIdx)
368 {
369 assert(0 <= phaseIdx && phaseIdx < numPhases);
370
371 const auto T = fluidState.temperature(phaseIdx);
372 const auto p = fluidState.pressure(phaseIdx);
373
374 if (phaseIdx == liquidPhaseIdx)
375 {
376 // assume pure brine
377 if (Policy::useBrineDensityAsLiquidMixtureDensity())
379
380 // assume one molecule of gas replaces one "brine" molecule
381 else
383 *(H2O::molarMass()*fluidState.moleFraction(liquidPhaseIdx, H2OIdx)
384 + NaCl::molarMass()*fluidState.moleFraction(liquidPhaseIdx, NaClIdx)
385 + Air::molarMass()*fluidState.moleFraction(liquidPhaseIdx, AirIdx));
386 }
387 else if (phaseIdx == phase1Idx)
388 {
389 // for the gas phase assume an ideal gas
390 if (Policy::useIdealGasDensity())
391 return IdealGas::density(fluidState.averageMolarMass(phase1Idx), T, p);
392
393 // if useComplexRelations = true, compute density. NaCl is assumed
394 // not to be present in gas phase, NaCl has only solid interfaces implemented
395 return H2O::gasDensity(T, fluidState.partialPressure(phase1Idx, H2OIdx))
396 + Air::gasDensity(T, fluidState.partialPressure(phase1Idx, AirIdx));
397 }
398 else
399 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
400 }
401
402 using Base::molarDensity;
412 template <class FluidState>
413 static Scalar molarDensity(const FluidState& fluidState, int phaseIdx)
414 {
415 if (phaseIdx == liquidPhaseIdx)
417 else if (phaseIdx == phase1Idx)
418 {
419 const Scalar T = fluidState.temperature(phaseIdx);
420
421 // for the gas phase assume an ideal gas
422 if (Policy::useIdealGasDensity())
423 return IdealGas::molarDensity(T, fluidState.pressure(phaseIdx));
424
425 // if useComplexRelations = true, compute density. NaCl is assumed
426 // not to be present in gas phase, NaCl has only solid interfaces implemented
427 return H2O::gasMolarDensity(T, fluidState.partialPressure(phase1Idx, H2OIdx))
428 + Air::gasMolarDensity(T, fluidState.partialPressure(phase1Idx, AirIdx));
429 }
430 else
431 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
432 }
433
434 using Base::viscosity;
445 template <class FluidState>
446 static Scalar viscosity(const FluidState& fluidState, int phaseIdx)
447 {
448 assert(0 <= phaseIdx && phaseIdx < numPhases);
449
450 if (phaseIdx == liquidPhaseIdx)
452 else
453 return Air::gasViscosity(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx));
454 }
455
478 template <class FluidState>
479 static Scalar fugacityCoefficient(const FluidState& fluidState, int phaseIdx, int compIdx)
480 {
481 assert(0 <= phaseIdx && phaseIdx < numPhases);
482 assert(0 <= compIdx && compIdx < numComponents);
483
484 Scalar T = fluidState.temperature(phaseIdx);
485 Scalar p = fluidState.pressure(phaseIdx);
486
487 if (phaseIdx == gasPhaseIdx)
488 return 1.0;
489
490 else if (phaseIdx == liquidPhaseIdx)
491 {
492 // TODO: Should we use the vapor pressure of the mixture (brine) here?
493 // The presence of NaCl lowers the vapor pressure, thus, we would
494 // expect the fugacity coefficient to be lower as well. However,
495 // with the fugacity coefficient being dependent on the salinity,
496 // the equation system for the phase equilibria becomes non-linear
497 // and our constraint solvers assume linear system of equations.
498 if (compIdx == H2OIdx)
499 return H2O::vaporPressure(T)/p;
500
501 else if (compIdx == AirIdx)
502 return BinaryCoeff::H2O_Air::henry(T)/p;
503
504 // we assume nacl always stays in the liquid phase
505 else
506 return 0.0;
507 }
508
509 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
510 }
511
513 template <class FluidState>
514 static Scalar diffusionCoefficient(const FluidState &fluidState,
515 int phaseIdx,
516 int compIdx)
517 {
518 DUNE_THROW(Dune::NotImplemented, "Diffusion coefficients");
519 }
520
532 template <class FluidState>
533 static Scalar binaryDiffusionCoefficient(const FluidState& fluidState,
534 int phaseIdx,
535 int compIIdx,
536 int compJIdx)
537 {
538 assert(0 <= phaseIdx && phaseIdx < numPhases);
539 assert(0 <= compIIdx && compIIdx < numComponents);
540 assert(0 <= compJIdx && compJIdx < numComponents);
541
542 const auto T = fluidState.temperature(phaseIdx);
543 const auto p = fluidState.pressure(phaseIdx);
544
545 if (compIIdx > compJIdx)
546 std::swap(compIIdx, compJIdx);
547
548 if (phaseIdx == liquidPhaseIdx)
549 {
550 if(compIIdx == H2OIdx && compJIdx == AirIdx)
551 return H2O_Air::liquidDiffCoeff(T, p);
552 else if (compIIdx == H2OIdx && compJIdx == NaClIdx)
554 else
555 DUNE_THROW(Dune::NotImplemented, "Binary diffusion coefficient of components "
556 << compIIdx << " and " << compJIdx
557 << " in phase " << phaseIdx);
558 }
559 else if (phaseIdx == gasPhaseIdx)
560 {
561 if (compIIdx == H2OIdx && compJIdx == AirIdx)
562 return H2O_Air::gasDiffCoeff(T, p);
563
564 // NaCl is expected to never be present in the gas phase. we need to
565 // return a diffusion coefficient that does not case numerical problems.
566 // We choose a very small value here.
567 else if (compIIdx == AirIdx && compJIdx == NaClIdx)
568 return 1e-12;
569
570 else
571 DUNE_THROW(Dune::NotImplemented, "Binary diffusion coefficient of components "
572 << compIIdx << " and " << compJIdx
573 << " in phase " << phaseIdx);
574 }
575
576 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
577 }
578
579 using Base::enthalpy;
599 template <class FluidState>
600 static Scalar enthalpy(const FluidState& fluidState, int phaseIdx)
601 {
602 assert(0 <= phaseIdx && phaseIdx < numPhases);
603
604 const Scalar T = fluidState.temperature(phaseIdx);
605 const Scalar p = fluidState.pressure(phaseIdx);
606
607 if (phaseIdx == liquidPhaseIdx)
609 else if (phaseIdx == gasPhaseIdx)
610 {
611 // This assumes NaCl not to be present in the gas phase
612 const Scalar XAir = fluidState.massFraction(gasPhaseIdx, AirIdx);
613 const Scalar XH2O = fluidState.massFraction(gasPhaseIdx, H2OIdx);
614 return XH2O*H2O::gasEnthalpy(T, p) + XAir*Air::gasEnthalpy(T, p);
615 }
616
617 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
618 }
619
626 template <class FluidState>
627 static Scalar componentEnthalpy(const FluidState& fluidState, int phaseIdx, int componentIdx)
628 {
629 const Scalar T = fluidState.temperature(gasPhaseIdx);
630 const Scalar p = fluidState.pressure(gasPhaseIdx);
631
632 if (phaseIdx == liquidPhaseIdx)
633 DUNE_THROW(Dune::NotImplemented, "The component enthalpies in the liquid phase are not implemented.");
634
635 else if (phaseIdx == gasPhaseIdx)
636 {
637 if (componentIdx == H2OIdx)
638 return H2O::gasEnthalpy(T, p);
639 else if (componentIdx == AirIdx)
640 return Air::gasEnthalpy(T, p);
641 else if (componentIdx == NaClIdx)
642 DUNE_THROW(Dune::InvalidStateException, "Implementation assumes NaCl not to be present in gas phase");
643 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << componentIdx);
644 }
645 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
646 }
647
660 template <class FluidState>
661 static Scalar thermalConductivity(const FluidState& fluidState, int phaseIdx)
662 {
663 if (phaseIdx == liquidPhaseIdx)
665 // This assumes NaCl not to be present in the gas phase
666 else if (phaseIdx == gasPhaseIdx)
667 return Air::gasThermalConductivity(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx)) * fluidState.massFraction(gasPhaseIdx, AirIdx)
668 + H2O::gasThermalConductivity(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx)) * fluidState.massFraction(gasPhaseIdx, H2OIdx);
669
670 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
671 }
672
683 using Base::heatCapacity;
684 template <class FluidState>
685 static Scalar heatCapacity(const FluidState &fluidState, int phaseIdx)
686 {
687 const Scalar T = fluidState.temperature(phaseIdx);
688 const Scalar p = fluidState.pressure(phaseIdx);
689
690 if (phaseIdx == liquidPhaseIdx)
692
693 // We assume NaCl not to be present in the gas phase here
694 else if (phaseIdx == gasPhaseIdx)
695 return Air::gasHeatCapacity(T, p)*fluidState.massFraction(gasPhaseIdx, AirIdx)
696 + H2O::gasHeatCapacity(T, p)*fluidState.massFraction(gasPhaseIdx, H2OIdx);
697
698 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
699 }
700};
701
702} // end namespace FluidSystems
703} // end namespace Dumux
704
705#endif
Some exceptions thrown in DuMux
A collection of input/output field names for common physical quantities.
Binary coefficients for water and air.
Material properties of pure water .
A simple class for the air fluid properties.
Material properties of pure salt .
Tabulates all thermodynamic properties of a given untabulated chemical species.
Adapter class for fluid states with different indices.
Relations valid for an ideal gas.
Definition: adapt.hh:29
std::string gaseousPhase() noexcept
I/O name of gaseous phase.
Definition: name.hh:123
std::string liquidPhase() noexcept
I/O name of liquid phase.
Definition: name.hh:119
Binary coefficients for water and air.
Definition: h2o_air.hh:37
static Scalar henry(Scalar temperature)
Henry coefficient for air in liquid water.
Definition: h2o_air.hh:48
static Scalar gasDiffCoeff(Scalar temperature, Scalar pressure)
Binary diffusion coefficient for molecular water and air.
Definition: h2o_air.hh:68
static Scalar liquidDiffCoeff(Scalar temperature, Scalar pressure)
Diffusion coefficient for molecular nitrogen in liquid water.
Definition: h2o_air.hh:101
A class for the air fluid properties.
Definition: air.hh:47
static Scalar gasDensity(Scalar temperature, Scalar pressure)
The density of Air at a given pressure and temperature.
Definition: air.hh:85
static constexpr Scalar molarMass()
The molar mass in of Air.
Definition: air.hh:62
static const Scalar gasHeatCapacity(Scalar temperature, Scalar pressure)
Specific isobaric heat capacity of pure air.
Definition: air.hh:313
static Scalar gasViscosity(Scalar temperature, Scalar pressure)
The dynamic viscosity of Air at a given pressure and temperature.
Definition: air.hh:193
static Scalar gasThermalConductivity(Scalar temperature, Scalar pressure)
Thermal conductivity of air.
Definition: air.hh:351
static constexpr bool gasIsIdeal()
Returns true, the gas phase is assumed to be ideal.
Definition: air.hh:109
static Scalar gasEnthalpy(Scalar temperature, Scalar pressure)
Specific enthalpy of Air with 273.15 as basis.
Definition: air.hh:276
static std::string name()
A human readable name for Air.
Definition: air.hh:54
static Scalar gasMolarDensity(Scalar temperature, Scalar pressure)
The molar density of air in , depending on pressure and temperature.
Definition: air.hh:97
A class for the NaCl properties.
Definition: nacl.hh:47
static std::string name()
A human readable name for the NaCl.
Definition: nacl.hh:52
static constexpr Scalar molarMass()
The molar mass of NaCl in .
Definition: nacl.hh:60
Tabulates all thermodynamic properties of a given untabulated chemical species.
Definition: tabulatedcomponent.hh:82
A central place for various physical constants occuring in some equations.
Definition: constants.hh:39
Adapter class for fluid states with different indices.
Definition: adapter.hh:44
Fluid system base class.
Definition: fluidsystems/base.hh:45
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
A compositional single phase fluid system consisting of two components, which are H2O and NaCl.
Definition: fluidsystems/brine.hh:47
static Scalar thermalConductivity(const FluidState &fluidState, int phaseIdx)
Thermal conductivity of a fluid phase .
Definition: fluidsystems/brine.hh:462
static Scalar viscosity(const FluidState &fluidState, int phaseIdx=liquidPhaseIdx)
Return the viscosity of the phase.
Definition: fluidsystems/brine.hh:292
static constexpr int H2OIdx
index of the water component
Definition: fluidsystems/brine.hh:62
static Scalar heatCapacity(const FluidState &fluidState, int phaseIdx)
Specific isobaric heat capacity of a fluid phase. .
Definition: fluidsystems/brine.hh:484
static Scalar enthalpy(const FluidState &fluidState, int phaseIdx)
Given a phase's composition, temperature and pressure, return its specific enthalpy .
Definition: fluidsystems/brine.hh:352
static Scalar vaporPressure(const FluidState &fluidState, int compIdx)
Vapor pressure of a component .
Definition: fluidsystems/brine.hh:323
static constexpr int NaClIdx
index of the NaCl component
Definition: fluidsystems/brine.hh:63
static Scalar molarDensity(const FluidState &fluidState, int phaseIdx=liquidPhaseIdx)
The molar density of the fluid phase in .
Definition: fluidsystems/brine.hh:397
static constexpr int liquidPhaseIdx
The one considered phase is liquid.
Definition: fluidsystems/brine.hh:60
static bool isCompressible(int phaseIdx=liquidPhaseIdx)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition: fluidsystems/brine.hh:125
static Scalar density(const FluidState &fluidState, int phaseIdx=liquidPhaseIdx)
Return the phase density [kg/m^3].
Definition: fluidsystems/brine.hh:230
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/brine.hh:429
Policy for the brine-air fluid system.
Definition: brineair.hh:56
static constexpr bool useIdealGasDensity()
Definition: brineair.hh:58
static constexpr bool useKelvinVaporPressure()
Definition: brineair.hh:59
static constexpr bool useBrineDensityAsLiquidMixtureDensity()
Definition: brineair.hh:57
A compositional two-phase fluid system with a liquid and a gaseous phase and , and (dissolved miner...
Definition: brineair.hh:75
static Scalar componentEnthalpy(const FluidState &fluidState, int phaseIdx, int componentIdx)
Returns the specific enthalpy of a component in a specific phase.
Definition: brineair.hh:627
static Scalar density(const FluidState &fluidState, int phaseIdx)
Given a phase's composition, temperature, pressure, and the partial pressures of all components,...
Definition: brineair.hh:367
static Scalar thermalConductivity(const FluidState &fluidState, int phaseIdx)
Thermal conductivity of a fluid phase .
Definition: brineair.hh:661
static constexpr int comp0Idx
Definition: brineair.hh:113
static constexpr bool isCompressible(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition: brineair.hh:206
static constexpr int liquidPhaseIdx
Definition: brineair.hh:101
static constexpr int AirIdx
Definition: brineair.hh:111
static constexpr int NaClIdx
Definition: brineair.hh:112
static void init()
Initialize the fluid system's static parameters generically.
Definition: brineair.hh:320
static constexpr int numPhases
Definition: brineair.hh:98
static Scalar heatCapacity(const FluidState &fluidState, int phaseIdx)
Definition: brineair.hh:685
static constexpr int numComponents
Definition: brineair.hh:99
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: brineair.hh:533
H2Otype H2O
export the involved components
Definition: brineair.hh:82
static Scalar viscosity(const FluidState &fluidState, int phaseIdx)
Calculate the dynamic viscosity of a fluid phase .
Definition: brineair.hh:446
static constexpr int comp2Idx
Definition: brineair.hh:115
static Scalar enthalpy(const FluidState &fluidState, int phaseIdx)
Given a phase's composition, temperature and pressure, return its specific enthalpy .
Definition: brineair.hh:600
static Scalar molarDensity(const FluidState &fluidState, int phaseIdx)
The molar density of a fluid phase in .
Definition: brineair.hh:413
static constexpr int getMainComponent(int phaseIdx)
Get the main component of a given phase if possible.
Definition: brineair.hh:234
static constexpr bool isGas(int phaseIdx)
Return whether a phase is gaseous.
Definition: brineair.hh:168
static Scalar fugacityCoefficient(const FluidState &fluidState, int phaseIdx, int compIdx)
Returns the fugacity coefficient of a component in a phase.
Definition: brineair.hh:479
static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp, Scalar pressMin, Scalar pressMax, unsigned nPress)
Initialize the fluid system's static parameters using problem specific temperature and pressure range...
Definition: brineair.hh:341
static Scalar molarMass(int compIdx)
Return the molar mass of a component in .
Definition: brineair.hh:266
static Scalar vaporPressure(const FluidState &fluidState, int compIdx)
Vapor pressure of a component .
Definition: brineair.hh:285
static std::string componentName(int compIdx)
Return the human readable name of a component.
Definition: brineair.hh:250
static std::string phaseName(int phaseIdx)
Return the human readable name of a fluid phase.
Definition: brineair.hh:147
static constexpr int comp1Idx
Definition: brineair.hh:114
static constexpr int phase1Idx
Definition: brineair.hh:105
static constexpr int phase0Idx
Definition: brineair.hh:104
static constexpr bool isMiscible()
Returns whether the fluids are miscible.
Definition: brineair.hh:161
static Scalar diffusionCoefficient(const FluidState &fluidState, int phaseIdx, int compIdx)
Definition: brineair.hh:514
static bool isIdealGas(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition: brineair.hh:221
static constexpr int gasPhaseIdx
Definition: brineair.hh:102
static bool isIdealMixture(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition: brineair.hh:188
Dumux::FluidSystems::Brine< Scalar, H2Otype > Brine
export the underlying brine fluid system for the liquid phase
Definition: brineair.hh:87
static constexpr int H2OIdx
Definition: brineair.hh:110
The a parameter cache which does nothing.
Definition: nullparametercache.hh:34
Relations valid for an ideal gas.
Definition: idealgas.hh:37
static constexpr Scalar density(Scalar avgMolarMass, Scalar temperature, Scalar pressure)
The density of the gas in , depending on pressure, temperature and average molar mass of the gas.
Definition: idealgas.hh:49
static constexpr Scalar molarDensity(Scalar temperature, Scalar pressure)
The molar density of the gas , depending on pressure and temperature.
Definition: idealgas.hh:70
Fluid system base class.
A fluid system for brine, i.e. H2O with dissolved NaCl.