32 template<
class MixtureType,
class BasicThermoType>
33 template<
class Mixture,
class Method,
class ... Args>
48 IOobject::groupName(psiName, this->
group()),
55 auto Yslicer = this->Yslicer();
59 auto composition = this->cellComposition(Yslicer, celli);
62 ((this->*mixture)(composition).*psiMethod)(
args[celli] ...);
72 this->patchFaceComposition(Yslicer,
patchi, patchFacei);
74 psiBf[
patchi][patchFacei] =
75 ((this->*mixture)(composition).*psiMethod)
86 template<
class MixtureType,
class BasicThermoType>
87 template<
class Mixture,
class Method,
class ... Args>
102 IOobject::groupName(psiName, this->
group()),
109 auto Yslicer = this->Yslicer();
113 auto composition = this->cellComposition(Yslicer, celli);
116 ((this->*mixture)(composition).*psiMethod)(
args[celli] ...);
123 template<
class MixtureType,
class BasicThermoType>
124 template<
class Mixture,
class Method,
class ... Args>
140 auto Yslicer = this->Yslicer();
144 auto composition = this->cellComposition(Yslicer,
cells[i]);
146 psi[i] = ((this->*mixture)(composition).*psiMethod)(
args[i] ...);
153 template<
class MixtureType,
class BasicThermoType>
154 template<
class Mixture,
class Method,
class ... Args>
170 auto Yslicer = this->Yslicer();
175 this->patchFaceComposition(Yslicer,
patchi, patchFacei);
178 ((this->*mixture)(composition).*psiMethod)(
args[patchFacei] ...);
185 template<
class MixtureType,
class BasicThermoType>
186 template<
class Mixture,
class Method,
class ... Args>
203 IOobject::groupName(psiName, this->
group()),
210 auto Yslicer = this->Yslicer(model, source);
214 auto composition = this->sourceCellComposition(Yslicer, celli);
217 ((this->*mixture)(composition).*psiMethod)(
args[celli] ...);
224 template<
class MixtureType,
class BasicThermoType>
225 template<
class Mixture,
class Method,
class ... Args>
240 auto Yslicer = this->Yslicer(model, source,
cells);
245 this->sourceCellComposition(Yslicer, i);
248 ((this->*mixture)(composition).*psiMethod)(
args[i] ...);
255 template<
class MixtureType,
class BasicThermoType>
267 template<
class MixtureType,
class BasicThermoType>
275 return psi.primitiveField();
279 template<
class MixtureType,
class BasicThermoType>
287 if (isA<gradientEnergyFvPatchScalarField>(hBf[
patchi]))
289 refCast<gradientEnergyFvPatchScalarField>(hBf[
patchi]).gradient()
290 = hBf[
patchi].fvPatchField::snGrad();
292 else if (isA<mixedEnergyFvPatchScalarField>(hBf[
patchi]))
294 refCast<mixedEnergyFvPatchScalarField>(hBf[
patchi]).refGrad()
295 = hBf[
patchi].fvPatchField::snGrad();
303 template<
class MixtureType,
class BasicThermoType>
307 const word& phaseName
311 MixtureType(properties()),
315 static_cast<const MixtureType&>(*this),
324 BasicThermoType::phasePropertyName
326 MixtureType::thermoType::heName(),
334 volScalarFieldProperty
338 &MixtureType::thermoMixture,
339 &MixtureType::thermoMixtureType::
he,
343 this->heBoundaryTypes(),
344 this->heBoundaryBaseTypes(),
345 this->heSourcesTypes()
352 BasicThermoType::phasePropertyName(
"Cp", phaseName),
364 BasicThermoType::phasePropertyName(
"Cv", phaseName),
378 template<
class MixtureType,
class BasicThermoType>
385 template<
class MixtureType,
class BasicThermoType>
389 return volScalarFieldProperty
393 &MixtureType::thermoMixture,
399 template<
class MixtureType,
class BasicThermoType>
406 return patchFieldProperty
408 &MixtureType::thermoMixture,
415 template<
class MixtureType,
class BasicThermoType>
419 if (MixtureType::thermoType::enthalpy())
430 template<
class MixtureType,
class BasicThermoType>
438 return volScalarFieldProperty
442 &MixtureType::thermoMixture,
450 template<
class MixtureType,
class BasicThermoType>
458 return volInternalScalarFieldProperty
462 &MixtureType::thermoMixture,
470 template<
class MixtureType,
class BasicThermoType>
478 return cellSetProperty
480 &MixtureType::thermoMixture,
483 cellSetScalarList(this->p_,
cells),
489 template<
class MixtureType,
class BasicThermoType>
497 return patchFieldProperty
499 &MixtureType::thermoMixture,
502 this->p_.boundaryField()[
patchi],
508 template<
class MixtureType,
class BasicThermoType>
517 return fieldSourceProperty
521 &MixtureType::thermoMixture,
525 this->p_.internalField(),
531 template<
class MixtureType,
class BasicThermoType>
541 return fieldSourceProperty
543 &MixtureType::thermoMixture,
548 cellSetScalarList(this->p_,
cells),
554 template<
class MixtureType,
class BasicThermoType>
558 return volScalarFieldProperty
562 &MixtureType::thermoMixture,
570 template<
class MixtureType,
class BasicThermoType>
578 return volScalarFieldProperty
582 &MixtureType::thermoMixture,
590 template<
class MixtureType,
class BasicThermoType>
598 return volInternalScalarFieldProperty
602 &MixtureType::thermoMixture,
610 template<
class MixtureType,
class BasicThermoType>
618 return cellSetProperty
620 &MixtureType::thermoMixture,
623 cellSetScalarList(this->p_,
cells),
629 template<
class MixtureType,
class BasicThermoType>
637 return patchFieldProperty
639 &MixtureType::thermoMixture,
642 this->p_.boundaryField()[
patchi],
648 template<
class MixtureType,
class BasicThermoType>
652 return volScalarFieldProperty
656 &MixtureType::thermoMixture,
664 template<
class MixtureType,
class BasicThermoType>
672 return volScalarFieldProperty
676 &MixtureType::thermoMixture,
684 template<
class MixtureType,
class BasicThermoType>
692 return volInternalScalarFieldProperty
696 &MixtureType::thermoMixture,
704 template<
class MixtureType,
class BasicThermoType>
712 return cellSetProperty
714 &MixtureType::thermoMixture,
717 cellSetScalarList(this->p_,
cells),
723 template<
class MixtureType,
class BasicThermoType>
731 return patchFieldProperty
733 &MixtureType::thermoMixture,
736 this->p_.boundaryField()[
patchi],
742 template<
class MixtureType,
class BasicThermoType>
750 return patchFieldProperty
752 &MixtureType::thermoMixture,
755 this->p_.boundaryField()[
patchi],
761 template<
class MixtureType,
class BasicThermoType>
769 return patchFieldProperty
771 &MixtureType::thermoMixture,
774 this->p_.boundaryField()[
patchi],
780 template<
class MixtureType,
class BasicThermoType>
788 if (MixtureType::thermoType::enthalpy())
799 template<
class MixtureType,
class BasicThermoType>
808 return volScalarFieldProperty
812 &MixtureType::thermoMixture,
813 &MixtureType::thermoMixtureType::The,
821 template<
class MixtureType,
class BasicThermoType>
830 return cellSetProperty
832 &MixtureType::thermoMixture,
833 &MixtureType::thermoMixtureType::The,
836 cellSetScalarList(this->p_,
cells),
842 template<
class MixtureType,
class BasicThermoType>
851 return patchFieldProperty
853 &MixtureType::thermoMixture,
854 &MixtureType::thermoMixtureType::The,
857 this->p_.boundaryField()[
patchi],
863 template<
class MixtureType,
class BasicThermoType>
scalar hs(const scalar p, const scalar T) const
scalar Cp(const scalar p, const scalar T) const
scalar ha(const scalar p, const scalar T) const
scalar Cv(const scalar p, const scalar T) const
#define forAll(list, i)
Loop across all elements in list.
BasicThermo(const fvMesh &, const word &phaseName)
Construct from mesh and phase name.
virtual const volScalarField & Cv() const
Heat capacity at constant volume [J/kg/K].
virtual tmp< volScalarField > W() const
Molecular weight [kg/kmol].
virtual const volScalarField & Cp() const
Heat capacity at constant pressure [J/kg/K].
tmp< volScalarField > volScalarFieldProperty(const word &psiName, const dimensionSet &psiDim, Mixture mixture, Method psiMethod, const Args &... args) const
Return a volScalarField of the given property.
tmp< scalarField > cellSetProperty(Mixture mixture, Method psiMethod, const labelList &cells, const Args &... args) const
Return a scalarField of the given property on a cell set.
static UIndirectList< scalar > cellSetScalarList(const volScalarField &psi, const labelUList &cells)
Return an indirect list of a field for the given set of cells.
virtual tmp< volScalarField > ha() const
Absolute enthalpy [J/kg/K].
virtual const volScalarField & he() const
Enthalpy/Internal energy [J/kg].
volScalarField he_
Energy field.
tmp< scalarField > patchFieldProperty(Mixture mixture, Method psiMethod, const label patchi, const Args &... args) const
Return a scalarField of the given property on a patch.
tmp< volScalarField::Internal > fieldSourceProperty(const word &psiName, const dimensionSet &psiDim, Mixture mixture, Method psiMethod, const fvSource &model, const volScalarField::Internal &source, const Args &... args) const
Return a scalarField of the given property for a source.
void heBoundaryCorrection(volScalarField &he)
Correct the enthalpy/internal energy field boundaries.
virtual const volScalarField & Cpv() const
Heat capacity at constant pressure/volume [J/kg/K].
virtual ~BasicThermo()
Destructor.
virtual tmp< volScalarField > hs() const
Sensible enthalpy [J/kg/K].
virtual tmp< volScalarField > The(const volScalarField &h, const volScalarField &p, const volScalarField &T0) const
Temperature from enthalpy/internal energy.
tmp< volScalarField::Internal > volInternalScalarFieldProperty(const word &psiName, const dimensionSet &psiDim, Mixture mixture, Method psiMethod, const Args &... args) const
Return a volScalarField::Internal of the given property.
virtual bool read()
Read thermophysical properties dictionary.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricBoundaryField class.
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
A List with indirect addressing.
Dimension set for the base types.
Mesh data needed to do the Finite Volume discretisation.
Base class for finite volume sources.
A base class for physical properties.
virtual bool read()
Read physicalProperties dictionary.
A class for managing temporary objects.
T & ref() const
Return non-const reference or generate a fatal error.
A class for handling words, derived from string.
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
volScalarField scalarField(fieldObject, mesh)
const volScalarField & psi
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
const char *const group
Group name for atomic constants.
const dimensionedScalar h
Planck constant.
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 dimEnergy
void T(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
const dimensionSet dimTemperature
const dimensionSet dimMoles
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
const dimensionSet dimMass
const scalarList W(::W(thermo))
Foam::argList args(argc, argv)