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-FileCopyrightText: 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::FluidSystems { | ||
26 | |||
27 | /*! | ||
28 | * \ingroup FluidSystems | ||
29 | * \brief A compositional single phase fluid system consisting of | ||
30 | * two components, which are H2O and NaCl. | ||
31 | */ | ||
32 | template< class Scalar, class H2OType = Components::TabulatedComponent<Dumux::Components::H2O<Scalar>> > | ||
33 | class Brine : public Base< Scalar, Brine<Scalar, H2OType>> | ||
34 | { | ||
35 | using ThisType = Brine<Scalar, H2OType>; | ||
36 | |||
37 | public: | ||
38 | //! export the involved components | ||
39 | using H2O = H2OType; | ||
40 | using NaCl = Dumux::Components::NaCl<Scalar>; | ||
41 | |||
42 | static const int numPhases = 1; //!< Number of phases in the fluid system | ||
43 | static const int numComponents = 2; //!< Number of components in the fluid system (H2O, NaCl) | ||
44 | |||
45 | static constexpr int phase0Idx = 0; //!< Index of the first (and only) phase | ||
46 | static constexpr int liquidPhaseIdx = phase0Idx; //!< The one considered phase is liquid | ||
47 | |||
48 | static constexpr int H2OIdx = 0; //!< index of the water component | ||
49 | static constexpr int NaClIdx = 1; //!< index of the NaCl component | ||
50 | static constexpr int comp0Idx = H2OIdx; //!< index of the first component | ||
51 | static constexpr int comp1Idx = NaClIdx; //!< index of the second component | ||
52 | |||
53 | /*! | ||
54 | * \brief Return the human readable name of the phase | ||
55 | * \param phaseIdx The index of the fluid phase to consider | ||
56 | */ | ||
57 | 39 | static const std::string phaseName(int phaseIdx = liquidPhaseIdx) | |
58 | { | ||
59 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 39 times.
|
39 | assert(phaseIdx == liquidPhaseIdx); |
60 | 39 | return IOName::liquidPhase(); | |
61 | } | ||
62 | |||
63 | /*! | ||
64 | * \brief Returns whether the fluids are miscible | ||
65 | * \note There is only one phase, so miscibility makes no sense here | ||
66 | */ | ||
67 | static constexpr bool isMiscible() | ||
68 | { | ||
69 | return false; | ||
70 | } | ||
71 | |||
72 | /*! | ||
73 | * \brief Return whether a phase is gaseous | ||
74 | * \param phaseIdx The index of the fluid phase to consider | ||
75 | */ | ||
76 | static constexpr bool isGas(int phaseIdx = liquidPhaseIdx) | ||
77 | { | ||
78 | assert(phaseIdx == liquidPhaseIdx); | ||
79 | return false; | ||
80 | } | ||
81 | |||
82 | /*! | ||
83 | * \brief Returns true if and only if a fluid phase is assumed to | ||
84 | * be an ideal mixture. | ||
85 | * | ||
86 | * We define an ideal mixture as a fluid phase where the fugacity | ||
87 | * coefficients of all components times the pressure of the phase | ||
88 | * are independent on the fluid composition. This assumption is true | ||
89 | * if Henry's law and Raoult's law apply. If you are unsure what | ||
90 | * this function should return, it is safe to return false. The | ||
91 | * only damage done will be (slightly) increased computation times | ||
92 | * in some cases. | ||
93 | * | ||
94 | * \param phaseIdx The index of the fluid phase to consider | ||
95 | */ | ||
96 | static bool isIdealMixture(int phaseIdx = liquidPhaseIdx) | ||
97 | { | ||
98 | assert(phaseIdx == liquidPhaseIdx); | ||
99 | return true; | ||
100 | } | ||
101 | |||
102 | /*! | ||
103 | * \brief Returns true if and only if a fluid phase is assumed to | ||
104 | * be compressible. | ||
105 | * | ||
106 | * Compressible means that the partial derivative of the density | ||
107 | * to the fluid pressure is always larger than zero. | ||
108 | * | ||
109 | * \param phaseIdx The index of the fluid phase to consider | ||
110 | */ | ||
111 | static bool isCompressible(int phaseIdx = liquidPhaseIdx) | ||
112 | { | ||
113 | assert(phaseIdx == liquidPhaseIdx); | ||
114 | return H2O::liquidIsCompressible(); | ||
115 | } | ||
116 | |||
117 | /*! | ||
118 | * \brief Returns true if and only if a fluid phase is assumed to | ||
119 | * be an ideal gas. | ||
120 | * | ||
121 | * \param phaseIdx The index of the fluid phase to consider | ||
122 | */ | ||
123 | |||
124 | static bool isIdealGas(int phaseIdx = liquidPhaseIdx) | ||
125 | { | ||
126 | assert(phaseIdx == liquidPhaseIdx); | ||
127 | return false; /*we're a liquid!*/ | ||
128 | } | ||
129 | |||
130 | /**************************************** | ||
131 | * Component related static parameters | ||
132 | ****************************************/ | ||
133 | |||
134 | /*! | ||
135 | * \brief Return the human readable name of a component | ||
136 | * \param compIdx The index of the component to consider | ||
137 | */ | ||
138 | 16 | static std::string componentName(int compIdx) | |
139 | { | ||
140 |
2/3✓ Branch 0 taken 7 times.
✓ Branch 1 taken 3 times.
✗ Branch 2 not taken.
|
16 | switch (compIdx) |
141 | { | ||
142 | 12 | case H2OIdx: return H2O::name(); | |
143 | 4 | case NaClIdx: return NaCl::name(); | |
144 | }; | ||
145 | ✗ | DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << compIdx); | |
146 | } | ||
147 | |||
148 | /*! | ||
149 | * \brief Return the molar mass of a component in \f$\mathrm{[kg/mol]}\f$. | ||
150 | * \param compIdx The index of the component to consider | ||
151 | */ | ||
152 | 16231339 | static Scalar molarMass(int compIdx) | |
153 | { | ||
154 |
2/3✓ Branch 0 taken 9023166 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 7119960 times.
|
16143152 | switch (compIdx) |
155 | { | ||
156 | case H2OIdx: return H2O::molarMass(); | ||
157 | 9023182 | case NaClIdx: return NaCl::molarMass(); | |
158 | }; | ||
159 | ✗ | DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << compIdx); | |
160 | } | ||
161 | |||
162 | /**************************************** | ||
163 | * thermodynamic relations | ||
164 | ****************************************/ | ||
165 | /*! | ||
166 | * \brief Initialize the fluid system's static parameters generically | ||
167 | * \note If a tabulated H2O component is used, we do our best to create | ||
168 | * tables that always work. | ||
169 | */ | ||
170 | 1 | static void init() | |
171 | { | ||
172 |
1/2✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
|
1 | init(/*tempMin=*/273.15, |
173 | /*tempMax=*/623.15, | ||
174 | /*numTemp=*/100, | ||
175 | /*pMin=*/-10., | ||
176 | /*pMax=*/20e6, | ||
177 | /*numP=*/200); | ||
178 | 1 | } | |
179 | |||
180 | /*! | ||
181 | * \brief Initialize the fluid system's static parameters using | ||
182 | * problem specific temperature and pressure ranges | ||
183 | * | ||
184 | * \param tempMin The minimum temperature used for tabulation of water [K] | ||
185 | * \param tempMax The maximum temperature used for tabulation of water [K] | ||
186 | * \param nTemp The number of ticks on the temperature axis of the table of water | ||
187 | * \param pressMin The minimum pressure used for tabulation of water [Pa] | ||
188 | * \param pressMax The maximum pressure used for tabulation of water [Pa] | ||
189 | * \param nPress The number of ticks on the pressure axis of the table of water | ||
190 | */ | ||
191 | 1 | static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp, | |
192 | Scalar pressMin, Scalar pressMax, unsigned nPress) | ||
193 | { | ||
194 | |||
195 | if (H2O::isTabulated) | ||
196 | { | ||
197 | 1 | std::cout << "Initializing tables for the H2O fluid properties (" | |
198 | 1 | << nTemp*nPress | |
199 | 1 | << " entries).\n"; | |
200 | |||
201 | 1 | H2O::init(tempMin, tempMax, nTemp, pressMin, pressMax, nPress); | |
202 | } | ||
203 | 1 | } | |
204 | |||
205 | using Base<Scalar, ThisType>::density; | ||
206 | /*! | ||
207 | * \brief Return the phase density [kg/m^3]. | ||
208 | * \note The density is computed as a function of the salt mass fraction, pressure and temperature. | ||
209 | * The used function is an empirical relationship fitted to experimental data. | ||
210 | * It is presented by Batzle and Wang, 1992 (DOI: 10.1190/1.1443207) \cite batzle1992, | ||
211 | * better description and comparison with other approaches in Adams and Bachu, 2002 | ||
212 | * (DOI: 10.1046/j.1468-8123.2002.00041.x) \cite adams2002. | ||
213 | * \param fluidState An arbitrary fluid state | ||
214 | * \param phaseIdx The index of the phase for which to compute the density (for compatibility, should be `liquidPhaseIdx`) | ||
215 | */ | ||
216 | template <class FluidState> | ||
217 | 5360378 | static Scalar density(const FluidState& fluidState, int phaseIdx = liquidPhaseIdx) | |
218 | { | ||
219 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 5360364 times.
|
5360378 | assert(phaseIdx == liquidPhaseIdx); |
220 | 5360378 | const Scalar temperature = fluidState.temperature(phaseIdx); | |
221 | 5360378 | const Scalar pressure = fluidState.pressure(phaseIdx); | |
222 | 5360378 | const Scalar xNaCl = fluidState.massFraction(phaseIdx, NaClIdx); | |
223 | |||
224 | using std::max; | ||
225 | 5360378 | const Scalar TempC = temperature - 273.15; | |
226 | 5360378 | const Scalar pMPa = pressure/1.0E6; | |
227 |
2/2✓ Branch 0 taken 5348320 times.
✓ Branch 1 taken 12044 times.
|
5360378 | const Scalar salinity = max(0.0, xNaCl); |
228 | |||
229 | 5360378 | const Scalar rhow = H2O::liquidDensity(temperature, pressure); | |
230 | 5360378 | const Scalar density = rhow + 1000*salinity*(0.668 + | |
231 | 5360378 | 0.44*salinity + | |
232 | 5360378 | 1.0E-6*(300*pMPa - | |
233 | 5360378 | 2400*pMPa*salinity + | |
234 | 5360378 | TempC*(80.0 + | |
235 | 5360378 | 3*TempC - | |
236 | 5360378 | 3300*salinity - | |
237 | 5360378 | 13*pMPa + | |
238 | 5360378 | 47*pMPa*salinity))); | |
239 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 5360364 times.
|
5360378 | assert(density > 0.0); |
240 | 5360378 | return density; | |
241 | } | ||
242 | |||
243 | using Base<Scalar, ThisType>::fugacityCoefficient; | ||
244 | //! \copydoc Base<Scalar,ThisType>::fugacityCoefficient(const FluidState&,int,int) | ||
245 | template <class FluidState> | ||
246 | 2 | static Scalar fugacityCoefficient(const FluidState &fluidState, | |
247 | int phaseIdx, | ||
248 | int compIdx) | ||
249 | { | ||
250 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
2 | assert(0 <= phaseIdx && phaseIdx < numPhases); |
251 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
2 | assert(0 <= compIdx && compIdx < numComponents); |
252 | |||
253 |
2/2✓ Branch 0 taken 1 times.
✓ Branch 1 taken 1 times.
|
2 | if (phaseIdx == compIdx) |
254 | // We could calculate the real fugacity coefficient of | ||
255 | // the component in the fluid. Probably that's not worth | ||
256 | // the effort, since the fugacity coefficient of the other | ||
257 | // component is infinite anyway... | ||
258 | 1 | return 1.0; | |
259 | return std::numeric_limits<Scalar>::infinity(); | ||
260 | } | ||
261 | |||
262 | using Base<Scalar, ThisType>::viscosity; | ||
263 | /*! | ||
264 | * \brief Return the viscosity of the phase. | ||
265 | * \note The viscosity is computed as a function of the salt mass fraction and temperature. | ||
266 | * The used function is an empirical relationship fitted to experimental data. | ||
267 | * It is presented by Batzle and Wang, 1992 (DOI: 10.1190/1.1443207) \cite batzle1992, | ||
268 | * better description and comparison with other approaches in Adams and Bachu, 2002 (DOI: 10.1046/j.1468-8123.2002.00041.x) \cite adams2002. | ||
269 | * However, the equation given in Adams and Bachu, 2002(DOI: 10.1046/j.1468-8123.2002.00041.x) \cite adams2002 | ||
270 | * is obviously wrong when compared to the original by Batzle and Wang, 1992 (DOI: 10.1190/1.1443207) \cite batzle1992. | ||
271 | * \param fluidState An arbitrary fluid state | ||
272 | * \param phaseIdx The index of the phase for which to compute the viscosity (for compatibility, should be `liquidPhaseIdx`) | ||
273 | */ | ||
274 | template <class FluidState> | ||
275 | 2680189 | static Scalar viscosity(const FluidState& fluidState, int phaseIdx = liquidPhaseIdx) | |
276 | { | ||
277 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2680182 times.
|
2680189 | assert(phaseIdx == liquidPhaseIdx); |
278 | 2680189 | const Scalar temperature = fluidState.temperature(phaseIdx); | |
279 | 2680189 | const Scalar xNaCl = fluidState.massFraction(phaseIdx, NaClIdx); | |
280 | |||
281 | using std::pow; | ||
282 | using Dune::power; | ||
283 | using std::exp; | ||
284 | using std::max; | ||
285 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2680182 times.
|
2680189 | const Scalar T = max(temperature, 275.0); |
286 |
2/2✓ Branch 0 taken 2674160 times.
✓ Branch 1 taken 6022 times.
|
2680189 | const Scalar salinity = max(0.0, xNaCl); |
287 | |||
288 | 2680189 | const Scalar T_C = T - 273.15; | |
289 | 2680189 | const Scalar A = ((0.42*power((pow(salinity, 0.8)-0.17), 2)) + 0.045)*pow(T_C, 0.8); | |
290 | 2680189 | const Scalar mu_brine = 0.1 + (0.333*salinity) + (1.65+(91.9*salinity*salinity*salinity))*exp(-A); // [cP] | |
291 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2680182 times.
|
2680189 | assert(mu_brine > 0.0); |
292 | 2680189 | return mu_brine/1000.0; // [Pa·s] | |
293 | } | ||
294 | |||
295 | /*! | ||
296 | * \brief Vapor pressure of a component \f$\mathrm{[Pa]}\f$. | ||
297 | * \note The vapor pressure of brine decreases with the mole fraction of water in the liquid phase. | ||
298 | * This is described by Raoult's law, see Thomas Fetzer's Dissertation Eq. 2.11. | ||
299 | * It is also the simplified version of the Kelvin equation, without the influence of the capillary pressure | ||
300 | * as we have one-phase flow. | ||
301 | * | ||
302 | * \param fluidState The fluid state | ||
303 | * \param compIdx The index of the component to consider | ||
304 | */ | ||
305 | template <class FluidState> | ||
306 | static Scalar vaporPressure(const FluidState& fluidState, int compIdx) | ||
307 | { | ||
308 | if (compIdx == H2OIdx) | ||
309 | { | ||
310 | const Scalar temperature = fluidState.temperature(H2OIdx); | ||
311 | // Raoult's law, see Thomas Fetzer's Dissertation Eq. 2.11. | ||
312 | return H2O::vaporPressure(temperature)*fluidState.moleFraction(phase0Idx, H2OIdx); | ||
313 | } | ||
314 | else if (compIdx == NaClIdx) | ||
315 | DUNE_THROW(Dune::NotImplemented, "NaCl::vaporPressure(t)"); | ||
316 | else | ||
317 | DUNE_THROW(Dune::NotImplemented, "Invalid component index " << compIdx); | ||
318 | } | ||
319 | |||
320 | using Base<Scalar, ThisType>::enthalpy; | ||
321 | /*! | ||
322 | * \brief Given a phase's composition, temperature and pressure, | ||
323 | * return its specific enthalpy \f$\mathrm{[J/kg]}\f$. | ||
324 | * | ||
325 | * \param fluidState The fluid state | ||
326 | * \param phaseIdx The index of the phase to consider | ||
327 | * | ||
328 | * Equations given in: | ||
329 | * - Palliser & McKibbin (1998) \cite palliser1998 <BR> | ||
330 | * - Michaelides (1981) \cite michaelides1981 <BR> | ||
331 | * - Daubert & Danner (1989) \cite daubert1989 | ||
332 | * | ||
333 | */ | ||
334 | template <class FluidState> | ||
335 | 85982 | static Scalar enthalpy(const FluidState& fluidState, int phaseIdx) | |
336 | { | ||
337 | //use private enthalpy function to recycle it for the heat capacity calculation | ||
338 | 85982 | return enthalpy_(fluidState.pressure(phaseIdx), | |
339 | fluidState.temperature(phaseIdx), | ||
340 | 85982 | fluidState.massFraction(phase0Idx, NaClIdx)); /*J/kg*/ | |
341 | } | ||
342 | |||
343 | /*! | ||
344 | * \brief Returns the specific enthalpy \f$\mathrm{[J/kg]}\f$ of a component in a specific phase | ||
345 | * \param fluidState An arbitrary fluid state | ||
346 | * \param phaseIdx The index of the fluid phase to consider | ||
347 | * \param componentIdx The index of the component to consider | ||
348 | * | ||
349 | */ | ||
350 | template <class FluidState> | ||
351 | static Scalar componentEnthalpy(const FluidState &fluidState, | ||
352 | int phaseIdx, | ||
353 | int componentIdx) | ||
354 | { | ||
355 | const Scalar T = fluidState.temperature(liquidPhaseIdx); | ||
356 | const Scalar p = fluidState.pressure(liquidPhaseIdx); | ||
357 | |||
358 | if (phaseIdx == liquidPhaseIdx) | ||
359 | { | ||
360 | if (componentIdx == H2OIdx) | ||
361 | return H2O::liquidEnthalpy(T, p); | ||
362 | else if (componentIdx == NaClIdx) | ||
363 | DUNE_THROW(Dune::NotImplemented, "The component enthalpy for NaCl is not implemented."); | ||
364 | DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << componentIdx); | ||
365 | } | ||
366 | DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx); | ||
367 | } | ||
368 | |||
369 | using Base<Scalar, ThisType>::molarDensity; | ||
370 | //! \copydoc Base<Scalar,ThisType>::molarDensity(const FluidState&,int) | ||
371 | template <class FluidState> | ||
372 | 5191558 | static Scalar molarDensity(const FluidState& fluidState, int phaseIdx = liquidPhaseIdx) | |
373 | { | ||
374 | 5191558 | return density(fluidState, phaseIdx)/fluidState.averageMolarMass(phaseIdx); | |
375 | } | ||
376 | |||
377 | using Base<Scalar, ThisType>::diffusionCoefficient; | ||
378 | //! \copydoc Base<Scalar,ThisType>::diffusionCoefficient(const FluidState&,int,int) | ||
379 | template <class FluidState> | ||
380 | 2 | static Scalar diffusionCoefficient(const FluidState& fluidState, int phaseIdx, int compIdx) | |
381 | { | ||
382 |
10/20✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 2 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✓ Branch 13 taken 2 times.
✗ Branch 14 not taken.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✓ Branch 19 taken 2 times.
✗ Branch 20 not taken.
✓ Branch 22 taken 2 times.
✗ Branch 23 not taken.
✓ Branch 25 taken 2 times.
✗ Branch 26 not taken.
✓ Branch 28 taken 2 times.
✗ Branch 29 not taken.
|
8 | DUNE_THROW(Dune::NotImplemented, "FluidSystems::Brine::diffusionCoefficient()"); |
383 | } | ||
384 | |||
385 | using Base<Scalar, ThisType>::binaryDiffusionCoefficient; | ||
386 | /*! | ||
387 | * \brief Given a phase's composition, temperature and pressure, | ||
388 | * return the binary diffusion coefficient \f$\mathrm{[m^2/s]}\f$ for components | ||
389 | * \f$\mathrm{i}\f$ and \f$\mathrm{j}\f$ in this phase. | ||
390 | * \param fluidState The fluid state | ||
391 | * \param phaseIdx Index of the fluid phase | ||
392 | * \param compIIdx Index of the component i | ||
393 | * \param compJIdx Index of the component j | ||
394 | * | ||
395 | * The implemented value for NaCl is for a molar concentration of 2.5984 mol/l and a temperature of 25°C, | ||
396 | * see Rard and Miller, 1979 (DOI: 10.1007/BF00648776) \cite Rard1979. | ||
397 | * 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. | ||
398 | * 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. | ||
399 | */ | ||
400 | template <class FluidState> | ||
401 | 2680207 | static Scalar binaryDiffusionCoefficient(const FluidState& fluidState, | |
402 | int phaseIdx, | ||
403 | int compIIdx, | ||
404 | int compJIdx) | ||
405 | { | ||
406 |
1/2✓ Branch 0 taken 2680191 times.
✗ Branch 1 not taken.
|
2680207 | if (phaseIdx == liquidPhaseIdx) |
407 | { | ||
408 |
2/2✓ Branch 0 taken 1 times.
✓ Branch 1 taken 2680190 times.
|
2680207 | if (compIIdx > compJIdx) |
409 | { | ||
410 | using std::swap; | ||
411 | 2 | swap(compIIdx, compJIdx); | |
412 | } | ||
413 | |||
414 |
2/2✓ Branch 0 taken 2680190 times.
✓ Branch 1 taken 1 times.
|
2680207 | if (compJIdx == NaClIdx) |
415 | 2680205 | return 1.54e-9; | |
416 | else | ||
417 |
15/30✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 1 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 1 times.
✗ Branch 11 not taken.
✓ Branch 13 taken 1 times.
✗ Branch 14 not taken.
✓ Branch 16 taken 1 times.
✗ Branch 17 not taken.
✓ Branch 19 taken 1 times.
✗ Branch 20 not taken.
✓ Branch 22 taken 1 times.
✗ Branch 23 not taken.
✓ Branch 25 taken 1 times.
✗ Branch 26 not taken.
✓ Branch 28 taken 1 times.
✗ Branch 29 not taken.
✓ Branch 31 taken 1 times.
✗ Branch 32 not taken.
✓ Branch 34 taken 1 times.
✗ Branch 35 not taken.
✓ Branch 37 taken 1 times.
✗ Branch 38 not taken.
✓ Branch 40 taken 1 times.
✗ Branch 41 not taken.
✓ Branch 43 taken 1 times.
✗ Branch 44 not taken.
|
10 | DUNE_THROW(Dune::NotImplemented, "Binary diffusion coefficient of components " |
418 | << compIIdx << " and " << compJIdx | ||
419 | << " in phase " << phaseIdx); | ||
420 | } | ||
421 | |||
422 | ✗ | DUNE_THROW(Dune::InvalidStateException, "Invalid phase index: " << phaseIdx); | |
423 | } | ||
424 | |||
425 | using Base<Scalar, ThisType>::thermalConductivity; | ||
426 | /*! | ||
427 | * \brief Thermal conductivity of a fluid phase \f$\mathrm{[W/(m K)]}\f$. | ||
428 | * \param fluidState An arbitrary fluid state | ||
429 | * \param phaseIdx The index of the fluid phase to consider | ||
430 | * | ||
431 | * 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 | ||
432 | */ | ||
433 | template <class FluidState> | ||
434 | 88194 | static Scalar thermalConductivity(const FluidState& fluidState, int phaseIdx) | |
435 | { | ||
436 |
1/2✓ Branch 0 taken 88187 times.
✗ Branch 1 not taken.
|
88194 | if (phaseIdx == liquidPhaseIdx) |
437 | { | ||
438 | 88194 | Scalar tempC = fluidState.temperature(phaseIdx)-273.15; | |
439 | 88194 | Scalar m = fluidState.moleFraction(phaseIdx, NaClIdx)/(molarMass(H2OIdx)*(1- fluidState.moleFraction(phaseIdx, NaClIdx))); // molality of NaCl | |
440 | 88194 | Scalar S = 5844.3 * m / (1000 + 58.443 *m); | |
441 | 88194 | 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; | |
442 | 88194 | return contribNaClFactor * H2O::liquidThermalConductivity(fluidState.temperature(phaseIdx), fluidState.pressure(phaseIdx)); | |
443 | } | ||
444 | ✗ | DUNE_THROW(Dune::InvalidStateException, "Invalid phase index: " << phaseIdx); | |
445 | } | ||
446 | |||
447 | using Base<Scalar, ThisType>::heatCapacity; | ||
448 | //! \copydoc Base<Scalar,ThisType>::heatCapacity(const FluidState&,int) | ||
449 | template <class FluidState> | ||
450 | 14 | static Scalar heatCapacity(const FluidState &fluidState, int phaseIdx) | |
451 | { | ||
452 |
1/2✓ Branch 0 taken 7 times.
✗ Branch 1 not taken.
|
14 | if (phaseIdx == liquidPhaseIdx){ |
453 | 14 | const Scalar eps = fluidState.temperature(phaseIdx)*1e-8; | |
454 | //calculate heat capacity from the difference in enthalpy with temperature at constant pressure. | ||
455 | 28 | return (enthalpy_(fluidState.pressure(phaseIdx), | |
456 | 14 | fluidState.temperature(phaseIdx) +eps , | |
457 | fluidState.massFraction(phase0Idx, NaClIdx)) | ||
458 | 14 | - enthalpy_(fluidState.pressure(phaseIdx), | |
459 | fluidState.temperature(phaseIdx), | ||
460 | 14 | fluidState.massFraction(phase0Idx, NaClIdx)))/eps; /*J/kg*/ | |
461 | } | ||
462 | ✗ | DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx); | |
463 | } | ||
464 | |||
465 | private: | ||
466 | /*! | ||
467 | * \brief Given a phase's composition, temperature and pressure, | ||
468 | * return its specific enthalpy \f$\mathrm{[J/kg]}\f$. | ||
469 | * | ||
470 | * Equations given in: | ||
471 | * - Palliser & McKibbin (1998) \cite palliser1998 <BR> | ||
472 | * - Michaelides (1981) \cite michaelides1981 <BR> | ||
473 | * - Daubert & Danner (1989) \cite daubert1989 | ||
474 | * | ||
475 | */ | ||
476 | 85989 | static Scalar enthalpy_(const Scalar p, const Scalar T, const Scalar xNaCl) | |
477 | { | ||
478 | /*Numerical coefficients from PALLISER*/ | ||
479 | static const Scalar f[] = { 2.63500E-1, 7.48368E-6, 1.44611E-6, -3.80860E-10 }; | ||
480 | |||
481 | /*Numerical coefficients from MICHAELIDES for the enthalpy of brine*/ | ||
482 | static const Scalar a[4][3] = { { +9633.6, -4080.0, +286.49 }, | ||
483 | { +166.58, +68.577, -4.6856 }, | ||
484 | { -0.90963, -0.36524, +0.249667E-1 }, | ||
485 | { +0.17965E-2, +0.71924E-3, -0.4900E-4 } }; | ||
486 | |||
487 | 85989 | const Scalar theta = T - 273.15; | |
488 | 85989 | const Scalar salSat = f[0] + f[1]*theta + f[2]*theta*theta + f[3]*theta*theta*theta; | |
489 | |||
490 | /*Regularization*/ | ||
491 | using std::min; | ||
492 | using std::max; | ||
493 |
3/4✗ Branch 0 not taken.
✓ Branch 1 taken 85968 times.
✓ Branch 2 taken 1360 times.
✓ Branch 3 taken 84608 times.
|
85989 | const Scalar salinity = min(max(xNaCl,0.0), salSat); |
494 | |||
495 | 85989 | const Scalar hw = H2O::liquidEnthalpy(T, p)/1E3; /* kJ/kg */ | |
496 | |||
497 | /*component enthalpy of soluted NaCl after DAUBERT and DANNER*/ | ||
498 | 85989 | const Scalar h_NaCl = (3.6710E4*T + 0.5*(6.2770E1)*T*T - ((6.6670E-2)/3)*T*T*T | |
499 | 85989 | + ((2.8000E-5)/4)*(T*T*T*T))/(58.44E3)- 2.045698e+02; /* U [kJ/kg] */ | |
500 | |||
501 | 85989 | const Scalar m = (1E3/58.44)*(salinity/(1-salinity)); | |
502 | |||
503 | using Dune::power; | ||
504 | 85989 | Scalar d_h = 0; | |
505 |
2/2✓ Branch 0 taken 343872 times.
✓ Branch 1 taken 85968 times.
|
429945 | for (int i = 0; i<=3; i++) |
506 |
2/2✓ Branch 0 taken 1031616 times.
✓ Branch 1 taken 343872 times.
|
1375824 | for (int j=0; j<=2; j++) |
507 | 3095604 | d_h = d_h + a[i][j] * power(theta, i) * power(m, j); | |
508 | |||
509 | /* heat of dissolution for halite according to Michaelides 1981 */ | ||
510 | 85989 | const Scalar delta_h = (4.184/(1E3 + (58.44 * m)))*d_h; | |
511 | |||
512 | /* Enthalpy of brine without any dissolved gas */ | ||
513 | 85989 | const Scalar h_ls1 =(1-salinity)*hw + salinity*h_NaCl + salinity*delta_h; /* kJ/kg */ | |
514 | 85989 | return h_ls1*1E3; /*J/kg*/ | |
515 | } | ||
516 | }; | ||
517 | |||
518 | } // end namespace Dumux::FluidSystems | ||
519 | |||
520 | #endif | ||
521 |