Go to the documentation of this file.
38 namespace incompressible
64 return scalar(1) - 0.3*
exp(-
sqr(Rt));
133 const word& propertiesName,
332 bound(epsilon_, epsilonMin_);
334 if (
type == typeName)
408 epsEqn.
ref().relax();
Defines the attributes of an object for which implicit objectRegistry management is supported,...
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
addToRunTimeSelectionTable(RASModel, kkLOmega, dictionary)
LienCubicKE(const geometricOneField &alpha, const geometricOneField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName)
Construct from components.
A class for handling words, derived from Foam::string.
RASModel< turbulenceModel > RASModel
dimensionedTensor skew(const dimensionedTensor &dt)
const dimensionedScalar G
Newtonian constant of gravitation.
A class for managing temporary objects.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
static const wallDist & New(const fvMesh &mesh, Args &&... args)
Get existing or create a new MeshObject.
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
dimensionedScalar Cgamma1_
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedScalar Cgamma4_
dimensionedSymmTensor innerSqr(const dimensionedSymmTensor &dt)
tmp< volScalarField > f2() const
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
dimensionedScalar Cbeta1_
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
dimensionedScalar sigmak_
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
bool readIfPresent(const dictionary &dict)
dimensionedScalar pow3(const dimensionedScalar &ds)
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
virtual bool read()
Re-read model coefficients if they have changed.
Bound the given scalar field if it has gone unbounded.
SolverPerformance< Type > solve(faMatrix< Type > &, Istream &)
Solve returning the solution statistics given convergence tolerance.
defineTypeNameAndDebug(kkLOmega, 0)
virtual void correctNonlinearStress(const volTensorField &gradU)
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Macros for easy insertion into run-time selection tables.
dimensionedScalar sigmaEps_
Generic dimensioned Type class.
const volScalarField & y_
Wall distance.
dimensionedScalar Cbeta2_
dimensionedScalar Cbeta3_
Base-class for all transport models used by the incompressible turbulence models.
virtual void correctNut()
void correctBoundaryConditions()
Correct boundary field.
dimensionedScalar Cgamma2_
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
tmp< volScalarField > fMu() const
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
tmp< volScalarField > E(const volScalarField &f2) const
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
volSymmTensorField nonlinearStress_
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)