Go to the documentation of this file.
31 template<
class BasePhaseModel>
34 const phaseSystem&
fluid,
35 const word& phaseName,
39 BasePhaseModel(
fluid, phaseName, index),
44 IOobject::groupName(
"K", this->
name()),
56 template<
class BasePhaseModel>
63 template<
class BasePhaseModel>
66 BasePhaseModel::correctKinematics();
71 template<
class BasePhaseModel>
76 this->thermo_->correct();
80 template<
class BasePhaseModel>
94 tmp<fvScalarMatrix> tEEqn
113 if (he.name() == this->thermo_->phasePropertyName(
"e"))
119 else if (this->thermo_->dpdt())
128 template<
class BasePhaseModel>
131 return !this->
thermo().incompressible();
135 template<
class BasePhaseModel>
virtual void correctThermo()
Correct the thermodynamics.
AnisothermalPhaseModel(const phaseSystem &fluid, const word &phaseName, const label index)
A class for managing temporary objects.
const dimensionSet dimVelocity
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
surfaceScalarField alphaPhi(IOobject("alphaPhi", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), phi *fvc::interpolate(alpha1))
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
autoPtr< compressible::turbulenceModel > turbulence
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
virtual ~AnisothermalPhaseModel()
Destructor.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
virtual tmp< fvScalarMatrix > heEqn()
Return the enthalpy equation.
const volScalarField & alphaEff
virtual void correctKinematics()
Correct the kinematics.
virtual bool compressible() const
Return true if the phase is compressible otherwise false.
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual const volScalarField & K() const
Return the phase kinetic energy.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Generic GeometricField class.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
word name(const complex &)
Return a string representation of a complex.