12#ifndef DUMUX_BRINE_FLUID_SYSTEM_HH
13#define DUMUX_BRINE_FLUID_SYSTEM_HH
15#include <dune/common/math.hh>
26namespace FluidSystems {
33template<
class Scalar,
class H2OType = Components::TabulatedComponent<Dumux::Components::H2O<Scalar>> >
34class Brine :
public Base< Scalar, Brine<Scalar, H2OType>>
115 return H2O::liquidIsCompressible();
143 case H2OIdx:
return H2O::name();
146 DUNE_THROW(Dune::InvalidStateException,
"Invalid component index " << compIdx);
157 case H2OIdx:
return H2O::molarMass();
160 DUNE_THROW(Dune::InvalidStateException,
"Invalid component index " << compIdx);
192 static void init(Scalar tempMin, Scalar tempMax,
unsigned nTemp,
193 Scalar pressMin, Scalar pressMax,
unsigned nPress)
196 if (H2O::isTabulated)
198 std::cout <<
"Initializing tables for the H2O fluid properties ("
202 H2O::init(tempMin, tempMax, nTemp, pressMin, pressMax, nPress);
217 template <
class Flu
idState>
221 const Scalar
temperature = fluidState.temperature(phaseIdx);
222 const Scalar
pressure = fluidState.pressure(phaseIdx);
223 const Scalar xNaCl = fluidState.massFraction(phaseIdx,
NaClIdx);
228 const Scalar salinity = max(0.0, xNaCl);
231 const Scalar
density = rhow + 1000*salinity*(0.668 +
246 template <
class Flu
idState>
251 assert(0 <= phaseIdx && phaseIdx <
numPhases);
254 if (phaseIdx == compIdx)
260 return std::numeric_limits<Scalar>::infinity();
275 template <
class Flu
idState>
279 const Scalar
temperature = fluidState.temperature(phaseIdx);
280 const Scalar xNaCl = fluidState.massFraction(phaseIdx,
NaClIdx);
287 const Scalar salinity = max(0.0, xNaCl);
289 const Scalar T_C = T - 273.15;
290 const Scalar A = ((0.42*power((pow(salinity, 0.8)-0.17), 2)) + 0.045)*pow(T_C, 0.8);
291 const Scalar mu_brine = 0.1 + (0.333*salinity) + (1.65+(91.9*salinity*salinity*salinity))*exp(-A);
292 assert(mu_brine > 0.0);
293 return mu_brine/1000.0;
306 template <
class Flu
idState>
316 DUNE_THROW(Dune::NotImplemented,
"NaCl::vaporPressure(t)");
318 DUNE_THROW(Dune::NotImplemented,
"Invalid component index " << compIdx);
335 template <
class Flu
idState>
336 static Scalar
enthalpy(
const FluidState& fluidState,
int phaseIdx)
339 return enthalpy_(fluidState.pressure(phaseIdx),
340 fluidState.temperature(phaseIdx),
351 template <
class Flu
idState>
361 if (componentIdx ==
H2OIdx)
362 return H2O::liquidEnthalpy(T, p);
363 else if (componentIdx ==
NaClIdx)
364 DUNE_THROW(Dune::NotImplemented,
"The component enthalpy for NaCl is not implemented.");
365 DUNE_THROW(Dune::InvalidStateException,
"Invalid component index " << componentIdx);
367 DUNE_THROW(Dune::InvalidStateException,
"Invalid phase index " << phaseIdx);
372 template <
class Flu
idState>
375 return density(fluidState, phaseIdx)/fluidState.averageMolarMass(phaseIdx);
380 template <
class Flu
idState>
383 DUNE_THROW(Dune::NotImplemented,
"FluidSystems::Brine::diffusionCoefficient()");
401 template <
class Flu
idState>
409 if (compIIdx > compJIdx)
412 swap(compIIdx, compJIdx);
418 DUNE_THROW(Dune::NotImplemented,
"Binary diffusion coefficient of components "
419 << compIIdx <<
" and " << compJIdx
420 <<
" in phase " << phaseIdx);
423 DUNE_THROW(Dune::InvalidStateException,
"Invalid phase index: " << phaseIdx);
434 template <
class Flu
idState>
439 Scalar tempC = fluidState.temperature(phaseIdx)-273.15;
441 Scalar S = 5844.3 * m / (1000 + 58.443 *m);
442 Scalar contribNaClFactor = 1.0 - (2.3434e-3 - 7.924e-6*tempC + 3.924e-8*tempC*tempC)*S + (1.06e-5 - 2.0e-8*tempC + 1.2e-10*tempC*tempC)*S*S;
443 return contribNaClFactor * H2O::liquidThermalConductivity(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx));
445 DUNE_THROW(Dune::InvalidStateException,
"Invalid phase index: " << phaseIdx);
450 template <
class Flu
idState>
454 const Scalar eps = fluidState.temperature(phaseIdx)*1e-8;
456 return (enthalpy_(fluidState.pressure(phaseIdx),
457 fluidState.temperature(phaseIdx) +eps ,
459 - enthalpy_(fluidState.pressure(phaseIdx),
460 fluidState.temperature(phaseIdx),
463 DUNE_THROW(Dune::InvalidStateException,
"Invalid phase index " << phaseIdx);
477 static Scalar enthalpy_(
const Scalar p,
const Scalar T,
const Scalar xNaCl)
480 static const Scalar f[] = { 2.63500E-1, 7.48368E-6, 1.44611E-6, -3.80860E-10 };
483 static const Scalar a[4][3] = { { +9633.6, -4080.0, +286.49 },
484 { +166.58, +68.577, -4.6856 },
485 { -0.90963, -0.36524, +0.249667E-1 },
486 { +0.17965E-2, +0.71924E-3, -0.4900E-4 } };
488 const Scalar theta = T - 273.15;
489 const Scalar salSat = f[0] + f[1]*theta + f[2]*theta*theta + f[3]*theta*theta*theta;
494 const Scalar salinity = min(max(xNaCl,0.0), salSat);
496 const Scalar hw = H2O::liquidEnthalpy(T, p)/1E3;
499 const Scalar h_NaCl = (3.6710E4*T + 0.5*(6.2770E1)*T*T - ((6.6670E-2)/3)*T*T*T
500 + ((2.8000E-5)/4)*(T*T*T*T))/(58.44E3)- 2.045698e+02;
502 const Scalar m = (1E3/58.44)*(salinity/(1-salinity));
506 for (
int i = 0; i<=3; i++)
507 for (
int j=0; j<=2; j++)
508 d_h = d_h + a[i][j] * power(theta, i) * power(m, j);
511 const Scalar delta_h = (4.184/(1E3 + (58.44 * m)))*d_h;
514 const Scalar h_ls1 =(1-salinity)*hw + salinity*h_NaCl + salinity*delta_h;
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
Fluid system base class.
Definition: fluidsystems/base.hh:33
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 Scalar componentEnthalpy(const FluidState &fluidState, int phaseIdx, int componentIdx)
Returns the specific enthalpy of a component in a specific phase.
Definition: fluidsystems/brine.hh:352
static bool isIdealGas(int phaseIdx=liquidPhaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal gas.
Definition: fluidsystems/brine.hh:125
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/brine.hh:381
static constexpr int phase0Idx
Index of the first (and only) phase.
Definition: fluidsystems/brine.hh:46
static constexpr bool isGas(int phaseIdx=liquidPhaseIdx)
Return whether a phase is gaseous.
Definition: fluidsystems/brine.hh:77
static constexpr int comp0Idx
index of the first component
Definition: fluidsystems/brine.hh:51
static const int numComponents
Number of components in the fluid system (H2O, NaCl)
Definition: fluidsystems/brine.hh:44
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: fluidsystems/brine.hh:192
static Scalar fugacityCoefficient(const FluidState &fluidState, int phaseIdx, int compIdx)
Calculate the fugacity coefficient of an individual component in a fluid phase.
Definition: fluidsystems/brine.hh:247
static constexpr int comp1Idx
index of the second component
Definition: fluidsystems/brine.hh:52
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 molarMass(int compIdx)
Return the molar mass of a component in .
Definition: fluidsystems/brine.hh:153
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 const std::string phaseName(int phaseIdx=liquidPhaseIdx)
Return the human readable name of the phase.
Definition: fluidsystems/brine.hh:58
static constexpr int NaClIdx
index of the NaCl component
Definition: fluidsystems/brine.hh:50
H2OType H2O
export the involved components
Definition: fluidsystems/brine.hh:40
static constexpr bool isMiscible()
Returns whether the fluids are miscible.
Definition: fluidsystems/brine.hh:68
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 std::string componentName(int compIdx)
Return the human readable name of a component.
Definition: fluidsystems/brine.hh:139
static Scalar density(const FluidState &fluidState, int phaseIdx=liquidPhaseIdx)
Return the phase density [kg/m^3].
Definition: fluidsystems/brine.hh:218
static bool isIdealMixture(int phaseIdx=liquidPhaseIdx)
Returns true if and only if a fluid phase is assumed to be an ideal mixture.
Definition: fluidsystems/brine.hh:97
static const int numPhases
Number of phases in the fluid system.
Definition: fluidsystems/brine.hh:43
static void init()
Initialize the fluid system's static parameters generically.
Definition: fluidsystems/brine.hh:171
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
A central place for various physical constants occurring in some equations.
Some exceptions thrown in DuMux
Material properties of pure water .
Material properties of pure salt .
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 liquidPhase() noexcept
I/O name of liquid phase.
Definition: name.hh:107
std::string pressure(int phaseIdx) noexcept
I/O name of pressure for multiphase systems.
Definition: name.hh:22
Tabulates all thermodynamic properties of a given untabulated chemical species.