Go to the documentation of this file.
32 template<
class Thermo,
class OtherThermo>
45 IOobject::groupName(
"YSolvent", pair.
name()),
46 pair.
phase1().mesh().time().timeName(),
53 if (k_.size() != this->speciesNames_.size())
56 <<
"Differing number of species and solubilities"
64 template<
class Thermo,
class OtherThermo>
70 YSolvent_ = scalar(1);
72 for (
const word& speciesName : this->speciesNames_)
74 YSolvent_ -= Yf(speciesName, Tf);
79 template<
class Thermo,
class OtherThermo>
83 const word& speciesName,
87 if (this->speciesNames_.found(speciesName))
89 const label index = this->speciesNames_[speciesName];
93 *this->otherThermo_.composition().Y(speciesName)
94 *this->otherThermo_.rhoThermo::rho()
95 /this->thermo_.rhoThermo::rho();
101 *this->thermo_.composition().Y(speciesName);
106 template<
class Thermo,
class OtherThermo>
110 const word& speciesName,
116 IOobject::groupName(
"YfPrime", this->pair_.name()),
117 this->pair_.phase1().mesh(),
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Base class for interface composition models, templated on the two thermodynamic models either side of...
A class for handling words, derived from Foam::string.
A class for managing temporary objects.
virtual void update(const volScalarField &Tf)
Update the composition.
Henry(const dictionary &dict, const phasePair &pair)
Construct from components.
virtual tmp< volScalarField > YfPrime(const word &speciesName, const volScalarField &Tf) const
The interface species fraction derivative w.r.t. temperature.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual tmp< volScalarField > Yf(const word &speciesName, const volScalarField &Tf) const
The interface species fraction.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
errorManipArg< error, int > exit(error &err, const int errNo=1)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
virtual word name() const
Pair name.
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
const phaseModel & phase1() const
const dimensionSet dimless
Dimensionless.