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().find(speciesName)
54 template<
class Thermo,
class OtherThermo>
55 template<
class ThermoType>
59 const word& speciesName,
60 const pureMixture<ThermoType>& globalThermo
63 return globalThermo.cellMixture(0);
69 template<
class Thermo,
class OtherThermo>
72 const dictionary&
dict,
76 interfaceCompositionModel(
dict, pair),
97 template<
class Thermo,
class OtherThermo>
101 const word& speciesName,
107 - thermo_.composition().Y()
109 thermo_.composition().species().find(speciesName)
114 template<
class Thermo,
class OtherThermo>
118 const word& speciesName
121 const typename Thermo::thermoType& localThermo =
147 localThermo.alphah(
p[celli],
T[celli])
148 /localThermo.rho(
p[celli],
T[celli]);
157 template<
class Thermo,
class OtherThermo>
161 const word& speciesName,
165 const typename Thermo::thermoType& localThermo =
171 const typename OtherThermo::thermoType& otherLocalThermo =
181 tmp<volScalarField> tmpL
196 localThermo.Ha(
p[celli], Tf[celli])
197 - otherLocalThermo.Ha(otherP[celli], Tf[celli]);
204 template<
class Thermo,
class OtherThermo>
213 for (
const word& speciesName : this->speciesNames_)
217 thermo_.rhoThermo::rho()
223 mDotL += rhoKDL*dY(speciesName, Tf);
224 mDotLPrime += rhoKDL*YfPrime(speciesName, Tf);
virtual tmp< volScalarField > D(const word &speciesName) const
Mass diffusivity.
const vector L(dict.get< vector >("L"))
const word dictName("faMeshDefinition")
const dimensionSet dimless
Dimensionless.
CGAL::Exact_predicates_exact_constructions_kernel K
#define forAll(list, i)
Loop across all elements in list.
virtual void addMDotL(const volScalarField &K, const volScalarField &Tf, volScalarField &mDotL, volScalarField &mDotLPrime) const
Add latent heat flow rate to total.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const pureMixture< ThermoType >::thermoType & getLocalThermo(const word &speciesName, const pureMixture< ThermoType > &globalThermo) const
Get a reference to the local thermo for a pure mixture.
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
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)
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
const Mesh & mesh() const noexcept
Return mesh.
virtual tmp< volScalarField > dY(const word &speciesName, const volScalarField &Tf) const
Mass fraction difference between the interface and the field.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const dimensionSet dimEnergy
ThermoType thermoType
The type of thermodynamics this mixture is instantiated for.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
ThermoType thermoType
The type of thermodynamics this mixture is instantiated for.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
const dimensionedScalar & D
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
A class for managing temporary objects.
InterfaceCompositionModel(const dictionary &dict, const phasePair &pair)
Construct from components.
const dimensionSet dimArea(sqr(dimLength))