GCC Code Coverage Report


Directory: ../../../builds/dumux-repositories/
File: /builds/dumux-repositories/dumux/dumux/material/fluidsystems/brine.hh
Date: 2024-09-21 20:52:54
Exec Total Coverage
Lines: 91 101 90.1%
Functions: 55 64 85.9%
Branches: 45 164 27.4%

Line Branch Exec Source
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 //
7 /*!
8 * \file
9 * \ingroup FluidSystems
10 * \brief A fluid system for brine, i.e. H2O with dissolved NaCl.
11 */
12 #ifndef DUMUX_BRINE_FLUID_SYSTEM_HH
13 #define DUMUX_BRINE_FLUID_SYSTEM_HH
14
15 #include <dune/common/math.hh>
16
17 #include <dumux/common/exceptions.hh>
18 #include <dumux/io/name.hh>
19 #include <dumux/material/fluidsystems/base.hh>
20 #include <dumux/material/constants.hh>
21 #include <dumux/material/components/h2o.hh>
22 #include <dumux/material/components/nacl.hh>
23 #include <dumux/material/components/tabulatedcomponent.hh>
24
25 namespace Dumux {
26 namespace FluidSystems {
27
28 /*!
29 * \ingroup FluidSystems
30 * \brief A compositional single phase fluid system consisting of
31 * two components, which are H2O and NaCl.
32 */
33 template< class Scalar, class H2OType = Components::TabulatedComponent<Dumux::Components::H2O<Scalar>> >
34 class Brine : public Base< Scalar, Brine<Scalar, H2OType>>
35 {
36 using ThisType = Brine<Scalar, H2OType>;
37
38 public:
39 //! export the involved components
40 using H2O = H2OType;
41 using NaCl = Dumux::Components::NaCl<Scalar>;
42
43 static const int numPhases = 1; //!< Number of phases in the fluid system
44 static const int numComponents = 2; //!< Number of components in the fluid system (H2O, NaCl)
45
46 static constexpr int phase0Idx = 0; //!< Index of the first (and only) phase
47 static constexpr int liquidPhaseIdx = phase0Idx; //!< The one considered phase is liquid
48
49 static constexpr int H2OIdx = 0; //!< index of the water component
50 static constexpr int NaClIdx = 1; //!< index of the NaCl component
51 static constexpr int comp0Idx = H2OIdx; //!< index of the first component
52 static constexpr int comp1Idx = NaClIdx; //!< index of the second component
53
54 /*!
55 * \brief Return the human readable name of the phase
56 * \param phaseIdx The index of the fluid phase to consider
57 */
58 39 static const std::string phaseName(int phaseIdx = liquidPhaseIdx)
59 {
60
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 39 times.
39 assert(phaseIdx == liquidPhaseIdx);
61 39 return IOName::liquidPhase();
62 }
63
64 /*!
65 * \brief Returns whether the fluids are miscible
66 * \note There is only one phase, so miscibility makes no sense here
67 */
68 static constexpr bool isMiscible()
69 {
70 return false;
71 }
72
73 /*!
74 * \brief Return whether a phase is gaseous
75 * \param phaseIdx The index of the fluid phase to consider
76 */
77 static constexpr bool isGas(int phaseIdx = liquidPhaseIdx)
78 {
79 assert(phaseIdx == liquidPhaseIdx);
80 return false;
81 }
82
83 /*!
84 * \brief Returns true if and only if a fluid phase is assumed to
85 * be an ideal mixture.
86 *
87 * We define an ideal mixture as a fluid phase where the fugacity
88 * coefficients of all components times the pressure of the phase
89 * are independent on the fluid composition. This assumption is true
90 * if Henry's law and Raoult's law apply. If you are unsure what
91 * this function should return, it is safe to return false. The
92 * only damage done will be (slightly) increased computation times
93 * in some cases.
94 *
95 * \param phaseIdx The index of the fluid phase to consider
96 */
97 static bool isIdealMixture(int phaseIdx = liquidPhaseIdx)
98 {
99 assert(phaseIdx == liquidPhaseIdx);
100 return true;
101 }
102
103 /*!
104 * \brief Returns true if and only if a fluid phase is assumed to
105 * be compressible.
106 *
107 * Compressible means that the partial derivative of the density
108 * to the fluid pressure is always larger than zero.
109 *
110 * \param phaseIdx The index of the fluid phase to consider
111 */
112 static bool isCompressible(int phaseIdx = liquidPhaseIdx)
113 {
114 assert(phaseIdx == liquidPhaseIdx);
115 return H2O::liquidIsCompressible();
116 }
117
118 /*!
119 * \brief Returns true if and only if a fluid phase is assumed to
120 * be an ideal gas.
121 *
122 * \param phaseIdx The index of the fluid phase to consider
123 */
124
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
135 /*!
136 * \brief Return the human readable name of a component
137 * \param compIdx The index of the component to consider
138 */
139 16 static std::string componentName(int compIdx)
140 {
141
2/3
✓ Branch 0 taken 7 times.
✓ Branch 1 taken 3 times.
✗ Branch 2 not taken.
16 switch (compIdx)
142 {
143 12 case H2OIdx: return H2O::name();
144 4 case NaClIdx: return NaCl::name();
145 };
146 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << compIdx);
147 }
148
149 /*!
150 * \brief Return the molar mass of a component in \f$\mathrm{[kg/mol]}\f$.
151 * \param compIdx The index of the component to consider
152 */
153 14699570 static Scalar molarMass(int compIdx)
154 {
155
2/3
✓ Branch 0 taken 8258210 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 6441327 times.
14699570 switch (compIdx)
156 {
157 case H2OIdx: return H2O::molarMass();
158 8258226 case NaClIdx: return NaCl::molarMass();
159 };
160 DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << compIdx);
161 }
162
163 /****************************************
164 * thermodynamic relations
165 ****************************************/
166 /*!
167 * \brief Initialize the fluid system's static parameters generically
168 * \note If a tabulated H2O component is used, we do our best to create
169 * tables that always work.
170 */
171 static void init()
172 {
173
1/2
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
2 init(/*tempMin=*/273.15,
174 /*tempMax=*/623.15,
175 /*numTemp=*/100,
176 /*pMin=*/-10.,
177 /*pMax=*/20e6,
178 /*numP=*/200);
179 }
180
181 /*!
182 * \brief Initialize the fluid system's static parameters using
183 * problem specific temperature and pressure ranges
184 *
185 * \param tempMin The minimum temperature used for tabulation of water [K]
186 * \param tempMax The maximum temperature used for tabulation of water [K]
187 * \param nTemp The number of ticks on the temperature axis of the table of water
188 * \param pressMin The minimum pressure used for tabulation of water [Pa]
189 * \param pressMax The maximum pressure used for tabulation of water [Pa]
190 * \param nPress The number of ticks on the pressure axis of the table of water
191 */
192 1 static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp,
193 Scalar pressMin, Scalar pressMax, unsigned nPress)
194 {
195
196 if (H2O::isTabulated)
197 {
198 1 std::cout << "Initializing tables for the H2O fluid properties ("
199 1 << nTemp*nPress
200 1 << " entries).\n";
201
202 1 H2O::init(tempMin, tempMax, nTemp, pressMin, pressMax, nPress);
203 }
204 1 }
205
206 using Base<Scalar, ThisType>::density;
207 /*!
208 * \brief Return the phase density [kg/m^3].
209 * \note The density is computed as a function of the salt mass fraction, pressure and temperature.
210 * The used function is an empirical relationship fitted to experimental data.
211 * It is presented by Batzle and Wang, 1992 (DOI: 10.1190/1.1443207) \cite batzle1992,
212 * better description and comparison with other approaches in Adams and Bachu, 2002
213 * (DOI: 10.1046/j.1468-8123.2002.00041.x) \cite adams2002.
214 * \param fluidState An arbitrary fluid state
215 * \param phaseIdx The index of the phase for which to compute the density (for compatibility, should be `liquidPhaseIdx`)
216 */
217 template <class FluidState>
218 4944622 static Scalar density(const FluidState& fluidState, int phaseIdx = liquidPhaseIdx)
219 {
220
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 4944608 times.
4944622 assert(phaseIdx == liquidPhaseIdx);
221 4944622 const Scalar temperature = fluidState.temperature(phaseIdx);
222 4944622 const Scalar pressure = fluidState.pressure(phaseIdx);
223 4944622 const Scalar xNaCl = fluidState.massFraction(phaseIdx, NaClIdx);
224
225 using std::max;
226 4944622 const Scalar TempC = temperature - 273.15;
227 4944622 const Scalar pMPa = pressure/1.0E6;
228
2/2
✓ Branch 0 taken 4932564 times.
✓ Branch 1 taken 12044 times.
4944622 const Scalar salinity = max(0.0, xNaCl);
229
230 4944622 const Scalar rhow = H2O::liquidDensity(temperature, pressure);
231 14833866 const Scalar density = rhow + 1000*salinity*(0.668 +
232 9889244 0.44*salinity +
233 14833866 1.0E-6*(300*pMPa -
234 9889244 2400*pMPa*salinity +
235 9889244 TempC*(80.0 +
236 9889244 3*TempC -
237 9889244 3300*salinity -
238 9889244 13*pMPa +
239 4944622 47*pMPa*salinity)));
240
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 4944608 times.
4944622 assert(density > 0.0);
241 4944622 return density;
242 }
243
244 using Base<Scalar, ThisType>::fugacityCoefficient;
245 //! \copydoc Base<Scalar,ThisType>::fugacityCoefficient(const FluidState&,int,int)
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
263 using Base<Scalar, ThisType>::viscosity;
264 /*!
265 * \brief Return the viscosity of the phase.
266 * \note The viscosity is computed as a function of the salt mass fraction and temperature.
267 * The used function is an empirical relationship fitted to experimental data.
268 * It is presented by Batzle and Wang, 1992 (DOI: 10.1190/1.1443207) \cite batzle1992,
269 * better description and comparison with other approaches in Adams and Bachu, 2002 (DOI: 10.1046/j.1468-8123.2002.00041.x) \cite adams2002.
270 * However, the equation given in Adams and Bachu, 2002(DOI: 10.1046/j.1468-8123.2002.00041.x) \cite adams2002
271 * is obviously wrong when compared to the original by Batzle and Wang, 1992 (DOI: 10.1190/1.1443207) \cite batzle1992.
272 * \param fluidState An arbitrary fluid state
273 * \param phaseIdx The index of the phase for which to compute the viscosity (for compatibility, should be `liquidPhaseIdx`)
274 */
275 template <class FluidState>
276 2472311 static Scalar viscosity(const FluidState& fluidState, int phaseIdx = liquidPhaseIdx)
277 {
278
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2472304 times.
2472311 assert(phaseIdx == liquidPhaseIdx);
279 2472311 const Scalar temperature = fluidState.temperature(phaseIdx);
280 2472311 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
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2472304 times.
2472311 const Scalar T = max(temperature, 275.0);
287
2/2
✓ Branch 0 taken 2466282 times.
✓ Branch 1 taken 6022 times.
2472311 const Scalar salinity = max(0.0, xNaCl);
288
289 2472311 const Scalar T_C = T - 273.15;
290 2472311 const Scalar A = ((0.42*power((pow(salinity, 0.8)-0.17), 2)) + 0.045)*pow(T_C, 0.8);
291 2472311 const Scalar mu_brine = 0.1 + (0.333*salinity) + (1.65+(91.9*salinity*salinity*salinity))*exp(-A); // [cP]
292
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 2472304 times.
2472311 assert(mu_brine > 0.0);
293 2472311 return mu_brine/1000.0; // [Pa·s]
294 }
295
296 /*!
297 * \brief Vapor pressure of a component \f$\mathrm{[Pa]}\f$.
298 * \note The vapor pressure of brine decreases with the mole fraction of water in the liquid phase.
299 * This is described by Raoult's law, see Thomas Fetzer's Dissertation Eq. 2.11.
300 * It is also the simplified version of the Kelvin equation, without the influence of the capillary pressure
301 * as we have one-phase flow.
302 *
303 * \param fluidState The fluid state
304 * \param compIdx The index of the component to consider
305 */
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
321 using Base<Scalar, ThisType>::enthalpy;
322 /*!
323 * \brief Given a phase's composition, temperature and pressure,
324 * return its specific enthalpy \f$\mathrm{[J/kg]}\f$.
325 *
326 * \param fluidState The fluid state
327 * \param phaseIdx The index of the phase to consider
328 *
329 * Equations given in:
330 * - Palliser & McKibbin (1998) \cite palliser1998 <BR>
331 * - Michaelides (1981) \cite michaelides1981 <BR>
332 * - Daubert & Danner (1989) \cite daubert1989
333 *
334 */
335 template <class FluidState>
336 2 static Scalar enthalpy(const FluidState& fluidState, int phaseIdx)
337 {
338 //use private enthalpy function to recycle it for the heat capacity calculation
339 2 return enthalpy_(fluidState.pressure(phaseIdx),
340 fluidState.temperature(phaseIdx),
341 2 fluidState.massFraction(phase0Idx, NaClIdx)); /*J/kg*/
342 }
343
344 /*!
345 * \brief Returns the specific enthalpy \f$\mathrm{[J/kg]}\f$ of a component in a specific phase
346 * \param fluidState An arbitrary fluid state
347 * \param phaseIdx The index of the fluid phase to consider
348 * \param componentIdx The index of the component to consider
349 *
350 */
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
370 using Base<Scalar, ThisType>::molarDensity;
371 //! \copydoc Base<Scalar,ThisType>::molarDensity(const FluidState&,int)
372 template <class FluidState>
373 1 static Scalar molarDensity(const FluidState& fluidState, int phaseIdx = liquidPhaseIdx)
374 {
375 4775798 return density(fluidState, phaseIdx)/fluidState.averageMolarMass(phaseIdx);
376 }
377
378 using Base<Scalar, ThisType>::diffusionCoefficient;
379 //! \copydoc Base<Scalar,ThisType>::diffusionCoefficient(const FluidState&,int,int)
380 template <class FluidState>
381 2 static Scalar diffusionCoefficient(const FluidState& fluidState, int phaseIdx, int compIdx)
382 {
383
7/16
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 11 taken 2 times.
✗ Branch 12 not taken.
✓ Branch 15 taken 2 times.
✗ Branch 16 not taken.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✓ Branch 21 taken 2 times.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✗ Branch 24 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
22 DUNE_THROW(Dune::NotImplemented, "FluidSystems::Brine::diffusionCoefficient()");
384 }
385
386 using Base<Scalar, ThisType>::binaryDiffusionCoefficient;
387 /*!
388 * \brief Given a phase's composition, temperature and pressure,
389 * return the binary diffusion coefficient \f$\mathrm{[m^2/s]}\f$ for components
390 * \f$\mathrm{i}\f$ and \f$\mathrm{j}\f$ in this phase.
391 * \param fluidState The fluid state
392 * \param phaseIdx Index of the fluid phase
393 * \param compIIdx Index of the component i
394 * \param compJIdx Index of the component j
395 *
396 * The implemented value for NaCl is for a molar concentration of 2.5984 mol/l and a temperature of 25°C,
397 * see Rard and Miller, 1979 (DOI: 10.1007/BF00648776) \cite Rard1979.
398 * Dependent on the salt concentration the coefficient can vary between 1.47e-9 m^2/s and 1.6e-9 m^2/s, see Rard and Miller, 1979.
399 * It also depends on temperature; values for different temperatures can e.g. found here: Alanis et al., 2000 (DOI: 10.1117/1.602422) \cite Alanis2000.
400 */
401 template <class FluidState>
402 2472329 static Scalar binaryDiffusionCoefficient(const FluidState& fluidState,
403 int phaseIdx,
404 int compIIdx,
405 int compJIdx)
406 {
407
1/2
✓ Branch 0 taken 2472313 times.
✗ Branch 1 not taken.
2472329 if (phaseIdx == liquidPhaseIdx)
408 {
409
2/2
✓ Branch 0 taken 1 times.
✓ Branch 1 taken 2472312 times.
2472329 if (compIIdx > compJIdx)
410 {
411 using std::swap;
412 2 swap(compIIdx, compJIdx);
413 }
414
415
2/2
✓ Branch 0 taken 2472312 times.
✓ Branch 1 taken 1 times.
2472329 if (compJIdx == NaClIdx)
416 2472327 return 1.54e-9;
417 else
418
9/20
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✓ Branch 11 taken 1 times.
✗ Branch 12 not taken.
✓ Branch 16 taken 1 times.
✗ Branch 17 not taken.
✓ Branch 20 taken 1 times.
✗ Branch 21 not taken.
✓ Branch 24 taken 1 times.
✗ Branch 25 not taken.
✓ Branch 27 taken 1 times.
✗ Branch 28 not taken.
✓ Branch 30 taken 1 times.
✗ Branch 31 not taken.
✓ Branch 32 taken 1 times.
✗ Branch 33 not taken.
✗ Branch 35 not taken.
✓ Branch 36 taken 1 times.
✗ Branch 38 not taken.
✗ Branch 39 not taken.
28 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
426 using Base<Scalar, ThisType>::thermalConductivity;
427 /*!
428 * \brief Thermal conductivity of a fluid phase \f$\mathrm{[W/(m K)]}\f$.
429 * \param fluidState An arbitrary fluid state
430 * \param phaseIdx The index of the fluid phase to consider
431 *
432 * The thermal conductivity of brine is implemented based on the contribution of NaCl (\f$\lambda_{brine}\f$/\f$\lambda_{H_2O}\f$) of \cite Yusufova1975 https://link.springer.com/content/pdf/10.1007/BF00867119.pdf, also discussed in \cite Ozbek1980 https://docecity.com/thermal-conductivity-of-aqueous-sodium-chloride-acs-publicat-5f10766acba00.html
433 */
434 template <class FluidState>
435 86330 static Scalar thermalConductivity(const FluidState& fluidState, int phaseIdx)
436 {
437
1/2
✓ Branch 0 taken 86323 times.
✗ Branch 1 not taken.
86330 if (phaseIdx == liquidPhaseIdx)
438 {
439 86330 Scalar tempC = fluidState.temperature(phaseIdx)-273.15;
440 172646 Scalar m = fluidState.moleFraction(phaseIdx, NaClIdx)/(molarMass(H2OIdx)*(1- fluidState.moleFraction(phaseIdx, NaClIdx))); // molality of NaCl
441 86330 Scalar S = 5844.3 * m / (1000 + 58.443 *m);
442 86330 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 258962 return contribNaClFactor * H2O::liquidThermalConductivity(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx));
444 }
445 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index: " << phaseIdx);
446 }
447
448 using Base<Scalar, ThisType>::heatCapacity;
449 //! \copydoc Base<Scalar,ThisType>::heatCapacity(const FluidState&,int)
450 template <class FluidState>
451 14 static Scalar heatCapacity(const FluidState &fluidState, int phaseIdx)
452 {
453
1/2
✓ Branch 0 taken 7 times.
✗ Branch 1 not taken.
14 if (phaseIdx == liquidPhaseIdx){
454 14 const Scalar eps = fluidState.temperature(phaseIdx)*1e-8;
455 //calculate heat capacity from the difference in enthalpy with temperature at constant pressure.
456 28 return (enthalpy_(fluidState.pressure(phaseIdx),
457 14 fluidState.temperature(phaseIdx) +eps ,
458 fluidState.massFraction(phase0Idx, NaClIdx))
459 14 - enthalpy_(fluidState.pressure(phaseIdx),
460 fluidState.temperature(phaseIdx),
461 14 fluidState.massFraction(phase0Idx, NaClIdx)))/eps; /*J/kg*/
462 }
463 DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
464 }
465
466 private:
467 /*!
468 * \brief Given a phase's composition, temperature and pressure,
469 * return its specific enthalpy \f$\mathrm{[J/kg]}\f$.
470 *
471 * Equations given in:
472 * - Palliser & McKibbin (1998) \cite palliser1998 <BR>
473 * - Michaelides (1981) \cite michaelides1981 <BR>
474 * - Daubert & Danner (1989) \cite daubert1989
475 *
476 */
477 86337 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 86337 const Scalar theta = T - 273.15;
489 86337 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
3/4
✗ Branch 0 not taken.
✓ Branch 1 taken 86316 times.
✓ Branch 2 taken 1520 times.
✓ Branch 3 taken 84796 times.
86337 const Scalar salinity = min(max(xNaCl,0.0), salSat);
495
496 86337 const Scalar hw = H2O::liquidEnthalpy(T, p)/1E3; /* kJ/kg */
497
498 /*component enthalpy of soluted NaCl after DAUBERT and DANNER*/
499 172674 const Scalar h_NaCl = (3.6710E4*T + 0.5*(6.2770E1)*T*T - ((6.6670E-2)/3)*T*T*T
500 86337 + ((2.8000E-5)/4)*(T*T*T*T))/(58.44E3)- 2.045698e+02; /* U [kJ/kg] */
501
502 86337 const Scalar m = (1E3/58.44)*(salinity/(1-salinity));
503
504 using Dune::power;
505 86337 Scalar d_h = 0;
506
2/2
✓ Branch 0 taken 345264 times.
✓ Branch 1 taken 86316 times.
431685 for (int i = 0; i<=3; i++)
507
2/2
✓ Branch 0 taken 1035792 times.
✓ Branch 1 taken 345264 times.
1381392 for (int j=0; j<=2; j++)
508 3108132 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 86337 const Scalar delta_h = (4.184/(1E3 + (58.44 * m)))*d_h;
512
513 /* Enthalpy of brine without any dissolved gas */
514 86337 const Scalar h_ls1 =(1-salinity)*hw + salinity*h_NaCl + salinity*delta_h; /* kJ/kg */
515 86337 return h_ls1*1E3; /*J/kg*/
516 }
517 };
518
519 } // end namespace FluidSystems
520 } // end namespace Dumux
521
522 #endif
523