39 template<
class BasicTurbulenceModel>
47 { shieldingMode::standard,
"standard" },
48 { shieldingMode::ZDES2020,
"ZDES2020" },
54 template<
class BasicTurbulenceModel>
60 const volScalarField r(this->r(this->nuEff(), magGradU, this->y_));
62 tmp<volScalarField> tfd = 1 -
tanh(
pow(Cd1_*r, Cd2_));
66 case shieldingMode::standard:
70 case shieldingMode::ZDES2020:
77 const auto& nuTilda = this->nuTilda_;
83 / (maxEps(magGradU, SMALL)*this->kappa_*this->y_)
90 *
sqrt(nuTilda/maxEps(
pow3(magGradU), SMALL))
98 *
pos(4*Cd4_/3.0 - GOmega)*
pos(GOmega - Cd4_)
106 (1.0 -
tanh(
pow(Cd1_*betaZDES_*r, Cd2_)))
107 / maxEps(fdStd, SMALL);
110 fdStd *= 1 - (1 - fdGnuTilda)*fRGOmega;
120 template<
class BasicTurbulenceModel>
138 Omega - fd(
mag(gradU))*
pos(lRAS - lLES)*(Omega - Ssigma)
144 SsigmaDES + fv2*this->nuTilda_/
sqr(this->kappa_*dTilda),
160 template<
class BasicTurbulenceModel>
174 lRAS - fd(
mag(gradU))*
max(lRAS - lLES, l0),
182 template<
class BasicTurbulenceModel>
185 const alphaField&
alpha,
191 const word& propertiesName,
208 shieldingModeNames.getOrDefault
212 shieldingMode::standard
277 if (
type == typeName)
279 this->printCoeffs(
type);
285 Info<<
"shielding function: standard DDES " 286 <<
"(Spalart et al., 2006)" 292 Info<<
"shielding function: ZDES mode 2 (Deck & Renard, 2020)" 299 <<
"Unrecognised 'shielding' option: " 307 Info<<
"fP2 term: active" <<
nl;
311 Info<<
"fP2 term: inactive" <<
nl;
319 template<
class BasicTurbulenceModel>
324 shieldingModeNames.readIfPresent
331 Cd1_.readIfPresent(this->coeffDict());
332 Cd2_.readIfPresent(this->coeffDict());
333 Cd3_.readIfPresent(this->coeffDict());
334 Cd4_.readIfPresent(this->coeffDict());
335 betaZDES_.readIfPresent(this->coeffDict());
336 usefP2_.readIfPresent(
"usefP2", this->coeffDict());
345 template<
class BasicTurbulenceModel>
dimensionedScalar tanh(const dimensionedScalar &ds)
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
errorManipArg< error, int > exit(error &err, const int errNo=1)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const volVectorField & n() const
Return reference to cached normal-to-wall field.
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
static const wallDist & New(const fvMesh &mesh, Args &&... args)
Get existing or create MeshObject registered with typeName.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
constexpr char nl
The newline '\n' character (0x0a)
dimensionedScalar sqrt(const dimensionedScalar &ds)
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, any/none. Also accepts 0/1 as a string and shortcuts t/f, y/n.
virtual bool read()
Re-read model coefficients if they have changed.
SpalartAllmarasDES DES turbulence model for incompressible and compressible flows.
GeometricField< vector, fvPatchField, volMesh > volVectorField
tmp< GeometricField< Type, fvPatchField, volMesh > > curl(const GeometricField< Type, fvPatchField, volMesh > &vf)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar pos(const dimensionedScalar &ds)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
dimensionedScalar exp(const dimensionedScalar &ds)
shieldingMode
Shielding modes.
A class for handling words, derived from Foam::string.
virtual tmp< volScalarField > dTilda(const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU) const
Return the length scale.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
static const Enum< shieldingMode > shieldingModeNames
Shielding mode names.
SpalartAllmaras DDES turbulence model for incompressible and compressible flows.
virtual tmp< volScalarField > Stilda(const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU, const volScalarField &dTilda) const
Return the production term.
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
virtual bool read()
Read from dictionary.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Base-class for all transport models used by the incompressible turbulence models. ...
messageStream Info
Information stream (stdout output on master, null elsewhere)
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef()
A class for managing temporary objects.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
virtual tmp< volScalarField > fd() const
Return the shielding function.
shieldingMode shielding_
Shielding mode.
static constexpr const zero Zero
Global zero (0)