41template<
class BasicTurbulenceModel>
46 this->
nut_.correctBoundaryConditions();
49 BasicTurbulenceModel::correctNut();
53template<
class BasicTurbulenceModel>
73template<
class BasicTurbulenceModel>
80 pow(k_, 1.5)/epsilon_,
98template<
class BasicTurbulenceModel>
99kEpsilonPhitF<BasicTurbulenceModel>::kEpsilonPhitF
101 const alphaField&
alpha,
107 const word& propertiesName,
254 IOobject::groupName(
"k", alphaRhoPhi.group()),
266 IOobject::groupName(
"epsilon", alphaRhoPhi.group()),
278 IOobject::groupName(
"phit", alphaRhoPhi.group()),
290 IOobject::groupName(
"f", alphaRhoPhi.group()),
325 this->printCoeffs(
type);
336 <<
"Non-zero values are required for the model constants:" <<
nl
347template<
class BasicTurbulenceModel>
352 includeNu_.readIfPresent(
"includeNu", this->coeffDict());
353 Cmu_.readIfPresent(this->coeffDict());
354 Ceps1a_.readIfPresent(this->coeffDict());
355 Ceps1b_.readIfPresent(this->coeffDict());
356 Ceps1c_.readIfPresent(this->coeffDict());
357 Ceps2_.readIfPresent(this->coeffDict());
358 Cf1_.readIfPresent(this->coeffDict());
359 Cf2_.readIfPresent(this->coeffDict());
360 CL_.readIfPresent(this->coeffDict());
361 Ceta_.readIfPresent(this->coeffDict());
362 CT_.readIfPresent(this->coeffDict());
363 sigmaK_.readIfPresent(this->coeffDict());
364 sigmaEps_.readIfPresent(this->coeffDict());
365 sigmaPhit_.readIfPresent(this->coeffDict());
374template<
class BasicTurbulenceModel>
377 if (!this->turbulence_)
383 const alphaField&
alpha = this->alpha_;
384 const rhoField&
rho = this->rho_;
390 eddyViscosity<RASModel<BasicTurbulenceModel>>
::correct();
409 Ceps1a_*(Ceps1b_ + Ceps1c_*
sqrt(1.0/phit_()))
414 epsilon_.boundaryFieldRef().updateCoeffs();
416 epsilon_.boundaryFieldRef().template evaluateCoupled<coupledFvPatch>();
420 tmp<fvScalarMatrix> epsEqn
436 epsEqn.ref().relax();
438 epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
441 bound(epsilon_, this->epsilonMin_);
446 tmp<fvScalarMatrix> kEqn
462 bound(k_, this->kMin_);
467 tmp<fvScalarMatrix> fEqn
473 (Cf1_ - 1.0)*(phit_() - 2.0/3.0)/T_()
475 +(Cf2_*(2.0/3.0)*
divU)
488 tmp<fvScalarMatrix> phitEqn
502 /(k_()*sigmaPhit_*phit_())
509 phitEqn.ref().relax();
513 bound(phit_, phitMin_);
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.
volScalarField epsilon_
Turbulent kinetic energy dissipation rate [m2/s3].
BasicTurbulenceModel::alphaField alphaField
BasicTurbulenceModel::rhoField rhoField
volScalarField k_
Turbulent kinetic energy [m2/s2].
dimensionedScalar sigmaEps_
dimensionedScalar sigmaPhit_
tmp< volScalarField > DphitEff() const
Return the effective diffusivity for phit (LUU:Eq. 17).
virtual void correct()
Solve the transport equations and correct the turbulent viscosity.
dimensionedScalar phitMin_
dimensionedScalar Ceps1b_
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon (LUU:Eq. 4).
volScalarField phit_
Normalised wall-normal fluctuating velocity scale [-].
dimensionedScalar Ceps1a_
volScalarField T_
Turbulent time scale [s].
tmp< volScalarField > Ts() const
Return the turbulent time scale, T.
virtual void correctNut()
Update nut with the latest available k, phit, and T.
volScalarField f_
Elliptic relaxation factor [1/s].
BasicTurbulenceModel::transportModel transportModel
dimensionedScalar sigmaK_
tmp< volScalarField > Ls() const
Return the turbulent length scale, L.
virtual bool read()
Re-read model coefficients if they have changed.
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k (LUU:Eq. 3).
dimensionedScalar Ceps1c_
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
const Type & value() const noexcept
Return const reference to value.
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.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
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< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
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.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
GeometricField< vector, fvPatchField, volMesh > volVectorField
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)
dimensionedScalar pow3(const dimensionedScalar &ds)
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
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")
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0).
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
dimensionedScalar pow025(const dimensionedScalar &ds)
constexpr char nl
The newline '\n' character (0x0a).