kOmegaSSTDDES.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) 2015 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 "kOmegaSSTDDES.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 1 - tanh(pow(Cd1_*this->r(this->nuEff(), magGradU), Cd2_));
48 }
49 
50 
51 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
52 
53 template<class BasicTurbulenceModel>
55 (
56  const volTensorField& gradU
57 ) const
58 {
61 
62  if (this->useSigma_)
63  {
64  volScalarField& S2 = tS2.ref();
65 
66  const volScalarField& k = this->k_;
67  const volScalarField& omega = this->omega_;
68 
69  const volScalarField CDkOmega
70  (
71  (2*this->alphaOmega2_)*(fvc::grad(k) & fvc::grad(omega))/omega
72  );
73 
74  const volScalarField F1(this->F1(CDkOmega));
75 
76  const volScalarField CDES(this->CDES(F1));
77  const volScalarField Ssigma(this->Ssigma(gradU));
78 
79  S2 -=
80  fd(mag(gradU))
81  *pos(this->lengthScaleRAS() - this->lengthScaleLES(CDES))
82  *(S2 - sqr(Ssigma));
83  }
84 
85  return tS2;
86 }
87 
88 
89 template<class BasicTurbulenceModel>
91 (
92  const volScalarField& magGradU,
93  const volScalarField& CDES
94 ) const
95 {
96  const volScalarField lRAS(this->lengthScaleRAS());
97  const volScalarField lLES(this->lengthScaleLES(CDES));
98  const dimensionedScalar l0(dimLength, Zero);
99 
100  return max
101  (
102  lRAS - fd(magGradU)*max(lRAS - lLES, l0),
103  dimensionedScalar("small", dimLength, SMALL)
104  );
105 }
106 
107 template<class BasicTurbulenceModel>
109 (
110  const volTensorField& gradU,
111  const volScalarField& S2
112 ) const
113 {
114  if (this->useSigma_)
115  {
116  return S2();
117  }
118 
119  return
121 }
122 
123 
124 template<class BasicTurbulenceModel>
126 (
127  const volScalarField::Internal& GbyNu0,
129  const volScalarField::Internal& S2
130 ) const
131 {
132  return GbyNu0; // Unlimited
133 }
134 
135 
136 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
137 
138 template<class BasicTurbulenceModel>
140 (
141  const alphaField& alpha,
142  const rhoField& rho,
143  const volVectorField& U,
144  const surfaceScalarField& alphaRhoPhi,
145  const surfaceScalarField& phi,
146  const transportModel& transport,
147  const word& propertiesName,
148  const word& type
149 )
150 :
151  kOmegaSSTDES<BasicTurbulenceModel>
152  (
153  alpha,
154  rho,
155  U,
156  alphaRhoPhi,
157  phi,
158  transport,
159  propertiesName,
160  type
161  ),
162 
163  Cd1_
164  (
165  this->useSigma_
166  ? dimensioned<scalar>::getOrAddToDict
167  (
168  "Cd1Sigma",
169  this->coeffDict_,
170  22
171  )
172  : dimensioned<scalar>::getOrAddToDict
173  (
174  "Cd1",
175  this->coeffDict_,
176  20
177  )
178  ),
179  Cd2_
180  (
181  dimensioned<scalar>::getOrAddToDict
182  (
183  "Cd2",
184  this->coeffDict_,
185  3
186  )
187  )
188 {
189  if (type == typeName)
190  {
191  this->printCoeffs(type);
192  }
193 }
194 
195 
196 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
197 
198 template<class BasicTurbulenceModel>
200 {
202  {
203  Cd1_.readIfPresent(this->coeffDict());
204  Cd2_.readIfPresent(this->coeffDict());
205 
206  return true;
207  }
208 
209  return false;
210 }
211 
212 
213 template<class BasicTurbulenceModel>
215 {
216  return fd(mag(fvc::grad(this->U_)));
217 }
218 
219 
220 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
221 
222 } // End namespace LESModels
223 } // End namespace Foam
224 
225 // ************************************************************************* //
k-omega-SST DDES turbulence model for incompressible and compressible flows.
Definition: kOmegaSSTDDES.H:64
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)
#define F1(B, C, D)
Definition: SHA1.C:149
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)
virtual tmp< volScalarField > dTilda(const volScalarField &magGradU, const volScalarField &CDES) const
Length scale.
Definition: kOmegaSSTDDES.C:84
Base class implementation of the k-omega-SST turbulence model for incompressible and compressible flo...
#define F2(B, C, D)
Definition: SHA1.C:150
label k
Boltzmann constant.
virtual tmp< volScalarField::Internal > GbyNu0(const volTensorField &gradU, const volScalarField &S2) const
Return (G/nu)_0.
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 > fd() const
Return the shielding function.
virtual tmp< volScalarField::Internal > GbyNu(const volScalarField::Internal &GbyNu0, const volScalarField::Internal &F2, const volScalarField::Internal &S2) const
Return G/nu.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
k-omega-SST DES turbulence model for incompressible and compressible flows.
Definition: kOmegaSSTDES.H:67
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.
Base-class for all transport models used by the incompressible turbulence models. ...
virtual bool read()
Re-read model coefficients if they have changed.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
virtual tmp< volScalarField > S2(const volTensorField &gradU) const
Return square of strain rate.
Definition: kOmegaSSTDDES.C:48
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Namespace for OpenFOAM.
virtual bool read()
Re-read model coefficients if they have changed.
Definition: kOmegaSSTDES.C:242
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127