59 #ifndef kEpsilonLopesdaCosta_H
60 #define kEpsilonLopesdaCosta_H
77 template<
class BasicMomentumTransportModel>
146 typedef typename BasicMomentumTransportModel::alphaField
alphaField;
147 typedef typename BasicMomentumTransportModel::rhoField
rhoField;
188 (this->
nut_/sigmak_ + this->nu())
198 (this->
nut_/sigmaEps_ + this->nu())
Graphite solid properties.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
Variant of the standard k-epsilon turbulence model with additional source terms to handle the changes...
BasicMomentumTransportModel::alphaField alphaField
volScalarField::Internal C1_
volScalarField::Internal C5_
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
virtual tmp< fvScalarMatrix > kSource(const volScalarField::Internal &magU, const volScalarField::Internal &magU3) const
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
virtual ~kEpsilonLopesdaCosta()
Destructor.
virtual tmp< volScalarField > omega() const
Return the turbulence specific dissipation rate.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
tmp< volScalarField > boundEpsilon()
Bound epsilon and return Cmu*sqr(k) for nut.
volScalarField::Internal betad_
volScalarField::Internal betap_
TypeName("kEpsilonLopesdaCosta")
Runtime type information.
volScalarField::Internal C4_
void setPorosityCoefficients()
kEpsilonLopesdaCosta(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity, const word &type=typeName)
Construct from components.
void operator=(const kEpsilonLopesdaCosta &)=delete
Disallow default bitwise assignment.
volScalarField::Internal C2_
void setCdAv(volScalarField::Internal &C, const porosityModels::powerLawLopesdaCosta &pm)
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
virtual void correctNut()
Correct the eddy-viscosity nut.
virtual tmp< fvScalarMatrix > epsilonSource(const volScalarField::Internal &magU, const volScalarField::Internal &magU3) const
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
void setPorosityCoefficient(volScalarField::Internal &C, const porosityModels::powerLawLopesdaCosta &pm)
volScalarField::Internal CdAv_
virtual bool read()
Re-read model coefficients if they have changed.
BasicMomentumTransportModel::rhoField rhoField
Eddy viscosity turbulence model base class.
Variant of the power law porosity model with spatially varying drag coefficient.
A class for managing temporary objects.
Abstract base class for all fluid physical properties.
A class for handling words, derived from string.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.