25#ifndef DUMUX_COMPOSITIONAL_FLASH_HH
26#define DUMUX_COMPOSITIONAL_FLASH_HH
28#include <dune/common/fvector.hh>
49template <
class Scalar,
class Flu
idSystem>
54 enum { numPhases = FluidSystem::numPhases,
55 numComponents = FluidSystem::numComponents
59 phase0Idx = FluidSystem::phase0Idx,
60 phase1Idx = FluidSystem::phase1Idx,
61 comp0Idx = FluidSystem::comp0Idx,
62 comp1Idx = FluidSystem::comp1Idx
91 template<
class Flu
idState>
101 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
102 assert(FluidSystem::isIdealMixture(phaseIdx));
109 fluidState.setPressure(phase0Idx, phasePressure[phase0Idx]);
110 fluidState.setPressure(phase1Idx, phasePressure[phase1Idx]);
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));
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));
125 Scalar K10 = fluidState.massFraction(phase1Idx, comp0Idx) / fluidState.massFraction(phase0Idx, comp0Idx);
126 Scalar K11 = (1. - fluidState.massFraction(phase1Idx, comp0Idx)) / (1. - fluidState.massFraction(phase0Idx, comp0Idx));
129 Scalar Nu0 = 1. + ((Z0 * (K10 - 1.)) + ((1. - Z0) * (K11 - 1.))) / ((K11 - 1.) * (K10 - 1.));
132 std::array<Scalar, 2> phaseMassFraction;
135 if (Nu0 > 0. && Nu0 < 1.)
136 phaseMassFraction[phase0Idx] = Nu0;
139 phaseMassFraction[phase0Idx] = 0.0;
140 fluidState.setMassFraction(phase1Idx,comp0Idx, Z0);
144 phaseMassFraction[phase0Idx] = 1.0;
145 fluidState.setMassFraction(phase0Idx, comp0Idx, Z0);
149 phaseMassFraction[phase1Idx] = 1.0 - phaseMassFraction[phase0Idx];
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);
182 phasePressure,
const int presentPhaseIdx,
const Scalar&
temperature)
186 fluidState.
setPressure(phase0Idx, phasePressure[phase0Idx]);
187 fluidState.
setPressure(phase1Idx, phasePressure[phase1Idx]);
193 fluidState.
setMoleFraction(presentPhaseIdx, comp0Idx, Z0 / FluidSystem::molarMass(comp0Idx)
194 / (Z0 / FluidSystem::molarMass(comp0Idx) + (1. - Z0) / FluidSystem::molarMass(comp1Idx)));
197 fluidState.
massFraction(presentPhaseIdx, comp0Idx) * FluidSystem::molarMass(comp0Idx)
198 + fluidState.
massFraction(presentPhaseIdx, comp1Idx) * FluidSystem::molarMass(comp1Idx));
223 template<
class Flu
idState>
233 for (
int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
234 assert(FluidSystem::isIdealMixture(phaseIdx));
241 fluidState.setPressure(phase0Idx, phasePressure[phase0Idx]);
242 fluidState.setPressure(phase1Idx, phasePressure[phase1Idx]);
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));
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));
266 fluidState.setSaturation(phase0Idx,
saturation);
267 fluidState.setSaturation(phase1Idx, 1.0-
saturation);
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