SpalartAllmarasDDES.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) 2015-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 "SpalartAllmarasDDES.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace LESModels
37 {
38 
39 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
40 
41 template<class BasicTurbulenceModel>
43 (
44  const volScalarField& magGradU
45 ) const
46 {
47  return
48  1
49  - tanh(pow(this->Cd1_*this->r(this->nuEff(), magGradU, this->y_), Cd2_));
50 }
51 
52 
53 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
54 
55 template<class BasicTurbulenceModel>
57 (
58  const volScalarField& chi,
59  const volScalarField& fv1,
60  const volTensorField& gradU,
61  const volScalarField& dTilda
62 ) const
63 {
64  if (this->useSigma_)
65  {
66  const volScalarField& lRAS(this->y_);
67  const volScalarField fv2(this->fv2(chi, fv1));
68  const volScalarField lLES(this->lengthScaleLES(chi, fv1));
69  const volScalarField Omega(this->Omega(gradU));
70  const volScalarField Ssigma(this->Ssigma(gradU));
71  const volScalarField SsigmaDES
72  (
73  Omega - fd(mag(gradU))*pos(lRAS - lLES)*(Omega - Ssigma)
74  );
75 
76  return
77  max
78  (
79  SsigmaDES + fv2*this->nuTilda_/sqr(this->kappa_*dTilda),
80  this->Cs_*SsigmaDES
81  );
82  }
83 
84  return
86  (
87  chi,
88  fv1,
89  gradU,
90  dTilda
91  );
92 }
93 
94 
95 template<class BasicTurbulenceModel>
97 (
98  const volScalarField& chi,
99  const volScalarField& fv1,
100  const volTensorField& gradU
101 ) const
102 {
103  const volScalarField& lRAS(this->y_);
104  const volScalarField lLES(this->lengthScaleLES(chi, fv1));
105  const dimensionedScalar l0(dimLength, Zero);
106 
107  return max
108  (
109  lRAS - fd(mag(gradU))*max(lRAS - lLES, l0),
110  dimensionedScalar("small", dimLength, SMALL)
111  );
112 }
113 
114 
115 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
116 
117 template<class BasicTurbulenceModel>
119 (
120  const alphaField& alpha,
121  const rhoField& rho,
122  const volVectorField& U,
123  const surfaceScalarField& alphaRhoPhi,
124  const surfaceScalarField& phi,
125  const transportModel& transport,
126  const word& propertiesName,
127  const word& type
128 )
129 :
130  SpalartAllmarasDES<BasicTurbulenceModel>
131  (
132  alpha,
133  rho,
134  U,
135  alphaRhoPhi,
136  phi,
137  transport,
138  propertiesName,
139  type
140  ),
141 
142  Cd1_
143  (
144  this->useSigma_
145  ? dimensioned<scalar>::getOrAddToDict
146  (
147  "Cd1Sigma",
148  this->coeffDict_,
149  10
150  )
151  : dimensioned<scalar>::getOrAddToDict
152  (
153  "Cd1",
154  this->coeffDict_,
155  8
156  )
157  ),
158  Cd2_
159  (
160  dimensioned<scalar>::getOrAddToDict
161  (
162  "Cd2",
163  this->coeffDict_,
164  3
165  )
166  )
167 {
168  if (type == typeName)
169  {
170  this->printCoeffs(type);
171  }
172 }
173 
174 
175 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
176 
177 template<class BasicTurbulenceModel>
179 {
181  {
182  Cd1_.readIfPresent(this->coeffDict());
183  Cd2_.readIfPresent(this->coeffDict());
184 
185  return true;
186  }
187 
188  return false;
189 }
190 
191 
192 template<class BasicTurbulenceModel>
194 {
195  return fd(mag(fvc::grad(this->U_)));
196 }
197 
198 
199 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
200 
201 } // End namespace LESModels
202 } // End namespace Foam
203 
204 // ************************************************************************* //
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
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
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 bool read()
Re-read model coefficients if they have changed.
SpalartAllmarasDES DES turbulence model for incompressible and compressible flows.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
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
A class for handling words, derived from Foam::string.
Definition: word.H:63
virtual tmp< volScalarField > dTilda(const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU) const
Return the length scale.
SpalartAllmaras DDES turbulence model for incompressible and compressible flows.
virtual tmp< volScalarField > Stilda(const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU, const volScalarField &dTilda) const
Return the production term.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
U
Definition: pEqn.H:72
virtual bool read()
Read from dictionary.
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. ...
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
virtual tmp< volScalarField > fd() const
Return the shielding function.
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127