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 FluidStates | ||
10 | * \brief Represents all relevant thermodynamic quantities of a | ||
11 | * multi-phase, multi-component fluid system assuming | ||
12 | * thermodynamic equilibrium. | ||
13 | */ | ||
14 | #ifndef DUMUX_COMPOSITIONAL_FLUID_STATE_HH | ||
15 | #define DUMUX_COMPOSITIONAL_FLUID_STATE_HH | ||
16 | |||
17 | #include <algorithm> | ||
18 | #include <cmath> | ||
19 | #include <type_traits> | ||
20 | #include <cassert> | ||
21 | #include <array> | ||
22 | |||
23 | #include <dune/common/exceptions.hh> | ||
24 | |||
25 | namespace Dumux { | ||
26 | |||
27 | /*! | ||
28 | * \ingroup FluidStates | ||
29 | * \brief Represents all relevant thermodynamic quantities of a | ||
30 | * multi-phase, multi-component fluid system assuming | ||
31 | * thermodynamic equilibrium. | ||
32 | */ | ||
33 | template <class ScalarType, class FluidSystem> | ||
34 | class CompositionalFluidState | ||
35 | { | ||
36 | public: | ||
37 | static constexpr int numPhases = FluidSystem::numPhases; | ||
38 | static constexpr int numComponents = FluidSystem::numComponents; | ||
39 | |||
40 | //! export the scalar type | ||
41 | using Scalar = ScalarType; | ||
42 | |||
43 | //! default constructor | ||
44 |
10/13✗ Branch 0 not taken.
✓ Branch 1 taken 201727 times.
✓ Branch 2 taken 17310 times.
✓ Branch 3 taken 22925 times.
✓ Branch 4 taken 24160 times.
✓ Branch 5 taken 46576 times.
✓ Branch 6 taken 185984 times.
✓ Branch 7 taken 83008 times.
✓ Branch 9 taken 1359115 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✓ Branch 12 taken 20328 times.
✓ Branch 8 taken 24 times.
|
36301225 | CompositionalFluidState() = default; |
45 | |||
46 | //! copy constructor from arbitrary fluid state | ||
47 | template <class FluidState, typename std::enable_if_t<!std::is_same<FluidState, CompositionalFluidState>::value, int> = 0> | ||
48 | explicit CompositionalFluidState(const FluidState &fs) | ||
49 | { assign(fs); } | ||
50 | |||
51 | // copy and move constructor / assignment operator | ||
52 | CompositionalFluidState(const CompositionalFluidState &fs) = default; | ||
53 | CompositionalFluidState(CompositionalFluidState &&fs) = default; | ||
54 | CompositionalFluidState& operator=(const CompositionalFluidState &fs) = default; | ||
55 | CompositionalFluidState& operator=(CompositionalFluidState &&fs) = default; | ||
56 | |||
57 | /***************************************************** | ||
58 | * Generic access to fluid properties (No assumptions | ||
59 | * on thermodynamic equilibrium required) | ||
60 | *****************************************************/ | ||
61 | /*! | ||
62 | * \brief Returns the index of the most wetting phase in the | ||
63 | * fluid-solid configuration (for porous medium systems). | ||
64 | */ | ||
65 |
1/2✓ Branch 1 taken 41660 times.
✗ Branch 2 not taken.
|
100901272 | int wettingPhase() const { return wPhaseIdx_; } |
66 | |||
67 | /*! | ||
68 | * \brief Returns the saturation \f$S_\alpha\f$ of a fluid phase \f$\alpha\f$ in \f$\mathrm{[-]}\f$. | ||
69 | * | ||
70 | * The saturation is defined as the pore space occupied by the fluid divided by the total pore space: | ||
71 | * \f[S_\alpha := \frac{\phi \mathcal{V}_\alpha}{\phi \mathcal{V}}\f] | ||
72 | * | ||
73 | * \param phaseIdx the index of the phase | ||
74 | */ | ||
75 | 3138605455 | Scalar saturation(int phaseIdx) const | |
76 |
22/24✓ Branch 0 taken 440689 times.
✓ Branch 1 taken 6825783 times.
✓ Branch 2 taken 48469778 times.
✓ Branch 3 taken 96112244 times.
✓ Branch 5 taken 445865 times.
✓ Branch 6 taken 190525 times.
✓ Branch 8 taken 399578 times.
✓ Branch 9 taken 120929 times.
✓ Branch 10 taken 159685 times.
✓ Branch 11 taken 127928 times.
✓ Branch 13 taken 75446013 times.
✓ Branch 14 taken 10258010 times.
✓ Branch 15 taken 651 times.
✓ Branch 16 taken 302023 times.
✓ Branch 18 taken 30553 times.
✓ Branch 19 taken 54025 times.
✓ Branch 21 taken 30553 times.
✓ Branch 22 taken 54025 times.
✓ Branch 4 taken 708295 times.
✓ Branch 12 taken 86633770 times.
✓ Branch 7 taken 949336 times.
✓ Branch 17 taken 124903 times.
✗ Branch 20 not taken.
✗ Branch 23 not taken.
|
2517912616 | { return saturation_[phaseIdx]; } |
77 | |||
78 | /*! | ||
79 | * \brief Returns the molar fraction \f$x^\kappa_\alpha\f$ of the component \f$\kappa\f$ in fluid phase \f$\alpha\f$ in \f$\mathrm{[-]}\f$. | ||
80 | * | ||
81 | * The molar fraction \f$x^\kappa_\alpha\f$ is defined as the ratio of the number of molecules | ||
82 | * of component \f$\kappa\f$ and the total number of molecules of the phase \f$\alpha\f$. | ||
83 | * | ||
84 | * \param phaseIdx the index of the phase | ||
85 | * \param compIdx the index of the component | ||
86 | */ | ||
87 | 10535490740 | Scalar moleFraction(int phaseIdx, int compIdx) const | |
88 |
13/21✓ Branch 0 taken 150202141 times.
✓ Branch 1 taken 22890971 times.
✓ Branch 2 taken 145721439 times.
✓ Branch 3 taken 42030072 times.
✓ Branch 4 taken 11950876 times.
✓ Branch 5 taken 1453869 times.
✓ Branch 7 taken 1206498 times.
✓ Branch 8 taken 2933131 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2595772 times.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 39 times.
✓ Branch 18 taken 42400 times.
✗ Branch 19 not taken.
✓ Branch 20 taken 479 times.
✓ Branch 6 taken 127929 times.
|
3579491083 | { return moleFraction_[phaseIdx][compIdx]; } |
89 | |||
90 | /*! | ||
91 | * \brief Returns the mass fraction \f$X^\kappa_\alpha\f$ of component \f$\kappa\f$ in fluid phase \f$\alpha\f$ in \f$\mathrm{[-]}\f$. | ||
92 | * | ||
93 | * The mass fraction \f$X^\kappa_\alpha\f$ is defined as the weight of all molecules of a | ||
94 | * component divided by the total mass of the fluid phase. It is related with the component's mole fraction by means of the relation | ||
95 | * | ||
96 | * \f[X^\kappa_\alpha = x^\kappa_\alpha \frac{M^\kappa}{\overline M_\alpha}\;,\f] | ||
97 | * | ||
98 | * where \f$M^\kappa\f$ is the molar mass of component \f$\kappa\f$ and | ||
99 | * \f$\overline M_\alpha\f$ is the mean molar mass of a molecule of phase | ||
100 | * \f$\alpha\f$. | ||
101 | * | ||
102 | * \param phaseIdx the index of the phase | ||
103 | * \param compIdx the index of the component | ||
104 | */ | ||
105 | 6217467920 | Scalar massFraction(int phaseIdx, int compIdx) const | |
106 | { | ||
107 | // calculate the mass fractions: | ||
108 | // for "mass" models this is just a back calculation | ||
109 |
7/9✓ Branch 2 taken 7951615 times.
✓ Branch 3 taken 316438171 times.
✓ Branch 5 taken 22372 times.
✗ Branch 6 not taken.
✓ Branch 0 taken 266466812 times.
✓ Branch 1 taken 277847878 times.
✓ Branch 4 taken 12577120 times.
✓ Branch 7 taken 5120 times.
✗ Branch 8 not taken.
|
2620148767 | return sumMoleFractions_[phaseIdx] |
110 |
8/10✓ Branch 3 taken 304613419 times.
✓ Branch 4 taken 8669435 times.
✓ Branch 6 taken 8058208 times.
✗ Branch 7 not taken.
✓ Branch 0 taken 266466812 times.
✓ Branch 1 taken 277847878 times.
✓ Branch 5 taken 14060060 times.
✓ Branch 9 taken 5120 times.
✗ Branch 10 not taken.
✓ Branch 2 taken 1588156 times.
|
5970809174 | * moleFraction(phaseIdx, compIdx) |
111 |
2/4✓ Branch 4 taken 4992 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 5120 times.
✗ Branch 7 not taken.
|
5973768770 | * FluidSystem::molarMass(compIdx) |
112 |
4/4✓ Branch 0 taken 16917 times.
✓ Branch 1 taken 168803 times.
✓ Branch 2 taken 325556 times.
✓ Branch 3 taken 12044 times.
|
1943295718 | / averageMolarMass_[phaseIdx]; |
113 | } | ||
114 | |||
115 | /*! | ||
116 | * \brief Returns the phase mass fraction, i.e. phase mass per total mass \f$\mathrm{[kg/kg]}\f$. | ||
117 | * \param phaseIdx the index of the phase | ||
118 | */ | ||
119 | Scalar phaseMassFraction(int phaseIdx) const | ||
120 | { | ||
121 | Scalar totalMass = 0.0; | ||
122 |
2/2✓ Branch 0 taken 20 times.
✓ Branch 1 taken 10 times.
|
30 | for (int pIdx = 0; pIdx < numPhases; ++pIdx) |
123 | 20 | totalMass += saturation(pIdx)*density(pIdx); | |
124 | |||
125 | 10 | return saturation(phaseIdx)*density(phaseIdx) / totalMass; | |
126 | } | ||
127 | |||
128 | /*! | ||
129 | * \brief The average molar mass \f$\overline M_\alpha\f$ of phase \f$\alpha\f$ in \f$\mathrm{[kg/mol]}\f$ | ||
130 | * | ||
131 | * The average molar mass is the mean mass of a mole of the | ||
132 | * fluid at current composition. It is defined as the sum of the | ||
133 | * component's molar masses weighted by the current mole fraction: | ||
134 | * \f[\mathrm{ \overline M_\alpha = \sum_\kappa M^\kappa x_\alpha^\kappa}\f] | ||
135 | */ | ||
136 | 92375929 | Scalar averageMolarMass(int phaseIdx) const | |
137 |
4/9✓ Branch 1 taken 56828 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 0 not taken.
|
91985600 | { return averageMolarMass_[phaseIdx]; } |
138 | |||
139 | /*! | ||
140 | * \brief The molar concentration \f$c^\kappa_\alpha\f$ of component \f$\kappa\f$ in fluid phase \f$\alpha\f$ in \f$\mathrm{[mol/m^3]}\f$ | ||
141 | * | ||
142 | * This quantity is usually called "molar concentration" or just | ||
143 | * "concentration", but there are many other (though less common) | ||
144 | * measures for concentration. | ||
145 | * | ||
146 | * http://en.wikipedia.org/wiki/Concentration | ||
147 | */ | ||
148 | 875226 | Scalar molarity(int phaseIdx, int compIdx) const | |
149 | 2673 | { return molarDensity(phaseIdx)*moleFraction(phaseIdx, compIdx); } | |
150 | |||
151 | /*! | ||
152 | * \brief The fugacity \f$f^\kappa_\alpha\f$ of component \f$\kappa\f$ | ||
153 | * in fluid phase \f$\alpha\f$ in \f$\mathrm{[Pa]}\f$ | ||
154 | * | ||
155 | * The fugacity is defined as: | ||
156 | * \f$f_\alpha^\kappa := \Phi^\kappa_\alpha x^\kappa_\alpha p_\alpha \;,\f$ | ||
157 | * where \f$\Phi^\kappa_\alpha\f$ is the fugacity coefficient \cite reid1987 . | ||
158 | * The physical meaning of fugacity becomes clear from the equation: | ||
159 | * \f[f_\alpha^\kappa = p_\alpha \exp\left\{\frac{\zeta^\kappa_\alpha}{R T_\alpha} \right\} \;,\f] | ||
160 | * where \f$\zeta^\kappa_\alpha\f$ represents the \f$\kappa\f$'s chemical | ||
161 | * potential in phase \f$\alpha\f$, \f$R\f$ stands for the ideal gas constant, | ||
162 | * and \f$T_\alpha\f$ for the absolute temperature of phase \f$\alpha\f$. Assuming thermal equilibrium, | ||
163 | * there is a one-to-one mapping between a component's chemical potential | ||
164 | * \f$\zeta^\kappa_\alpha\f$ and its fugacity \f$f^\kappa_\alpha\f$. In this | ||
165 | * case chemical equilibrium can thus be expressed by: | ||
166 | * \f[f^\kappa := f^\kappa_\alpha = f^\kappa_\beta\quad\forall \alpha, \beta\f] | ||
167 | */ | ||
168 | 51090861 | Scalar fugacity(int phaseIdx, int compIdx) const | |
169 | 51047431 | { return fugacityCoefficient(phaseIdx, compIdx)*moleFraction(phaseIdx, compIdx)*pressure(phaseIdx); } | |
170 | |||
171 | /*! | ||
172 | * \brief The fugacity coefficient \f$\Phi^\kappa_\alpha\f$ of component \f$\kappa\f$ in fluid phase \f$\alpha\f$ in \f$\mathrm{[-]}\f$ | ||
173 | */ | ||
174 | 74937962 | Scalar fugacityCoefficient(int phaseIdx, int compIdx) const | |
175 |
2/6✗ Branch 0 not taken.
✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 5 not taken.
✓ Branch 6 taken 1182423 times.
|
125362471 | { return fugacityCoefficient_[phaseIdx][compIdx]; } |
176 | |||
177 | /*! | ||
178 | * \brief The molar volume \f$v_{mol,\alpha}\f$ of a fluid phase \f$\alpha\f$ in \f$\mathrm{[m^3/mol]}\f$ | ||
179 | * | ||
180 | * This quantity is the inverse of the molar density. | ||
181 | */ | ||
182 | Scalar molarVolume(int phaseIdx) const | ||
183 |
4/8✓ 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.
|
4 | { return 1.0/molarDensity(phaseIdx); } |
184 | |||
185 | /*! | ||
186 | * \brief The mass density \f$\rho_\alpha\f$ of the fluid phase | ||
187 | * \f$\alpha\f$ in \f$\mathrm{[kg/m^3]}\f$ | ||
188 | */ | ||
189 | 9451847790 | Scalar density(int phaseIdx) const | |
190 |
30/50✓ Branch 10 taken 58968056 times.
✓ Branch 11 taken 55981282 times.
✓ Branch 13 taken 2767555 times.
✗ Branch 14 not taken.
✓ Branch 16 taken 13625 times.
✓ Branch 17 taken 20650 times.
✓ Branch 19 taken 14075 times.
✗ Branch 20 not taken.
✓ Branch 22 taken 7025 times.
✗ Branch 23 not taken.
✓ Branch 25 taken 25 times.
✗ Branch 26 not taken.
✓ Branch 28 taken 25 times.
✗ Branch 29 not taken.
✓ Branch 31 taken 25 times.
✗ Branch 32 not taken.
✓ Branch 33 taken 25 times.
✗ Branch 34 not taken.
✓ Branch 35 taken 25 times.
✗ Branch 36 not taken.
✗ Branch 37 not taken.
✗ Branch 38 not taken.
✗ Branch 39 not taken.
✗ Branch 40 not taken.
✗ Branch 41 not taken.
✗ Branch 42 not taken.
✗ Branch 43 not taken.
✗ Branch 44 not taken.
✓ Branch 45 taken 25 times.
✗ Branch 46 not taken.
✓ Branch 47 taken 25 times.
✗ Branch 48 not taken.
✓ Branch 49 taken 25 times.
✗ Branch 50 not taken.
✓ Branch 51 taken 25 times.
✗ Branch 52 not taken.
✓ Branch 2 taken 78511068 times.
✓ Branch 3 taken 14303092 times.
✓ Branch 4 taken 18528739 times.
✓ Branch 5 taken 50231695 times.
✓ Branch 8 taken 60938241 times.
✓ Branch 9 taken 63373607 times.
✓ Branch 12 taken 11461310 times.
✓ Branch 6 taken 4348592 times.
✓ Branch 15 taken 14450 times.
✓ Branch 7 taken 5637626 times.
✓ Branch 18 taken 14400 times.
✓ Branch 24 taken 6200 times.
✓ Branch 1 taken 54721888 times.
✓ Branch 0 taken 1885144 times.
|
8801908181 | { return density_[phaseIdx]; } |
191 | |||
192 | /*! | ||
193 | * \brief The molar density \f$\rho_\alpha\f$ of the fluid phase | ||
194 | * \f$\alpha\f$ in \f$\mathrm{[mol/m^3]}\f$ | ||
195 | */ | ||
196 | 2095651040 | Scalar molarDensity(int phaseIdx) const | |
197 |
9/11✓ Branch 1 taken 75461916 times.
✓ Branch 2 taken 34606 times.
✓ Branch 4 taken 129330 times.
✓ Branch 5 taken 16157538 times.
✓ Branch 7 taken 1 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 1 times.
✗ Branch 11 not taken.
✓ Branch 3 taken 135044219 times.
✓ Branch 6 taken 3053753 times.
✓ Branch 0 taken 270887 times.
|
2097616086 | { return molarDensity_[phaseIdx]; } |
198 | |||
199 | /*! | ||
200 | * \brief The absolute temperature\f$T_\alpha\f$ of a fluid phase \f$\alpha\f$ in \f$\mathrm{[K]}\f$ | ||
201 | */ | ||
202 |
2/4✓ Branch 0 taken 1382068 times.
✓ Branch 1 taken 1382068 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
2447972041 | Scalar temperature(int phaseIdx) const |
203 |
31/32✓ Branch 0 taken 132998141 times.
✓ Branch 1 taken 152809017 times.
✓ Branch 3 taken 188984641 times.
✓ Branch 4 taken 135495740 times.
✓ Branch 5 taken 98838062 times.
✓ Branch 6 taken 174253372 times.
✓ Branch 7 taken 121387615 times.
✓ Branch 8 taken 93242589 times.
✓ Branch 9 taken 161340380 times.
✓ Branch 10 taken 96114826 times.
✓ Branch 11 taken 60879562 times.
✓ Branch 12 taken 77745390 times.
✓ Branch 13 taken 104841214 times.
✓ Branch 14 taken 26696460 times.
✓ Branch 15 taken 19095078 times.
✓ Branch 16 taken 39894657 times.
✓ Branch 17 taken 3191433 times.
✓ Branch 18 taken 29176989 times.
✓ Branch 19 taken 20819222 times.
✓ Branch 20 taken 2532627 times.
✓ Branch 22 taken 296074 times.
✓ Branch 23 taken 987888 times.
✓ Branch 24 taken 3532476 times.
✓ Branch 25 taken 8864489 times.
✓ Branch 26 taken 1184296 times.
✓ Branch 27 taken 4755174 times.
✓ Branch 28 taken 251840 times.
✓ Branch 29 taken 251840 times.
✗ Branch 30 not taken.
✓ Branch 31 taken 4266828 times.
✓ Branch 2 taken 95536692 times.
✓ Branch 21 taken 6001794 times.
|
2832106816 | { return temperature_[phaseIdx]; } |
204 | |||
205 | /*! | ||
206 | * \brief The pressure \f$p_\alpha\f$ of a fluid phase \f$\alpha\f$ in \f$\mathrm{[Pa]}\f$ | ||
207 | */ | ||
208 | 4737520492 | Scalar pressure(int phaseIdx) const | |
209 |
32/34✓ Branch 1 taken 170154308 times.
✓ Branch 2 taken 127789666 times.
✓ Branch 3 taken 190836411 times.
✓ Branch 4 taken 113288547 times.
✓ Branch 5 taken 85911882 times.
✓ Branch 6 taken 127293551 times.
✓ Branch 7 taken 109766946 times.
✓ Branch 8 taken 98659821 times.
✓ Branch 10 taken 28738782 times.
✓ Branch 11 taken 38374081 times.
✓ Branch 12 taken 74690316 times.
✓ Branch 13 taken 42337193 times.
✓ Branch 14 taken 29956232 times.
✓ Branch 15 taken 1752328 times.
✓ Branch 16 taken 12537314 times.
✓ Branch 17 taken 16986686 times.
✓ Branch 18 taken 23613217 times.
✓ Branch 19 taken 365506 times.
✓ Branch 20 taken 1200891 times.
✓ Branch 21 taken 1090128 times.
✓ Branch 22 taken 163300 times.
✓ Branch 23 taken 7410166 times.
✓ Branch 24 taken 8254351 times.
✓ Branch 25 taken 6018171 times.
✓ Branch 26 taken 251840 times.
✓ Branch 27 taken 251840 times.
✗ Branch 28 not taken.
✓ Branch 29 taken 48356 times.
✓ Branch 30 taken 119459 times.
✓ Branch 31 taken 6461 times.
✗ Branch 32 not taken.
✓ Branch 33 taken 4266828 times.
✓ Branch 9 taken 66070630 times.
✓ Branch 0 taken 82204540 times.
|
4309367073 | { return pressure_[phaseIdx]; } |
210 | |||
211 | /*! | ||
212 | * \brief The partial pressure of a component in a phase \f$\mathrm{[Pa]}\f$ | ||
213 | * \todo is this needed? | ||
214 | */ | ||
215 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 5191606 times.
|
81067418 | Scalar partialPressure(int phaseIdx, int compIdx) const |
216 | { | ||
217 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 5191606 times.
|
45725315 | assert(FluidSystem::isGas(phaseIdx)); |
218 | 40533771 | return moleFraction(phaseIdx, compIdx) * pressure(phaseIdx); | |
219 | } | ||
220 | |||
221 | /*! | ||
222 | * \brief The specific enthalpy \f$h_\alpha\f$ of a fluid phase \f$\alpha\f$ in \f$\mathrm{[J/kg]}\f$ | ||
223 | */ | ||
224 | 426389552 | Scalar enthalpy(int phaseIdx) const | |
225 |
7/10✓ Branch 1 taken 2959 times.
✓ Branch 2 taken 781 times.
✓ Branch 4 taken 48357 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 1 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 1 times.
✗ Branch 11 not taken.
✓ Branch 0 taken 14551 times.
✓ Branch 3 taken 1475 times.
|
424236056 | { return enthalpy_[phaseIdx]; } |
226 | |||
227 | /*! | ||
228 | * \brief The specific internal energy \f$u_\alpha\f$ of a fluid phase \f$\alpha\f$ in \f$\mathrm{[J/kg]}\f$ | ||
229 | * | ||
230 | * The specific internal energy is defined by the relation: | ||
231 | * | ||
232 | * \f[u_\alpha = h_\alpha - \frac{p_\alpha}{\rho_\alpha}\f] | ||
233 | */ | ||
234 | 574082044 | Scalar internalEnergy(int phaseIdx) const | |
235 |
4/8✓ 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.
|
573890084 | { return enthalpy_[phaseIdx] - pressure(phaseIdx)/density(phaseIdx); } |
236 | |||
237 | /*! | ||
238 | * \brief The dynamic viscosity \f$\mu_\alpha\f$ of fluid phase \f$\alpha\f$ in \f$\mathrm{[Pa s]}\f$ | ||
239 | */ | ||
240 | 3464163002 | Scalar viscosity(int phaseIdx) const | |
241 |
16/34✓ Branch 2 taken 6644644 times.
✓ Branch 3 taken 15498810 times.
✓ Branch 4 taken 8724342 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 6332783 times.
✓ Branch 7 taken 1 times.
✓ Branch 8 taken 6731811 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 4295040 times.
✗ Branch 11 not taken.
✓ Branch 12 taken 25 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 25 times.
✗ Branch 15 not taken.
✓ Branch 16 taken 25 times.
✗ Branch 17 not taken.
✓ Branch 18 taken 25 times.
✗ Branch 19 not taken.
✓ Branch 20 taken 25 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✗ Branch 27 not taken.
✗ Branch 28 not taken.
✗ Branch 29 not taken.
✓ Branch 30 taken 25 times.
✗ Branch 31 not taken.
✓ Branch 32 taken 25 times.
✗ Branch 33 not taken.
✓ Branch 0 taken 3568529 times.
✓ Branch 1 taken 25268679 times.
|
935924047 | { return viscosity_[phaseIdx]; } |
242 | |||
243 | |||
244 | /***************************************************** | ||
245 | * Access to fluid properties which only make sense | ||
246 | * if assuming thermodynamic equilibrium | ||
247 | *****************************************************/ | ||
248 | |||
249 | /*! | ||
250 | * \brief The temperature within the domain \f$\mathrm{[K]}\f$ | ||
251 | */ | ||
252 | 164018253 | Scalar temperature() const | |
253 |
3/6✓ Branch 0 taken 9858644 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 34736 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 36288 times.
✗ Branch 5 not taken.
|
135705760 | { return temperature_[0]; } |
254 | |||
255 | /*! | ||
256 | * \brief The fugacity of a component \f$\mathrm{[Pa]}\f$ | ||
257 | * | ||
258 | * This assumes chemical equilibrium. | ||
259 | */ | ||
260 | 38400 | Scalar fugacity(int compIdx) const | |
261 | 38400 | { return fugacity(0, compIdx); } | |
262 | |||
263 | |||
264 | /***************************************************** | ||
265 | * Setter methods. Note that these are not part of the | ||
266 | * generic FluidState interface but specific for each | ||
267 | * implementation... | ||
268 | *****************************************************/ | ||
269 | |||
270 | /*! | ||
271 | * \brief Retrieve all parameters from an arbitrary fluid | ||
272 | * state. | ||
273 | * | ||
274 | * \note If the other fluid state object is inconsistent with the | ||
275 | * thermodynamic equilibrium, the result of this method is | ||
276 | * undefined. | ||
277 | */ | ||
278 | template <class FluidState> | ||
279 | 1 | void assign(const FluidState &fs) | |
280 | { | ||
281 |
2/2✓ Branch 0 taken 3 times.
✓ Branch 1 taken 1 times.
|
4 | for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) |
282 | { | ||
283 | 3 | averageMolarMass_[phaseIdx] = 0; | |
284 | 3 | sumMoleFractions_[phaseIdx] = 0; | |
285 | 3 | temperature_[phaseIdx] = fs.temperature(); | |
286 |
2/2✓ Branch 0 taken 21 times.
✓ Branch 1 taken 3 times.
|
24 | for (int compIdx = 0; compIdx < numComponents; ++compIdx) |
287 | { | ||
288 | 21 | moleFraction_[phaseIdx][compIdx] = fs.moleFraction(phaseIdx, compIdx); | |
289 | 21 | fugacityCoefficient_[phaseIdx][compIdx] = fs.fugacityCoefficient(phaseIdx, compIdx); | |
290 | 21 | averageMolarMass_[phaseIdx] += moleFraction_[phaseIdx][compIdx]*FluidSystem::molarMass(compIdx); | |
291 | 21 | sumMoleFractions_[phaseIdx] += moleFraction_[phaseIdx][compIdx]; | |
292 | } | ||
293 | 3 | pressure_[phaseIdx] = fs.pressure(phaseIdx); | |
294 | 3 | saturation_[phaseIdx] = fs.saturation(phaseIdx); | |
295 | 3 | density_[phaseIdx] = fs.density(phaseIdx); | |
296 | 3 | molarDensity_[phaseIdx] = fs.molarDensity(phaseIdx); | |
297 | 3 | enthalpy_[phaseIdx] = fs.enthalpy(phaseIdx); | |
298 | 3 | viscosity_[phaseIdx] = fs.viscosity(phaseIdx); | |
299 | } | ||
300 | 1 | wPhaseIdx_ = fs.wettingPhase(); | |
301 | 1 | } | |
302 | |||
303 | /*! | ||
304 | * \brief Set the temperature \f$\mathrm{[K]}\f$ of all phases. | ||
305 | */ | ||
306 | 25963306 | void setTemperature(Scalar value) | |
307 | { | ||
308 |
78/78✓ Branch 0 taken 11571700 times.
✓ Branch 1 taken 3858352 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 5 times.
✓ Branch 4 taken 3 times.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 1 times.
✓ Branch 7 taken 1 times.
✓ Branch 9 taken 3 times.
✓ Branch 10 taken 1 times.
✓ Branch 11 taken 1 times.
✓ Branch 12 taken 1 times.
✓ Branch 13 taken 1 times.
✓ Branch 14 taken 1 times.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 1 times.
✓ Branch 17 taken 2 times.
✓ Branch 18 taken 1 times.
✓ Branch 19 taken 2 times.
✓ Branch 20 taken 1 times.
✓ Branch 21 taken 2 times.
✓ Branch 22 taken 1 times.
✓ Branch 23 taken 2 times.
✓ Branch 24 taken 1 times.
✓ Branch 25 taken 2 times.
✓ Branch 26 taken 1 times.
✓ Branch 28 taken 3 times.
✓ Branch 29 taken 1 times.
✓ Branch 31 taken 3 times.
✓ Branch 32 taken 1 times.
✓ Branch 34 taken 3 times.
✓ Branch 35 taken 1 times.
✓ Branch 36 taken 1 times.
✓ Branch 37 taken 1 times.
✓ Branch 38 taken 2 times.
✓ Branch 39 taken 1 times.
✓ Branch 40 taken 2 times.
✓ Branch 41 taken 1 times.
✓ Branch 42 taken 2 times.
✓ Branch 43 taken 1 times.
✓ Branch 44 taken 2 times.
✓ Branch 45 taken 1 times.
✓ Branch 46 taken 2 times.
✓ Branch 47 taken 1 times.
✓ Branch 48 taken 2 times.
✓ Branch 49 taken 1 times.
✓ Branch 50 taken 2 times.
✓ Branch 51 taken 1 times.
✓ Branch 52 taken 2 times.
✓ Branch 53 taken 1 times.
✓ Branch 54 taken 2 times.
✓ Branch 55 taken 1 times.
✓ Branch 56 taken 2 times.
✓ Branch 57 taken 1 times.
✓ Branch 58 taken 2 times.
✓ Branch 59 taken 1 times.
✓ Branch 60 taken 2 times.
✓ Branch 61 taken 1 times.
✓ Branch 62 taken 2 times.
✓ Branch 63 taken 1 times.
✓ Branch 64 taken 2 times.
✓ Branch 65 taken 1 times.
✓ Branch 66 taken 2 times.
✓ Branch 67 taken 1 times.
✓ Branch 68 taken 2 times.
✓ Branch 69 taken 1 times.
✓ Branch 70 taken 1 times.
✓ Branch 71 taken 1 times.
✓ Branch 73 taken 3 times.
✓ Branch 74 taken 1 times.
✓ Branch 75 taken 2 times.
✓ Branch 76 taken 1 times.
✓ Branch 77 taken 2 times.
✓ Branch 78 taken 1 times.
✓ Branch 79 taken 2 times.
✓ Branch 80 taken 1 times.
✓ Branch 81 taken 2 times.
✓ Branch 82 taken 1 times.
|
59445088 | for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) |
309 |
1/2✓ Branch 1 taken 17 times.
✗ Branch 2 not taken.
|
55631238 | temperature_[phaseIdx] = value; |
310 | } | ||
311 | |||
312 | /*! | ||
313 | * \brief Set the temperature \f$\mathrm{[K]}\f$ of a specific phase. | ||
314 | * This is not implemented in this fluidstate. | ||
315 | */ | ||
316 | 271754777 | void setTemperature(const int phaseIdx, const Scalar value) | |
317 |
1/2✓ Branch 0 taken 329394 times.
✗ Branch 1 not taken.
|
271754777 | { temperature_[phaseIdx] = value; } |
318 | |||
319 | /*! | ||
320 | * \brief Set the fluid pressure of a phase \f$\mathrm{[Pa]}\f$ | ||
321 | */ | ||
322 | 199050970 | void setPressure(int phaseIdx, Scalar value) | |
323 |
2/3✓ Branch 0 taken 20279268 times.
✓ Branch 1 taken 76201401 times.
✗ Branch 2 not taken.
|
154546999 | { pressure_[phaseIdx] = value; } |
324 | |||
325 | /*! | ||
326 | * \brief Set the saturation of a phase \f$\mathrm{[-]}\f$ | ||
327 | */ | ||
328 | 197813977 | void setSaturation(int phaseIdx, Scalar value) | |
329 |
2/2✓ Branch 0 taken 1207452 times.
✓ Branch 1 taken 3268398 times.
|
197813969 | { saturation_[phaseIdx] = value; } |
330 | |||
331 | /*! | ||
332 | * \brief Set the mole fraction of a component in a phase \f$\mathrm{[-]}\f$ | ||
333 | * and update the average molar mass \f$\mathrm{[kg/mol]}\f$ according | ||
334 | * to the current composition of the phase | ||
335 | */ | ||
336 | 593352574 | void setMoleFraction(int phaseIdx, int compIdx, Scalar value) | |
337 | { | ||
338 | 593397119 | moleFraction_[phaseIdx][compIdx] = value; | |
339 | |||
340 | // re-calculate the mean molar mass | ||
341 | 593397119 | sumMoleFractions_[phaseIdx] = 0.0; | |
342 | 593397119 | averageMolarMass_[phaseIdx] = 0.0; | |
343 |
26/26✓ Branch 0 taken 1323809985 times.
✓ Branch 1 taken 466503550 times.
✓ Branch 2 taken 183177738 times.
✓ Branch 3 taken 63913634 times.
✓ Branch 4 taken 50077293 times.
✓ Branch 5 taken 19442503 times.
✓ Branch 6 taken 42697697 times.
✓ Branch 7 taken 20818086 times.
✓ Branch 8 taken 24581637 times.
✓ Branch 9 taken 12290816 times.
✓ Branch 10 taken 9 times.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 9 times.
✓ Branch 13 taken 2 times.
✓ Branch 14 taken 9 times.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 24 times.
✓ Branch 17 taken 12 times.
✓ Branch 18 taken 2 times.
✓ Branch 19 taken 1 times.
✓ Branch 20 taken 2 times.
✓ Branch 21 taken 1 times.
✓ Branch 22 taken 2 times.
✓ Branch 23 taken 1 times.
✓ Branch 24 taken 2 times.
✓ Branch 25 taken 1 times.
|
2226218212 | for (int compJIdx = 0; compJIdx < numComponents; ++compJIdx) |
344 | { | ||
345 |
6/10✓ Branch 0 taken 98154810 times.
✓ Branch 1 taken 98154820 times.
✓ Branch 2 taken 54415452 times.
✓ Branch 3 taken 108830904 times.
✓ Branch 4 taken 9069242 times.
✓ Branch 5 taken 18138484 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
|
1645112512 | sumMoleFractions_[phaseIdx] += moleFraction_[phaseIdx][compJIdx]; |
346 |
6/10✓ Branch 0 taken 98154810 times.
✓ Branch 1 taken 98154820 times.
✓ Branch 2 taken 54415452 times.
✓ Branch 3 taken 108830904 times.
✓ Branch 4 taken 9069242 times.
✓ Branch 5 taken 18138484 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
|
2039822843 | averageMolarMass_[phaseIdx] += moleFraction_[phaseIdx][compJIdx]*FluidSystem::molarMass(compJIdx); |
347 | } | ||
348 | 521473781 | } | |
349 | |||
350 | /*! | ||
351 | * \brief Set the mass fraction of a component in a phase \f$\mathrm{[-]}\f$ | ||
352 | * and update the average molar mass \f$\mathrm{[kg/mol]}\f$ according | ||
353 | * to the current composition of the phase | ||
354 | */ | ||
355 | 68034467 | void setMassFraction(int phaseIdx, int compIdx, Scalar value) | |
356 | { | ||
357 | if (numComponents != 2) | ||
358 |
2/2✓ Branch 0 taken 1950 times.
✓ Branch 1 taken 1969 times.
|
68034467 | DUNE_THROW(Dune::NotImplemented, "This currently only works for 2 components."); |
359 | else | ||
360 | { | ||
361 | // calculate average molar mass of the gas phase | ||
362 | 68034467 | Scalar M1 = FluidSystem::molarMass(compIdx); | |
363 |
2/2✓ Branch 0 taken 1965 times.
✓ Branch 1 taken 1950 times.
|
68034467 | Scalar M2 = FluidSystem::molarMass(1-compIdx); |
364 | 68034467 | Scalar X2 = 1.0-value; | |
365 | 68034467 | Scalar avgMolarMass = M1*M2/(M2 + X2*(M1 - M2)); | |
366 | |||
367 | 68034467 | moleFraction_[phaseIdx][compIdx] = value * avgMolarMass / M1; | |
368 | 68034467 | moleFraction_[phaseIdx][1-compIdx] = 1.0-moleFraction_[phaseIdx][compIdx]; | |
369 | |||
370 | // re-calculate the mean molar mass | ||
371 | 68034467 | sumMoleFractions_[phaseIdx] = 0.0; | |
372 | 68034467 | averageMolarMass_[phaseIdx] = 0.0; | |
373 |
2/2✓ Branch 0 taken 136068934 times.
✓ Branch 1 taken 68034467 times.
|
204103401 | for (int compJIdx = 0; compJIdx < numComponents; ++compJIdx) { |
374 |
2/2✓ Branch 0 taken 3915 times.
✓ Branch 1 taken 3915 times.
|
136068934 | sumMoleFractions_[phaseIdx] += moleFraction_[phaseIdx][compJIdx]; |
375 |
2/2✓ Branch 0 taken 3915 times.
✓ Branch 1 taken 3915 times.
|
136076764 | averageMolarMass_[phaseIdx] += moleFraction_[phaseIdx][compJIdx]*FluidSystem::molarMass(compJIdx); |
376 | } | ||
377 | } | ||
378 | 68034467 | } | |
379 | |||
380 | /*! | ||
381 | * \brief Set the fugacity coefficient \f$\Phi^\kappa_\alpha\f$ of component \f$\kappa\f$ | ||
382 | * in fluid phase \f$\alpha\f$ in \f$\mathrm{[-]}\f$ | ||
383 | */ | ||
384 | 221254550 | void setFugacityCoefficient(int phaseIdx, int compIdx, Scalar value) | |
385 | 221254550 | { fugacityCoefficient_[phaseIdx][compIdx] = value; } | |
386 | |||
387 | /*! | ||
388 | * \brief Set the density of a phase \f$\mathrm{[kg / m^3]}\f$ | ||
389 | */ | ||
390 | 298762337 | void setDensity(int phaseIdx, Scalar value) | |
391 |
2/2✓ Branch 0 taken 23268 times.
✓ Branch 1 taken 11634 times.
|
297918667 | { density_[phaseIdx] = value; } |
392 | |||
393 | /*! | ||
394 | * \brief Set the molar density of a phase \f$\mathrm{[mol / m^3]}\f$ | ||
395 | */ | ||
396 | 298762337 | void setMolarDensity(int phaseIdx, Scalar value) | |
397 |
2/2✓ Branch 0 taken 23268 times.
✓ Branch 1 taken 11634 times.
|
279665717 | { molarDensity_[phaseIdx] = value; } |
398 | |||
399 | /*! | ||
400 | * \brief Set the specific enthalpy of a phase\f$\mathrm{[J/kg]}\f$ | ||
401 | */ | ||
402 | 315765238 | void setEnthalpy(int phaseIdx, Scalar value) | |
403 | 311043077 | { enthalpy_[phaseIdx] = value; } | |
404 | |||
405 | /*! | ||
406 | * \brief Set the dynamic viscosity of a phase \f$\mathrm{[Pa s]}\f$ | ||
407 | */ | ||
408 | 312690349 | void setViscosity(int phaseIdx, Scalar value) | |
409 | 312029998 | { viscosity_[phaseIdx] = value; } | |
410 | |||
411 | /*! | ||
412 | * \brief Set the index of the most wetting phase | ||
413 | */ | ||
414 | 102092394 | void setWettingPhase(int phaseIdx) | |
415 |
2/2✓ Branch 0 taken 12291843 times.
✓ Branch 1 taken 87380262 times.
|
102092394 | { wPhaseIdx_ = phaseIdx; } |
416 | |||
417 | protected: | ||
418 | //! zero-initialize all data members with braces syntax | ||
419 | std::array<std::array<Scalar, numComponents>, numPhases> moleFraction_ = {}; | ||
420 | std::array<std::array<Scalar, numComponents>, numPhases> fugacityCoefficient_ = {}; | ||
421 | std::array<Scalar, numPhases> averageMolarMass_ = {}; | ||
422 | std::array<Scalar, numPhases> sumMoleFractions_ = {}; | ||
423 | std::array<Scalar, numPhases> pressure_ = {}; | ||
424 | std::array<Scalar, numPhases> saturation_ = {}; | ||
425 | std::array<Scalar, numPhases> density_ = {}; | ||
426 | std::array<Scalar, numPhases> molarDensity_ = {}; | ||
427 | std::array<Scalar, numPhases> enthalpy_ = {}; | ||
428 | std::array<Scalar, numPhases> viscosity_ = {}; | ||
429 | std::array<Scalar, numPhases> temperature_ = {}; | ||
430 | |||
431 | int wPhaseIdx_{0}; | ||
432 | }; | ||
433 | |||
434 | } // end namespace Dumux | ||
435 | |||
436 | #endif | ||
437 |