Flash calculation routines for compositional sequential models. More...
#include <dumux/material/constraintsolvers/compositionalflash.hh>
Flash calculation routines for compositional sequential models.
Routines for isothermal and isobaric 2p2c and 1p2c flash, assuming an ideal mixture. The flash assumes that the fugacities of a component \( \kappa \) in each phase are the same. The fugacity is defined as:
\( f^\kappa = \Phi^\kappa_\alpha(T_\alpha, p_\alpha) p_\alpha x^\kappa_\alpha\; \),
where \( Phi^\kappa_\alpha \) is a fixed fugacity coefficient independent of the phase's composition (ideal mixture), \( T_\alpha \) is a fixed temperature, and \( p_\alpha \) a fixed phase pressure. From the equality of fugacities, the mole (and mass) fractions \( x^\kappa_\alpha \) in equilibrium are calculated.
Public Types | |
using | ComponentVector = Dune::FieldVector< Scalar, numComponents > |
using | PhaseVector = Dune::FieldVector< Scalar, numPhases > |
Static Public Member Functions | |||||||||
Concentration flash for a given feed fraction | |||||||||
2p2c Flash for constant p & T if concentration (feed mass fraction) is given. This flash uses the Rachford-Rice equation: Rachford Jr, H. H., & Rice, J. D. (1952). Procedure for use of electronic digital computers in calculating flash vaporization hydrocarbon equilibrium. Journal of Petroleum Technology, 4(10), 19-3. Routine goes as follows:
| |||||||||
template<class FluidState > | |||||||||
static void | concentrationFlash2p2c (FluidState &fluidState, const Scalar Z0, const PhaseVector &phasePressure, const Scalar temperature) | ||||||||
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. More... | |||||||||
Saturation flash for a given saturation (e.g. at boundary) | |||||||||
2p2c flash for constant p & T if the saturation instead of concentration (feed mass fraction) is known. Routine goes as follows:
| |||||||||
template<class FluidState > | |||||||||
static void | saturationFlash2p2c (FluidState &fluidState, const Scalar saturation, const PhaseVector &phasePressure, const Scalar temperature) | ||||||||
using Dumux::CompositionalFlash< Scalar, FluidSystem >::ComponentVector = Dune::FieldVector<Scalar, numComponents> |
using Dumux::CompositionalFlash< Scalar, FluidSystem >::PhaseVector = Dune::FieldVector<Scalar, numPhases> |
|
inlinestatic |
The simplest possible update routine for 1p2c "flash" calculations.
Routine goes as follows:
fluidState | The sequential fluid state |
Z0 | Feed mass fraction: Mass of first component per total mass \(\mathrm{[-]}\) |
phasePressure | Vector holding the pressure \(\mathrm{[Pa]}\) |
presentPhaseIdx | Subdomain Index = Indication which phase is present |
temperature | Temperature \(\mathrm{[K]}\) |
|
inlinestatic |
|
inlinestatic |