32 template<
class BasicThermo,
class MixtureType>
36 class PatchFaceMixture,
45 CellMixture cellMixture,
46 PatchFaceMixture patchFaceMixture,
55 IOobject::groupName(psiName, this->
group()),
65 psi[celli] = ((this->*cellMixture)(celli).*psiMethod)(args[celli] ...);
77 ((this->*patchFaceMixture)(patchi, facei).*psiMethod)
79 args.boundaryField()[
patchi][facei] ...
88 template<
class BasicThermo,
class MixtureType>
89 template<
class CellMixture,
class Method,
class ... Args>
93 CellMixture cellMixture,
108 ((this->*cellMixture)(cells[celli]).*psiMethod)(args[celli] ...);
115 template<
class BasicThermo,
class MixtureType>
116 template<
class PatchFaceMixture,
class Method,
class ... Args>
120 PatchFaceMixture patchFaceMixture,
135 ((this->*patchFaceMixture)(patchi, facei).*psiMethod)
145 template<
class BasicThermo,
class MixtureType>
157 template<
class BasicThermo,
class MixtureType>
169 template<
class BasicThermo,
class MixtureType>
177 if (isA<gradientEnergyFvPatchScalarField>(hBf[patchi]))
179 refCast<gradientEnergyFvPatchScalarField>(hBf[
patchi]).gradient()
180 = hBf[
patchi].fvPatchField::snGrad();
182 else if (isA<mixedEnergyFvPatchScalarField>(hBf[patchi]))
184 refCast<mixedEnergyFvPatchScalarField>(hBf[
patchi]).refGrad()
185 = hBf[
patchi].fvPatchField::snGrad();
193 template<
class BasicThermo,
class MixtureType>
197 const word& phaseName
200 BasicThermo(mesh, phaseName),
201 MixtureType(*
this, mesh, phaseName),
207 BasicThermo::phasePropertyName
209 MixtureType::thermoType::heName(),
217 volScalarFieldProperty
221 &MixtureType::cellThermoMixture,
222 &MixtureType::patchFaceThermoMixture,
223 &MixtureType::thermoMixtureType::HE,
227 this->heBoundaryTypes(),
228 this->heBoundaryBaseTypes()
235 BasicThermo::phasePropertyName(
"Cp", phaseName),
247 BasicThermo::phasePropertyName(
"Cv", phaseName),
255 heBoundaryCorrection(he_);
261 template<
class BasicThermo,
class MixtureType>
268 template<
class BasicThermo,
class MixtureType>
275 return volScalarFieldProperty
279 &MixtureType::cellThermoMixture,
280 &MixtureType::patchFaceThermoMixture,
281 &MixtureType::thermoMixtureType::HE,
288 template<
class BasicThermo,
class MixtureType>
295 return cellSetProperty
297 &MixtureType::cellThermoMixture,
298 &MixtureType::thermoMixtureType::HE,
300 cellSetScalarList(this->p_, cells),
306 template<
class BasicThermo,
class MixtureType>
313 return patchFieldProperty
315 &MixtureType::patchFaceThermoMixture,
316 &MixtureType::thermoMixtureType::HE,
318 this->p_.boundaryField()[
patchi],
324 template<
class BasicThermo,
class MixtureType>
328 return volScalarFieldProperty
332 &MixtureType::cellThermoMixture,
333 &MixtureType::patchFaceThermoMixture,
341 template<
class BasicThermo,
class MixtureType>
348 return volScalarFieldProperty
352 &MixtureType::cellThermoMixture,
353 &MixtureType::patchFaceThermoMixture,
361 template<
class BasicThermo,
class MixtureType>
368 return cellSetProperty
370 &MixtureType::cellThermoMixture,
373 cellSetScalarList(this->p_, cells),
379 template<
class BasicThermo,
class MixtureType>
386 return patchFieldProperty
388 &MixtureType::patchFaceThermoMixture,
391 this->p_.boundaryField()[
patchi],
397 template<
class BasicThermo,
class MixtureType>
401 return volScalarFieldProperty
405 &MixtureType::cellThermoMixture,
406 &MixtureType::patchFaceThermoMixture,
414 template<
class BasicThermo,
class MixtureType>
421 return volScalarFieldProperty
425 &MixtureType::cellThermoMixture,
426 &MixtureType::patchFaceThermoMixture,
434 template<
class BasicThermo,
class MixtureType>
441 return cellSetProperty
443 &MixtureType::cellThermoMixture,
446 cellSetScalarList(this->p_, cells),
452 template<
class BasicThermo,
class MixtureType>
459 return patchFieldProperty
461 &MixtureType::patchFaceThermoMixture,
464 this->p_.boundaryField()[
patchi],
470 template<
class BasicThermo,
class MixtureType>
474 return volScalarFieldProperty
478 &MixtureType::cellThermoMixture,
479 &MixtureType::patchFaceThermoMixture,
480 &MixtureType::thermoMixtureType::Hf
485 template<
class BasicThermo,
class MixtureType>
492 return patchFieldProperty
494 &MixtureType::patchFaceThermoMixture,
497 this->p_.boundaryField()[
patchi],
503 template<
class BasicThermo,
class MixtureType>
511 return patchFieldProperty
513 &MixtureType::patchFaceThermoMixture,
516 this->p_.boundaryField()[
patchi],
522 template<
class BasicThermo,
class MixtureType>
529 return patchFieldProperty
531 &MixtureType::patchFaceThermoMixture,
532 &MixtureType::thermoMixtureType::gamma,
534 this->p_.boundaryField()[
patchi],
540 template<
class BasicThermo,
class MixtureType>
548 template<
class BasicThermo,
class MixtureType>
555 if (MixtureType::thermoType::enthalpy())
557 return Cp(T, patchi);
561 return Cv(T, patchi);
566 template<
class BasicThermo,
class MixtureType>
570 if (MixtureType::thermoType::enthalpy())
581 template<
class BasicThermo,
class MixtureType>
589 return volScalarFieldProperty
593 &MixtureType::cellThermoMixture,
594 &MixtureType::patchFaceThermoMixture,
595 &MixtureType::thermoMixtureType::THE,
603 template<
class BasicThermo,
class MixtureType>
611 return cellSetProperty
613 &MixtureType::cellThermoMixture,
614 &MixtureType::thermoMixtureType::THE,
617 cellSetScalarList(this->p_, cells),
623 template<
class BasicThermo,
class MixtureType>
631 return patchFieldProperty
633 &MixtureType::patchFaceThermoMixture,
634 &MixtureType::thermoMixtureType::THE,
637 this->p_.boundaryField()[
patchi],
643 template<
class BasicThermo,
class MixtureType>
647 return volScalarFieldProperty
651 &MixtureType::cellThermoMixture,
652 &MixtureType::patchFaceThermoMixture,
658 template<
class BasicThermo,
class MixtureType>
664 return patchFieldProperty
666 &MixtureType::patchFaceThermoMixture,
673 template<
class BasicThermo,
class MixtureType>
677 if (MixtureType::thermoType::enthalpy())
688 template<
class BasicThermo,
class MixtureType>
692 if (MixtureType::thermoType::enthalpy())
694 return this->kappa_.boundaryField()[
patchi]/Cp_.boundaryField()[
patchi];
698 return this->kappa_.boundaryField()[
patchi]/Cv_.boundaryField()[
patchi];
703 template<
class BasicThermo,
class MixtureType>
714 template<
class BasicThermo,
class MixtureType>
723 this->kappa_.boundaryField()[
patchi]
728 template<
class BasicThermo,
class MixtureType>
735 if (MixtureType::thermoType::enthalpy())
744 (this->kappa_ + Cp_*alphat)/Cv_
750 template<
class BasicThermo,
class MixtureType>
758 if (MixtureType::thermoType::enthalpy())
761 this->kappa_.boundaryField()[
patchi]/Cp_.boundaryField()[
patchi]
768 this->kappa_.boundaryField()[
patchi]
769 + Cp_.boundaryField()[
patchi]*alphat
770 )/Cv_.boundaryField()[
patchi];
775 template<
class BasicThermo,
class MixtureType>
const char *const group
Group name for atomic constants.
scalar Cv(const scalar p, const scalar T) const
#define forAll(list, i)
Loop across all elements in list.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
virtual tmp< volScalarField > gamma() const
Gamma = Cp/Cv [].
scalar Hs(const scalar p, const scalar T) const
T & ref() const
Return non-const reference or generate a fatal error.
virtual ~heThermo()
Destructor.
void size(const label)
Override size to be inconsistent with allocated storage.
virtual tmp< volScalarField > Cpv() const
Heat capacity at constant pressure/volume [J/kg/K].
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
const Time & time() const
Return the top-level database.
virtual volScalarField & he()
Enthalpy/Internal energy [J/kg].
virtual tmp< volScalarField > hc() const
Enthalpy of formation [J/kg].
tmp< scalarField > patchFieldProperty(PatchFaceMixture patchFaceMixture, Method psiMethod, const label patchi, const Args &... args) const
Return a scalarField of the given property on a patch.
Dimension set for the base types.
static word timeName(const scalar, const int precision=curPrecision_)
Return time name of given scalar time.
A class for handling words, derived from string.
static UIndirectList< scalar > cellSetScalarList(const volScalarField &psi, const labelList &cells)
Return an indirect list of a field for the given set of cells.
const volScalarField & psi
virtual tmp< volScalarField > THE(const volScalarField &h, const volScalarField &p, const volScalarField &T0) const
Temperature from enthalpy/internal energy.
virtual tmp< volScalarField > alphahe() const
Thermal diffusivity of energy of mixture [kg/m/s].
virtual tmp< volScalarField > hs() const
Sensible enthalpy [J/kg/K].
virtual tmp< volScalarField > ha() const
Absolute enthalpy [J/kg/K].
volScalarField scalarField(fieldObject, mesh)
virtual tmp< volScalarField > alphaEff(const volScalarField &alphat) const
Effective turbulent thermal diffusivity of energy.
virtual const volScalarField & Cv() const
Heat capacity at constant volume [J/kg/K].
tmp< scalarField > cellSetProperty(CellMixture cellMixture, Method psiMethod, const labelList &cells, const Args &... args) const
Return a scalarField of the given property on a cell set.
const dimensionSet dimEnergy
const dimensionSet dimMass
heThermo(const fvMesh &, const word &phaseName)
Construct from mesh.
virtual bool read()
Read thermophysical properties dictionary.
virtual tmp< volScalarField > kappaEff(const volScalarField &) const
Effective thermal turbulent conductivity of mixture [W/m/K].
const dimensionSet dimMoles
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual const volScalarField & Cp() const
Heat capacity at constant pressure [J/kg/K].
const scalarList W(::W(thermo))
Mesh data needed to do the Finite Volume discretisation.
A List with indirect addressing.
scalar Cp(const scalar p, const scalar T) const
scalar Ha(const scalar p, const scalar T) const
A class for managing temporary objects.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
void heBoundaryCorrection(volScalarField &he)
Correct the enthalpy/internal energy field boundaries.
tmp< volScalarField > volScalarFieldProperty(const word &psiName, const dimensionSet &psiDim, CellMixture cellMixture, PatchFaceMixture patchFaceMixture, Method psiMethod, const Args &... args) const
Return a volScalarField of the given property.
const dimensionSet dimTemperature
virtual tmp< volScalarField > W() const
Molecular weight [kg/kmol].