43template<
class BasicTurbulenceModel>
59 lEps =
y*ClStar*(1.0 -
exp(-1.0*
Rey/(2*ClStar)));
66 this->
nut_ *= lambdaEps;
69 const scalar Amu = 70.0;
74 this->
nut_.correctBoundaryConditions();
77 BasicTurbulenceModel::correctNut();
81template<
class BasicTurbulenceModel>
86 if (twoLayerTreatment_)
89 const alphaField&
alpha = this->alpha_;
90 const rhoField&
rho = this->rho_;
98 *(
sqrt(k_())/lEps - epsilon_()/k_()),
111template<
class BasicTurbulenceModel>
124template<
class BasicTurbulenceModel>
125kEpsilon<BasicTurbulenceModel>::kEpsilon
127 const alphaField&
alpha,
133 const word& propertiesName,
257 Info<<
"Two-layer wall treatment activated" <<
endl;
259 lEpsPtr_ = std::make_unique<volScalarField::Internal>
300template<
class BasicTurbulenceModel>
305 Cmu_.readIfPresent(this->coeffDict());
306 C1_.readIfPresent(this->coeffDict());
307 C2_.readIfPresent(this->coeffDict());
308 C3_.readIfPresent(this->coeffDict());
309 sigmak_.readIfPresent(this->coeffDict());
310 sigmaEps_.readIfPresent(this->coeffDict());
314 ReyStar_.readIfPresent(this->coeffDict());
315 ReyFactor_.readIfPresent(this->coeffDict());
324template<
class BasicTurbulenceModel>
327 if (!this->turbulence_)
333 const alphaField&
alpha = this->alpha_;
334 const rhoField&
rho = this->rho_;
341 eddyViscosity<RASModel<BasicTurbulenceModel>>
::correct();
358 epsilon_.boundaryFieldRef().updateCoeffs();
360 epsilon_.boundaryFieldRef().template evaluateCoupled<coupledFvPatch>();
363 tmp<fvScalarMatrix> epsEqn
376 epsEqn.ref().relax();
378 epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
381 bound(epsilon_, this->epsilonMin_);
384 tmp<fvScalarMatrix> kEqn
401 bound(k_, this->kMin_);
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
Bound the given scalar field if it has gone unbounded.
DimensionedField< scalar, volMesh > Internal
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ MUST_READ
Reading required.
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
static word group(const word &name)
Return group (extension part of name).
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
static FOAM_NO_DANGLING_REFERENCE const wallDist & New(const fvMesh &mesh, Args &&... args)
Templated abstract base class for RAS turbulence models.
virtual void printCoeffs(const word &type)
dimensionedScalar epsilonMin_
BasicTurbulenceModel::alphaField alphaField
virtual tmp< fvScalarMatrix > epsilonSource() const
dimensionedScalar sigmak_
BasicTurbulenceModel::rhoField rhoField
std::unique_ptr< volScalarField > lambdaEpsPtr_
dimensionedScalar sigmaEps_
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
dimensionedScalar ReyStar_
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
std::unique_ptr< volScalarField::Internal > lEpsPtr_
dimensionedScalar ReyFactor_
virtual void correctNut()
virtual tmp< fvScalarMatrix > kSource() const
BasicTurbulenceModel::transportModel transportModel
virtual bool read()
Re-read model coefficients if they have changed.
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
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.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
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.
const surfaceScalarField & alphaRhoPhi() const
Access function to phase flux field.
const volScalarField & y() const noexcept
Return reference to cached distance-to-wall field.
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.
dimensionedScalar exp(const dimensionedScalar &ds)
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
bool read(const char *buf, int32_t &val)
Same as readInt32.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
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)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.
dimensionedScalar tanh(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
dimensionedScalar atanh(const dimensionedScalar &ds)
static constexpr const zero Zero
Global zero (0).
const dimensionSet dimVolume(pow3(dimLength))
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.