SpalartAllmarasDES.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) 2016-2022 OpenCFD Ltd.
11 -------------------------------------------------------------------------------
12 License
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 "SpalartAllmarasDES.H"
31 #include "fvOptions.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace LESModels
38 {
39 
40 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
41 
42 template<class BasicTurbulenceModel>
44 (
45  const volScalarField& chi,
46  const volScalarField& fv1
47 ) const
48 {
49  auto tpsi = tmp<volScalarField>::New
50  (
51  IOobject
52  (
53  IOobject::scopedName(type(), "psi"),
54  this->time().timeName(),
55  this->mesh(),
58  ),
59  this->mesh(),
61  );
62 
63  if (lowReCorrection_)
64  {
65  auto& psi = tpsi.ref();
66 
67  const volScalarField fv2(this->fv2(chi, fv1));
68  const volScalarField ft2(this->ft2(chi));
69 
70  psi =
71  sqrt
72  (
73  min
74  (
75  scalar(100),
76  (1 - this->Cb1_/(this->Cw1_*sqr(this->kappa_)*fwStar_)
77  *(ft2 + (1 - ft2)*fv2))
78  /max(SMALL, (fv1*max(scalar(1e-10), 1 - ft2)))
79  )
80  );
81  }
82 
83  return tpsi;
84 }
85 
86 
87 template<class BasicTurbulenceModel>
89 (
90  const volScalarField& chi,
91  const volScalarField& fv1
92 ) const
93 {
94  return psi(chi, fv1)*CDES_*this->delta();
95 }
96 
97 
98 template<class BasicTurbulenceModel>
100 (
101  const volScalarField& chi,
102  const volScalarField& fv1,
103  const volTensorField& gradU,
104  const volScalarField& dTilda
105 ) const
106 {
107  if (useSigma_)
108  {
109  const volScalarField& lRAS = this->y_;
110  const volScalarField fv2(this->fv2(chi, fv1));
111  const volScalarField lLES(this->lengthScaleLES(chi, fv1));
112  const volScalarField Omega(this->Omega(gradU));
113  const volScalarField Ssigma(this->Ssigma(gradU));
114  const volScalarField SsigmaDES
115  (
116  pos(dTilda - lRAS)*Omega + (scalar(1) - pos(dTilda - lRAS))*Ssigma
117  );
118 
119  return
120  max
121  (
122  SsigmaDES + fv2*this->nuTilda_/sqr(this->kappa_*dTilda),
123  this->Cs_*SsigmaDES
124  );
125  }
126 
127  return
129  (
130  chi,
131  fv1,
132  gradU,
133  dTilda
134  );
135 }
136 
137 
138 template<class BasicTurbulenceModel>
140 (
141  const volScalarField& chi,
142  const volScalarField& fv1,
143  const volTensorField& gradU
144 ) const
145 {
146  // Initialise with LES length scale
147  tmp<volScalarField> tdTilda = lengthScaleLES(chi, fv1);
148 
149  // Take min vs wall distance
150  min(tdTilda.ref(), tdTilda(), this->y_);
151  return tdTilda;
152 }
153 
154 
155 template<class BasicTurbulenceModel>
157 {
158  // Correct the turbulence viscosity
160 
161  // Correct the turbulence thermal diffusivity
162  BasicTurbulenceModel::correctNut();
163 }
164 
165 
166 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
167 
168 template<class BasicTurbulenceModel>
170 (
171  const alphaField& alpha,
172  const rhoField& rho,
173  const volVectorField& U,
174  const surfaceScalarField& alphaRhoPhi,
175  const surfaceScalarField& phi,
176  const transportModel& transport,
177  const word& propertiesName,
178  const word& type
179 )
180 :
181  SpalartAllmarasBase<DESModel<BasicTurbulenceModel>>
182  (
183  type,
184  alpha,
185  rho,
186  U,
187  alphaRhoPhi,
188  phi,
189  transport,
190  propertiesName
191  ),
192 
193  useSigma_
194  (
195  Switch::getOrAddToDict
196  (
197  "useSigma",
198  this->coeffDict_,
199  false
200  )
201  ),
202  CDES_
203  (
204  dimensioned<scalar>::getOrAddToDict
205  (
206  "CDES",
207  this->coeffDict_,
208  0.65
209  )
210  ),
211  lowReCorrection_
212  (
213  Switch::getOrAddToDict
214  (
215  "lowReCorrection",
216  this->coeffDict_,
217  true
218  )
219  ),
220  fwStar_
221  (
222  dimensioned<scalar>::getOrAddToDict
223  (
224  "fwStar",
225  this->coeffDict_,
226  0.424
227  )
228  )
229 {
230  // Note: Ctrans coeff is model specific; for S-A = 67.7
231  this->Ctrans_ =
232  dimensioned<scalar>::getOrAddToDict("Ctrans", this->coeffDict_, 67.7);
233 
234  if (type == typeName)
235  {
237  << "This model is not recommended and will be deprecated in future "
238  << "releases. Please consider using the DDES/IDDES versions instead"
239  << endl;
240 
241  this->printCoeffs(type);
242  }
243 }
244 
245 
246 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
247 
248 template<class BasicTurbulenceModel>
250 {
251  if (SpalartAllmarasBase<DESModel<BasicTurbulenceModel>>::read())
252  {
253  useSigma_.readIfPresent("useSigma", this->coeffDict());
254  CDES_.readIfPresent(this->coeffDict());
255  lowReCorrection_.readIfPresent("lowReCorrection", this->coeffDict());
256  fwStar_.readIfPresent(this->coeffDict());
257 
258  return true;
259  }
261  return false;
262 }
263 
264 
265 template<class BasicTurbulenceModel>
268 {
269  const volScalarField chi(this->chi());
270  const volScalarField fv1(this->fv1(chi));
271 
273  (
274  IOobject
275  (
276  IOobject::scopedName("DES", "LESRegion"),
277  this->mesh_.time().timeName(),
278  this->mesh_,
281  ),
282  neg(dTilda(chi, fv1, fvc::grad(this->U_)) - this->y_)
283  );
284 }
285 
286 
287 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
288 
289 } // End namespace LESModels
290 } // End namespace Foam
291 
292 // ************************************************************************* //
scalar delta
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition: tmpI.H:235
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)
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
virtual tmp< volScalarField > Stilda(const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU, const volScalarField &dTilda) const
Return the production term.
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.
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
dimensionedScalar neg(const dimensionedScalar &ds)
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:127
Base class to handle various characteristics for SpalartAllmaras based LES/DES turbulence models for ...
virtual tmp< volScalarField > LESRegion() const
Return the LES field indicator.
word timeName
Definition: getTimeIndex.H:3
dimensionedScalar pos(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: POSIX.C:799
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition: IOobjectI.H:40
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
A class for handling words, derived from Foam::string.
Definition: word.H:63
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
virtual tmp< volScalarField > psi(const volScalarField &chi, const volScalarField &fv1) const
Return the low Reynolds number correction function.
static dimensioned< Type > getOrAddToDict(const word &name, dictionary &dict, const dimensionSet &dims=dimless, const Type &deflt=Type(Zero))
Construct dimensioned from dictionary, with default value.
U
Definition: pEqn.H:72
#define WarningInFunction
Report a warning using Foam::Warning.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Base-class for all transport models used by the incompressible turbulence models. ...
Nothing to be read.
virtual tmp< volScalarField > lengthScaleLES(const volScalarField &chi, const volScalarField &fv1) const
Return the LES length scale.
const volScalarField & psi
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
virtual tmp< volScalarField > dTilda(const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU) const
Return the length scale.
Templated abstract base class for DES models.
Definition: DESModel.H:51
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127