42template<
class BasicTurbulenceModel>
65 auto&
psi = tpsi.ref();
76 (1 - this->Cb1_/(this->Cw1_*
sqr(this->kappa_)*
fwStar_)
77 *(ft2 + (1 - ft2)*fv2))
78 /
max(SMALL, (fv1*
max(scalar(1
e-10), 1 - ft2)))
87template<
class BasicTurbulenceModel>
98template<
class BasicTurbulenceModel>
116 pos(dTilda - lRAS)*Omega + (scalar(1) -
pos(dTilda - lRAS))*Ssigma
122 SsigmaDES + fv2*this->nuTilda_/
sqr(this->kappa_*dTilda),
138template<
class BasicTurbulenceModel>
150 min(tdTilda.
ref(), tdTilda(), this->y_);
155template<
class BasicTurbulenceModel>
162 BasicTurbulenceModel::correctNut();
168template<
class BasicTurbulenceModel>
169SpalartAllmarasDES<BasicTurbulenceModel>::SpalartAllmarasDES
171 const alphaField&
alpha,
177 const word& propertiesName,
237 <<
"This model is not recommended and will be deprecated in future "
238 <<
"releases. Please consider using the DDES/IDDES versions instead"
241 this->printCoeffs(
type);
248template<
class BasicTurbulenceModel>
253 useSigma_.readIfPresent(
"useSigma", this->coeffDict());
254 CDES_.readIfPresent(this->coeffDict());
255 lowReCorrection_.readIfPresent(
"lowReCorrection", this->coeffDict());
256 fwStar_.readIfPresent(this->coeffDict());
265template<
class BasicTurbulenceModel>
277 this->mesh_.time().timeName(),
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing 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.
Templated abstract base class for DES models.
BasicTurbulenceModel::alphaField alphaField
Switch useSigma_
Switch to activate grey-area enhanced sigma-(D)DES.
BasicTurbulenceModel::rhoField rhoField
dimensionedScalar CDES_
DES coefficient.
virtual tmp< volScalarField > LESRegion() const
Return the LES field indicator.
virtual tmp< volScalarField > psi(const volScalarField &chi, const volScalarField &fv1) const
Return the low Reynolds number correction function.
dimensionedScalar fwStar_
virtual void correctNut()
virtual tmp< volScalarField > lengthScaleLES(const volScalarField &chi, const volScalarField &fv1) const
Return the LES length scale.
BasicTurbulenceModel::transportModel transportModel
virtual tmp< volScalarField > Stilda(const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU, const volScalarField &dTilda) const
Return the production term.
Switch lowReCorrection_
Flag for low Reynolds number correction.
virtual tmp< volScalarField > dTilda(const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU) const
Return the length scale.
virtual bool read()
Re-read model coefficients if they have changed.
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Generic dimensioned Type class.
static dimensioned< Type > getOrAddToDict(const word &name, dictionary &dict, const dimensionSet &dims=dimless, const Type &deflt=Type(Zero))
Construct dimensioned from dictionary, with default value.
A class for managing temporary objects.
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.
A class for handling words, derived from Foam::string.
const volScalarField & psi
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for LES SGS models.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
dimensionedScalar pos(const dimensionedScalar &ds)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
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)
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.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
GeometricField< tensor, fvPatchField, volMesh > volTensorField
dimensionedScalar sqrt(const dimensionedScalar &ds)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
dimensionedScalar neg(const dimensionedScalar &ds)
static constexpr const zero Zero
Global zero (0).
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.