version 3.8
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// 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_BRINE_AIR_FLUID_SYSTEM_HH
13#define DUMUX_BRINE_AIR_FLUID_SYSTEM_HH
14
15#include <array>
16#include <cassert>
17#include <iomanip>
18
28
30
31#include <dumux/io/name.hh>
32
33#include "brine.hh"
34
35namespace Dumux {
36namespace FluidSystems {
37
42template<bool fastButSimplifiedRelations = false>
44{
45 static constexpr bool useBrineDensityAsLiquidMixtureDensity() { return fastButSimplifiedRelations;}
46 static constexpr bool useIdealGasDensity() { return fastButSimplifiedRelations; }
47 static constexpr bool useKelvinVaporPressure() { return false; }
48};
49
58template <class Scalar,
60 class Policy = BrineAirDefaultPolicy<>>
62: public Base<Scalar, BrineAir<Scalar, H2Otype, Policy>>
63{
66
67public:
69 using H2O = H2Otype;
72
75
78
81
82 /****************************************
83 * Fluid phase related static parameters
84 ****************************************/
85 static constexpr int numPhases = 2; // one liquid and one gas phase
86 static constexpr int numComponents = 3; // H2O, Air, NaCl
87
88 static constexpr int liquidPhaseIdx = 0; // index of the liquid phase
89 static constexpr int gasPhaseIdx = 1; // index of the gas phase
90
91 static constexpr int phase0Idx = liquidPhaseIdx; // index of the first phase
92 static constexpr int phase1Idx = gasPhaseIdx; // index of the second phase
93
94 // export component indices to indicate the main component
95 // of the corresponding phase at atmospheric pressure 1 bar
96 // and room temperature 20°C:
97 static constexpr int H2OIdx = 0;
98 static constexpr int AirIdx = 1;
99 static constexpr int NaClIdx = 2;
100 static constexpr int comp0Idx = H2OIdx;
101 static constexpr int comp1Idx = AirIdx;
102 static constexpr int comp2Idx = NaClIdx;
103
104private:
105 struct BrineAdapterPolicy
106 {
107 using FluidSystem = Brine;
108
109 static constexpr int phaseIdx(int brinePhaseIdx) { return liquidPhaseIdx; }
110 static constexpr int compIdx(int brineCompIdx)
111 {
112 switch (brineCompIdx)
113 {
114 case Brine::H2OIdx: return H2OIdx;
115 case Brine::NaClIdx: return NaClIdx;
116 default: return 0; // this will never be reached, only needed to suppress compiler warning
117 }
118 }
119 };
120
121 template<class FluidState>
123
124public:
125
126 /****************************************
127 * phase related static parameters
128 ****************************************/
129
134 static std::string phaseName(int phaseIdx)
135 {
136 assert(0 <= phaseIdx && phaseIdx < numPhases);
137 switch (phaseIdx)
138 {
139 case liquidPhaseIdx: return IOName::liquidPhase();
140 case gasPhaseIdx: return IOName::gaseousPhase();
141 }
142 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
143 }
144
148 static constexpr bool isMiscible()
149 { return true; }
150
155 static constexpr bool isGas(int phaseIdx)
156 {
157 assert(0 <= phaseIdx && phaseIdx < numPhases);
158 return phaseIdx == gasPhaseIdx;
159 }
160
175 static bool isIdealMixture(int phaseIdx)
176 {
177 assert(0 <= phaseIdx && phaseIdx < numPhases);
178 // we assume Henry's and Raoult's laws for the water phase and
179 // and no interaction between gas molecules of different
180 // components, so all phases are ideal mixtures!
181 return true;
182 }
183
193 static constexpr bool isCompressible(int phaseIdx)
194 {
195 assert(0 <= phaseIdx && phaseIdx < numPhases);
196 // ideal gases are always compressible
197 if (phaseIdx == gasPhaseIdx)
198 return true;
199 // let brine decide for the liquid phase...
201 }
202
208 static bool isIdealGas(int phaseIdx)
209 {
210 assert(0 <= phaseIdx && phaseIdx < numPhases);
211 // let the fluids decide
212 if (phaseIdx == gasPhaseIdx)
213 return H2O::gasIsIdeal() && Air::gasIsIdeal();
214 return false; // not a gas
215 }
216
221 static constexpr int getMainComponent(int phaseIdx)
222 {
223 assert(0 <= phaseIdx && phaseIdx < numPhases);
224 if (phaseIdx == liquidPhaseIdx)
225 return H2OIdx;
226 else
227 return AirIdx;
228 }
229
230 /****************************************
231 * Component related static parameters
232 ****************************************/
237 static std::string componentName(int compIdx)
238 {
239 assert(0 <= compIdx && compIdx < numComponents);
240 switch (compIdx)
241 {
242 case H2OIdx: return H2O::name();
243 case AirIdx: return Air::name();
244 case NaClIdx: return NaCl::name();
245 }
246 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << compIdx);
247 }
248
253 static Scalar molarMass(int compIdx)
254 {
255 assert(0 <= compIdx && compIdx < numComponents);
256 switch (compIdx)
257 {
258 case H2OIdx: return H2O::molarMass();
259 case AirIdx: return Air::molarMass();
260 case NaClIdx: return NaCl::molarMass();
261 }
262 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << compIdx);
263 }
264
271 template <class FluidState>
272 static Scalar vaporPressure(const FluidState& fluidState, int compIdx)
273 {
274 // The vapor pressure of the water is affected by the
275 // salinity, thus, we forward to the interface of Brine here
276 if (compIdx == H2OIdx)
277 {
278 if (!Policy::useKelvinVaporPressure())
280 else
281 {
282 //additionally taking the influence of the capillary pressure on the vapor pressure, described by Kelvin's equation, into account.
283 const auto t = fluidState.temperature(liquidPhaseIdx);
284 const auto pc = (fluidState.wettingPhase() == (int) liquidPhaseIdx)
285 ? fluidState.pressure(gasPhaseIdx)-fluidState.pressure(liquidPhaseIdx)
286 : fluidState.pressure(liquidPhaseIdx)-fluidState.pressure(gasPhaseIdx);
287 //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
289 *exp( -pc / (Brine::molarDensity(fluidState, liquidPhaseIdx) * (Dumux::Constants<Scalar>::R*t)));
290 }
291 }
292 else if (compIdx == NaClIdx)
293 DUNE_THROW(Dune::NotImplemented, "NaCl::vaporPressure(t)");
294 else
295 DUNE_THROW(Dune::NotImplemented, "Invalid component index " << compIdx);
296 }
297
298 /****************************************
299 * thermodynamic relations
300 ****************************************/
307 static void init()
308 {
309 init(/*tempMin=*/273.15,
310 /*tempMax=*/800.0,
311 /*numTemptempSteps=*/200,
312 /*startPressure=*/-10,
313 /*endPressure=*/20e6,
314 /*pressureSteps=*/200);
315 }
316
328 static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp,
329 Scalar pressMin, Scalar pressMax, unsigned nPress)
330 {
331 std::cout << "The brine-air fluid system was configured with the following policy:\n";
332 std::cout << " - use brine density as liquid mixture density: " << std::boolalpha << Policy::useBrineDensityAsLiquidMixtureDensity() << "\n";
333 std::cout << " - use ideal gas density: " << std::boolalpha << Policy::useIdealGasDensity() << std::endl;
334
335 if (H2O::isTabulated)
336 H2O::init(tempMin, tempMax, nTemp, pressMin, pressMax, nPress);
337 }
338
353 template <class FluidState>
354 static Scalar density(const FluidState &fluidState, int phaseIdx)
355 {
356 assert(0 <= phaseIdx && phaseIdx < numPhases);
357
358 const auto T = fluidState.temperature(phaseIdx);
359 const auto p = fluidState.pressure(phaseIdx);
360
361 if (phaseIdx == liquidPhaseIdx)
362 {
363 // assume pure brine
364 if (Policy::useBrineDensityAsLiquidMixtureDensity())
366
367 // assume one molecule of gas replaces one "brine" molecule
368 else
370 *(H2O::molarMass()*fluidState.moleFraction(liquidPhaseIdx, H2OIdx)
371 + NaCl::molarMass()*fluidState.moleFraction(liquidPhaseIdx, NaClIdx)
372 + Air::molarMass()*fluidState.moleFraction(liquidPhaseIdx, AirIdx));
373 }
374 else if (phaseIdx == phase1Idx)
375 {
376 // for the gas phase assume an ideal gas
377 if (Policy::useIdealGasDensity())
378 return IdealGas::density(fluidState.averageMolarMass(phase1Idx), T, p);
379
380 // if useComplexRelations = true, compute density. NaCl is assumed
381 // not to be present in gas phase, NaCl has only solid interfaces implemented
382 return H2O::gasDensity(T, fluidState.partialPressure(phase1Idx, H2OIdx))
383 + Air::gasDensity(T, fluidState.partialPressure(phase1Idx, AirIdx));
384 }
385 else
386 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
387 }
388
391 template <class FluidState>
392 static Scalar molarDensity(const FluidState& fluidState, int phaseIdx)
393 {
394 if (phaseIdx == liquidPhaseIdx)
396 else if (phaseIdx == phase1Idx)
397 {
398 const Scalar T = fluidState.temperature(phaseIdx);
399
400 // for the gas phase assume an ideal gas
401 if (Policy::useIdealGasDensity())
402 return IdealGas::molarDensity(T, fluidState.pressure(phaseIdx));
403
404 // if useComplexRelations = true, compute density. NaCl is assumed
405 // not to be present in gas phase, NaCl has only solid interfaces implemented
406 return H2O::gasMolarDensity(T, fluidState.partialPressure(phase1Idx, H2OIdx))
407 + Air::gasMolarDensity(T, fluidState.partialPressure(phase1Idx, AirIdx));
408 }
409 else
410 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
411 }
412
424 template <class FluidState>
425 static Scalar viscosity(const FluidState& fluidState, int phaseIdx)
426 {
427 assert(0 <= phaseIdx && phaseIdx < numPhases);
428
429 if (phaseIdx == liquidPhaseIdx)
431 else
432 return Air::gasViscosity(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx));
433 }
434
457 template <class FluidState>
458 static Scalar fugacityCoefficient(const FluidState& fluidState, int phaseIdx, int compIdx)
459 {
460 assert(0 <= phaseIdx && phaseIdx < numPhases);
461 assert(0 <= compIdx && compIdx < numComponents);
462
463 Scalar T = fluidState.temperature(phaseIdx);
464 Scalar p = fluidState.pressure(phaseIdx);
465
466 if (phaseIdx == gasPhaseIdx)
467 return 1.0;
468
469 else if (phaseIdx == liquidPhaseIdx)
470 {
471 // TODO: Should we use the vapor pressure of the mixture (brine) here?
472 // The presence of NaCl lowers the vapor pressure, thus, we would
473 // expect the fugacity coefficient to be lower as well. However,
474 // with the fugacity coefficient being dependent on the salinity,
475 // the equation system for the phase equilibria becomes non-linear
476 // and our constraint solvers assume linear system of equations.
477 if (compIdx == H2OIdx)
478 return H2O::vaporPressure(T)/p;
479
480 else if (compIdx == AirIdx)
481 return BinaryCoeff::H2O_Air::henry(T)/p;
482
483 // we assume nacl always stays in the liquid phase
484 else
485 return 0.0;
486 }
487
488 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
489 }
490
493 template <class FluidState>
494 static Scalar diffusionCoefficient(const FluidState &fluidState,
495 int phaseIdx,
496 int compIdx)
497 {
498 DUNE_THROW(Dune::NotImplemented, "Diffusion coefficients");
499 }
500
503 template <class FluidState>
504 static Scalar binaryDiffusionCoefficient(const FluidState& fluidState,
505 int phaseIdx,
506 int compIIdx,
507 int compJIdx)
508 {
509 assert(0 <= phaseIdx && phaseIdx < numPhases);
510 assert(0 <= compIIdx && compIIdx < numComponents);
511 assert(0 <= compJIdx && compJIdx < numComponents);
512
513 const auto T = fluidState.temperature(phaseIdx);
514 const auto p = fluidState.pressure(phaseIdx);
515
516 if (compIIdx > compJIdx)
517 std::swap(compIIdx, compJIdx);
518
519 if (phaseIdx == liquidPhaseIdx)
520 {
521 if(compIIdx == H2OIdx && compJIdx == AirIdx)
522 return H2O_Air::liquidDiffCoeff(T, p);
523 else if (compIIdx == H2OIdx && compJIdx == NaClIdx)
525 else
526 DUNE_THROW(Dune::NotImplemented, "Binary diffusion coefficient of components "
527 << compIIdx << " and " << compJIdx
528 << " in phase " << phaseIdx);
529 }
530 else if (phaseIdx == gasPhaseIdx)
531 {
532 if (compIIdx == H2OIdx && compJIdx == AirIdx)
533 return H2O_Air::gasDiffCoeff(T, p);
534
535 // NaCl is expected to never be present in the gas phase. we need to
536 // return a diffusion coefficient that does not case numerical problems.
537 // We choose a very small value here.
538 else if (compIIdx == AirIdx && compJIdx == NaClIdx)
539 return 1e-12;
540
541 else
542 DUNE_THROW(Dune::NotImplemented, "Binary diffusion coefficient of components "
543 << compIIdx << " and " << compJIdx
544 << " in phase " << phaseIdx);
545 }
546
547 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
548 }
549
570 template <class FluidState>
571 static Scalar enthalpy(const FluidState& fluidState, int phaseIdx)
572 {
573 assert(0 <= phaseIdx && phaseIdx < numPhases);
574
575 const Scalar T = fluidState.temperature(phaseIdx);
576 const Scalar p = fluidState.pressure(phaseIdx);
577
578 if (phaseIdx == liquidPhaseIdx)
580 else if (phaseIdx == gasPhaseIdx)
581 {
582 // This assumes NaCl not to be present in the gas phase
583 const Scalar XAir = fluidState.massFraction(gasPhaseIdx, AirIdx);
584 const Scalar XH2O = fluidState.massFraction(gasPhaseIdx, H2OIdx);
585 return XH2O*H2O::gasEnthalpy(T, p) + XAir*Air::gasEnthalpy(T, p);
586 }
587
588 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
589 }
590
597 template <class FluidState>
598 static Scalar componentEnthalpy(const FluidState& fluidState, int phaseIdx, int componentIdx)
599 {
600 const Scalar T = fluidState.temperature(gasPhaseIdx);
601 const Scalar p = fluidState.pressure(gasPhaseIdx);
602
603 if (phaseIdx == liquidPhaseIdx)
604 DUNE_THROW(Dune::NotImplemented, "The component enthalpies in the liquid phase are not implemented.");
605
606 else if (phaseIdx == gasPhaseIdx)
607 {
608 if (componentIdx == H2OIdx)
609 return H2O::gasEnthalpy(T, p);
610 else if (componentIdx == AirIdx)
611 return Air::gasEnthalpy(T, p);
612 else if (componentIdx == NaClIdx)
613 DUNE_THROW(Dune::InvalidStateException, "Implementation assumes NaCl not to be present in gas phase");
614 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << componentIdx);
615 }
616 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
617 }
618
631 template <class FluidState>
632 static Scalar thermalConductivity(const FluidState& fluidState, int phaseIdx)
633 {
634 if (phaseIdx == liquidPhaseIdx)
636 // This assumes NaCl not to be present in the gas phase
637 else if (phaseIdx == gasPhaseIdx)
638 return Air::gasThermalConductivity(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx)) * fluidState.massFraction(gasPhaseIdx, AirIdx)
639 + H2O::gasThermalConductivity(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx)) * fluidState.massFraction(gasPhaseIdx, H2OIdx);
640
641 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
642 }
643
655 template <class FluidState>
656 static Scalar heatCapacity(const FluidState &fluidState, int phaseIdx)
657 {
658 const Scalar T = fluidState.temperature(phaseIdx);
659 const Scalar p = fluidState.pressure(phaseIdx);
660
661 if (phaseIdx == liquidPhaseIdx)
663
664 // We assume NaCl not to be present in the gas phase here
665 else if (phaseIdx == gasPhaseIdx)
666 return Air::gasHeatCapacity(T, p)*fluidState.massFraction(gasPhaseIdx, AirIdx)
667 + H2O::gasHeatCapacity(T, p)*fluidState.massFraction(gasPhaseIdx, H2OIdx);
668
669 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
670 }
671};
672
673} // end namespace FluidSystems
674} // end namespace Dumux
675
676#endif
Adapter class for fluid states with different indices.
A simple class for the air fluid properties.
Binary coefficients for water and air.
Definition: h2o_air.hh:25
static Scalar henry(Scalar temperature)
Henry coefficient for air in liquid water.
Definition: h2o_air.hh:36
static Scalar gasDiffCoeff(Scalar temperature, Scalar pressure)
Binary diffusion coefficient for molecular water and air.
Definition: h2o_air.hh:56
static Scalar liquidDiffCoeff(Scalar temperature, Scalar pressure)
Diffusion coefficient for molecular nitrogen in liquid water.
Definition: h2o_air.hh:89
A class for the air fluid properties.
Definition: air.hh:35
static Scalar gasDensity(Scalar temperature, Scalar pressure)
The density of Air at a given pressure and temperature.
Definition: air.hh:73
static constexpr Scalar molarMass()
The molar mass in of Air.
Definition: air.hh:50
static const Scalar gasHeatCapacity(Scalar temperature, Scalar pressure)
Specific isobaric heat capacity of pure air.
Definition: air.hh:301
static Scalar gasViscosity(Scalar temperature, Scalar pressure)
The dynamic viscosity of Air at a given pressure and temperature.
Definition: air.hh:181
static Scalar gasThermalConductivity(Scalar temperature, Scalar pressure)
Thermal conductivity of air.
Definition: air.hh:339
static constexpr bool gasIsIdeal()
Returns true, the gas phase is assumed to be ideal.
Definition: air.hh:97
static Scalar gasEnthalpy(Scalar temperature, Scalar pressure)
Specific enthalpy of Air with 273.15 as basis.
Definition: air.hh:264
static std::string name()
A human readable name for Air.
Definition: air.hh:42
static Scalar gasMolarDensity(Scalar temperature, Scalar pressure)
The molar density of air in , depending on pressure and temperature.
Definition: air.hh:85
A class for the NaCl properties.
Definition: nacl.hh:35
static std::string name()
A human readable name for the NaCl.
Definition: nacl.hh:40
static constexpr Scalar molarMass()
The molar mass of NaCl in .
Definition: nacl.hh:48
Tabulates all thermodynamic properties of a given component.
Definition: tabulatedcomponent.hh:672
A central place for various physical constants occurring in some equations.
Definition: constants.hh:27
Adapter class for fluid states with different indices.
Definition: adapter.hh:32
Fluid system base class.
Definition: fluidsystems/base.hh:33
A compositional two-phase fluid system with a liquid and a gaseous phase and , and (dissolved miner...
Definition: brineair.hh:63
static Scalar componentEnthalpy(const FluidState &fluidState, int phaseIdx, int componentIdx)
Returns the specific enthalpy of a component in a specific phase.
Definition: brineair.hh:598
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:354
static Scalar thermalConductivity(const FluidState &fluidState, int phaseIdx)
Thermal conductivity of a fluid phase .
Definition: brineair.hh:632
static constexpr int comp0Idx
Definition: brineair.hh:100
static constexpr bool isCompressible(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition: brineair.hh:193
static constexpr int liquidPhaseIdx
Definition: brineair.hh:88
static constexpr int AirIdx
Definition: brineair.hh:98
static constexpr int NaClIdx
Definition: brineair.hh:99
static void init()
Initialize the fluid system's static parameters generically.
Definition: brineair.hh:307
static constexpr int numPhases
Definition: brineair.hh:85
static Scalar heatCapacity(const FluidState &fluidState, int phaseIdx)
Specific isobaric heat capacity of a fluid phase. .
Definition: brineair.hh:656
static constexpr int numComponents
Definition: brineair.hh:86
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:504
H2Otype H2O
export the involved components
Definition: brineair.hh:69
static Scalar viscosity(const FluidState &fluidState, int phaseIdx)
Calculate the dynamic viscosity of a fluid phase .
Definition: brineair.hh:425
static constexpr int comp2Idx
Definition: brineair.hh:102
static Scalar enthalpy(const FluidState &fluidState, int phaseIdx)
Given a phase's composition, temperature and pressure, return its specific enthalpy .
Definition: brineair.hh:571
static Scalar molarDensity(const FluidState &fluidState, int phaseIdx)
Calculate the molar density of a fluid phase.
Definition: brineair.hh:392
static constexpr int getMainComponent(int phaseIdx)
Get the main component of a given phase if possible.
Definition: brineair.hh:221
static constexpr bool isGas(int phaseIdx)
Return whether a phase is gaseous.
Definition: brineair.hh:155
static Scalar fugacityCoefficient(const FluidState &fluidState, int phaseIdx, int compIdx)
Returns the fugacity coefficient of a component in a phase.
Definition: brineair.hh:458
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:328
static Scalar molarMass(int compIdx)
Return the molar mass of a component in .
Definition: brineair.hh:253
static Scalar vaporPressure(const FluidState &fluidState, int compIdx)
Vapor pressure of a component .
Definition: brineair.hh:272
static std::string componentName(int compIdx)
Return the human readable name of a component.
Definition: brineair.hh:237
static std::string phaseName(int phaseIdx)
Return the human readable name of a fluid phase.
Definition: brineair.hh:134
static constexpr int comp1Idx
Definition: brineair.hh:101
static constexpr int phase1Idx
Definition: brineair.hh:92
static constexpr int phase0Idx
Definition: brineair.hh:91
static constexpr bool isMiscible()
Returns whether the fluids are miscible.
Definition: brineair.hh:148
static Scalar diffusionCoefficient(const FluidState &fluidState, int phaseIdx, int compIdx)
Calculate the binary molecular diffusion coefficient for a component in a fluid phase .
Definition: brineair.hh:494
static bool isIdealGas(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition: brineair.hh:208
static constexpr int gasPhaseIdx
Definition: brineair.hh:89
static bool isIdealMixture(int phaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition: brineair.hh:175
Dumux::FluidSystems::Brine< Scalar, H2Otype > Brine
export the underlying brine fluid system for the liquid phase
Definition: brineair.hh:74
static constexpr int H2OIdx
Definition: brineair.hh:97
A compositional single phase fluid system consisting of two components, which are H2O and NaCl.
Definition: fluidsystems/brine.hh:35
static Scalar thermalConductivity(const FluidState &fluidState, int phaseIdx)
Thermal conductivity of a fluid phase .
Definition: fluidsystems/brine.hh:435
static Scalar viscosity(const FluidState &fluidState, int phaseIdx=liquidPhaseIdx)
Return the viscosity of the phase.
Definition: fluidsystems/brine.hh:276
static constexpr int H2OIdx
index of the water component
Definition: fluidsystems/brine.hh:49
static Scalar heatCapacity(const FluidState &fluidState, int phaseIdx)
Specific isobaric heat capacity of a fluid phase .
Definition: fluidsystems/brine.hh:451
static Scalar enthalpy(const FluidState &fluidState, int phaseIdx)
Given a phase's composition, temperature and pressure, return its specific enthalpy .
Definition: fluidsystems/brine.hh:336
static Scalar vaporPressure(const FluidState &fluidState, int compIdx)
Vapor pressure of a component .
Definition: fluidsystems/brine.hh:307
static constexpr int NaClIdx
index of the NaCl component
Definition: fluidsystems/brine.hh:50
static Scalar molarDensity(const FluidState &fluidState, int phaseIdx=liquidPhaseIdx)
Calculate the molar density of a fluid phase.
Definition: fluidsystems/brine.hh:373
static constexpr int liquidPhaseIdx
The one considered phase is liquid.
Definition: fluidsystems/brine.hh:47
static bool isCompressible(int phaseIdx=liquidPhaseIdx)
Returns true if and only if a fluid phase is assumed to be compressible.
Definition: fluidsystems/brine.hh:112
static Scalar density(const FluidState &fluidState, int phaseIdx=liquidPhaseIdx)
Return the phase density [kg/m^3].
Definition: fluidsystems/brine.hh:218
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:402
Relations valid for an ideal gas.
Definition: idealgas.hh:25
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:37
static constexpr Scalar molarDensity(Scalar temperature, Scalar pressure)
The molar density of the gas , depending on pressure and temperature.
Definition: idealgas.hh:58
The a parameter cache which does nothing.
Definition: nullparametercache.hh:22
Some exceptions thrown in DuMux
Fluid system base class.
A fluid system for brine, i.e. H2O with dissolved NaCl.
Material properties of pure water .
Binary coefficients for water and air.
Relations valid for an ideal gas.
Material properties of pure salt .
A collection of input/output field names for common physical quantities.
std::string gaseousPhase() noexcept
I/O name of gaseous phase.
Definition: name.hh:111
std::string liquidPhase() noexcept
I/O name of liquid phase.
Definition: name.hh:107
Definition: adapt.hh:17
Policy for the brine-air fluid system.
Definition: brineair.hh:44
static constexpr bool useIdealGasDensity()
Definition: brineair.hh:46
static constexpr bool useKelvinVaporPressure()
Definition: brineair.hh:47
static constexpr bool useBrineDensityAsLiquidMixtureDensity()
Definition: brineair.hh:45
Tabulates all thermodynamic properties of a given untabulated chemical species.