SpalartAllmarasDES.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) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2015-2022 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::SpalartAllmarasDES
29 
30 Group
31  grpDESTurbulence
32 
33 Description
34  SpalartAllmarasDES DES turbulence model for incompressible and
35  compressible flows.
36 
37  Reference:
38  \verbatim
39  Spalart, P. R., Jou, W. H., Strelets, M., & Allmaras, S. R. (1997).
40  Comments on the feasibility of LES for wings, and on a hybrid
41  RANS/LES approach.
42  Advances in DNS/LES, 1, 4-8.
43  \endverbatim
44 
45  Including the low Reynolds number correction documented in
46  \verbatim
47  Spalart, P. R., Deck, S., Shur, M. L., Squires,
48  K. D., Strelets, M. K., & Travin, A. (2006).
49  A new version of detached-eddy simulation,
50  resistant to ambiguous grid densities.
51  Theoretical and computational fluid dynamics, 20(3), 181-195.
52  DOI:10.1007/s00162-006-0015-0
53  \endverbatim
54 
55 Note
56  The default value of the DES constant implemented was calibrated for
57  OpenFOAM using decaying isotropic turbulence and agrees with the value
58  suggested in the reference publication.
59 
60 SourceFiles
61  SpalartAllmarasDES.C
62 
63 \*---------------------------------------------------------------------------*/
64 
65 #ifndef Foam_SpalartAllmarasDES_H
66 #define Foam_SpalartAllmarasDES_H
67 
68 #include "SpalartAllmarasBase.H"
69 #include "DESModel.H"
70 
71 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
72 
73 namespace Foam
74 {
75 namespace LESModels
76 {
77 
78 /*---------------------------------------------------------------------------*\
79  Class SpalartAllmarasDES Declaration
80 \*---------------------------------------------------------------------------*/
81 
82 template<class BasicTurbulenceModel>
84 :
85  public SpalartAllmarasBase<DESModel<BasicTurbulenceModel>>
86 {
87  // Private Member Functions
88 
89  //- No copy construct
90  SpalartAllmarasDES(const SpalartAllmarasDES&) = delete;
91 
92  //- No copy assignment
93  void operator=(const SpalartAllmarasDES&) = delete;
94 
95 
96 protected:
97 
98  // Protected Data
99 
100  //- Switch to activate grey-area enhanced sigma-(D)DES
103  // Model constants
104 
105  //- DES coefficient
107 
108  //- Flag for low Reynolds number correction
111 
112 
113  // Protected Member Functions
115  //- Return the low Reynolds number correction function
116  virtual tmp<volScalarField> psi
117  (
118  const volScalarField& chi,
119  const volScalarField& fv1
120  ) const;
121 
122  //- Return the LES length scale
124  (
125  const volScalarField& chi,
126  const volScalarField& fv1
127  ) const;
128 
129  //- Return the production term
131  (
132  const volScalarField& chi,
133  const volScalarField& fv1,
134  const volTensorField& gradU,
135  const volScalarField& dTilda
136  ) const;
137 
138  //- Return the length scale
140  (
141  const volScalarField& chi,
142  const volScalarField& fv1,
143  const volTensorField& gradU
144  ) const;
145 
146  virtual void correctNut();
147 
148 public:
149 
150  typedef typename BasicTurbulenceModel::alphaField alphaField;
151  typedef typename BasicTurbulenceModel::rhoField rhoField;
152  typedef typename BasicTurbulenceModel::transportModel transportModel;
153 
154 
155  //- Runtime type information
156  TypeName("SpalartAllmarasDES");
157 
158 
159  // Constructors
160 
161  //- Construct from components
163  (
165  const rhoField& rho,
166  const volVectorField& U,
167  const surfaceScalarField& alphaRhoPhi,
168  const surfaceScalarField& phi,
169  const transportModel& transport,
170  const word& propertiesName = turbulenceModel::propertiesName,
171  const word& type = typeName
172  );
173 
174 
175  //- Destructor
176  virtual ~SpalartAllmarasDES() = default;
177 
178 
179  // Member Functions
180 
181  //- Re-read model coefficients if they have changed
182  virtual bool read();
183 
184  //- Return the LES field indicator
185  virtual tmp<volScalarField> LESRegion() const;
186 };
187 
188 
189 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
190 
191 } // End namespace LESModels
192 } // End namespace Foam
193 
194 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
195 
196 #ifdef NoRepository
197  #include "SpalartAllmarasDES.C"
198 #endif
199 
200 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
201 
202 #endif
203 
204 // ************************************************************************* //
TypeName("SpalartAllmarasDES")
Runtime type information.
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 ~SpalartAllmarasDES()=default
Destructor.
virtual bool read()
Re-read model coefficients if they have changed.
dimensionedScalar CDES_
DES coefficient.
SpalartAllmarasDES DES turbulence model for incompressible and compressible flows.
Switch useSigma_
Switch to activate grey-area enhanced sigma-(D)DES.
BasicTurbulenceModel::transportModel transportModel
Base class to handle various characteristics for SpalartAllmaras based LES/DES turbulence models for ...
virtual tmp< volScalarField > LESRegion() const
Return the LES field indicator.
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 const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:63
virtual tmp< volScalarField > psi(const volScalarField &chi, const volScalarField &fv1) const
Return the low Reynolds number correction function.
BasicTurbulenceModel::rhoField rhoField
U
Definition: pEqn.H:72
virtual tmp< volScalarField > lengthScaleLES(const volScalarField &chi, const volScalarField &fv1) const
Return the LES length scale.
BasicTurbulenceModel::alphaField alphaField
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Switch lowReCorrection_
Flag for low Reynolds number correction.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
virtual tmp< volScalarField > dTilda(const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU) const
Return the length scale.
Namespace for OpenFOAM.