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-2018 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 BasicTurbulenceModel>
38 tmp<volScalarField> SpalartAllmarasIDDES<BasicTurbulenceModel>::alpha() const
39 {
40  return max
41  (
42  0.25 - this->y_/static_cast<const volScalarField&>(IDDESDelta_.hmax()),
43  scalar(-5)
44  );
45 }
46 
47 
48 template<class BasicTurbulenceModel>
49 tmp<volScalarField> SpalartAllmarasIDDES<BasicTurbulenceModel>::ft
50 (
51  const volScalarField& magGradU
52 ) const
53 {
54  return tanh(pow3(sqr(ct_)*rd(this->nut_, magGradU)));
55 }
56 
57 
58 template<class BasicTurbulenceModel>
59 tmp<volScalarField> SpalartAllmarasIDDES<BasicTurbulenceModel>::fl
60 (
61  const volScalarField& magGradU
62 ) const
63 {
64  return tanh(pow(sqr(cl_)*rd(this->nu(), magGradU), 10));
65 }
66 
67 
68 template<class BasicTurbulenceModel>
69 tmp<volScalarField> SpalartAllmarasIDDES<BasicTurbulenceModel>::rd
70 (
71  const volScalarField& nur,
72  const volScalarField& magGradU
73 ) const
74 {
75  return min
76  (
77  nur
78  /(
79  max
80  (
81  magGradU,
82  dimensionedScalar(magGradU.dimensions(), small)
83  )*sqr(this->kappa_*this->y_)
84  ),
85  scalar(10)
86  );
87 }
88 
89 
90 template<class BasicTurbulenceModel>
91 tmp<volScalarField> SpalartAllmarasIDDES<BasicTurbulenceModel>::fd
92 (
93  const volScalarField& magGradU
94 ) const
95 {
96  return 1 - tanh(pow3(8*rd(this->nuEff(), magGradU)));
97 }
98 
99 
100 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
101 
102 template<class BasicTurbulenceModel>
104 (
105  const volScalarField& chi,
106  const volScalarField& fv1,
107  const volTensorField& gradU
108 ) const
109 {
110  const volScalarField alpha(this->alpha());
111  const volScalarField expTerm(exp(sqr(alpha)));
112  const volScalarField magGradU(mag(gradU));
113 
114  tmp<volScalarField> fHill =
115  2*(pos0(alpha)*pow(expTerm, -11.09) + neg(alpha)*pow(expTerm, -9.0));
116 
117  tmp<volScalarField> fStep = min(2*pow(expTerm, -9.0), scalar(1));
118  const volScalarField fHyb(max(1 - fd(magGradU), fStep));
119  tmp<volScalarField> fAmp = 1 - max(ft(magGradU), fl(magGradU));
120  tmp<volScalarField> fRestore = max(fHill - 1, scalar(0))*fAmp;
121 
122  // IGNORING ft2 terms
123  const volScalarField Psi
124  (
125  sqrt
126  (
127  min
128  (
129  scalar(100),
130  (
131  1
132  - this->Cb1_*this->fv2(chi, fv1)
133  /(this->Cw1_*sqr(this->kappa_)*fwStar_)
134  )/max(small, fv1)
135  )
136  )
137  );
138 
139  return max
140  (
142  fHyb*(1 + fRestore*Psi)*this->y_
143  + (1 - fHyb)*this->CDES_*Psi*this->delta()
144  );
145 }
146 
147 
148 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
149 
150 template<class BasicTurbulenceModel>
152 (
153  const alphaField& alpha,
154  const rhoField& rho,
155  const volVectorField& U,
156  const surfaceScalarField& alphaRhoPhi,
157  const surfaceScalarField& phi,
158  const transportModel& transport,
159  const word& propertiesName,
160  const word& type
161 )
162 :
164  (
165  alpha,
166  rho,
167  U,
168  alphaRhoPhi,
169  phi,
170  transport,
171  propertiesName
172  ),
173  fwStar_
174  (
176  (
177  "fwStar",
178  this->coeffDict_,
179  0.424
180  )
181  ),
182  cl_
183  (
185  (
186  "cl",
187  this->coeffDict_,
188  3.55
189  )
190  ),
191  ct_
192  (
194  (
195  "ct",
196  this->coeffDict_,
197  1.63
198  )
199  ),
200  IDDESDelta_(refCast<IDDESDelta>(this->delta_()))
201 {}
202 
203 
204 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
205 
206 template<class BasicTurbulenceModel>
208 {
210  {
211  fwStar_.readIfPresent(this->coeffDict());
212  cl_.readIfPresent(this->coeffDict());
213  ct_.readIfPresent(this->coeffDict());
214 
215  return true;
216  }
217  else
218  {
219  return false;
220  }
221 }
222 
223 
224 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
225 
226 } // End namespace LESModels
227 } // End namespace Foam
228 
229 // ************************************************************************* //
scalar delta
dimensionedScalar tanh(const dimensionedScalar &ds)
surfaceScalarField & phi
BasicTurbulenceModel::transportModel transportModel
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual tmp< volScalarField > dTilda(const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU) const
Length scale.
dimensionedScalar sqrt(const dimensionedScalar &ds)
SpalartAllmarasDES DES turbulence model for incompressible and compressible flows.
Generic dimensioned Type class.
virtual bool read()
Read model coefficients if they have changed.
dimensionedScalar neg(const dimensionedScalar &ds)
BasicTurbulenceModel::alphaField alphaField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dimensionedScalar exp(const dimensionedScalar &ds)
SpalartAllmarasIDDES(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName)
Construct from components.
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)
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.
BasicTurbulenceModel::rhoField rhoField
dimensioned< scalar > mag(const dimensioned< Type > &)
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
volScalarField & nu
Namespace for OpenFOAM.