SpalartAllmarasIDDES.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2020 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "SpalartAllmarasIDDES.H"
27 
28 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32 namespace LESModels
33 {
34 
35 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
36 
37 template<class BasicMomentumTransportModel>
38 tmp<volScalarField::Internal>
39 SpalartAllmarasIDDES<BasicMomentumTransportModel>::IDDESalpha() const
40 {
42  (
43  modelName("alpha"),
44  max(0.25 - this->y_()/IDDESDelta_.hmax(), scalar(-5))
45  );
46 }
47 
48 
49 template<class BasicMomentumTransportModel>
50 tmp<volScalarField::Internal>
51 SpalartAllmarasIDDES<BasicMomentumTransportModel>::ft
52 (
53  const volScalarField::Internal& magGradU
54 ) const
55 {
57  (
58  modelName("ft"),
59  tanh(pow3(sqr(ct_)*rd(this->nut_, magGradU)))
60  );
61 }
62 
63 
64 template<class BasicMomentumTransportModel>
65 tmp<volScalarField::Internal>
66 SpalartAllmarasIDDES<BasicMomentumTransportModel>::fl
67 (
68  const volScalarField::Internal& magGradU
69 ) const
70 {
72  (
73  modelName("fl"),
74  tanh(pow(sqr(cl_)*rd(this->nu(), magGradU), 10))
75  );
76 }
77 
78 
79 template<class BasicMomentumTransportModel>
80 tmp<volScalarField::Internal>
81 SpalartAllmarasIDDES<BasicMomentumTransportModel>::rd
82 (
83  const volScalarField::Internal& nur,
84  const volScalarField::Internal& magGradU
85 ) const
86 {
88  (
89  modelName("rd"),
90  min
91  (
92  nur
93  /(
94  max
95  (
96  magGradU,
97  dimensionedScalar(magGradU.dimensions(), small)
98  )*sqr(this->kappa_*this->y_())
99  ),
100  scalar(10)
101  )
102  );
103 }
104 
105 
106 template<class BasicMomentumTransportModel>
107 tmp<volScalarField::Internal>
108 SpalartAllmarasIDDES<BasicMomentumTransportModel>::fd
109 (
110  const volScalarField::Internal& magGradU
111 ) const
112 {
114  (
115  modelName("fd"),
116  1 - tanh(pow3(8*rd(this->nuEff(), magGradU)))
117  );
118 }
119 
120 
121 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
122 
123 template<class BasicMomentumTransportModel>
124 tmp<volScalarField::Internal>
126 (
127  const volScalarField::Internal& chi,
128  const volScalarField::Internal& fv1,
129  const volTensorField::Internal& gradU
130 ) const
131 {
132  const volScalarField::Internal alpha(IDDESalpha());
133 
134  const volScalarField::Internal expTerm
135  (
136  modelName("expTerm"),
137  exp(sqr(alpha))
138  );
139 
140  const volScalarField::Internal magGradU(modelName("magGradU"), mag(gradU));
141 
143  (
145  (
146  modelName("fHill"),
147  2*(pos0(alpha)*pow(expTerm, -11.09) + neg(alpha)*pow(expTerm, -9.0))
148  )
149  );
150 
152  (
154  (
155  modelName("fStep"),
156  min(2*pow(expTerm, -9.0), scalar(1))
157  )
158  );
159 
160  const volScalarField::Internal fHyb
161  (
162  modelName("fHyb"),
163  max(1 - fd(magGradU), fStep)
164  );
165 
167  (
169  (
170  modelName("fAmp"),
171  1 - max(ft(magGradU), fl(magGradU))
172  )
173  );
174 
176  (
178  (
179  modelName("fRestore"),
180  max(fHill - 1, scalar(0))*fAmp
181  )
182  );
183 
184  // IGNORING ft2 terms
185  const volScalarField::Internal Psi
186  (
187  modelName("Psi"),
188  sqrt
189  (
190  min
191  (
192  scalar(100),
193  (
194  1
195  - this->Cb1_*this->fv2(chi, fv1)
196  /(this->Cw1_*sqr(this->kappa_)*fwStar_)
197  )/max(small, fv1)
198  )
199  )
200  );
201 
203  (
204  modelName("dTilda"),
205  max
206  (
208  fHyb*(1 + fRestore*Psi)*this->y_
209  + (1 - fHyb)*this->CDES_*Psi*this->delta()
210  )
211  );
212 }
213 
214 
215 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
216 
217 template<class BasicMomentumTransportModel>
219 (
220  const alphaField& alpha,
221  const rhoField& rho,
222  const volVectorField& U,
223  const surfaceScalarField& alphaRhoPhi,
224  const surfaceScalarField& phi,
225  const transportModel& transport,
226  const word& type
227 )
228 :
230  (
231  alpha,
232  rho,
233  U,
234  alphaRhoPhi,
235  phi,
236  transport
237  ),
238  fwStar_
239  (
241  (
242  "fwStar",
243  this->coeffDict_,
244  0.424
245  )
246  ),
247  cl_
248  (
250  (
251  "cl",
252  this->coeffDict_,
253  3.55
254  )
255  ),
256  ct_
257  (
259  (
260  "ct",
261  this->coeffDict_,
262  1.63
263  )
264  ),
265  IDDESDelta_(refCast<IDDESDelta>(this->delta_()))
266 {}
267 
268 
269 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
270 
271 template<class BasicMomentumTransportModel>
273 {
275  {
276  fwStar_.readIfPresent(this->coeffDict());
277  cl_.readIfPresent(this->coeffDict());
278  ct_.readIfPresent(this->coeffDict());
279 
280  return true;
281  }
282  else
283  {
284  return false;
285  }
286 }
287 
288 
289 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
290 
291 } // End namespace LESModels
292 } // End namespace Foam
293 
294 // ************************************************************************* //
scalar delta
dimensionedScalar tanh(const dimensionedScalar &ds)
BasicMomentumTransportModel::rhoField rhoField
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &)
Return a temporary field constructed from name, mesh.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
SpalartAllmarasDES DES turbulence model for incompressible and compressible flows.
Generic dimensioned Type class.
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
virtual bool read()
Read model coefficients if they have changed.
dimensionedScalar neg(const dimensionedScalar &ds)
phi
Definition: pEqn.H:104
dimensionedScalar exp(const dimensionedScalar &ds)
A class for handling words, derived from string.
Definition: word.H:59
dimensionedScalar pos0(const dimensionedScalar &ds)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
BasicMomentumTransportModel::alphaField alphaField
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.
BasicMomentumTransportModel::transportModel transportModel
SpalartAllmarasIDDES(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &type=typeName)
Construct from components.
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual tmp< volScalarField::Internal > dTilda(const volScalarField::Internal &chi, const volScalarField::Internal &fv1, const volTensorField::Internal &gradU) const
Length scale.
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.