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-2020 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_)*rd(this->nut_, magGradU)));
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_)*rd(this->nu(), magGradU), 10));
78}
79
80
81template<class BasicTurbulenceModel>
82tmp<volScalarField> SpalartAllmarasIDDES<BasicTurbulenceModel>::rd
83(
84 const volScalarField& nur,
85 const volScalarField& magGradU
86) const
87{
88 tmp<volScalarField> tr
89 (
90 min
91 (
92 nur
93 /(
94 max
95 (
96 magGradU,
97 dimensionedScalar("SMALL", magGradU.dimensions(), SMALL)
98 )
99 *sqr(this->kappa_*this->y_)
100 ),
101 scalar(10)
102 )
103 );
104 tr.ref().boundaryFieldRef() == 0.0;
105
106 return tr;
107}
108
109
110template<class BasicTurbulenceModel>
111tmp<volScalarField> SpalartAllmarasIDDES<BasicTurbulenceModel>::fdt
112(
113 const volScalarField& magGradU
114) const
115{
116 return 1 - tanh(pow(Cdt1_*rd(this->nut_, magGradU), Cdt2_));
117}
118
119
120// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
121
122template<class BasicTurbulenceModel>
124(
125 const volScalarField& chi,
126 const volScalarField& fv1,
127 const volTensorField& gradU
128) const
129{
130 const volScalarField magGradU(mag(gradU));
131 const volScalarField psi(this->psi(chi, fv1));
132
133 const volScalarField& lRAS(this->y_);
134 const volScalarField lLES(psi*this->CDES_*this->delta());
135
136 const volScalarField alpha(this->alpha());
137 const volScalarField expTerm(exp(sqr(alpha)));
138
139 tmp<volScalarField> fB = min(2*pow(expTerm, -9.0), scalar(1));
141 2*(pos0(alpha)*pow(expTerm, -11.09) + neg(alpha)*pow(expTerm, -9.0));
142 tmp<volScalarField> fe2 = 1 - max(ft(magGradU), fl(magGradU));
143 tmp<volScalarField> fe = max(fe1 - 1, scalar(0))*psi*fe2;
144
145 const volScalarField fdTilda(max(1 - fdt(magGradU), fB));
146
147 // Simplified formulation from Gritskevich et al. paper (2011) where fe = 0
148 // return max
149 // (
150 // fdTilda*lRAS + (1 - fdTilda)*lLES,
151 // dimensionedScalar("SMALL", dimLength, SMALL)
152 // );
153
154 // Original formulation from Shur et al. paper (2008)
155 return max
156 (
157 fdTilda*(1 + fe)*lRAS + (1 - fdTilda)*lLES,
158 dimensionedScalar("SMALL", dimLength, SMALL)
159 );
160}
161
162
163// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
164
165template<class BasicTurbulenceModel>
167(
168 const alphaField& alpha,
169 const rhoField& rho,
170 const volVectorField& U,
171 const surfaceScalarField& alphaRhoPhi,
172 const surfaceScalarField& phi,
173 const transportModel& transport,
174 const word& propertiesName,
175 const word& type
176)
177:
178 SpalartAllmarasDES<BasicTurbulenceModel>
179 (
180 alpha,
181 rho,
182 U,
183 alphaRhoPhi,
184 phi,
185 transport,
186 propertiesName,
187 type
188 ),
189
190 Cdt1_
191 (
192 dimensioned<scalar>::getOrAddToDict
193 (
194 "Cdt1",
195 this->coeffDict_,
196 8
197 )
198 ),
199 Cdt2_
200 (
201 dimensioned<scalar>::getOrAddToDict
202 (
203 "Cdt2",
204 this->coeffDict_,
205 3
206 )
207 ),
208 Cl_
209 (
210 dimensioned<scalar>::getOrAddToDict
211 (
212 "Cl",
213 this->coeffDict_,
214 3.55
215 )
216 ),
217 Ct_
218 (
219 dimensioned<scalar>::getOrAddToDict
220 (
221 "Ct",
222 this->coeffDict_,
223 1.63
224 )
225 ),
226 IDDESDelta_(setDelta())
227{
228 if (type == typeName)
229 {
230 this->printCoeffs(type);
231 }
232}
233
234
235// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
236
237template<class BasicTurbulenceModel>
239{
241 {
242 Cdt1_.readIfPresent(this->coeffDict());
243 Cdt2_.readIfPresent(this->coeffDict());
244 Cl_.readIfPresent(this->coeffDict());
245 Ct_.readIfPresent(this->coeffDict());
246
247 return true;
248 }
249
250 return false;
251}
252
253
254// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
255
256} // End namespace LESModels
257} // End namespace Foam
258
259// ************************************************************************* //
scalar delta
surfaceScalarField & phi
IDDESDelta used by the IDDES (improved low Re Spalart-Allmaras DES model) The min and max delta are c...
SpalartAllmarasDES DES turbulence model for incompressible and compressible flows.
SpalartAllmaras IDDES turbulence model for incompressible and compressible flows.
BasicTurbulenceModel::alphaField alphaField
BasicTurbulenceModel::rhoField rhoField
BasicTurbulenceModel::transportModel transportModel
virtual tmp< volScalarField > dTilda(const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU) const
Length scale.
virtual bool read()
Re-read model coefficients if they have changed.
Generic dimensioned Type class.
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
const volScalarField & psi
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
dimensionedScalar pos0(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
dimensionedScalar pow3(const dimensionedScalar &ds)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
dimensionedScalar tanh(const dimensionedScalar &ds)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
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:33
dimensionedScalar neg(const dimensionedScalar &ds)
error FatalError
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
volScalarField & nu
volScalarField & alpha
static const char *const typeName
The type name used in ensight case files.