kOmegaSSTDES.H
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) 2016-2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Class
28  Foam::LESModels::kOmegaSSTDES
29 
30 Group
31  grpDESTurbulence
32 
33 Description
34  k-omega-SST DES turbulence model for incompressible and compressible flows.
35 
36  Reference:
37  \verbatim
38  Strelets, M. (2001).
39  Detached Eddy Simulation of Massively Separated Flows.
40  39th AIAA Aerospace Sciences Meeting and Exhibit, Reno, NV.
41  \endverbatim
42 
43 Note
44  The default values of the DES constants implemented are code-specific
45  values calibrated for OpenFOAM using decaying isotropic turbulence, and
46  hence deviate slightly from the values suggested in the reference
47  publication.
48 
49 SourceFiles
50  kOmegaSSTDES.C
51 
52 \*---------------------------------------------------------------------------*/
53 
54 #ifndef Foam_kOmegaSSTDES_H
55 #define Foam_kOmegaSSTDES_H
56 
57 #include "DESModel.H"
58 #include "kOmegaSSTBase.H"
59 
60 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
61 
62 namespace Foam
63 {
64 namespace LESModels
65 {
66 
67 /*---------------------------------------------------------------------------*\
68  Class kOmegaSSTDES Declaration
69 \*---------------------------------------------------------------------------*/
70 
71 template<class BasicTurbulenceModel>
72 class kOmegaSSTDES
73 :
74  public kOmegaSSTBase<DESModel<BasicTurbulenceModel>>
75 {
76  // Private Member Functions
77 
78  //- No copy construct
79  kOmegaSSTDES(const kOmegaSSTDES&) = delete;
80 
81  //- No copy assignment
82  void operator=(const kOmegaSSTDES&) = delete;
83 
84 
85 protected:
86 
87  // Protected Data
88 
89  //- Switch to activate grey-area enhanced sigma-(D)DES
91 
92  // Model coefficients
93 
95 
96  //- DES coefficients
99 
101  // Protected Member Functions
102 
103  //- Blending for CDES parameter
104  virtual tmp<volScalarField> CDES(const volScalarField& F1) const
105  {
106  return this->blend(F1, CDESkom_, CDESkeps_);
107  }
108 
109  virtual void correctNut(const volScalarField& S2);
110  virtual void correctNut();
111 
113  (
114  const volScalarField& nur,
115  const volScalarField& magGradU
116  ) const;
117 
118  //- Return square of strain rate
119  virtual tmp<volScalarField> S2
120  (
121  const volTensorField& gradU
122  ) const;
123 
124  //- Return length scale
126  (
127  const volScalarField& magGradU,
128  const volScalarField& CDES
129  ) const;
130 
131  //- Return epsilon/k
133  (
134  const volScalarField& F1,
135  const volTensorField& gradU
136  ) const;
137 
138  //- Return (G/nu)_0
140  (
141  const volTensorField& gradU,
142  const volScalarField& S2
143  ) const;
144 
145  //- Return G/nu
147  (
151  ) const;
152 
153 
154 public:
155 
156  typedef typename BasicTurbulenceModel::alphaField alphaField;
157  typedef typename BasicTurbulenceModel::rhoField rhoField;
158  typedef typename BasicTurbulenceModel::transportModel transportModel;
159 
160 
161  //- Runtime type information
162  TypeName("kOmegaSSTDES");
163 
164 
165  // Constructors
166 
167  //- Construct from components
169  (
170  const alphaField& alpha,
171  const rhoField& rho,
173  const surfaceScalarField& alphaRhoPhi,
174  const surfaceScalarField& phi,
175  const transportModel& transport,
176  const word& propertiesName = turbulenceModel::propertiesName,
177  const word& type = typeName
178  );
179 
180 
181  //- Destructor
182  virtual ~kOmegaSSTDES() = default;
183 
184 
185  // Member Functions
186 
187  //- Re-read model coefficients if they have changed
188  virtual bool read();
189 
190  //- RAS length scale
191  virtual tmp<volScalarField> lengthScaleRAS() const;
192 
193  //- LES length scale
195  (
196  const volScalarField& CDES
197  ) const;
198 
199  //- Return the LES field indicator
200  virtual tmp<volScalarField> LESRegion() const;
201 };
202 
203 
204 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
205 
206 } // End namespace LESModels
207 } // End namespace Foam
208 
209 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
210 #ifdef NoRepository
211  #include "kOmegaSSTDES.C"
212 #endif
213 
214 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
215 #endif
216 
217 // ************************************************************************* //
tmp< volScalarField > blend(const volScalarField &F1, const dimensionedScalar &psi1, const dimensionedScalar &psi2) const
Return the blended field.
Switch useSigma_
Switch to activate grey-area enhanced sigma-(D)DES.
Definition: kOmegaSSTDES.H:91
dimensionedScalar CDESkom_
DES coefficients.
Definition: kOmegaSSTDES.H:100
TypeName("kOmegaSSTDES")
Runtime type information.
virtual tmp< volScalarField::Internal > GbyNu0(const volTensorField &gradU, const volScalarField &S2) const
Return (G/nu)_0.
Definition: kOmegaSSTDES.C:132
virtual tmp< volScalarField > CDES(const volScalarField &F1) const
Blending for CDES parameter.
Definition: kOmegaSSTDES.H:109
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
Base class implementation of the k-omega-SST turbulence model for incompressible and compressible flo...
BasicTurbulenceModel::rhoField rhoField
Definition: kOmegaSSTDES.H:172
virtual tmp< volScalarField > dTilda(const volScalarField &magGradU, const volScalarField &CDES) const
Return length scale.
Definition: kOmegaSSTDES.C:109
BasicTurbulenceModel::transportModel transportModel
Definition: kOmegaSSTDES.H:173
dimensionedScalar kappa_
Definition: kOmegaSSTDES.H:95
virtual tmp< volScalarField::Internal > epsilonByk(const volScalarField &F1, const volTensorField &gradU) const
Return epsilon/k.
Definition: kOmegaSSTDES.C:120
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
virtual tmp< volScalarField::Internal > GbyNu(const volScalarField::Internal &GbyNu0, const volScalarField::Internal &F2, const volScalarField::Internal &S2) const
Return G/nu.
Definition: kOmegaSSTDES.C:149
virtual tmp< volScalarField > LESRegion() const
Return the LES field indicator.
Definition: kOmegaSSTDES.C:281
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:63
tmp< volScalarField > r(const volScalarField &nur, const volScalarField &magGradU) const
Definition: kOmegaSSTDES.C:54
virtual tmp< volScalarField > lengthScaleRAS() const
RAS length scale.
Definition: kOmegaSSTDES.C:260
virtual tmp< volScalarField > F2() const
Definition: kOmegaSSTBase.C:66
virtual ~kOmegaSSTDES()=default
Destructor.
dimensionedScalar CDESkeps_
Definition: kOmegaSSTDES.H:101
k-omega-SST DES turbulence model for incompressible and compressible flows.
Definition: kOmegaSSTDES.H:67
U
Definition: pEqn.H:72
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 > F1(const volScalarField &CDkOmega) const
Definition: kOmegaSSTBase.C:37
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
virtual tmp< volScalarField > lengthScaleLES(const volScalarField &CDES) const
LES length scale.
Definition: kOmegaSSTDES.C:272
virtual tmp< volScalarField > S2(const volTensorField &gradU) const
Return square of strain rate.
Definition: kOmegaSSTDES.C:72
Namespace for OpenFOAM.
virtual bool read()
Re-read model coefficients if they have changed.
Definition: kOmegaSSTDES.C:242
BasicTurbulenceModel::alphaField alphaField
Definition: kOmegaSSTDES.H:171