31 template<
class MixtureType>
36 const dimensionSet& psiDim,
37 scalar (MixtureType::thermoType::*psiMethod)
47 const typename MixtureType::thermoType& thermo =
48 this->specieThermo(speciei);
50 tmp<volScalarField> tPsi
54 IOobject::groupName(psiName, T.group()),
64 psi[celli] = (thermo.*psiMethod)(p[celli], T[celli]);
67 volScalarField::Boundary& psiBf = psi.boundaryFieldRef();
77 ppsi[facei] = (thermo.*psiMethod)(pp[facei], pT[facei]);
85 template<
class MixtureType>
88 scalar (MixtureType::thermoType::*psiMethod)
98 const typename MixtureType::thermoType& thermo =
99 this->specieThermo(speciei);
107 psi[facei] = (thermo.*psiMethod)(p[facei], T[facei]);
116 template<
class MixtureType>
121 const word& phaseName
124 MixtureType(thermoDict, mesh, phaseName)
130 template<
class MixtureType>
133 return this->specieThermo(speciei).W();
137 template<
class MixtureType>
140 return this->specieThermo(speciei).Hf();
144 template<
class MixtureType>
152 return this->specieThermo(speciei).rho(p, T);
156 template<
class MixtureType>
164 return volScalarFieldProperty
176 template<
class MixtureType>
184 return this->specieThermo(speciei).Cp(p, T);
188 template<
class MixtureType>
196 return volScalarFieldProperty
208 template<
class MixtureType>
216 return this->specieThermo(speciei).HE(p, T);
220 template<
class MixtureType>
228 return fieldProperty(&MixtureType::thermoType::HE, speciei, p, T);
232 template<
class MixtureType>
240 return volScalarFieldProperty
244 &MixtureType::thermoType::HE,
252 template<
class MixtureType>
260 return this->specieThermo(speciei).Hs(p, T);
264 template<
class MixtureType>
276 template<
class MixtureType>
284 return volScalarFieldProperty
296 template<
class MixtureType>
304 return this->specieThermo(speciei).Ha(p, T);
308 template<
class MixtureType>
320 template<
class MixtureType>
328 return volScalarFieldProperty
340 template<
class MixtureType>
348 return this->specieThermo(speciei).mu(p, T);
352 template<
class MixtureType>
360 return volScalarFieldProperty
372 template<
class MixtureType>
380 return this->specieThermo(speciei).kappa(p, T);
384 template<
class MixtureType>
392 return volScalarFieldProperty
#define forAll(list, i)
Loop across all elements in list.
virtual scalar Hf(const label speciei) const
Enthalpy of formation [J/kg].
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
const dimensionSet dimPower
scalar Hs(const scalar p, const scalar T) const
virtual scalar HE(const label speciei, const scalar p, const scalar T) const
Enthalpy/Internal energy [J/kg].
A list of keyword definitions, which are a keyword followed by any number of values (e...
virtual scalar Hs(const label speciei, const scalar p, const scalar T) const
Sensible enthalpy [J/kg].
virtual scalar Wi(const label speciei) const
Molecular weight of the given specie [kg/kmol].
const dimensionSet dimLength
GeometricField< scalar, fvPatchField, volMesh > volScalarField
SpecieMixture(const dictionary &, const fvMesh &, const word &phaseName)
Construct from dictionary, mesh and phase name.
const dimensionSet dimTime
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
fvPatchField< scalar > fvPatchScalarField
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
A class for handling words, derived from string.
const dimensionSet dimDensity
virtual scalar Ha(const label speciei, const scalar p, const scalar T) const
Absolute enthalpy [J/kg].
volScalarField scalarField(fieldObject, mesh)
const dimensionedScalar mu
Atomic mass unit.
const dimensionSet dimEnergy
const dimensionSet dimMass
Internal & ref()
Return a reference to the dimensioned internal field.
virtual scalar rho(const label speciei, const scalar p, const scalar T) const
Density [kg/m^3].
virtual scalar mu(const label speciei, const scalar p, const scalar T) const
Dynamic viscosity [kg/m/s].
virtual tmp< volScalarField > kappa() const =0
Thermal diffusivity for temperature of mixture [W/m/K].
Mesh data needed to do the Finite Volume discretisation.
scalar Cp(const scalar p, const scalar T) const
scalar Ha(const scalar p, const scalar T) const
A class for managing temporary objects.
virtual tmp< volScalarField > Cp() const =0
Heat capacity at constant pressure for patch [J/kg/K].
const dimensionSet dimTemperature