Loading...
Searching...
No Matches
SpalartAllmarasDDES.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-2016 OpenFOAM Foundation
9 Copyright (C) 2022 Upstream CFD GmbH
10 Copyright (C) 2015-2022 OpenCFD Ltd.
11-------------------------------------------------------------------------------
12License
13 This file is part of OpenFOAM.
14
15 OpenFOAM is free software: you can redistribute it and/or modify it
16 under the terms of the GNU General Public License as published by
17 the Free Software Foundation, either version 3 of the License, or
18 (at your option) any later version.
19
20 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23 for more details.
24
25 You should have received a copy of the GNU General Public License
26 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27
28\*---------------------------------------------------------------------------*/
29
30#include "SpalartAllmarasDDES.H"
31
32// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33
34namespace Foam
35{
36namespace LESModels
37{
38
39template<class BasicTurbulenceModel>
40const Foam::Enum
41<
42 typename Foam::LESModels::
43 SpalartAllmarasDDES<BasicTurbulenceModel>::shieldingMode
44>
46({
47 { shieldingMode::standard, "standard" },
48 { shieldingMode::ZDES2020, "ZDES2020" },
49});
50
51
52// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
53
54template<class BasicTurbulenceModel>
56(
57 const volScalarField& magGradU
58) const
59{
60 const volScalarField r(this->r(this->nuEff(), magGradU, this->y_));
61
62 tmp<volScalarField> tfd = 1 - tanh(pow(Cd1_*r, Cd2_));
63
64 switch (shielding_)
65 {
66 case shieldingMode::standard:
67 {
68 return tfd;
69 }
70 case shieldingMode::ZDES2020:
71 {
72 auto maxEps = [](const volScalarField& fld, const scalar eps){
73 return max(fld, dimensionedScalar(fld.dimensions(), eps));
74 };
75
76 volScalarField& fdStd = tfd.ref();
77 const auto& nuTilda = this->nuTilda_;
78 const volVectorField& n = wallDist::New(this->mesh_).n();
79
80 const volScalarField GnuTilda
81 (
82 Cd3_*maxEps(fvc::grad(nuTilda) & n, Zero)
83 / (maxEps(magGradU, SMALL)*this->kappa_*this->y_)
84 );
85
86 volScalarField fdGnuTilda(1 - tanh(pow(Cd1_*GnuTilda, Cd2_)));
87 const volScalarField GOmega
88 (
89 - (fvc::grad(mag(fvc::curl(this->U_))) & n)
90 * sqrt(nuTilda/maxEps(pow3(magGradU), SMALL))
91 );
92 const volScalarField alpha((7.0/6.0*Cd4_ - GOmega)/(Cd4_/6.0));
93 const volScalarField fRGOmega
94 (
95 pos(Cd4_ - GOmega)
96 + 1.0
97 /(1 + exp(min(-6*alpha/max(1 - sqr(alpha), SMALL), scalar(50))))
98 *pos(4*Cd4_/3.0 - GOmega)*pos(GOmega - Cd4_)
99 );
100
101 // Use more conservative fP2-function in case switch is true;
102 // otherwise use simplified formulation
103 if (usefP2_)
104 {
105 fdGnuTilda *=
106 (1.0 - tanh(pow(Cd1_*betaZDES_*r, Cd2_)))
107 / maxEps(fdStd, SMALL);
108 }
109
110 fdStd *= 1 - (1 - fdGnuTilda)*fRGOmega;
111 }
112 }
113
114 return tfd;
115}
116
117
118// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
119
120template<class BasicTurbulenceModel>
122(
123 const volScalarField& chi,
124 const volScalarField& fv1,
125 const volTensorField& gradU,
126 const volScalarField& dTilda
127) const
128{
129 if (this->useSigma_)
130 {
131 const volScalarField& lRAS(this->y_);
132 const volScalarField fv2(this->fv2(chi, fv1));
133 const volScalarField lLES(this->lengthScaleLES(chi, fv1));
134 const volScalarField Omega(this->Omega(gradU));
135 const volScalarField Ssigma(this->Ssigma(gradU));
136 const volScalarField SsigmaDES
137 (
138 Omega - fd(mag(gradU))*pos(lRAS - lLES)*(Omega - Ssigma)
139 );
140
141 return
142 max
143 (
144 SsigmaDES + fv2*this->nuTilda_/sqr(this->kappa_*dTilda),
145 this->Cs_*SsigmaDES
146 );
147 }
148
149 return
151 (
152 chi,
153 fv1,
154 gradU,
155 dTilda
156 );
157}
158
159
160template<class BasicTurbulenceModel>
162(
163 const volScalarField& chi,
164 const volScalarField& fv1,
165 const volTensorField& gradU
166) const
167{
168 const volScalarField& lRAS(this->y_);
169 const volScalarField lLES(this->lengthScaleLES(chi, fv1));
171
172 return max
173 (
174 lRAS - fd(mag(gradU))*max(lRAS - lLES, l0),
175 dimensionedScalar("small", dimLength, SMALL)
176 );
177}
178
179
180// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
181
182template<class BasicTurbulenceModel>
183SpalartAllmarasDDES<BasicTurbulenceModel>::SpalartAllmarasDDES
184(
185 const alphaField& alpha,
186 const rhoField& rho,
187 const volVectorField& U,
188 const surfaceScalarField& alphaRhoPhi,
189 const surfaceScalarField& phi,
190 const transportModel& transport,
191 const word& propertiesName,
192 const word& type
193)
194:
195 SpalartAllmarasDES<BasicTurbulenceModel>
196 (
197 alpha,
198 rho,
199 U,
200 alphaRhoPhi,
201 phi,
202 transport,
203 propertiesName,
204 type
205 ),
206 shielding_
207 (
208 shieldingModeNames.getOrDefault
209 (
210 "shielding",
211 this->coeffDict_,
212 shieldingMode::standard
213 )
214 ),
215 Cd1_
216 (
217 this->useSigma_
218 ? dimensioned<scalar>::getOrAddToDict
219 (
220 "Cd1Sigma",
221 this->coeffDict_,
222 10
223 )
224 : dimensioned<scalar>::getOrAddToDict
225 (
226 "Cd1",
227 this->coeffDict_,
228 8
229 )
230 ),
231 Cd2_
232 (
233 dimensioned<scalar>::getOrAddToDict
234 (
235 "Cd2",
236 this->coeffDict_,
237 3
238 )
239 ),
240 Cd3_
241 (
242 dimensioned<scalar>::getOrAddToDict
243 (
244 "Cd3",
245 this->coeffDict_,
246 25
247 )
248 ),
249 Cd4_
250 (
251 dimensioned<scalar>::getOrAddToDict
252 (
253 "Cd4",
254 this->coeffDict_,
255 0.03
256 )
257 ),
258 betaZDES_
259 (
260 dimensioned<scalar>::getOrAddToDict
261 (
262 "betaZDES",
263 this->coeffDict_,
264 2.5
265 )
266 ),
267 usefP2_
268 (
269 Switch::getOrAddToDict
270 (
271 "usefP2",
272 this->coeffDict_,
273 false
274 )
275 )
276{
277 if (type == typeName)
278 {
279 this->printCoeffs(type);
280
281 switch (shielding_)
282 {
284 {
285 Info<< "shielding function: standard DDES "
286 << "(Spalart et al., 2006)"
287 << nl;
288 break;
289 }
291 {
292 Info<< "shielding function: ZDES mode 2 (Deck & Renard, 2020)"
293 << nl;
294 break;
295 }
296 default:
297 {
299 << "Unrecognised 'shielding' option: "
301 << exit(FatalError);
302 }
303 }
304
305 if (usefP2_)
306 {
307 Info<< "fP2 term: active" << nl;
308 }
309 else
310 {
311 Info<< "fP2 term: inactive" << nl;
312 }
314}
315
316
317// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
318
319template<class BasicTurbulenceModel>
321{
323 {
324 shieldingModeNames.readIfPresent
325 (
326 "shielding",
327 this->coeffDict(),
328 shielding_
329 );
330
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());
337
338 return true;
340
341 return false;
342}
343
344
345template<class BasicTurbulenceModel>
347{
348 return fd(mag(fvc::grad(this->U_)));
349}
350
351
352// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
353
354} // End namespace LESModels
355} // End namespace Foam
356
357// ************************************************************************* //
label n
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.
BasicTurbulenceModel::transportModel transportModel
static const Enum< shieldingMode > shieldingModeNames
Shielding mode names.
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,...
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.
const volVectorField & n() const
Return reference to cached normal-to-wall field.
Definition wallDist.C:168
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
#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
tmp< GeometricField< Type, fvPatchField, volMesh > > curl(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition fvcCurl.C:40
Namespace for OpenFOAM.
dimensionedScalar pos(const dimensionedScalar &ds)
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
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.
Definition hashSets.C:26
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
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
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
volScalarField & alpha