39template<
class BasicTurbulenceModel>
42 typename Foam::LESModels::
43 SpalartAllmarasDDES<BasicTurbulenceModel>::shieldingMode
47 { shieldingMode::standard,
"standard" },
48 { shieldingMode::ZDES2020,
"ZDES2020" },
54template<
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;
120template<
class BasicTurbulenceModel>
138 Omega - fd(
mag(gradU))*
pos(lRAS - lLES)*(Omega - Ssigma)
144 SsigmaDES + fv2*this->nuTilda_/
sqr(this->kappa_*dTilda),
160template<
class BasicTurbulenceModel>
174 lRAS - fd(
mag(gradU))*
max(lRAS - lLES, l0),
182template<
class BasicTurbulenceModel>
183SpalartAllmarasDDES<BasicTurbulenceModel>::SpalartAllmarasDDES
185 const alphaField&
alpha,
191 const word& propertiesName,
208 shieldingModeNames.getOrDefault
212 shieldingMode::standard
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;
319template<
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());
345template<
class BasicTurbulenceModel>
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef().
BasicTurbulenceModel::alphaField alphaField
BasicTurbulenceModel::rhoField rhoField
virtual tmp< volScalarField > fd() const
Return the shielding function.
dimensionedScalar betaZDES_
BasicTurbulenceModel::transportModel transportModel
static const Enum< shieldingMode > shieldingModeNames
Shielding mode names.
shieldingMode
Shielding modes.
virtual tmp< volScalarField > Stilda(const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU, const volScalarField &dTilda) const
Return the production term.
shieldingMode shielding_
Shielding mode.
virtual tmp< volScalarField > dTilda(const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU) const
Return the length scale.
virtual bool read()
Read from dictionary.
SpalartAllmarasDES DES turbulence model for incompressible and compressible flows.
Switch useSigma_
Switch to activate grey-area enhanced sigma-(D)DES.
virtual tmp< volScalarField > lengthScaleLES(const volScalarField &chi, const volScalarField &fv1) const
Return the LES length scale.
virtual bool read()
Re-read model coefficients if they have changed.
static FOAM_NO_DANGLING_REFERENCE const wallDist & New(const fvMesh &mesh, Args &&... args)
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Generic dimensioned Type class.
A class for managing temporary objects.
Base-class for all transport models used by the incompressible turbulence models.
const volVectorField & n() const
Return reference to cached normal-to-wall field.
A class for handling words, derived from Foam::string.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Namespace for LES SGS models.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< GeometricField< Type, fvPatchField, volMesh > > curl(const GeometricField< Type, fvPatchField, volMesh > &vf)
dimensionedScalar pos(const dimensionedScalar &ds)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
dimensionedScalar exp(const dimensionedScalar &ds)
GeometricField< vector, fvPatchField, volMesh > volVectorField
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
dimensionedScalar pow3(const dimensionedScalar &ds)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar tanh(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
GeometricField< tensor, fvPatchField, volMesh > volTensorField
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
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)
constexpr char nl
The newline '\n' character (0x0a).