48 #ifndef kineticTheoryModel_H 49 #define kineticTheoryModel_H 53 #include "phaseCompressibleTurbulenceModel.H" 55 #include "phaseModel.H" 56 #include "dragModel.H" 57 #include "viscosityModel.H" 58 #include "conductivityModel.H" 59 #include "radialModel.H" 60 #include "granularPressureModel.H" 61 #include "frictionalStressModel.H" 75 class kineticTheoryModel
79 RASModel<EddyDiffusivity<phaseCompressibleTurbulenceModel>>
86 const phaseModel& phase_;
92 autoPtr<kineticTheoryModels::viscosityModel> viscosityModel_;
95 autoPtr<kineticTheoryModels::conductivityModel> conductivityModel_;
98 autoPtr<kineticTheoryModels::radialModel> radialModel_;
101 autoPtr<kineticTheoryModels::granularPressureModel>
102 granularPressureModel_;
105 autoPtr<kineticTheoryModels::frictionalStressModel>
106 frictionalStressModel_;
154 kineticTheoryModel(
const kineticTheoryModel&) =
delete;
157 void operator=(
const kineticTheoryModel&) =
delete;
178 const word&
type = typeName
192 virtual tmp<volScalarField>
nuEff()
const 198 virtual tmp<scalarField>
nuEff(
const label patchi)
const 200 return this->
nut(patchi);
204 virtual tmp<volScalarField>
k()
const;
207 virtual tmp<volScalarField>
epsilon()
const;
210 virtual tmp<volScalarField>
omega()
const;
213 virtual tmp<volSymmTensorField>
R()
const;
217 virtual tmp<volScalarField>
pPrime()
const;
221 virtual tmp<surfaceScalarField>
pPrimef()
const;
224 virtual tmp<volSymmTensorField>
devRhoReff()
const;
virtual tmp< volScalarField > nut() const
Return the turbulence viscosity.
const volVectorField & U() const
Access function to velocity field.
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.
TypeName("kineticTheory")
Runtime type information.
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
GeometricField< vector, fvPatchField, volMesh > volVectorField
virtual bool read()
Re-read model coefficients if they have changed.
const transportModel & transport() const
Access function to incompressible transport model.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
virtual tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U) const
Return the source term for the momentum equation.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
virtual void correct()
Solve the kinetic theory equations and correct the viscosity.
virtual ~kineticTheoryModel()
Destructor.
static const word propertiesName
Default name of the turbulence properties dictionary.
virtual tmp< surfaceScalarField > pPrimef() const
Return the face-phase-pressure'.
const surfaceScalarField & alphaRhoPhi() const
Access function to phase flux field.
virtual tmp< volScalarField > omega() const
Return the specific dissipation rate.
const volScalarField & rho() const
Return the density field.
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
virtual tmp< volSymmTensorField > devRhoReff() const
Return the effective stress tensor.
virtual tmp< volScalarField > pPrime() const
Return the phase-pressure'.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual tmp< surfaceScalarField > phi() const
Return the volumetric flux field.
A class for managing temporary objects.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
virtual tmp< volScalarField > alpha() const
Return the laminar thermal diffusivity for enthalpy [kg/m/s].
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.