Go to the documentation of this file.
29 #include "phaseModel.H"
30 #include "phasePair.H"
38 template<
class Thermo,
class OtherThermo>
39 template<
class ThermoType>
43 const word& speciesName,
50 globalThermo.species()
58 template<
class Thermo,
class OtherThermo>
59 template<
class ThermoType>
63 const word& speciesName,
71 template<
class Thermo,
class OtherThermo>
72 template<
class ThermoType>
76 const word& speciesName,
87 mesh.time().timeName(),
94 zeroGradientFvPatchScalarField::typeName
105 template<
class Thermo,
class OtherThermo>
106 template<
class ThermoType>
110 const word& speciesName,
121 mesh.time().timeName(),
128 zeroGradientFvPatchScalarField::typeName
133 template<
class Thermo,
class OtherThermo>
134 template<
class ThermoType>
148 mesh.time().timeName(),
160 zeroGradientFvPatchScalarField::typeName
165 template<
class Thermo,
class OtherThermo>
166 template<
class ThermoType>
173 return refCast<const basicSpecieMixture>(
mixture).W();
179 template<
class Thermo,
class OtherThermo>
189 pair.
from().mesh().lookupObject<Thermo>
200 pair.
to().mesh().lookupObject<OtherThermo>
215 template<
class Thermo,
class OtherThermo>
219 const word& speciesName
222 const typename Thermo::thermoType& fromThermo =
237 IOobject::groupName(
"D", pair_.name()),
245 auto&
D = tmpD.ref();
250 fromThermo.alphah(
p[cellI],
T[cellI])
251 /fromThermo.rho(
p[cellI],
T[cellI]);
255 D.correctBoundaryConditions();
261 template<
class Thermo,
class OtherThermo>
265 const word& speciesName,
269 const typename Thermo::thermoType& fromThermo =
270 getLocalThermo(speciesName, fromThermo_);
271 const typename OtherThermo::thermoType& toThermo =
272 getLocalThermo(speciesName, toThermo_);
280 IOobject::groupName(
"L", pair_.name()),
286 zeroGradientFvPatchScalarField::typeName
289 auto&
L = tmpL.ref();
294 L[cellI] = fromThermo.Hc() - toThermo.Hc();
297 L.correctBoundaryConditions();
303 template<
class Thermo,
class OtherThermo>
307 const word& speciesName,
316 template<
class Thermo,
class OtherThermo>
320 const word& speciesName,
Defines the attributes of an object for which implicit objectRegistry management is supported,...
const vector L(dict.get< vector >("L"))
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
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.
static constexpr const zero Zero
Global zero (0)
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")
virtual const phaseModel & to() const
To phase.
const dimensionSet dimMoles(0, 0, 0, 0, 1, 0, 0)
const ThermoType & cellMixture(const label) const
ThermoType thermoType
The type of thermodynamics this mixture is instantiated for.
Generic base class for interface models. Mass transfer models are interface models between two thermo...
#define forAll(list, i)
Loop across all elements in list.
virtual tmp< volScalarField > Yf(const word &speciesName, const volScalarField &Tf) const
Reference mass fraction for species based models.
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
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))
Foam::multiComponentMixture.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
tmp< volScalarField > MwMixture(const pureMixture< ThermoType > &thermo) const
Return moleculas weight of the mixture for pureMixture [Kg/mol].
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Mesh data needed to do the Finite Volume discretisation.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
tmp< volScalarField > getSpecieMassFraction(const word &speciesName, const pureMixture< ThermoType > &thermo) const
Return mass fraction for a pureMixture equal to one.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
const word & name() const
const ThermoType & getLocalThermo(const label speciei) const
Return thermo based on index.
virtual const phaseModel & from() const
From phase.
const dimensionedScalar e
Elementary charge.
const dimensionedScalar & D
ThermoType thermoType
The type of thermodynamics this mixture is instantiated for.
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.