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 -------------------------------------------------------------------------------
11 License
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 
29 #include "SpalartAllmarasIDDES.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace LESModels
36 {
37 
38 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
39 
40 template<class BasicTurbulenceModel>
41 const 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 
54 template<class BasicTurbulenceModel>
55 tmp<volScalarField> SpalartAllmarasIDDES<BasicTurbulenceModel>::alpha() const
56 {
57  return max(0.25 - this->y_/IDDESDelta_.hmax(), scalar(-5));
58 }
59 
60 
61 template<class BasicTurbulenceModel>
62 tmp<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 
71 template<class BasicTurbulenceModel>
72 tmp<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 
81 template<class BasicTurbulenceModel>
82 tmp<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 
93 template<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  {
115  tmp<volScalarField> fe1 =
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 
140 template<class BasicTurbulenceModel>
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);
216  }
217 }
218 
219 
220 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
221 
222 template<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;
233  }
234 
235  return false;
236 }
237 
238 
239 template<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
dimensionedScalar tanh(const dimensionedScalar &ds)
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual tmp< volScalarField > dTilda(const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU) const
Return the length scale.
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.
Definition: Switch.H:77
virtual bool read()
Re-read model coefficients if they have changed.
SpalartAllmarasDES DES turbulence model for incompressible and compressible flows.
Generic dimensioned Type class.
virtual bool read()
Re-read model coefficients if they have changed.
SpalartAllmaras IDDES turbulence model for incompressible and compressible flows. ...
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
dimensionedScalar exp(const dimensionedScalar &ds)
A class for handling words, derived from Foam::string.
Definition: word.H:63
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
dimensionedScalar pos0(const dimensionedScalar &ds)
virtual tmp< volScalarField > fd() const
Return the shielding function.
dimensioned< Type > lerp(const dimensioned< Type > &a, const dimensioned< Type > &b, const scalar t)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
U
Definition: pEqn.H:72
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Base-class for all transport models used by the incompressible turbulence models. ...
const volScalarField & psi
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
volScalarField & nu
Namespace for OpenFOAM.