Go to the documentation of this file.
38 #ifndef InterfaceCompositionModel_H
39 #define InterfaceCompositionModel_H
41 #include "interfaceCompositionModel.H"
57 template<
class Thermo,
class OtherThermo>
79 template<
class ThermoType>
83 const word& speciesName,
88 template<
class ThermoType>
92 const word& speciesName,
118 const word& speciesName,
125 const word& speciesName
131 const word& speciesName,
153 #define makeInterfaceCompositionType(Type, Thermo, Comp, Mix, Phys, OtherThermo, OtherComp, OtherMix, OtherPhys)\
155 typedef Thermo<Comp, SpecieMixture<Mix<Phys>>> \
156 Type##Thermo##Comp##Mix##Phys; \
158 typedef OtherThermo<OtherComp, OtherMix<OtherPhys>> \
159 Type##Other##OtherThermo##OtherComp##OtherMix##OtherPhys; \
161 addInterfaceCompositionToRunTimeSelectionTable \
164 Type##Thermo##Comp##Mix##Phys, \
165 Type##Other##OtherThermo##OtherComp##OtherMix##OtherPhys \
169 #define makeSpecieInterfaceCompositionType(Type, Thermo, Comp, Mix, Phys, OtherThermo, OtherComp, OtherMix, OtherPhys)\
171 typedef Thermo<Comp, SpecieMixture<Mix<Phys>>> \
172 Type##Thermo##Comp##Mix##Phys; \
174 typedef OtherThermo<OtherComp, SpecieMixture<OtherMix<OtherPhys>>> \
175 Type##Other##OtherThermo##OtherComp##OtherMix##OtherPhys; \
177 addInterfaceCompositionToRunTimeSelectionTable \
180 Type##Thermo##Comp##Mix##Phys, \
181 Type##Other##OtherThermo##OtherComp##OtherMix##OtherPhys \
185 #define addInterfaceCompositionToRunTimeSelectionTable(Type, Thermo, OtherThermo)\
187 typedef Type<Thermo, OtherThermo> \
188 Type##Thermo##OtherThermo; \
190 defineTemplateTypeNameAndDebugWithName \
192 Type##Thermo##OtherThermo, \
194 word(Type##Thermo##OtherThermo::typeName_()) + "<" \
195 + word(Thermo::typeName) + "," \
196 + word(OtherThermo::typeName) + ">" \
201 addToRunTimeSelectionTable \
203 interfaceCompositionModel, \
204 Type##Thermo##OtherThermo, \
211 #include "InterfaceCompositionModel.C"
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.
virtual tmp< volScalarField > L(const word &speciesName, const volScalarField &Tf) const
Latent heat.
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 dimensionedScalar Le_
Lewis number.
~InterfaceCompositionModel()=default
Destructor.
const Thermo & thermo_
Thermo.
ThermoType thermoType
The type of thermodynamics this mixture is instantiated for.
Generic base class for interface composition models. These models describe the composition in phase 1...
const phasePair & pair() const
Return pair.
CGAL::Exact_predicates_exact_constructions_kernel K
InterfaceCompositionModel(const dictionary &dict, const phasePair &pair)
Construct from components.
const OtherThermo & otherThermo_
Other Thermo.
Foam::multiComponentMixture.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
virtual void addMDotL(const volScalarField &K, const volScalarField &Tf, volScalarField &mDotL, volScalarField &mDotLPrime) const
Add latent heat flow rate to total.
ThermoType thermoType
The type of thermodynamics this mixture is instantiated for.
virtual tmp< volScalarField > D(const word &speciesName) const
Mass diffusivity.
virtual tmp< volScalarField > dY(const word &speciesName, const volScalarField &Tf) const
Mass fraction difference between the interface and the field.