Loading...
Searching...
No Matches
SpalartAllmarasIDDES.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2015 OpenFOAM Foundation
9 Copyright (C) 2015-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
30
31// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32
33namespace Foam
34{
35namespace LESModels
36{
37
38// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
39
40template<class BasicTurbulenceModel>
41const IDDESDelta& SpalartAllmarasIDDES<BasicTurbulenceModel>::setDelta() const
42{
43 if (!isA<IDDESDelta>(this->delta_()))
44 {
46 << "The delta function must be set to a " << IDDESDelta::typeName
47 << " -based model" << exit(FatalError);
48 }
49
50 return refCast<const IDDESDelta>(this->delta_());
51}
52
53
54template<class BasicTurbulenceModel>
55tmp<volScalarField> SpalartAllmarasIDDES<BasicTurbulenceModel>::alpha() const
56{
57 return max(0.25 - this->y_/IDDESDelta_.hmax(), scalar(-5));
58}
59
60
61template<class BasicTurbulenceModel>
62tmp<volScalarField> SpalartAllmarasIDDES<BasicTurbulenceModel>::ft
63(
64 const volScalarField& magGradU
65) const
66{
67 return tanh(pow3(sqr(Ct_)*this->r(this->nut_, magGradU, this->y_)));
68}
69
70
71template<class BasicTurbulenceModel>
72tmp<volScalarField> SpalartAllmarasIDDES<BasicTurbulenceModel>::fl
73(
74 const volScalarField& magGradU
75) const
76{
77 return tanh(pow(sqr(Cl_)*this->r(this->nu(), magGradU, this->y_), 10));
78}
79
80
81template<class BasicTurbulenceModel>
82tmp<volScalarField> SpalartAllmarasIDDES<BasicTurbulenceModel>::fdt
83(
84 const volScalarField& magGradU
85) const
86{
87 return 1 - tanh(pow(Cdt1_*this->r(this->nut_, magGradU, this->y_), Cdt2_));
88}
89
90
91// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
92
93template<class BasicTurbulenceModel>
95(
96 const volScalarField& chi,
97 const volScalarField& fv1,
98 const volTensorField& gradU
99) const
100{
101 const volScalarField magGradU(mag(gradU));
102 const volScalarField psi(this->psi(chi, fv1));
103
104 const volScalarField& lRAS = this->y_;
105 const volScalarField lLES(psi*this->CDES_*this->delta());
106
107 const volScalarField alpha(this->alpha());
108 const volScalarField expTerm(exp(sqr(alpha)));
109
110 tmp<volScalarField> fB = min(2*pow(expTerm, -9.0), scalar(1));
111 const volScalarField fdTilda(max(1 - fdt(magGradU), fB));
112
113 if (fe_)
114 {
116 2*lerp(pow(expTerm, -9.), pow(expTerm, -11.09), pos0(alpha));
117 tmp<volScalarField> fe2 = 1 - max(ft(magGradU), fl(magGradU));
118 tmp<volScalarField> fe = max(fe1 - 1, scalar(0))*psi*fe2;
119
120 // Original formulation from Shur et al. paper (2008)
121 return max
122 (
123 fdTilda*(1 + fe)*lRAS + (1 - fdTilda)*lLES,
124 dimensionedScalar("SMALL", dimLength, SMALL)
125 );
126 }
127
128
129 // Simplified formulation from Gritskevich et al. paper (2011) where fe = 0
130 return max
131 (
132 lerp(lLES, lRAS, fdTilda),
133 dimensionedScalar("SMALL", dimLength, SMALL)
134 );
135}
136
137
138// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
139
140template<class BasicTurbulenceModel>
141SpalartAllmarasIDDES<BasicTurbulenceModel>::SpalartAllmarasIDDES
142(
143 const alphaField& alpha,
144 const rhoField& rho,
145 const volVectorField& U,
146 const surfaceScalarField& alphaRhoPhi,
147 const surfaceScalarField& phi,
148 const transportModel& transport,
149 const word& propertiesName,
150 const word& type
151)
152:
153 SpalartAllmarasDES<BasicTurbulenceModel>
154 (
155 alpha,
156 rho,
157 U,
158 alphaRhoPhi,
159 phi,
160 transport,
161 propertiesName,
162 type
163 ),
164
165 Cdt1_
166 (
167 dimensioned<scalar>::getOrAddToDict
168 (
169 "Cdt1",
170 this->coeffDict_,
171 8
172 )
173 ),
174 Cdt2_
175 (
176 dimensioned<scalar>::getOrAddToDict
177 (
178 "Cdt2",
179 this->coeffDict_,
180 3
181 )
182 ),
183 Cl_
184 (
185 dimensioned<scalar>::getOrAddToDict
186 (
187 "Cl",
188 this->coeffDict_,
189 3.55
190 )
191 ),
192 Ct_
193 (
194 dimensioned<scalar>::getOrAddToDict
195 (
196 "Ct",
197 this->coeffDict_,
198 1.63
199 )
200 ),
201 fe_
202 (
203 Switch::getOrAddToDict
204 (
205 "fe",
206 this->coeffDict_,
207 true
208 )
209 ),
210
211 IDDESDelta_(setDelta())
212{
213 if (type == typeName)
214 {
215 this->printCoeffs(type);
217}
218
219
220// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
221
222template<class BasicTurbulenceModel>
224{
226 {
227 Cdt1_.readIfPresent(this->coeffDict());
228 Cdt2_.readIfPresent(this->coeffDict());
229 Cl_.readIfPresent(this->coeffDict());
230 Ct_.readIfPresent(this->coeffDict());
231
232 return true;
234
235 return false;
236}
237
238
239template<class BasicTurbulenceModel>
241{
242 const volScalarField alpha(this->alpha());
243 const volScalarField expTerm(exp(sqr(alpha)));
244
245 tmp<volScalarField> fB = min(2*pow(expTerm, -9.0), scalar(1));
246 return max(1 - fdt(mag(fvc::grad(this->U_))), fB);
247}
248
249
250// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
251
252} // End namespace LESModels
253} // End namespace Foam
254
255// ************************************************************************* //
scalar delta
SpalartAllmarasDES DES turbulence model for incompressible and compressible flows.
dimensionedScalar CDES_
DES coefficient.
virtual tmp< volScalarField > psi(const volScalarField &chi, const volScalarField &fv1) const
Return the low Reynolds number correction function.
virtual bool read()
Re-read model coefficients if they have changed.
BasicTurbulenceModel::alphaField alphaField
BasicTurbulenceModel::rhoField rhoField
virtual tmp< volScalarField > fd() const
Return the shielding function.
BasicTurbulenceModel::transportModel transportModel
const IDDESDelta & IDDESDelta_
IDDES delta.
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,...
Definition Switch.H:81
Generic dimensioned Type class.
A class for managing temporary objects.
Definition tmp.H:75
Base-class for all transport models used by the incompressible turbulence models.
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
const volScalarField & psi
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
Namespace for LES SGS models.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
dimensionedScalar exp(const dimensionedScalar &ds)
GeometricField< vector, fvPatchField, volMesh > volVectorField
dimensionedScalar pos0(const dimensionedScalar &ds)
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)
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
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
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.
Definition hashSets.C:26
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)
Definition errorManip.H:125
dimensioned< Type > lerp(const dimensioned< Type > &a, const dimensioned< Type > &b, const scalar t)
volScalarField & nu
volScalarField & alpha