Go to the documentation of this file.
42 template<
class BasicTurbulenceModel>
45 this->nut_ = k_/omega_;
46 this->nut_.correctBoundaryConditions();
49 BasicTurbulenceModel::correctNut();
55 template<
class BasicTurbulenceModel>
64 const word& propertiesName,
131 this->runTime_.timeName(),
143 this->runTime_.timeName(),
151 bound(k_, this->kMin_);
152 bound(omega_, this->omegaMin_);
154 if (
type == typeName)
156 this->printCoeffs(
type);
163 template<
class BasicTurbulenceModel>
168 Cmu_.readIfPresent(this->coeffDict());
169 beta_.readIfPresent(this->coeffDict());
170 gamma_.readIfPresent(this->coeffDict());
171 alphaK_.readIfPresent(this->coeffDict());
172 alphaOmega_.readIfPresent(this->coeffDict());
181 template<
class BasicTurbulenceModel>
184 if (!this->turbulence_)
213 omega_.boundaryFieldRef().updateCoeffs();
228 omegaEqn.
ref().relax();
230 omegaEqn.
ref().boundaryManipulate(omega_.boundaryFieldRef());
233 bound(omega_, this->omegaMin_);
253 bound(k_, this->kMin_);
virtual bool read()
Read RASProperties dictionary.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
BasicTurbulenceModel::alphaField alphaField
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
A class for handling words, derived from Foam::string.
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
void clear() const noexcept
virtual void correctNut()
const dimensionedScalar G
Newtonian constant of gravitation.
A class for managing temporary objects.
Templated abstract base class for RAS turbulence models.
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
BasicTurbulenceModel::rhoField rhoField
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.
bool read(const char *buf, int32_t &val)
Same as readInt32.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
BasicTurbulenceModel::transportModel transportModel
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Bound the given scalar field if it has gone unbounded.
SolverPerformance< Type > solve(faMatrix< Type > &, Istream &)
Solve returning the solution statistics given convergence tolerance.
kOmega(const alphaField &alpha, const rhoField &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.
Generic dimensioned Type class.
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Eddy viscosity turbulence model base class.
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)