33 namespace interfaceCompositionModels
64 saturatedName_(species()[0]),
65 saturatedIndex_(
thermo().species()[saturatedName_]),
71 <<
"saturated model is suitable for one species only."
94 const word& speciesName,
98 if (saturatedName_ == speciesName)
100 return wRatioByP()*saturationModel_->pSat(Tf);
104 const label speciesIndex =
thermo().species()[speciesName];
107 thermo().Y()[speciesIndex]
108 *(scalar(1) - wRatioByP()*saturationModel_->pSat(Tf))
109 /
max(scalar(1) -
thermo().Y()[saturatedIndex_], small);
117 const word& speciesName,
121 if (saturatedName_ == speciesName)
123 return wRatioByP()*saturationModel_->pSatPrime(Tf);
127 const label speciesIndex =
thermo().species()[speciesName];
130 -
thermo().Y()[speciesIndex]
131 *wRatioByP()*saturationModel_->pSatPrime(Tf)
132 /
max(scalar(1) -
thermo().
Y()[saturatedIndex_], small);
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
virtual tmp< volScalarField > W() const =0
Molecular weight [kg/kmol].
A list of keyword definitions, which are a keyword followed by any number of values (e....
virtual const volScalarField & p() const =0
Pressure [Pa].
Generic base class for interface composition models. These models describe the composition in phase 1...
const rhoFluidMulticomponentThermo & thermo() const
Return the thermo.
const hashedWordList & species() const
Return the transferring species names.
Model which uses a saturation pressure model for a single species to calculate the interface composit...
virtual void update(const volScalarField &Tf)
Update the composition.
tmp< volScalarField > wRatioByP() const
Constant of proportionality between partial pressure and mass.
label saturatedIndex_
saturated species index
virtual ~saturated()
Destructor.
virtual tmp< volScalarField > Yf(const word &speciesName, const volScalarField &Tf) const
The interface species fraction.
virtual tmp< volScalarField > YfPrime(const word &speciesName, const volScalarField &Tf) const
The interface species fraction derivative w.r.t. temperature.
saturated(const dictionary &dict, const phaseInterface &interface)
Construct from a dictionary and an interface.
virtual dimensionedScalar Wi(const label speciei) const =0
Molecular weight [kg/kmol].
Class to represent an interface between phases. Derivations can further specify the configuration of ...
Model to describe the dependence of saturation pressure on temperature, and vice versa.
A class for managing temporary objects.
A class for handling words, derived from string.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
addToRunTimeSelectionTable(interfaceCompositionModel, Henry, dictionary)
defineTypeNameAndDebug(Henry, 0)
errorManipArg< error, int > exit(error &err, const int errNo=1)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
PtrList< volScalarField > & Y
fluidMulticomponentThermo & thermo