3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
compositionalflash.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
25#ifndef DUMUX_COMPOSITIONAL_FLASH_HH
26#define DUMUX_COMPOSITIONAL_FLASH_HH
27
28#include <dune/common/fvector.hh>
29
31
32namespace Dumux {
33
49template <class Scalar, class FluidSystem>
51{
53
54 enum { numPhases = FluidSystem::numPhases,
55 numComponents = FluidSystem::numComponents
56 };
57
58 enum{
59 phase0Idx = FluidSystem::phase0Idx,
60 phase1Idx = FluidSystem::phase1Idx,
61 comp0Idx = FluidSystem::comp0Idx,
62 comp1Idx = FluidSystem::comp1Idx
63 };
64
65public:
66 using ComponentVector = Dune::FieldVector<Scalar, numComponents>;
67 using PhaseVector = Dune::FieldVector<Scalar, numPhases>;
91 template<class FluidState>
92 static void concentrationFlash2p2c(FluidState &fluidState,
93 const Scalar &Z0,
94 const PhaseVector &phasePressure,
95 const Scalar &porosity,
96 const Scalar &temperature)
97 {
98#ifndef NDEBUG
99 // this solver can only handle fluid systems which
100 // assume ideal mixtures of all fluids.
101 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
102 assert(FluidSystem::isIdealMixture(phaseIdx));
103
104 }
105#endif
106
107 // set the temperature, pressure
108 fluidState.setTemperature(temperature);
109 fluidState.setPressure(phase0Idx, phasePressure[phase0Idx]);
110 fluidState.setPressure(phase1Idx, phasePressure[phase1Idx]);
111
112 // mole equilibrium ratios k for in case first phase is reference phase
113 Scalar k10 = FluidSystem::fugacityCoefficient(fluidState, phase0Idx, comp0Idx) * fluidState.pressure(phase0Idx)
114 / (FluidSystem::fugacityCoefficient(fluidState, phase1Idx, comp0Idx) * fluidState.pressure(phase1Idx));
115 Scalar k11 = FluidSystem::fugacityCoefficient(fluidState, phase0Idx, comp1Idx) * fluidState.pressure(phase0Idx)
116 / (FluidSystem::fugacityCoefficient(fluidState, phase1Idx, comp1Idx) * fluidState.pressure(phase1Idx));
117
118 // get mole fraction from equilibrium constants
119 fluidState.setMoleFraction(phase0Idx, comp0Idx, ((1. - k11) / (k10 - k11)));
120 fluidState.setMoleFraction(phase1Idx, comp0Idx, (fluidState.moleFraction(phase0Idx,comp0Idx) * k10));
121 fluidState.setMoleFraction(phase0Idx, comp1Idx, 1.0 - fluidState.moleFraction(phase0Idx,comp0Idx));
122 fluidState.setMoleFraction(phase1Idx, comp1Idx, 1.0 - fluidState.moleFraction(phase1Idx,comp0Idx));
123
124 // mass equilibrium ratios K for in case first phase is reference phase
125 Scalar K10 = fluidState.massFraction(phase1Idx, comp0Idx) / fluidState.massFraction(phase0Idx, comp0Idx);
126 Scalar K11 = (1. - fluidState.massFraction(phase1Idx, comp0Idx)) / (1. - fluidState.massFraction(phase0Idx, comp0Idx));
127
128 // phase mass fraction Nu (ratio of phase mass to total phase mass) of first phase
129 Scalar Nu0 = 1. + ((Z0 * (K10 - 1.)) + ((1. - Z0) * (K11 - 1.))) / ((K11 - 1.) * (K10 - 1.));
130
131 // an array of the phase mass fractions from which we will compute the saturations
132 std::array<Scalar, 2> phaseMassFraction;
133
134 // check phase presence
135 if (Nu0 > 0. && Nu0 < 1.) // two phases present
136 phaseMassFraction[phase0Idx] = Nu0;
137 else if (Nu0 <= 0.) // only second phase present
138 {
139 phaseMassFraction[phase0Idx] = 0.0; // no first phase
140 fluidState.setMassFraction(phase1Idx,comp0Idx, Z0); // assign complete mass dissolved into second phase
141 }
142 else // only first phase present
143 {
144 phaseMassFraction[phase0Idx] = 1.0; // no second phase
145 fluidState.setMassFraction(phase0Idx, comp0Idx, Z0); // assign complete mass dissolved into first phase
146 }
147
148 // complete phase mass fractions
149 phaseMassFraction[phase1Idx] = 1.0 - phaseMassFraction[phase0Idx];
150
151 // get densities with correct composition
152 fluidState.setDensity(phase0Idx, FluidSystem::density(fluidState, phase0Idx));
153 fluidState.setDensity(phase1Idx, FluidSystem::density(fluidState, phase1Idx));
154 fluidState.setMolarDensity(phase0Idx, FluidSystem::molarDensity(fluidState, phase0Idx));
155 fluidState.setMolarDensity(phase1Idx, FluidSystem::molarDensity(fluidState, phase1Idx));
156
157 fluidState.setViscosity(phase0Idx, FluidSystem::viscosity(fluidState, phase0Idx));
158 fluidState.setViscosity(phase1Idx, FluidSystem::viscosity(fluidState, phase1Idx));
159
160 Scalar sw = phaseMassFraction[phase0Idx] / fluidState.density(phase0Idx);
161 sw /= (phaseMassFraction[phase0Idx] / fluidState.density(phase0Idx)
162 + phaseMassFraction[phase1Idx] / fluidState.density(phase1Idx));
163 fluidState.setSaturation(phase0Idx, sw);
164 fluidState.setSaturation(phase1Idx, 1.0-sw);
165 }
166
181 static void concentrationFlash1p2c(FluidState1p2c& fluidState, const Scalar& Z0,const Dune::FieldVector<Scalar,numPhases>
182 phasePressure,const int presentPhaseIdx, const Scalar& temperature)
183 {
184 // set the temperature, pressure
185 fluidState.setTemperature(temperature);
186 fluidState.setPressure(phase0Idx, phasePressure[phase0Idx]);
187 fluidState.setPressure(phase1Idx, phasePressure[phase1Idx]);
188
189 fluidState.setPresentPhaseIdx(presentPhaseIdx);
190 fluidState.setMassFraction(presentPhaseIdx,comp0Idx, Z0);
191
192 // transform mass to mole fractions
193 fluidState.setMoleFraction(presentPhaseIdx, comp0Idx, Z0 / FluidSystem::molarMass(comp0Idx)
194 / (Z0 / FluidSystem::molarMass(comp0Idx) + (1. - Z0) / FluidSystem::molarMass(comp1Idx)));
195
196 fluidState.setAverageMolarMass(presentPhaseIdx,
197 fluidState.massFraction(presentPhaseIdx, comp0Idx) * FluidSystem::molarMass(comp0Idx)
198 + fluidState.massFraction(presentPhaseIdx, comp1Idx) * FluidSystem::molarMass(comp1Idx));
199
200 fluidState.setDensity(presentPhaseIdx, FluidSystem::density(fluidState, presentPhaseIdx));
201 fluidState.setMolarDensity(presentPhaseIdx, FluidSystem::molarDensity(fluidState, presentPhaseIdx));
202
203 fluidState.setViscosity(presentPhaseIdx, FluidSystem::viscosity(fluidState, presentPhaseIdx));
204 }
206
223 template<class FluidState>
224 static void saturationFlash2p2c(FluidState &fluidState,
225 const Scalar &saturation,
226 const PhaseVector &phasePressure,
227 const Scalar &porosity,
228 const Scalar &temperature)
229 {
230#ifndef NDEBUG
231 // this solver can only handle fluid systems which
232 // assume ideal mixtures of all fluids.
233 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
234 assert(FluidSystem::isIdealMixture(phaseIdx));
235
236 }
237#endif
238
239 // set the temperature, pressure
240 fluidState.setTemperature(temperature);
241 fluidState.setPressure(phase0Idx, phasePressure[phase0Idx]);
242 fluidState.setPressure(phase1Idx, phasePressure[phase1Idx]);
243
244 // mole equilibrium ratios k for in case first phase is reference phase
245 Scalar k10 = FluidSystem::fugacityCoefficient(fluidState, phase0Idx, comp0Idx) * fluidState.pressure(phase0Idx)
246 / (FluidSystem::fugacityCoefficient(fluidState, phase1Idx, comp0Idx) * fluidState.pressure(phase1Idx));
247 Scalar k11 = FluidSystem::fugacityCoefficient(fluidState, phase0Idx, comp1Idx) * fluidState.pressure(phase0Idx)
248 / (FluidSystem::fugacityCoefficient(fluidState, phase1Idx, comp1Idx) * fluidState.pressure(phase1Idx));
249
250 // get mole fraction from equilibrium constants
251 fluidState.setMoleFraction(phase0Idx,comp0Idx, ((1. - k11) / (k10 - k11)));
252 fluidState.setMoleFraction(phase1Idx,comp0Idx, (fluidState.moleFraction(phase0Idx,comp0Idx) * k10));
253 fluidState.setMoleFraction(phase0Idx, comp1Idx, 1.0 - fluidState.moleFraction(phase0Idx,comp0Idx));
254 fluidState.setMoleFraction(phase1Idx, comp1Idx, 1.0 - fluidState.moleFraction(phase1Idx,comp0Idx));
255
256 // get densities with correct composition
257 fluidState.setDensity(phase0Idx, FluidSystem::density(fluidState, phase0Idx));
258 fluidState.setDensity(phase1Idx, FluidSystem::density(fluidState, phase1Idx));
259 fluidState.setMolarDensity(phase0Idx, FluidSystem::molarDensity(fluidState, phase0Idx));
260 fluidState.setMolarDensity(phase1Idx, FluidSystem::molarDensity(fluidState, phase1Idx));
261
262 fluidState.setViscosity(phase0Idx, FluidSystem::viscosity(fluidState, phase0Idx));
263 fluidState.setViscosity(phase1Idx, FluidSystem::viscosity(fluidState, phase1Idx));
264
265 // set saturation
266 fluidState.setSaturation(phase0Idx, saturation);
267 fluidState.setSaturation(phase1Idx, 1.0-saturation);
268 }
270};
271
272} // end namespace Dumux
273
274#endif
Calculates phase state for a single phase but two-component state.
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
std::string temperature() noexcept
I/O name of temperature for equilibrium models.
Definition: name.hh:51
std::string saturation(int phaseIdx) noexcept
I/O name of saturation for multiphase systems.
Definition: name.hh:43
std::string viscosity(int phaseIdx) noexcept
I/O name of viscosity for multiphase systems.
Definition: name.hh:74
std::string molarDensity(int phaseIdx) noexcept
I/O name of molar density for multiphase systems.
Definition: name.hh:83
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:65
std::string porosity() noexcept
I/O name of porosity.
Definition: name.hh:139
Flash calculation routines for compositional sequential models.
Definition: compositionalflash.hh:51
static void saturationFlash2p2c(FluidState &fluidState, const Scalar &saturation, const PhaseVector &phasePressure, const Scalar &porosity, const Scalar &temperature)
2p2c flash for constant p & T if the saturation instead of concentration (feed mass fraction) is know...
Definition: compositionalflash.hh:224
Dune::FieldVector< Scalar, numComponents > ComponentVector
Definition: compositionalflash.hh:66
static void concentrationFlash1p2c(FluidState1p2c &fluidState, const Scalar &Z0, const Dune::FieldVector< Scalar, numPhases > phasePressure, const int presentPhaseIdx, const Scalar &temperature)
The simplest possible update routine for 1p2c "flash" calculations.
Definition: compositionalflash.hh:181
Dune::FieldVector< Scalar, numPhases > PhaseVector
Definition: compositionalflash.hh:67
static void concentrationFlash2p2c(FluidState &fluidState, const Scalar &Z0, const PhaseVector &phasePressure, const Scalar &porosity, const Scalar &temperature)
2p2c Flash for constant p & T if concentration (feed mass fraction) is given.
Definition: compositionalflash.hh:92
Container for compositional variables in a 1p2c situation.
Definition: pseudo1p2c.hh:44
void setAverageMolarMass(int phaseIdx, Scalar value)
Set the average molar mass of a fluid phase [kg/mol].
Definition: pseudo1p2c.hh:276
void setMassFraction(int phaseIdx, int compIdx, Scalar value)
Sets the mass fraction of a component in a phase.
Definition: pseudo1p2c.hh:216
void setPresentPhaseIdx(int phaseIdx)
Sets the phase Index that is present in this fluidState.
Definition: pseudo1p2c.hh:257
void setPressure(int phaseIdx, Scalar value)
Sets the phase pressure .
Definition: pseudo1p2c.hh:282
void setDensity(int phaseIdx, Scalar value)
Sets the density of a phase .
Definition: pseudo1p2c.hh:235
void setViscosity(int phaseIdx, Scalar value)
Sets the viscosity of a phase .
Definition: pseudo1p2c.hh:203
void setMolarDensity(int phaseIdx, Scalar value)
Set the molar density of a phase .
Definition: pseudo1p2c.hh:247
void setTemperature(Scalar value)
Sets the temperature.
Definition: pseudo1p2c.hh:265
Scalar massFraction(int phaseIdx, int compIdx) const
Returns the mass fraction of component in fluid phase in .
Definition: pseudo1p2c.hh:121
void setMoleFraction(int phaseIdx, int compIdx, Scalar value)
Sets the molar fraction of a component in a fluid phase.
Definition: pseudo1p2c.hh:226