42template<
class BasicTurbulenceModel>
46 this->
nut_.correctBoundaryConditions();
49 BasicTurbulenceModel::correctNut();
55template<
class BasicTurbulenceModel>
58 const alphaField&
alpha,
64 const word& propertiesName,
130 IOobject::groupName(
"k", alphaRhoPhi.group()),
142 IOobject::groupName(
"omega", alphaRhoPhi.group()),
156 this->printCoeffs(
type);
163template<
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());
181template<
class BasicTurbulenceModel>
184 if (!this->turbulence_)
190 const alphaField&
alpha = this->alpha_;
191 const rhoField&
rho = this->rho_;
197 eddyViscosity<RASModel<BasicTurbulenceModel>>
::correct();
213 omega_.boundaryFieldRef().updateCoeffs();
215 omega_.boundaryFieldRef().template evaluateCoupled<coupledFvPatch>();
218 tmp<fvScalarMatrix> omegaEqn
230 omegaEqn.ref().relax();
232 omegaEqn.ref().boundaryManipulate(omega_.boundaryFieldRef());
235 bound(omega_, this->omegaMin_);
239 tmp<fvScalarMatrix> kEqn
255 bound(k_, this->kMin_);
Bound the given scalar field if it has gone unbounded.
DimensionedField< scalar, volMesh > Internal
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Templated abstract base class for RAS turbulence models.
BasicTurbulenceModel::alphaField alphaField
BasicTurbulenceModel::rhoField rhoField
tmp< volScalarField > DomegaEff() const
Return the effective diffusivity for omega.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
dimensionedScalar alphaOmega_
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.
virtual void correctNut()
BasicTurbulenceModel::transportModel transportModel
dimensionedScalar alphaK_
virtual bool read()
Read RASProperties dictionary.
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
Generic dimensioned Type class.
Eddy viscosity turbulence model base class.
eddyViscosity(const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName)
virtual tmp< volScalarField > nut() const
void correct(GeometricField< Type, PatchField, GeoMesh > &field)
Apply correction to field.
Finite-volume options, which is an IOdictionary of values and a fv::optionList.
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present otherwise lookup and return.
A class for managing temporary objects.
void clear() const noexcept
If object pointer points to valid object: delete object and set pointer to nullptr.
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Base-class for all transport models used by the incompressible turbulence models.
A class for handling words, derived from Foam::string.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
zeroField SuSp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
GeometricField< vector, fvPatchField, volMesh > volVectorField
bool read(const char *buf, int32_t &val)
Same as readInt32.
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")