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 | * \copybrief Dumux::FluidSystems::H2ON2 | ||
11 | */ | ||
12 | #ifndef DUMUX_H2O_N2_FLUID_SYSTEM_HH | ||
13 | #define DUMUX_H2O_N2_FLUID_SYSTEM_HH | ||
14 | |||
15 | #include <cassert> | ||
16 | #include <iomanip> | ||
17 | |||
18 | #include <dumux/common/exceptions.hh> | ||
19 | |||
20 | #include <dumux/material/idealgas.hh> | ||
21 | |||
22 | #include <dumux/material/components/n2.hh> | ||
23 | #include <dumux/material/components/h2o.hh> | ||
24 | #include <dumux/material/components/tabulatedcomponent.hh> | ||
25 | #include <dumux/material/binarycoefficients/h2o_n2.hh> | ||
26 | |||
27 | #include <dumux/io/name.hh> | ||
28 | |||
29 | #include "base.hh" | ||
30 | |||
31 | namespace Dumux::FluidSystems { | ||
32 | |||
33 | /*! | ||
34 | * \ingroup FluidSystems | ||
35 | * \brief Policy for the H2O-N2 fluid system | ||
36 | */ | ||
37 | template<bool fastButSimplifiedRelations = false> | ||
38 | struct H2ON2DefaultPolicy | ||
39 | { | ||
40 | static constexpr bool useH2ODensityAsLiquidMixtureDensity() { return fastButSimplifiedRelations; } | ||
41 | static constexpr bool useIdealGasDensity() { return fastButSimplifiedRelations; } | ||
42 | static constexpr bool useN2ViscosityAsGasMixtureViscosity() { return fastButSimplifiedRelations; } | ||
43 | static constexpr bool useN2HeatConductivityAsGasMixtureHeatConductivity() { return fastButSimplifiedRelations; } | ||
44 | static constexpr bool useIdealGasHeatCapacities() { return fastButSimplifiedRelations; } | ||
45 | }; | ||
46 | |||
47 | /*! | ||
48 | * \ingroup FluidSystems | ||
49 | * | ||
50 | * \brief A two-phase fluid system with two components water \f$(\mathrm{H_2O})\f$ | ||
51 | * Nitrogen \f$(\mathrm{N_2})\f$ for non-equilibrium models. | ||
52 | */ | ||
53 | template <class Scalar, class Policy = H2ON2DefaultPolicy<>> | ||
54 | class H2ON2 | ||
55 | : public Base<Scalar, H2ON2<Scalar, Policy> > | ||
56 | { | ||
57 | using ThisType = H2ON2<Scalar, Policy>; | ||
58 | |||
59 | // convenience aliases using declarations | ||
60 | using IdealGas = Dumux::IdealGas<Scalar>; | ||
61 | using TabulatedH2O = Components::TabulatedComponent<Dumux::Components::H2O<Scalar> >; | ||
62 | using SimpleN2 = Dumux::Components::N2<Scalar>; | ||
63 | |||
64 | public: | ||
65 | using H2O = TabulatedH2O; //!< The component for pure water | ||
66 | using N2 = SimpleN2; //!< The component for pure nitrogen | ||
67 | |||
68 | static constexpr int numPhases = 2; //!< Number of phases in the fluid system | ||
69 | static constexpr int numComponents = 2; //!< Number of components in the fluid system | ||
70 | |||
71 | static constexpr int liquidPhaseIdx = 0; //!< index of the liquid phase | ||
72 | static constexpr int gasPhaseIdx = 1; //!< index of the gas phase | ||
73 | static constexpr int phase0Idx = liquidPhaseIdx; //!< index of the first phase | ||
74 | static constexpr int phase1Idx = gasPhaseIdx; //!< index of the second phase | ||
75 | |||
76 | static constexpr int H2OIdx = 0; | ||
77 | static constexpr int N2Idx = 1; | ||
78 | static constexpr int comp0Idx = H2OIdx; //!< index of the first component | ||
79 | static constexpr int comp1Idx = N2Idx; //!< index of the second component | ||
80 | static constexpr int liquidCompIdx = H2OIdx; //!< index of the liquid component | ||
81 | static constexpr int gasCompIdx = N2Idx; //!< index of the gas component | ||
82 | |||
83 | /**************************************** | ||
84 | * Fluid phase related static parameters | ||
85 | ****************************************/ | ||
86 | /*! | ||
87 | * \brief Return the human readable name of a fluid phase | ||
88 | * | ||
89 | * \param phaseIdx The index of the fluid phase to consider | ||
90 | */ | ||
91 | 936 | static std::string phaseName(int phaseIdx) | |
92 | { | ||
93 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 482 times.
|
490 | assert(0 <= phaseIdx && phaseIdx < numPhases); |
94 |
2/2✓ Branch 0 taken 274 times.
✓ Branch 1 taken 208 times.
|
490 | switch (phaseIdx) |
95 | { | ||
96 | 278 | case liquidPhaseIdx: return IOName::liquidPhase(); | |
97 | 212 | case gasPhaseIdx: return IOName::gaseousPhase(); | |
98 | } | ||
99 | DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx); | ||
100 | } | ||
101 | |||
102 | /*! | ||
103 | * \brief Returns whether the fluids are miscible | ||
104 | */ | ||
105 | static constexpr bool isMiscible() | ||
106 | { return true; } | ||
107 | |||
108 | /*! | ||
109 | * \brief Return whether a phase is gaseous | ||
110 | * | ||
111 | * \param phaseIdx The index of the fluid phase to consider | ||
112 | */ | ||
113 | 99336 | static constexpr bool isGas(int phaseIdx) | |
114 | { | ||
115 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 99332 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 4 times.
|
99336 | assert(0 <= phaseIdx && phaseIdx < numPhases); |
116 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 99332 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 4 times.
|
99336 | return phaseIdx == gasPhaseIdx; |
117 | } | ||
118 | |||
119 | /*! | ||
120 | * \brief Returns true if and only if a fluid phase is assumed to | ||
121 | * be an ideal mixture. | ||
122 | * | ||
123 | * We define an ideal mixture as a fluid phase where the fugacity | ||
124 | * coefficients of all components times the pressure of the phase | ||
125 | * are independent on the fluid composition. This assumption is true | ||
126 | * if Henry's law and Raoult's law apply. If you are unsure what | ||
127 | * this function should return, it is safe to return false. The | ||
128 | * only damage done will be (slightly) increased computation times | ||
129 | * in some cases. | ||
130 | * | ||
131 | * \param phaseIdx The index of the fluid phase to consider | ||
132 | */ | ||
133 | 15071263 | static bool isIdealMixture(int phaseIdx) | |
134 | { | ||
135 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 15070937 times.
|
15070937 | assert(0 <= phaseIdx && phaseIdx < numPhases); |
136 | // we assume Henry's and Raoult's laws for the water phase and | ||
137 | // and no interaction between gas molecules of different | ||
138 | // components, so all phases are ideal mixtures! | ||
139 | return true; | ||
140 | } | ||
141 | |||
142 | /*! | ||
143 | * \brief Returns true if and only if a fluid phase is assumed to | ||
144 | * be compressible. | ||
145 | * | ||
146 | * Compressible means that the partial derivative of the density | ||
147 | * to the fluid pressure is always larger than zero. | ||
148 | * | ||
149 | * \param phaseIdx The index of the fluid phase to consider | ||
150 | */ | ||
151 | static constexpr bool isCompressible(int phaseIdx) | ||
152 | { | ||
153 | assert(0 <= phaseIdx && phaseIdx < numPhases); | ||
154 | // gases are always compressible | ||
155 | if (phaseIdx == gasPhaseIdx) | ||
156 | return true; | ||
157 | // the water component decides for the liquid phase... | ||
158 | return H2O::liquidIsCompressible(); | ||
159 | } | ||
160 | |||
161 | /*! | ||
162 | * \brief Returns true if and only if a fluid phase is assumed to | ||
163 | * be an ideal gas. | ||
164 | * | ||
165 | * \param phaseIdx The index of the fluid phase to consider | ||
166 | */ | ||
167 | static bool isIdealGas(int phaseIdx) | ||
168 | { | ||
169 | assert(0 <= phaseIdx && phaseIdx < numPhases); | ||
170 | |||
171 | if (phaseIdx == gasPhaseIdx) | ||
172 | // let the components decide | ||
173 | return H2O::gasIsIdeal() && N2::gasIsIdeal(); | ||
174 | return false; // not a gas | ||
175 | } | ||
176 | |||
177 | /**************************************** | ||
178 | * Component related static parameters | ||
179 | ****************************************/ | ||
180 | /*! | ||
181 | * \brief Return the human readable name of a component | ||
182 | * | ||
183 | * \param compIdx The index of the component to consider | ||
184 | */ | ||
185 | 380 | static std::string componentName(int compIdx) | |
186 | { | ||
187 |
2/3✓ Branch 0 taken 187 times.
✓ Branch 1 taken 185 times.
✗ Branch 2 not taken.
|
380 | switch (compIdx) |
188 | { | ||
189 | 191 | case H2OIdx: return H2O::name(); | |
190 | 189 | case N2Idx: return N2::name(); | |
191 | } | ||
192 | |||
193 | ✗ | DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << compIdx); | |
194 | } | ||
195 | |||
196 | /*! | ||
197 | * \brief Return the molar mass of a component in \f$\mathrm{[kg/mol]}\f$. | ||
198 | * | ||
199 | * \param compIdx The index of the component to consider | ||
200 | */ | ||
201 | 1118345125 | static Scalar molarMass(int compIdx) | |
202 | { | ||
203 | static const Scalar M[] = { | ||
204 | H2O::molarMass(), | ||
205 | N2::molarMass(), | ||
206 | }; | ||
207 | |||
208 |
5/10✗ Branch 0 not taken.
✓ Branch 1 taken 531617108 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 28982222 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 4 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 12 not taken.
✓ Branch 13 taken 7 times.
|
560599343 | assert(0 <= compIdx && compIdx < numComponents); |
209 |
7/11✓ Branch 0 taken 21451240 times.
✓ Branch 1 taken 2 times.
✓ Branch 3 taken 2 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✓ Branch 2 taken 41348792 times.
|
585064540 | return M[compIdx]; |
210 | } | ||
211 | |||
212 | /*! | ||
213 | * \brief Critical temperature of a component \f$\mathrm{[K]}\f$. | ||
214 | * | ||
215 | * \param compIdx The index of the component to consider | ||
216 | */ | ||
217 | static Scalar criticalTemperature(int compIdx) | ||
218 | { | ||
219 | static const Scalar Tcrit[] = { | ||
220 | H2O::criticalTemperature(), | ||
221 | N2::criticalTemperature() | ||
222 | }; | ||
223 | |||
224 | assert(0 <= compIdx && compIdx < numComponents); | ||
225 | return Tcrit[compIdx]; | ||
226 | } | ||
227 | |||
228 | /*! | ||
229 | * \brief Critical pressure of a component \f$\mathrm{[Pa]}\f$. | ||
230 | * | ||
231 | * \param compIdx The index of the component to consider | ||
232 | */ | ||
233 | static Scalar criticalPressure(int compIdx) | ||
234 | { | ||
235 | static const Scalar pcrit[] = { | ||
236 | H2O::criticalPressure(), | ||
237 | N2::criticalPressure() | ||
238 | }; | ||
239 | |||
240 | assert(0 <= compIdx && compIdx < numComponents); | ||
241 | return pcrit[compIdx]; | ||
242 | } | ||
243 | |||
244 | /*! | ||
245 | * \brief Vapor pressure including the Kelvin equation in \f$\mathrm{[Pa]}\f$ | ||
246 | * | ||
247 | * Calculate the decreased vapor pressure due to capillarity | ||
248 | * | ||
249 | * \param fluidState An arbitrary fluid state | ||
250 | * \param phaseIdx The index of the fluid phase to consider | ||
251 | * \param compIdx The index of the component to consider | ||
252 | */ | ||
253 | template <class FluidState> | ||
254 | static Scalar kelvinVaporPressure(const FluidState &fluidState, | ||
255 | const int phaseIdx, | ||
256 | const int compIdx) | ||
257 | { | ||
258 | assert(compIdx == H2OIdx && phaseIdx == liquidPhaseIdx); | ||
259 | |||
260 | using std::exp; | ||
261 | return fugacityCoefficient(fluidState, phaseIdx, compIdx) | ||
262 | * fluidState.pressure(phaseIdx) | ||
263 | * exp(-(fluidState.pressure(gasPhaseIdx)-fluidState.pressure(liquidPhaseIdx)) | ||
264 | / density(fluidState, phaseIdx) | ||
265 | / (Dumux::Constants<Scalar>::R / molarMass(compIdx)) | ||
266 | / fluidState.temperature()); | ||
267 | } | ||
268 | |||
269 | /*! | ||
270 | * \brief Molar volume of a component at the critical point \f$\mathrm{[m^3/mol]}\f$. | ||
271 | * | ||
272 | * \param compIdx The index of the component to consider | ||
273 | */ | ||
274 | static Scalar criticalMolarVolume(int compIdx) | ||
275 | { | ||
276 | DUNE_THROW(Dune::NotImplemented, | ||
277 | "H2ON2FluidSystem::criticalMolarVolume()"); | ||
278 | } | ||
279 | |||
280 | /*! | ||
281 | * \brief The acentric factor of a component \f$\mathrm{[-]}\f$. | ||
282 | * | ||
283 | * \param compIdx The index of the component to consider | ||
284 | */ | ||
285 | static Scalar acentricFactor(int compIdx) | ||
286 | { | ||
287 | static const Scalar accFac[] = { | ||
288 | H2O::acentricFactor(), | ||
289 | N2::acentricFactor() | ||
290 | }; | ||
291 | |||
292 | assert(0 <= compIdx && compIdx < numComponents); | ||
293 | return accFac[compIdx]; | ||
294 | } | ||
295 | |||
296 | /**************************************** | ||
297 | * thermodynamic relations | ||
298 | ****************************************/ | ||
299 | |||
300 | /*! | ||
301 | * \brief Initialize the fluid system's static parameters generically | ||
302 | * | ||
303 | * If a tabulated H2O component is used, we do our best to create | ||
304 | * tables that always work. | ||
305 | */ | ||
306 | 55 | static void init() | |
307 | { | ||
308 |
1/2✓ Branch 1 taken 51 times.
✗ Branch 2 not taken.
|
55 | init(/*tempMin=*/273.15, |
309 | /*tempMax=*/623.15, | ||
310 | /*numTemp=*/100, | ||
311 | /*pMin=*/0.0, | ||
312 | /*pMax=*/20e6, | ||
313 | /*numP=*/200); | ||
314 | 11 | } | |
315 | |||
316 | /*! | ||
317 | * \brief Initialize the fluid system's static parameters using | ||
318 | * problem specific temperature and pressure ranges | ||
319 | * | ||
320 | * \param tempMin The minimum temperature used for tabulation of water \f$\mathrm{[K]}\f$ | ||
321 | * \param tempMax The maximum temperature used for tabulation of water \f$\mathrm{[K]}\f$ | ||
322 | * \param nTemp The number of ticks on the temperature axis of the table of water | ||
323 | * \param pressMin The minimum pressure used for tabulation of water \f$\mathrm{[Pa]}\f$ | ||
324 | * \param pressMax The maximum pressure used for tabulation of water \f$\mathrm{[Pa]}\f$ | ||
325 | * \param nPress The number of ticks on the pressure axis of the table of water | ||
326 | */ | ||
327 | 81 | static void init(Scalar tempMin, Scalar tempMax, unsigned nTemp, | |
328 | Scalar pressMin, Scalar pressMax, unsigned nPress) | ||
329 | { | ||
330 | 81 | std::cout << "The H2O-N2 fluid system was configured with the following policy:\n"; | |
331 | 81 | std::cout << " - use H2O density as liquid mixture density: " << std::boolalpha << Policy::useH2ODensityAsLiquidMixtureDensity() << "\n"; | |
332 | 81 | std::cout << " - use ideal gas density: " << std::boolalpha << Policy::useIdealGasDensity() << "\n"; | |
333 | 81 | std::cout << " - use N2 viscosity as gas mixture viscosity: " << std::boolalpha << Policy::useN2ViscosityAsGasMixtureViscosity() << "\n"; | |
334 | 81 | std::cout << " - use N2 heat conductivity as gas mixture heat conductivity: " << std::boolalpha << Policy::useN2HeatConductivityAsGasMixtureHeatConductivity() << "\n"; | |
335 | 81 | std::cout << " - use ideal gas heat capacities: " << std::boolalpha << Policy::useIdealGasHeatCapacities() << std::endl; | |
336 | |||
337 | if constexpr (H2O::isTabulated) | ||
338 | 81 | H2O::init(tempMin, tempMax, nTemp, pressMin, pressMax, nPress); | |
339 | 81 | } | |
340 | |||
341 | using Base<Scalar, ThisType>::density; | ||
342 | /*! | ||
343 | * \brief Given a phase's composition, temperature, pressure, and | ||
344 | * the partial pressures of all components, return its | ||
345 | * density \f$\mathrm{[kg/m^3]}\f$. | ||
346 | * | ||
347 | * If Policy::useH2ODensityAsLiquidMixtureDensity() == false, we apply Eq. (7) | ||
348 | * in Class et al. (2002a) \cite A3:class:2002b <BR> | ||
349 | * for the liquid density. | ||
350 | * | ||
351 | * \param fluidState An arbitrary fluid state | ||
352 | * \param phaseIdx The index of the fluid phase to consider | ||
353 | */ | ||
354 | template <class FluidState> | ||
355 | 55007568 | static Scalar density(const FluidState &fluidState, | |
356 | int phaseIdx) | ||
357 | { | ||
358 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 37542842 times.
|
37542860 | assert(0 <= phaseIdx && phaseIdx < numPhases); |
359 | |||
360 |
2/2✓ Branch 0 taken 18656363 times.
✓ Branch 1 taken 18886471 times.
|
55007568 | Scalar T = fluidState.temperature(phaseIdx); |
361 | 55007568 | Scalar p = fluidState.pressure(phaseIdx); | |
362 | |||
363 | // liquid phase | ||
364 |
2/2✓ Branch 0 taken 18656367 times.
✓ Branch 1 taken 18886475 times.
|
37542860 | if (phaseIdx == liquidPhaseIdx) { |
365 | if (Policy::useH2ODensityAsLiquidMixtureDensity()) | ||
366 | // assume pure water | ||
367 | 24508456 | return H2O::liquidDensity(T, p); | |
368 | else | ||
369 | { | ||
370 | // See: Eq. (7) in Class et al. (2002a) | ||
371 | // This assumes each gas molecule displaces exactly one | ||
372 | // molecule in the liquid. | ||
373 | 11612628 | return H2O::liquidMolarDensity(T, p) | |
374 | 11612494 | * (H2O::molarMass()*fluidState.moleFraction(liquidPhaseIdx, H2OIdx) | |
375 | 11612628 | + N2::molarMass()*fluidState.moleFraction(liquidPhaseIdx, N2Idx)); | |
376 | } | ||
377 | } | ||
378 | |||
379 | // gas phase | ||
380 | using std::max; | ||
381 | if (Policy::useIdealGasDensity()) | ||
382 | // for the gas phase assume an ideal gas | ||
383 | { | ||
384 | 7249032 | const Scalar averageMolarMass = fluidState.averageMolarMass(gasPhaseIdx); | |
385 | 7249032 | return IdealGas::density(averageMolarMass, T, p); | |
386 | } | ||
387 | |||
388 | // assume ideal mixture: steam and nitrogen don't "see" each other | ||
389 | 11637460 | Scalar rho_gH2O = H2O::gasDensity(T, fluidState.partialPressure(gasPhaseIdx, H2OIdx)); | |
390 | 11637460 | Scalar rho_gN2 = N2::gasDensity(T, fluidState.partialPressure(gasPhaseIdx, N2Idx)); | |
391 | 11637460 | return (rho_gH2O + rho_gN2); | |
392 | } | ||
393 | |||
394 | using Base<Scalar, ThisType>::molarDensity; | ||
395 | //! \copydoc Dumux::FluidSystems::Base<Scalar,ThisType>::molarDensity(const FluidState&,int) | ||
396 | template <class FluidState> | ||
397 | 50390042 | static Scalar molarDensity(const FluidState &fluidState, int phaseIdx) | |
398 | { | ||
399 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 32925316 times.
|
32925334 | assert(0 <= phaseIdx && phaseIdx < numPhases); |
400 | |||
401 |
2/2✓ Branch 0 taken 16347600 times.
✓ Branch 1 taken 16577708 times.
|
50390042 | Scalar T = fluidState.temperature(phaseIdx); |
402 | 50390042 | Scalar p = fluidState.pressure(phaseIdx); | |
403 | |||
404 | // liquid phase | ||
405 |
2/2✓ Branch 0 taken 16347604 times.
✓ Branch 1 taken 16577712 times.
|
32925334 | if (phaseIdx == liquidPhaseIdx) |
406 | { | ||
407 | // assume pure water or that each gas molecule displaces exactly one | ||
408 | // molecule in the liquid. | ||
409 | 33812321 | return H2O::liquidMolarDensity(T, p); | |
410 | } | ||
411 | |||
412 | // gas phase | ||
413 | using std::max; | ||
414 | if (Policy::useIdealGasDensity()) | ||
415 | // for the gas phase assume an ideal gas | ||
416 | { | ||
417 | 4940269 | return IdealGas::molarDensity(T, p); | |
418 | } | ||
419 | |||
420 | // assume ideal mixture: steam and nitrogen don't "see" each other | ||
421 | 11637460 | Scalar rho_gH2O = H2O::gasMolarDensity(T, fluidState.partialPressure(gasPhaseIdx, H2OIdx)); | |
422 | 11637460 | Scalar rho_gN2 = N2::gasMolarDensity(T, fluidState.partialPressure(gasPhaseIdx, N2Idx)); | |
423 | 11637460 | return rho_gH2O + rho_gN2; | |
424 | } | ||
425 | |||
426 | using Base<Scalar, ThisType>::viscosity; | ||
427 | /*! | ||
428 | * \brief Calculate the dynamic viscosity of a fluid phase \f$\mathrm{[Pa*s]}\f$ | ||
429 | * | ||
430 | * Compositional effects in the gas phase are accounted by the Wilke method. | ||
431 | * See Reid et al. (1987) \cite reid1987 <BR> | ||
432 | * 4th edition, McGraw-Hill, 1987, 407-410 | ||
433 | * 5th edition, McGraw-Hill, 20001, p. 9.21/22 | ||
434 | * | ||
435 | * \param fluidState An arbitrary fluid state | ||
436 | * \param phaseIdx The index of the fluid phase to consider | ||
437 | * \note Compositional effects for a liquid mixture have to be implemented. | ||
438 | */ | ||
439 | template <class FluidState> | ||
440 | 51686242 | static Scalar viscosity(const FluidState &fluidState, | |
441 | int phaseIdx) | ||
442 | { | ||
443 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 34266073 times.
|
34266081 | assert(0 <= phaseIdx && phaseIdx < numPhases); |
444 | |||
445 |
2/2✓ Branch 0 taken 17659382 times.
✓ Branch 1 taken 16606683 times.
|
51686244 | Scalar T = fluidState.temperature(phaseIdx); |
446 | 51686252 | Scalar p = fluidState.pressure(phaseIdx); | |
447 | |||
448 | // liquid phase | ||
449 |
2/2✓ Branch 0 taken 17659386 times.
✓ Branch 1 taken 16606687 times.
|
34266081 | if (phaseIdx == liquidPhaseIdx) { |
450 | // assume pure water for the liquid phase | ||
451 | 35079553 | return H2O::liquidViscosity(T, p); | |
452 | } | ||
453 | |||
454 | // gas phase | ||
455 | if (Policy::useN2ViscosityAsGasMixtureViscosity()) | ||
456 | { | ||
457 | // assume pure nitrogen for the gas phase | ||
458 | 4969373 | return N2::gasViscosity(T, p); | |
459 | } | ||
460 | else | ||
461 | { | ||
462 | // Wilke method (Reid et al.): | ||
463 | 11637326 | Scalar muResult = 0; | |
464 | 23274652 | const Scalar mu[numComponents] = { | |
465 | 11637326 | h2oGasViscosityInMixture(T, p), | |
466 | 11637326 | N2::gasViscosity(T, p) | |
467 | }; | ||
468 | |||
469 | 11637326 | Scalar sumx = 0.0; | |
470 | using std::max; | ||
471 |
2/2✓ Branch 0 taken 23274648 times.
✓ Branch 1 taken 11637324 times.
|
34911978 | for (int compIdx = 0; compIdx < numComponents; ++compIdx) |
472 | 23274652 | sumx += fluidState.moleFraction(phaseIdx, compIdx); | |
473 |
1/2✓ Branch 0 taken 11637324 times.
✗ Branch 1 not taken.
|
11637326 | sumx = max(1e-10, sumx); |
474 | |||
475 |
2/2✓ Branch 0 taken 23274648 times.
✓ Branch 1 taken 11637324 times.
|
34911978 | for (int i = 0; i < numComponents; ++i) { |
476 | Scalar divisor = 0; | ||
477 | // using std::sqrt; | ||
478 | // using std::pow; | ||
479 |
2/2✓ Branch 0 taken 46549296 times.
✓ Branch 1 taken 23274648 times.
|
69823956 | for (int j = 0; j < numComponents; ++j) { |
480 | 46549304 | Scalar phiIJ = 1 + sqrt(mu[i]/mu[j]) * pow(molarMass(j)/molarMass(i), 1/4.0); | |
481 | 46549304 | phiIJ *= phiIJ; | |
482 | 46549304 | phiIJ /= sqrt(8*(1 + molarMass(i)/molarMass(j))); | |
483 | 46549304 | divisor += fluidState.moleFraction(phaseIdx, j)/sumx * phiIJ; | |
484 | } | ||
485 | 23274652 | muResult += fluidState.moleFraction(phaseIdx, i)/sumx * mu[i] / divisor; | |
486 | } | ||
487 | |||
488 | return muResult; | ||
489 | } | ||
490 | } | ||
491 | |||
492 | using Base<Scalar, ThisType>::fugacityCoefficient; | ||
493 | //! \copydoc Dumux::FluidSystems::Base<Scalar,ThisType>::fugacityCoefficient(const FluidState&,int,int) | ||
494 | template <class FluidState> | ||
495 | 70164444 | static Scalar fugacityCoefficient(const FluidState &fluidState, | |
496 | int phaseIdx, | ||
497 | int compIdx) | ||
498 | { | ||
499 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 70164444 times.
|
70164460 | assert(0 <= phaseIdx && phaseIdx < numPhases); |
500 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 70164444 times.
|
70164460 | assert(0 <= compIdx && compIdx < numComponents); |
501 | |||
502 |
2/2✓ Branch 0 taken 35082218 times.
✓ Branch 1 taken 35082210 times.
|
70164468 | Scalar T = fluidState.temperature(phaseIdx); |
503 | 70164444 | Scalar p = fluidState.pressure(phaseIdx); | |
504 | |||
505 | // liquid phase | ||
506 |
2/2✓ Branch 0 taken 35082226 times.
✓ Branch 1 taken 35082218 times.
|
70164460 | if (phaseIdx == liquidPhaseIdx) { |
507 |
2/2✓ Branch 0 taken 17541109 times.
✓ Branch 1 taken 17541117 times.
|
35082234 | if (compIdx == H2OIdx) |
508 | 17541121 | return H2O::vaporPressure(T)/p; | |
509 | 17541121 | return BinaryCoeff::H2O_N2::henry(T)/p; | |
510 | } | ||
511 | |||
512 | // for the gas phase, assume an ideal gas when it comes to | ||
513 | // fugacity (-> fugacity == partial pressure) | ||
514 | return 1.0; | ||
515 | } | ||
516 | |||
517 | using Base<Scalar, ThisType>::diffusionCoefficient; | ||
518 | //! \copydoc Dumux::FluidSystems::Base<Scalar,ThisType>::diffusionCoefficient(const FluidState&,int,int) | ||
519 | template <class FluidState> | ||
520 | 32 | static Scalar diffusionCoefficient(const FluidState &fluidState, | |
521 | int phaseIdx, | ||
522 | int compIdx) | ||
523 | { | ||
524 |
10/20✓ Branch 1 taken 16 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 16 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 16 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 16 times.
✗ Branch 11 not taken.
✓ Branch 13 taken 16 times.
✗ Branch 14 not taken.
✓ Branch 16 taken 16 times.
✗ Branch 17 not taken.
✓ Branch 19 taken 16 times.
✗ Branch 20 not taken.
✓ Branch 22 taken 16 times.
✗ Branch 23 not taken.
✓ Branch 25 taken 16 times.
✗ Branch 26 not taken.
✓ Branch 28 taken 16 times.
✗ Branch 29 not taken.
|
128 | DUNE_THROW(Dune::NotImplemented, "Diffusion coefficients"); |
525 | } | ||
526 | |||
527 | using Base<Scalar, ThisType>::binaryDiffusionCoefficient; | ||
528 | //! \copydoc Dumux::FluidSystems::Base<Scalar,ThisType>::binaryDiffusionCoefficient(const FluidState&,int,int,int) | ||
529 | template <class FluidState> | ||
530 | 55561647 | static Scalar binaryDiffusionCoefficient(const FluidState &fluidState, | |
531 | int phaseIdx, | ||
532 | int compIIdx, | ||
533 | int compJIdx) | ||
534 | |||
535 | { | ||
536 |
2/2✓ Branch 0 taken 18076700 times.
✓ Branch 1 taken 37484915 times.
|
55561647 | if (compIIdx > compJIdx) |
537 | { | ||
538 | using std::swap; | ||
539 | 18076708 | swap(compIIdx, compJIdx); | |
540 | } | ||
541 | |||
542 |
2/2✓ Branch 0 taken 36499151 times.
✓ Branch 1 taken 19062432 times.
|
55561647 | const Scalar T = fluidState.temperature(phaseIdx); |
543 | 55561647 | const Scalar p = fluidState.pressure(phaseIdx); | |
544 | |||
545 |
6/6✓ Branch 0 taken 36499167 times.
✓ Branch 1 taken 19062448 times.
✓ Branch 2 taken 36499163 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 36499159 times.
✓ Branch 5 taken 4 times.
|
55561647 | if (phaseIdx == liquidPhaseIdx && compIIdx == H2OIdx && compJIdx == N2Idx) |
546 | 36499167 | return BinaryCoeff::H2O_N2::liquidDiffCoeff(T, p); | |
547 | |||
548 |
6/6✓ Branch 0 taken 19062448 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 19062444 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 19062440 times.
✓ Branch 5 taken 4 times.
|
19062480 | else if (phaseIdx == gasPhaseIdx && compIIdx == H2OIdx && compJIdx == N2Idx) |
549 | 19062448 | return BinaryCoeff::H2O_N2::gasDiffCoeff(T, p); | |
550 | |||
551 | else | ||
552 |
16/32✓ Branch 1 taken 16 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 16 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 16 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 16 times.
✗ Branch 11 not taken.
✓ Branch 13 taken 16 times.
✗ Branch 14 not taken.
✓ Branch 16 taken 16 times.
✗ Branch 17 not taken.
✓ Branch 19 taken 16 times.
✗ Branch 20 not taken.
✓ Branch 22 taken 16 times.
✗ Branch 23 not taken.
✓ Branch 25 taken 16 times.
✗ Branch 26 not taken.
✓ Branch 28 taken 16 times.
✗ Branch 29 not taken.
✓ Branch 31 taken 16 times.
✗ Branch 32 not taken.
✓ Branch 34 taken 16 times.
✗ Branch 35 not taken.
✓ Branch 37 taken 16 times.
✗ Branch 38 not taken.
✓ Branch 40 taken 16 times.
✗ Branch 41 not taken.
✓ Branch 43 taken 16 times.
✗ Branch 44 not taken.
✓ Branch 46 taken 16 times.
✗ Branch 47 not taken.
|
128 | DUNE_THROW(Dune::InvalidStateException, |
553 | "Binary diffusion coefficient of components " | ||
554 | << compIIdx << " and " << compJIdx | ||
555 | << " in phase " << phaseIdx << " is unavailable!\n"); | ||
556 | } | ||
557 | |||
558 | using Base<Scalar, ThisType>::enthalpy; | ||
559 | /*! | ||
560 | * \brief Given a phase's composition, temperature, pressure and | ||
561 | * density, calculate its specific enthalpy \f$\mathrm{[J/kg]}\f$. | ||
562 | * | ||
563 | * \note This fluid system neglects the contribution of | ||
564 | * gas-molecules in the liquid phase. This contribution is | ||
565 | * probably not big. Somebody would have to find out the | ||
566 | * enthalpy of solution for this system. ... | ||
567 | * | ||
568 | * \param fluidState An arbitrary fluid state | ||
569 | * \param phaseIdx The index of the fluid phase to consider | ||
570 | */ | ||
571 | template <class FluidState> | ||
572 |
2/2✓ Branch 0 taken 27330178 times.
✓ Branch 1 taken 16147205 times.
|
43477399 | static Scalar enthalpy(const FluidState &fluidState, |
573 | int phaseIdx) | ||
574 | { | ||
575 |
2/2✓ Branch 0 taken 27330178 times.
✓ Branch 1 taken 16147205 times.
|
43477399 | const Scalar T = fluidState.temperature(phaseIdx); |
576 | 43477399 | const Scalar p = fluidState.pressure(phaseIdx); | |
577 | |||
578 | // liquid phase | ||
579 |
2/2✓ Branch 0 taken 27330182 times.
✓ Branch 1 taken 16147209 times.
|
43477399 | if (phaseIdx == liquidPhaseIdx) { |
580 | 27330186 | return H2O::liquidEnthalpy(T, p); | |
581 | } | ||
582 | // gas phase | ||
583 | else { | ||
584 | // assume ideal mixture: which means | ||
585 | // that the total specific enthalpy is the sum of the | ||
586 | // "partial specific enthalpies" of the components. | ||
587 | 16147213 | Scalar hH2O = | |
588 | 16147213 | fluidState.massFraction(gasPhaseIdx, H2OIdx) | |
589 | 16147213 | * H2O::gasEnthalpy(T, p); | |
590 | 16147213 | Scalar hN2 = | |
591 | 16147213 | fluidState.massFraction(gasPhaseIdx, N2Idx) | |
592 | 16147213 | * N2::gasEnthalpy(T, p); | |
593 | 16147213 | return hH2O + hN2; | |
594 | } | ||
595 | } | ||
596 | |||
597 | /*! | ||
598 | * \brief Returns the specific enthalpy \f$\mathrm{[J/kg]}\f$ of a component in the specified phase | ||
599 | * \param fluidState The fluid state | ||
600 | * \param phaseIdx The index of the phase | ||
601 | * \param componentIdx The index of the component | ||
602 | */ | ||
603 | template <class FluidState> | ||
604 | static Scalar componentEnthalpy(const FluidState &fluidState, | ||
605 | int phaseIdx, | ||
606 | int componentIdx) | ||
607 | { | ||
608 | const Scalar T = fluidState.temperature(phaseIdx); | ||
609 | const Scalar p = fluidState.pressure(phaseIdx); | ||
610 | |||
611 | if (phaseIdx == liquidPhaseIdx) | ||
612 | { | ||
613 | if (componentIdx == H2OIdx) | ||
614 | return H2O::liquidEnthalpy(T, p); | ||
615 | else if (componentIdx == N2Idx) | ||
616 | DUNE_THROW(Dune::NotImplemented, "Component enthalpy of nitrogen in liquid phase"); | ||
617 | else | ||
618 | DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << componentIdx); | ||
619 | } | ||
620 | else if (phaseIdx == gasPhaseIdx) | ||
621 | { | ||
622 | if (componentIdx == H2OIdx) | ||
623 | return H2O::gasEnthalpy(T, p); | ||
624 | else if (componentIdx == N2Idx) | ||
625 | return N2::gasEnthalpy(T, p); | ||
626 | else | ||
627 | DUNE_THROW(Dune::InvalidStateException, "Invalid component index " << componentIdx); | ||
628 | } | ||
629 | else | ||
630 | DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx); | ||
631 | } | ||
632 | |||
633 | using Base<Scalar, ThisType>::thermalConductivity; | ||
634 | /*! | ||
635 | * \brief Thermal conductivity of a fluid phase \f$\mathrm{[W/(m K)]}\f$. | ||
636 | * | ||
637 | * Use the conductivity of air and water as a first approximation. | ||
638 | * | ||
639 | * \param fluidState An arbitrary fluid state | ||
640 | * \param phaseIdx The index of the fluid phase to consider | ||
641 | */ | ||
642 | template <class FluidState> | ||
643 | 50071224 | static Scalar thermalConductivity(const FluidState &fluidState, | |
644 | const int phaseIdx) | ||
645 | { | ||
646 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 35736736 times.
|
35736744 | assert(0 <= phaseIdx && phaseIdx < numPhases); |
647 | |||
648 |
2/2✓ Branch 0 taken 17855948 times.
✓ Branch 1 taken 17880780 times.
|
51148548 | const Scalar temperature = fluidState.temperature(phaseIdx) ; |
649 | 50071226 | const Scalar pressure = fluidState.pressure(phaseIdx); | |
650 |
2/2✓ Branch 0 taken 17855952 times.
✓ Branch 1 taken 17880784 times.
|
35736744 | if (phaseIdx == liquidPhaseIdx) |
651 | { | ||
652 | 33256690 | return H2O::liquidThermalConductivity(temperature, pressure); | |
653 | } | ||
654 | else | ||
655 | { | ||
656 | 17353198 | Scalar lambdaPureN2 = N2::gasThermalConductivity(temperature, pressure); | |
657 | if (!Policy::useN2HeatConductivityAsGasMixtureHeatConductivity()) | ||
658 | { | ||
659 | 11637326 | Scalar xN2 = fluidState.moleFraction(phaseIdx, N2Idx); | |
660 | 11637326 | Scalar xH2O = fluidState.moleFraction(phaseIdx, H2OIdx); | |
661 | 11637326 | Scalar lambdaN2 = xN2 * lambdaPureN2; | |
662 | 11637326 | Scalar partialPressure = pressure * xH2O; | |
663 | 11637326 | Scalar lambdaH2O = xH2O * H2O::gasThermalConductivity(temperature, partialPressure); | |
664 | 11637326 | return lambdaN2 + lambdaH2O; | |
665 | } | ||
666 | else | ||
667 | 6243462 | return lambdaPureN2; | |
668 | } | ||
669 | } | ||
670 | |||
671 | using Base<Scalar, ThisType>::heatCapacity; | ||
672 | //! \copydoc Dumux::FluidSystems::Base<Scalar,ThisType>::heatCapacity(const FluidState&,int) | ||
673 | template <class FluidState> | ||
674 | 2964311 | static Scalar heatCapacity(const FluidState &fluidState, | |
675 | int phaseIdx) | ||
676 | { | ||
677 |
2/2✓ Branch 0 taken 786240 times.
✓ Branch 1 taken 786240 times.
|
1572488 | if (phaseIdx == liquidPhaseIdx) { |
678 | 2178067 | return H2O::liquidHeatCapacity(fluidState.temperature(phaseIdx), | |
679 | 786244 | fluidState.pressure(phaseIdx)); | |
680 | } | ||
681 | |||
682 | // for the gas phase, assume ideal mixture | ||
683 | Scalar c_pN2; | ||
684 | Scalar c_pH2O; | ||
685 | // let the water and nitrogen components do things their own way | ||
686 | if (!Policy::useIdealGasHeatCapacities()) { | ||
687 | 2 | c_pN2 = N2::gasHeatCapacity(fluidState.temperature(phaseIdx), | |
688 | 2 | fluidState.pressure(phaseIdx) | |
689 | 2 | * fluidState.moleFraction(phaseIdx, N2Idx)); | |
690 | |||
691 | 2 | c_pH2O = H2O::gasHeatCapacity(fluidState.temperature(phaseIdx), | |
692 | 2 | fluidState.pressure(phaseIdx) | |
693 | 2 | * fluidState.moleFraction(phaseIdx, H2OIdx)); | |
694 | } | ||
695 | else { | ||
696 | // assume an ideal gas for both components. See: | ||
697 | // http://en.wikipedia.org/wiki/Heat_capacity | ||
698 | 786242 | Scalar c_vN2molar = Constants<Scalar>::R*2.39; | |
699 | 786242 | Scalar c_pN2molar = Constants<Scalar>::R + c_vN2molar; | |
700 | |||
701 | 786242 | Scalar c_vH2Omolar = Constants<Scalar>::R*3.37; // <- correct?? | |
702 | 786242 | Scalar c_pH2Omolar = Constants<Scalar>::R + c_vH2Omolar; | |
703 | |||
704 | 786242 | c_pN2 = c_pN2molar/molarMass(N2Idx); | |
705 | 786242 | c_pH2O = c_pH2Omolar/molarMass(H2OIdx); | |
706 | } | ||
707 | |||
708 | // mangle both components together | ||
709 | 786244 | return c_pH2O*fluidState.massFraction(gasPhaseIdx, H2OIdx) | |
710 | 786244 | + c_pN2*fluidState.massFraction(gasPhaseIdx, N2Idx); | |
711 | } | ||
712 | }; | ||
713 | |||
714 | } // end namespace Dumux::FluidSystems | ||
715 | |||
716 | #endif | ||
717 |