14#ifndef DUMUX_COMPOSITIONAL_FLUID_STATE_HH
15#define DUMUX_COMPOSITIONAL_FLUID_STATE_HH
23#include <dune/common/exceptions.hh>
33template <
class ScalarType,
class Flu
idSystem>
37 static constexpr int numPhases = FluidSystem::numPhases;
47 template <class FluidState, typename std::enable_if_t<!std::is_same<FluidState, CompositionalFluidState>::value,
int> = 0>
111 * FluidSystem::molarMass(compIdx)
122 for (
int pIdx = 0; pIdx <
numPhases; ++pIdx)
217 assert(FluidSystem::isGas(phaseIdx));
278 template <
class Flu
idState>
281 for (
int phaseIdx = 0; phaseIdx <
numPhases; ++phaseIdx)
288 moleFraction_[phaseIdx][compIdx] = fs.moleFraction(phaseIdx, compIdx);
293 pressure_[phaseIdx] = fs.pressure(phaseIdx);
295 density_[phaseIdx] = fs.density(phaseIdx);
297 enthalpy_[phaseIdx] = fs.enthalpy(phaseIdx);
298 viscosity_[phaseIdx] = fs.viscosity(phaseIdx);
308 for (
int phaseIdx = 0; phaseIdx <
numPhases; ++phaseIdx)
358 DUNE_THROW(Dune::NotImplemented,
"This currently only works for 2 components.");
362 Scalar M1 = FluidSystem::molarMass(compIdx);
363 Scalar M2 = FluidSystem::molarMass(1-compIdx);
365 Scalar avgMolarMass = M1*M2/(M2 + X2*(M1 - M2));
367 moleFraction_[phaseIdx][compIdx] = value * avgMolarMass / M1;
373 for (
int compJIdx = 0; compJIdx <
numComponents; ++compJIdx) {
385 template <
class Flu
idState>
389 assert(phaseIdx == FluidSystem::nPhaseIdx);
390 assert(compIdx == FluidSystem::wCompIdx);
392 assert(FluidSystem::isGas(phaseIdx));
394 Scalar moleFraction = value * FluidSystem::vaporPressure(fluidState, FluidSystem::wCompIdx)
395 / fluidState.pressure(phaseIdx);
396 fluidState.setMoleFraction(phaseIdx, FluidSystem::wCompIdx,
moleFraction);
397 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:35
CompositionalFluidState & operator=(CompositionalFluidState &&fs)=default
int wPhaseIdx_
Definition: compositional.hh:451
void setMolarDensity(int phaseIdx, Scalar value)
Set the molar density of a phase .
Definition: compositional.hh:416
CompositionalFluidState(CompositionalFluidState &&fs)=default
Scalar pressure(int phaseIdx) const
The pressure of a fluid phase in .
Definition: compositional.hh:208
Scalar moleFraction(int phaseIdx, int compIdx) const
Returns the molar fraction of the component in fluid phase in .
Definition: compositional.hh:87
Scalar viscosity(int phaseIdx) const
The dynamic viscosity of fluid phase in .
Definition: compositional.hh:240
Scalar temperature(int phaseIdx) const
The absolute temperature of a fluid phase in .
Definition: compositional.hh:202
CompositionalFluidState & operator=(const CompositionalFluidState &fs)=default
Scalar internalEnergy(int phaseIdx) const
The specific internal energy of a fluid phase in .
Definition: compositional.hh:234
Scalar molarDensity(int phaseIdx) const
The molar density of the fluid phase in .
Definition: compositional.hh:196
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:386
void setSaturation(int phaseIdx, Scalar value)
Set the saturation of a phase .
Definition: compositional.hh:328
Scalar saturation(int phaseIdx) const
Returns the saturation of a fluid phase in .
Definition: compositional.hh:75
int wettingPhase() const
Returns the index of the most wetting phase in the fluid-solid configuration (for porous medium syste...
Definition: compositional.hh:65
void setFugacityCoefficient(int phaseIdx, int compIdx, Scalar value)
Set the fugacity coefficient of component in fluid phase in .
Definition: compositional.hh:404
void assign(const FluidState &fs)
Retrieve all parameters from an arbitrary fluid state.
Definition: compositional.hh:279
Scalar enthalpy(int phaseIdx) const
The specific enthalpy of a fluid phase in .
Definition: compositional.hh:224
Scalar density(int phaseIdx) const
The mass density of the fluid phase in .
Definition: compositional.hh:189
void setWettingPhase(int phaseIdx)
Set the index of the most wetting phase.
Definition: compositional.hh:434
Scalar partialPressure(int phaseIdx, int compIdx) const
The partial pressure of a component in a phase .
Definition: compositional.hh:215
std::array< Scalar, numPhases > pressure_
Definition: compositional.hh:443
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:336
void setTemperature(Scalar value)
Set the temperature of all phases.
Definition: compositional.hh:306
std::array< std::array< Scalar, numComponents >, numPhases > fugacityCoefficient_
Definition: compositional.hh:440
Scalar molarity(int phaseIdx, int compIdx) const
The molar concentration of component in fluid phase in .
Definition: compositional.hh:148
Scalar phaseMassFraction(int phaseIdx) const
Returns the phase mass fraction, i.e. phase mass per total mass .
Definition: compositional.hh:119
static constexpr int numPhases
Definition: compositional.hh:37
Scalar averageMolarMass(int phaseIdx) const
The average molar mass of phase in .
Definition: compositional.hh:136
CompositionalFluidState()=default
default constructor
std::array< Scalar, numPhases > averageMolarMass_
Definition: compositional.hh:441
Scalar molarVolume(int phaseIdx) const
The molar volume of a fluid phase in .
Definition: compositional.hh:182
CompositionalFluidState(const FluidState &fs)
copy constructor from arbitrary fluid state
Definition: compositional.hh:48
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:355
std::array< Scalar, numPhases > saturation_
Definition: compositional.hh:444
Scalar temperature() const
The temperature within the domain .
Definition: compositional.hh:252
std::array< Scalar, numPhases > sumMoleFractions_
Definition: compositional.hh:442
std::array< Scalar, numPhases > density_
Definition: compositional.hh:445
std::array< Scalar, numPhases > enthalpy_
Definition: compositional.hh:447
static constexpr int numComponents
Definition: compositional.hh:38
Scalar fugacity(int compIdx) const
The fugacity of a component .
Definition: compositional.hh:260
std::array< Scalar, numPhases > temperature_
Definition: compositional.hh:449
Scalar massFraction(int phaseIdx, int compIdx) const
Returns the mass fraction of component in fluid phase in .
Definition: compositional.hh:105
std::array< Scalar, numPhases > molarDensity_
Definition: compositional.hh:446
Scalar fugacityCoefficient(int phaseIdx, int compIdx) const
The fugacity coefficient of component in fluid phase in .
Definition: compositional.hh:174
CompositionalFluidState(const CompositionalFluidState &fs)=default
std::array< Scalar, numPhases > viscosity_
Definition: compositional.hh:448
ScalarType Scalar
export the scalar type
Definition: compositional.hh:41
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:316
void setDensity(int phaseIdx, Scalar value)
Set the density of a phase .
Definition: compositional.hh:410
void setViscosity(int phaseIdx, Scalar value)
Set the dynamic viscosity of a phase .
Definition: compositional.hh:428
Scalar fugacity(int phaseIdx, int compIdx) const
The fugacity of component in fluid phase in .
Definition: compositional.hh:168
void setEnthalpy(int phaseIdx, Scalar value)
Set the specific enthalpy of a phase .
Definition: compositional.hh:422
void setPressure(int phaseIdx, Scalar value)
Set the fluid pressure of a phase .
Definition: compositional.hh:322
std::array< std::array< Scalar, numComponents >, numPhases > moleFraction_
zero-initialize all data members with braces syntax
Definition: compositional.hh:439