version 3.10-dev
fluidsystems/brine.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_FLUID_SYSTEM_HH
13#define DUMUX_BRINE_FLUID_SYSTEM_HH
14
15#include <dune/common/math.hh>
16
18#include <dumux/io/name.hh>
24
25namespace Dumux {
26namespace FluidSystems {
27
33template< class Scalar, class H2OType = Components::TabulatedComponent<Dumux::Components::H2O<Scalar>> >
34class Brine : public Base< Scalar, Brine<Scalar, H2OType>>
35{
37
38public:
40 using H2O = H2OType;
42
43 static const int numPhases = 1;
44 static const int numComponents = 2;
45
46 static constexpr int phase0Idx = 0;
47 static constexpr int liquidPhaseIdx = phase0Idx;
48
49 static constexpr int H2OIdx = 0;
50 static constexpr int NaClIdx = 1;
51 static constexpr int comp0Idx = H2OIdx;
52 static constexpr int comp1Idx = NaClIdx;
53
58 static const std::string phaseName(int phaseIdx = liquidPhaseIdx)
59 {
60 assert(phaseIdx == liquidPhaseIdx);
61 return IOName::liquidPhase();
62 }
63
68 static constexpr bool isMiscible()
69 {
70 return false;
71 }
72
77 static constexpr bool isGas(int phaseIdx = liquidPhaseIdx)
78 {
79 assert(phaseIdx == liquidPhaseIdx);
80 return false;
81 }
82
97 static bool isIdealMixture(int phaseIdx = liquidPhaseIdx)
98 {
99 assert(phaseIdx == liquidPhaseIdx);
100 return true;
101 }
102
112 static bool isCompressible(int phaseIdx = liquidPhaseIdx)
113 {
114 assert(phaseIdx == liquidPhaseIdx);
115 return H2O::liquidIsCompressible();
116 }
117
125 static bool isIdealGas(int phaseIdx = liquidPhaseIdx)
126 {
127 assert(phaseIdx == liquidPhaseIdx);
128 return false; /*we're a liquid!*/
129 }
130
131 /****************************************
132 * Component related static parameters
133 ****************************************/
134
139 static std::string componentName(int compIdx)
140 {
141 switch (compIdx)
142 {
143 case H2OIdx: return H2O::name();
144 case NaClIdx: return NaCl::name();
145 };
146 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << compIdx);
147 }
148
153 static Scalar molarMass(int compIdx)
154 {
155 switch (compIdx)
156 {
157 case H2OIdx: return H2O::molarMass();
158 case NaClIdx: return NaCl::molarMass();
159 };
160 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << compIdx);
161 }
162
163 /****************************************
164 * thermodynamic relations
165 ****************************************/
171 static void init()
172 {
173 init(/*tempMin=*/273.15,
174 /*tempMax=*/623.15,
175 /*numTemp=*/100,
176 /*pMin=*/-10.,
177 /*pMax=*/20e6,
178 /*numP=*/200);
179 }
180
192 static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp,
193 Scalar pressMin, Scalar pressMax, unsigned nPress)
194 {
195
196 if (H2O::isTabulated)
197 {
198 std::cout << "Initializing tables for the H2O fluid properties ("
199 << nTemp*nPress
200 << " entries).\n";
201
202 H2O::init(tempMin, tempMax, nTemp, pressMin, pressMax, nPress);
203 }
204 }
205
217 template <class FluidState>
218 static Scalar density(const FluidState& fluidState, int phaseIdx = liquidPhaseIdx)
219 {
220 assert(phaseIdx == liquidPhaseIdx);
221 const Scalar temperature = fluidState.temperature(phaseIdx);
222 const Scalar pressure = fluidState.pressure(phaseIdx);
223 const Scalar xNaCl = fluidState.massFraction(phaseIdx, NaClIdx);
224
225 using std::max;
226 const Scalar TempC = temperature - 273.15;
227 const Scalar pMPa = pressure/1.0E6;
228 const Scalar salinity = max(0.0, xNaCl);
229
230 const Scalar rhow = H2O::liquidDensity(temperature, pressure);
231 const Scalar density = rhow + 1000*salinity*(0.668 +
232 0.44*salinity +
233 1.0E-6*(300*pMPa -
234 2400*pMPa*salinity +
235 TempC*(80.0 +
236 3*TempC -
237 3300*salinity -
238 13*pMPa +
239 47*pMPa*salinity)));
240 assert(density > 0.0);
241 return density;
242 }
243
246 template <class FluidState>
247 static Scalar fugacityCoefficient(const FluidState &fluidState,
248 int phaseIdx,
249 int compIdx)
250 {
251 assert(0 <= phaseIdx && phaseIdx < numPhases);
252 assert(0 <= compIdx && compIdx < numComponents);
253
254 if (phaseIdx == compIdx)
255 // We could calculate the real fugacity coefficient of
256 // the component in the fluid. Probably that's not worth
257 // the effort, since the fugacity coefficient of the other
258 // component is infinite anyway...
259 return 1.0;
260 return std::numeric_limits<Scalar>::infinity();
261 }
262
275 template <class FluidState>
276 static Scalar viscosity(const FluidState& fluidState, int phaseIdx = liquidPhaseIdx)
277 {
278 assert(phaseIdx == liquidPhaseIdx);
279 const Scalar temperature = fluidState.temperature(phaseIdx);
280 const Scalar xNaCl = fluidState.massFraction(phaseIdx, NaClIdx);
281
282 using std::pow;
283 using Dune::power;
284 using std::exp;
285 using std::max;
286 const Scalar T = max(temperature, 275.0);
287 const Scalar salinity = max(0.0, xNaCl);
288
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); // [cP]
292 assert(mu_brine > 0.0);
293 return mu_brine/1000.0; // [Pa·s]
294 }
295
306 template <class FluidState>
307 static Scalar vaporPressure(const FluidState& fluidState, int compIdx)
308 {
309 if (compIdx == H2OIdx)
310 {
311 const Scalar temperature = fluidState.temperature(H2OIdx);
312 // Raoult's law, see Thomas Fetzer's Dissertation Eq. 2.11.
313 return H2O::vaporPressure(temperature)*fluidState.moleFraction(phase0Idx, H2OIdx);
314 }
315 else if (compIdx == NaClIdx)
316 DUNE_THROW(Dune::NotImplemented, "NaCl::vaporPressure(t)");
317 else
318 DUNE_THROW(Dune::NotImplemented, "Invalid component index " << compIdx);
319 }
320
335 template <class FluidState>
336 static Scalar enthalpy(const FluidState& fluidState, int phaseIdx)
337 {
338 //use private enthalpy function to recycle it for the heat capacity calculation
339 return enthalpy_(fluidState.pressure(phaseIdx),
340 fluidState.temperature(phaseIdx),
341 fluidState.massFraction(phase0Idx, NaClIdx)); /*J/kg*/
342 }
343
351 template <class FluidState>
352 static Scalar componentEnthalpy(const FluidState &fluidState,
353 int phaseIdx,
354 int componentIdx)
355 {
356 const Scalar T = fluidState.temperature(liquidPhaseIdx);
357 const Scalar p = fluidState.pressure(liquidPhaseIdx);
358
359 if (phaseIdx == liquidPhaseIdx)
360 {
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);
366 }
367 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
368 }
369
372 template <class FluidState>
373 static Scalar molarDensity(const FluidState& fluidState, int phaseIdx = liquidPhaseIdx)
374 {
375 return density(fluidState, phaseIdx)/fluidState.averageMolarMass(phaseIdx);
376 }
377
380 template <class FluidState>
381 static Scalar diffusionCoefficient(const FluidState& fluidState, int phaseIdx, int compIdx)
382 {
383 DUNE_THROW(Dune::NotImplemented, "FluidSystems::Brine::diffusionCoefficient()");
384 }
385
401 template <class FluidState>
402 static Scalar binaryDiffusionCoefficient(const FluidState& fluidState,
403 int phaseIdx,
404 int compIIdx,
405 int compJIdx)
406 {
407 if (phaseIdx == liquidPhaseIdx)
408 {
409 if (compIIdx > compJIdx)
410 {
411 using std::swap;
412 swap(compIIdx, compJIdx);
413 }
414
415 if (compJIdx == NaClIdx)
416 return 1.54e-9;
417 else
418 DUNE_THROW(Dune::NotImplemented, "Binary diffusion coefficient of components "
419 << compIIdx << " and " << compJIdx
420 << " in phase " << phaseIdx);
421 }
422
423 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index: " << phaseIdx);
424 }
425
434 template <class FluidState>
435 static Scalar thermalConductivity(const FluidState& fluidState, int phaseIdx)
436 {
437 if (phaseIdx == liquidPhaseIdx)
438 {
439 Scalar tempC = fluidState.temperature(phaseIdx)-273.15;
440 Scalar m = fluidState.moleFraction(phaseIdx, NaClIdx)/(molarMass(H2OIdx)*(1- fluidState.moleFraction(phaseIdx, NaClIdx))); // molality of NaCl
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));
444 }
445 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index: " << phaseIdx);
446 }
447
450 template <class FluidState>
451 static Scalar heatCapacity(const FluidState &fluidState, int phaseIdx)
452 {
453 if (phaseIdx == liquidPhaseIdx){
454 const Scalar eps = fluidState.temperature(phaseIdx)*1e-8;
455 //calculate heat capacity from the difference in enthalpy with temperature at constant pressure.
456 return (enthalpy_(fluidState.pressure(phaseIdx),
457 fluidState.temperature(phaseIdx) +eps ,
458 fluidState.massFraction(phase0Idx, NaClIdx))
459 - enthalpy_(fluidState.pressure(phaseIdx),
460 fluidState.temperature(phaseIdx),
461 fluidState.massFraction(phase0Idx, NaClIdx)))/eps; /*J/kg*/
462 }
463 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
464 }
465
466private:
477 static Scalar enthalpy_(const Scalar p, const Scalar T, const Scalar xNaCl)
478 {
479 /*Numerical coefficients from PALLISER*/
480 static const Scalar f[] = { 2.63500E-1, 7.48368E-6, 1.44611E-6, -3.80860E-10 };
481
482 /*Numerical coefficients from MICHAELIDES for the enthalpy of brine*/
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 } };
487
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;
490
491 /*Regularization*/
492 using std::min;
493 using std::max;
494 const Scalar salinity = min(max(xNaCl,0.0), salSat);
495
496 const Scalar hw = H2O::liquidEnthalpy(T, p)/1E3; /* kJ/kg */
497
498 /*component enthalpy of soluted NaCl after DAUBERT and DANNER*/
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; /* U [kJ/kg] */
501
502 const Scalar m = (1E3/58.44)*(salinity/(1-salinity));
503
504 using Dune::power;
505 Scalar d_h = 0;
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);
509
510 /* heat of dissolution for halite according to Michaelides 1981 */
511 const Scalar delta_h = (4.184/(1E3 + (58.44 * m)))*d_h;
512
513 /* Enthalpy of brine without any dissolved gas */
514 const Scalar h_ls1 =(1-salinity)*hw + salinity*h_NaCl + salinity*delta_h; /* kJ/kg */
515 return h_ls1*1E3; /*J/kg*/
516 }
517};
518
519} // end namespace FluidSystems
520} // end namespace Dumux
521
522#endif
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
Fluid system base class.
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
Definition: adapt.hh:17
Tabulates all thermodynamic properties of a given untabulated chemical species.