42template<
class BasicTurbulenceModel>
45 this->nut_ = Cmu_*
sqr(k_)/epsilon_;
46 this->nut_.correctBoundaryConditions();
49 BasicTurbulenceModel::correctNut();
53template<
class BasicTurbulenceModel>
61 dimVolume*this->rho_.dimensions()*k_.dimensions()
68template<
class BasicTurbulenceModel>
76 dimVolume*this->rho_.dimensions()*epsilon_.dimensions()
85template<
class BasicTurbulenceModel>
94 const word& propertiesName,
169 IOobject::groupName(
"k", alphaRhoPhi.group()),
181 IOobject::groupName(
"epsilon", alphaRhoPhi.group()),
190 bound(k_, this->kMin_);
191 bound(epsilon_, this->epsilonMin_);
193 if (type == typeName)
195 this->printCoeffs(type);
202template<
class BasicTurbulenceModel>
207 Cmu_.readIfPresent(this->coeffDict());
208 C1_.readIfPresent(this->coeffDict());
209 C2_.readIfPresent(this->coeffDict());
210 C3_.readIfPresent(this->coeffDict());
211 sigmak_.readIfPresent(this->coeffDict());
212 sigmaEps_.readIfPresent(this->coeffDict());
221template<
class BasicTurbulenceModel>
224 if (!this->turbulence_)
242 fvc::div(fvc::absolute(this->
phi(), U))().v()
248 this->type() +
":GbyNu",
255 epsilon_.boundaryFieldRef().updateCoeffs();
261 + fvm::div(alphaRhoPhi, epsilon_)
262 - fvm::laplacian(
alpha*
rho*DepsilonEff(), epsilon_)
265 - fvm::SuSp(((2.0/3.0)*C1_ - C3_)*
alpha()*
rho()*
divU, epsilon_)
266 - fvm::Sp(C2_*
alpha()*
rho()*epsilon_()/k_(), epsilon_)
271 epsEqn.
ref().relax();
273 epsEqn.
ref().boundaryManipulate(epsilon_.boundaryFieldRef());
276 bound(epsilon_, this->epsilonMin_);
282 + fvm::div(alphaRhoPhi, k_)
283 - fvm::laplacian(
alpha*
rho*DkEff(), k_)
287 - fvm::Sp(
alpha()*
rho()*epsilon_()/k_(), k_)
296 bound(k_, this->kMin_);
Bound the given scalar field if it has gone unbounded.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Templated abstract base class for RAS turbulence models.
Standard k-epsilon turbulence model for incompressible and compressible flows including rapid distort...
BasicTurbulenceModel::alphaField alphaField
virtual tmp< fvScalarMatrix > epsilonSource() const
BasicTurbulenceModel::rhoField rhoField
virtual void correctNut()
virtual tmp< fvScalarMatrix > kSource() const
BasicTurbulenceModel::transportModel transportModel
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Generic dimensioned Type class.
Eddy viscosity turbulence model base class.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
A class for managing temporary objects.
void clear() const noexcept
A class for handling words, derived from Foam::string.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
const dimensionSet dimVolume(pow3(dimLength))