26#ifndef DUMUX_COMPOSITIONAL_FLUID_STATE_HH
27#define DUMUX_COMPOSITIONAL_FLUID_STATE_HH
35#include <dune/common/exceptions.hh>
45template <
class ScalarType,
class Flu
idSystem>
49 static constexpr int numPhases = FluidSystem::numPhases;
59 template <class FluidState, typename std::enable_if_t<!std::is_same<FluidState, CompositionalFluidState>::value,
int> = 0>
123 * FluidSystem::molarMass(compIdx)
134 for (
int pIdx = 0; pIdx <
numPhases; ++pIdx)
229 assert(FluidSystem::isGas(phaseIdx));
290 template <
class Flu
idState>
293 for (
int phaseIdx = 0; phaseIdx <
numPhases; ++phaseIdx)
300 moleFraction_[phaseIdx][compIdx] = fs.moleFraction(phaseIdx, compIdx);
305 pressure_[phaseIdx] = fs.pressure(phaseIdx);
307 density_[phaseIdx] = fs.density(phaseIdx);
309 enthalpy_[phaseIdx] = fs.enthalpy(phaseIdx);
310 viscosity_[phaseIdx] = fs.viscosity(phaseIdx);
320 for (
int phaseIdx = 0; phaseIdx <
numPhases; ++phaseIdx)
370 DUNE_THROW(Dune::NotImplemented,
"This currently only works for 2 components.");
374 Scalar M1 = FluidSystem::molarMass(compIdx);
375 Scalar M2 = FluidSystem::molarMass(1-compIdx);
377 Scalar avgMolarMass = M1*M2/(M2 + X2*(M1 - M2));
379 moleFraction_[phaseIdx][compIdx] = value * avgMolarMass / M1;
385 for (
int compJIdx = 0; compJIdx <
numComponents; ++compJIdx) {
397 template <
class Flu
idState>
401 assert(phaseIdx == FluidSystem::nPhaseIdx);
402 assert(compIdx == FluidSystem::wCompIdx);
404 assert(FluidSystem::isGas(phaseIdx));
406 Scalar moleFraction = value * FluidSystem::vaporPressure(fluidState, FluidSystem::wCompIdx)
407 / fluidState.pressure(phaseIdx);
408 fluidState.setMoleFraction(phaseIdx, FluidSystem::wCompIdx,
moleFraction);
409 fluidState.setMoleFraction(phaseIdx, FluidSystem::nCompIdx, 1.0-
moleFraction);
Represents all relevant thermodynamic quantities of a multi-phase, multi-component fluid system assum...
Definition: compositional.hh:47
CompositionalFluidState & operator=(CompositionalFluidState &&fs)=default
int wPhaseIdx_
Definition: compositional.hh:463
void setMolarDensity(int phaseIdx, Scalar value)
Set the molar density of a phase .
Definition: compositional.hh:428
CompositionalFluidState(CompositionalFluidState &&fs)=default
Scalar pressure(int phaseIdx) const
The pressure of a fluid phase in .
Definition: compositional.hh:220
Scalar moleFraction(int phaseIdx, int compIdx) const
Returns the molar fraction of the component in fluid phase in .
Definition: compositional.hh:99
Scalar viscosity(int phaseIdx) const
The dynamic viscosity of fluid phase in .
Definition: compositional.hh:252
Scalar temperature(int phaseIdx) const
The absolute temperature of a fluid phase in .
Definition: compositional.hh:214
CompositionalFluidState & operator=(const CompositionalFluidState &fs)=default
Scalar internalEnergy(int phaseIdx) const
The specific internal energy of a fluid phase in .
Definition: compositional.hh:246
Scalar molarDensity(int phaseIdx) const
The molar density of the fluid phase in .
Definition: compositional.hh:208
void setRelativeHumidity(FluidState &fluidState, int phaseIdx, int compIdx, Scalar value)
Set the relative humidity of a component in a phase and update the average molar mass according to ...
Definition: compositional.hh:398
void setSaturation(int phaseIdx, Scalar value)
Set the saturation of a phase .
Definition: compositional.hh:340
Scalar saturation(int phaseIdx) const
Returns the saturation of a fluid phase in .
Definition: compositional.hh:87
int wettingPhase() const
Returns the index of the most wetting phase in the fluid-solid configuration (for porous medium syste...
Definition: compositional.hh:77
void setFugacityCoefficient(int phaseIdx, int compIdx, Scalar value)
Set the fugacity coefficient of component in fluid phase in .
Definition: compositional.hh:416
void assign(const FluidState &fs)
Retrieve all parameters from an arbitrary fluid state.
Definition: compositional.hh:291
Scalar enthalpy(int phaseIdx) const
The specific enthalpy of a fluid phase in .
Definition: compositional.hh:236
Scalar density(int phaseIdx) const
The mass density of the fluid phase in .
Definition: compositional.hh:201
void setWettingPhase(int phaseIdx)
Set the index of the most wetting phase.
Definition: compositional.hh:446
Scalar partialPressure(int phaseIdx, int compIdx) const
The partial pressure of a component in a phase .
Definition: compositional.hh:227
std::array< Scalar, numPhases > pressure_
Definition: compositional.hh:455
void setMoleFraction(int phaseIdx, int compIdx, Scalar value)
Set the mole fraction of a component in a phase and update the average molar mass according to the ...
Definition: compositional.hh:348
void setTemperature(Scalar value)
Set the temperature of all phases.
Definition: compositional.hh:318
std::array< std::array< Scalar, numComponents >, numPhases > fugacityCoefficient_
Definition: compositional.hh:452
Scalar molarity(int phaseIdx, int compIdx) const
The molar concentration of component in fluid phase in .
Definition: compositional.hh:160
Scalar phaseMassFraction(int phaseIdx) const
Returns the phase mass fraction, i.e. phase mass per total mass .
Definition: compositional.hh:131
static constexpr int numPhases
Definition: compositional.hh:49
Scalar averageMolarMass(int phaseIdx) const
The average molar mass of phase in .
Definition: compositional.hh:148
CompositionalFluidState()=default
default constructor
std::array< Scalar, numPhases > averageMolarMass_
Definition: compositional.hh:453
Scalar molarVolume(int phaseIdx) const
The molar volume of a fluid phase in .
Definition: compositional.hh:194
CompositionalFluidState(const FluidState &fs)
copy constructor from arbitrary fluid state
Definition: compositional.hh:60
void setMassFraction(int phaseIdx, int compIdx, Scalar value)
Set the mass fraction of a component in a phase and update the average molar mass according to the ...
Definition: compositional.hh:367
std::array< Scalar, numPhases > saturation_
Definition: compositional.hh:456
Scalar temperature() const
The temperature within the domain .
Definition: compositional.hh:264
std::array< Scalar, numPhases > sumMoleFractions_
Definition: compositional.hh:454
std::array< Scalar, numPhases > density_
Definition: compositional.hh:457
std::array< Scalar, numPhases > enthalpy_
Definition: compositional.hh:459
static constexpr int numComponents
Definition: compositional.hh:50
Scalar fugacity(int compIdx) const
The fugacity of a component .
Definition: compositional.hh:272
std::array< Scalar, numPhases > temperature_
Definition: compositional.hh:461
Scalar massFraction(int phaseIdx, int compIdx) const
Returns the mass fraction of component in fluid phase in .
Definition: compositional.hh:117
std::array< Scalar, numPhases > molarDensity_
Definition: compositional.hh:458
Scalar fugacityCoefficient(int phaseIdx, int compIdx) const
The fugacity coefficient of component in fluid phase in .
Definition: compositional.hh:186
CompositionalFluidState(const CompositionalFluidState &fs)=default
std::array< Scalar, numPhases > viscosity_
Definition: compositional.hh:460
ScalarType Scalar
export the scalar type
Definition: compositional.hh:53
void setTemperature(const int phaseIdx, const Scalar value)
Set the temperature of a specific phase. This is not implemented in this fluidstate.
Definition: compositional.hh:328
void setDensity(int phaseIdx, Scalar value)
Set the density of a phase .
Definition: compositional.hh:422
void setViscosity(int phaseIdx, Scalar value)
Set the dynamic viscosity of a phase .
Definition: compositional.hh:440
Scalar fugacity(int phaseIdx, int compIdx) const
The fugacity of component in fluid phase in .
Definition: compositional.hh:180
void setEnthalpy(int phaseIdx, Scalar value)
Set the specific enthalpy of a phase .
Definition: compositional.hh:434
void setPressure(int phaseIdx, Scalar value)
Set the fluid pressure of a phase .
Definition: compositional.hh:334
std::array< std::array< Scalar, numComponents >, numPhases > moleFraction_
zero-initialize all data members with braces syntax
Definition: compositional.hh:451