Go to the documentation of this file.
32 template<
class TransportModel>
43 const word& propertiesName
67 template<
class TransportModel>
76 const word& propertiesName
104 template<
class TransportModel>
112 IOobject::groupName(
"pPrime", this->alphaRhoPhi_.group()),
113 this->runTime_.timeName(),
124 template<
class TransportModel>
132 IOobject::groupName(
"pPrimef", this->alphaRhoPhi_.group()),
133 this->runTime_.timeName(),
144 template<
class TransportModel>
152 template<
class TransportModel>
159 return divDevRhoReff(
U);
163 template<
class TransportModel>
173 template<
class TransportModel>
182 return divDevReff(
U);
virtual tmp< volSymmTensorField > devRhoReff() const
Return the effective stress tensor.
const dimensionSet dimPressure
Defines the attributes of an object for which implicit objectRegistry management is supported,...
A class for handling words, derived from Foam::string.
static autoPtr< PhaseIncompressibleTurbulenceModel > New(const alphaField &alpha, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const TransportModel &transportModel, const word &propertiesName=turbulenceModel::propertiesName)
Return a reference to the selected turbulence model.
A class for managing temporary objects.
static constexpr const zero Zero
Global zero (0)
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
virtual tmp< fvVectorMatrix > divDevReff(volVectorField &U) const
Return the source term for the momentum equation.
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
virtual tmp< volScalarField > pPrime() const
Return the phase-pressure'.
virtual tmp< volSymmTensorField > devReff() const
Return the effective stress tensor.
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Templated abstract base class for turbulence models.
virtual tmp< surfaceScalarField > pPrimef() const
Return the face-phase-pressure'.
Templated abstract base class for multiphase incompressible turbulence models.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
PhaseIncompressibleTurbulenceModel(const word &type, const alphaField &alpha, const geometricOneField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const TransportModel &transportModel, const word &propertiesName)
Construct.
virtual tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U) const
Return the source term for the momentum equation.
TransportModel transportModel
Abstract base class for turbulence models (RAS, LES and laminar).