Go to the documentation of this file.
29 #include "phaseModel.H"
30 #include "phasePair.H"
37 template<
class Thermo,
class OtherThermo>
38 template<
class ThermoType>
42 const word& speciesName,
43 const multiComponentMixture<ThermoType>& globalThermo
47 globalThermo.getLocalThermo
49 globalThermo.species()
57 template<
class Thermo,
class OtherThermo>
58 template<
class ThermoType>
62 const word& speciesName,
63 const pureMixture<ThermoType>& globalThermo
66 return globalThermo.cellMixture(0);
72 template<
class Thermo,
class OtherThermo>
75 const dictionary&
dict,
79 interfaceCompositionModel(
dict, pair),
100 template<
class Thermo,
class OtherThermo>
104 const word& speciesName,
110 - thermo_.composition().Y()
112 thermo_.composition().species()[speciesName]
117 template<
class Thermo,
class OtherThermo>
121 const word& speciesName
124 const typename Thermo::thermoType& localThermo =
135 tmp<volScalarField> tmpD
150 localThermo.alphah(
p[celli],
T[celli])
151 /localThermo.rho(
p[celli],
T[celli]);
160 template<
class Thermo,
class OtherThermo>
164 const word& speciesName,
168 const typename Thermo::thermoType& localThermo =
174 const typename OtherThermo::thermoType& otherLocalThermo =
184 tmp<volScalarField> tmpL
199 localThermo.Ha(
p[celli], Tf[celli])
200 - otherLocalThermo.Ha(otherP[celli], Tf[celli]);
207 template<
class Thermo,
class OtherThermo>
216 for (
const word& speciesName : this->speciesNames_)
220 thermo_.rhoThermo::rho()
226 mDotL += rhoKDL*dY(speciesName, Tf);
227 mDotLPrime += rhoKDL*YfPrime(speciesName, Tf);
const vector L(dict.get< vector >("L"))
A class for handling words, derived from Foam::string.
virtual tmp< volScalarField > L(const word &speciesName, const volScalarField &Tf) const
Latent heat (to - from)(thermo - otherThermo)
A class for managing temporary objects.
const dimensionSet dimEnergy
const pureMixture< ThermoType >::thermoType & getLocalThermo(const word &speciesName, const pureMixture< ThermoType > &globalThermo) const
Get a reference to the local thermo for a pure mixture.
const word dictName("faMeshDefinition")
ThermoType thermoType
The type of thermodynamics this mixture is instantiated for.
#define forAll(list, i)
Loop across all elements in list.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &ds, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return tmp field from name, mesh, dimensions and patch type.
CGAL::Exact_predicates_exact_constructions_kernel K
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
InterfaceCompositionModel(const dictionary &dict, const phasePair &pair)
Construct from components.
const dimensionSet dimArea(sqr(dimLength))
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
virtual void addMDotL(const volScalarField &K, const volScalarField &Tf, volScalarField &mDotL, volScalarField &mDotLPrime) const
Add latent heat flow rate to total.
const dimensionedScalar & D
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
ThermoType thermoType
The type of thermodynamics this mixture is instantiated for.
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
virtual tmp< volScalarField > D(const word &speciesName) const
Mass diffusivity of the local thermo.
virtual tmp< volScalarField > dY(const word &speciesName, const volScalarField &Tf) const
Mass fraction difference between the interface and the field.
const dimensionSet dimless
Dimensionless.