SpalartAllmarasBase.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::SpalartAllmarasBase
29 
30 Group
31  grpDESTurbulence
32 
33 Description
34  Base class to handle various characteristics for \c SpalartAllmaras based
35  LES/DES turbulence models for incompressible and compressible flows.
36 
37  References:
38  \verbatim
39  Standard model:
40  Spalart, P.R., & Allmaras, S.R. (1994).
41  A one-equation turbulence model for aerodynamic flows.
42  La Recherche Aerospatiale, 1, 5-21.
43 
44  Standard model:
45  Spalart, P. R., Jou, W. H., Strelets, M., & Allmaras, S. R. (1997).
46  Comments on the feasibility of LES for wings, and on a hybrid
47  RANS/LES approach.
48  Advances in DNS/LES, 1, 4-8.
49 
50  Estimation expression for k and epsilon (tag:B), Eq. 4.50:
51  Bourgoin, A. (2019).
52  Bathymetry induced turbulence modelling the
53  Alderney Race site: regional approach with TELEMAC-LES.
54  Normandie Université.
55 
56  Estimation expressions for omega (tag:P):
57  Pope, S. B. (2000).
58  Turbulent flows.
59  Cambridge, UK: Cambridge Univ. Press
60  DOI:10.1017/CBO9780511840531
61  \endverbatim
62 
63 SourceFiles
64  SpalartAllmarasBase.C
65 
66 \*---------------------------------------------------------------------------*/
67 
68 #ifndef Foam_SpalartAllmarasBase_H
69 #define Foam_SpalartAllmarasBase_H
70 
71 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
72 
73 namespace Foam
74 {
75 
76 /*---------------------------------------------------------------------------*\
77  Class SpalartAllmarasBase Declaration
78 \*---------------------------------------------------------------------------*/
79 
80 template<class BasicEddyViscosityModel>
82 :
83  public BasicEddyViscosityModel
84 {
85  // Private Member Functions
86 
87  //- No copy construct
89 
90  //- No copy assignment
91  void operator=(const SpalartAllmarasBase&) = delete;
92 
93 
94 protected:
95 
96  // Protected Data
97 
98  // Model constants
99 
115 
116 
117  // Fields
118 
119  //- Modified kinematic viscosity [m^2/s]
122  //- Wall distance
123  // Note: different to wall distance in parent RASModel
124  // which is for near-wall cells only
125  const volScalarField& y_;
126 
127 
128  // Protected Member Functions
130  tmp<volScalarField> chi() const;
131 
133 
135  (
136  const volScalarField& chi,
137  const volScalarField& fv1
138  ) const;
139 
141 
142  tmp<volScalarField> Omega(const volTensorField& gradU) const;
143 
145  (
146  const volScalarField& nur,
147  const volScalarField& Stilda,
148  const volScalarField& dTilda
149  ) const;
150 
152  (
153  const volScalarField& Stilda,
154  const volScalarField& dTilda
155  ) const;
156 
158  (
159  const volScalarField& chi,
160  const volScalarField& fv1,
161  const volTensorField& gradU,
162  const volScalarField& dTilda
163  ) const;
164 
165  //- Length scale
167  (
168  const volScalarField& chi,
169  const volScalarField& fv1,
170  const volTensorField& gradU
171  ) const = 0;
172 
173  void correctNut(const volScalarField& fv1);
174  virtual void correctNut();
175 
176 
177 public:
178 
179  typedef typename BasicEddyViscosityModel::alphaField alphaField;
180  typedef typename BasicEddyViscosityModel::rhoField rhoField;
181  typedef typename BasicEddyViscosityModel::transportModel transportModel;
182 
183 
184  // Constructors
186  //- Construct from components
188  (
189  const word& type,
190  const alphaField& alpha,
191  const rhoField& rho,
192  const volVectorField& U,
193  const surfaceScalarField& alphaRhoPhi,
194  const surfaceScalarField& phi,
195  const transportModel& transport,
196  const word& propertiesName = turbulenceModel::propertiesName
197  );
198 
199 
200  //- Destructor
201  virtual ~SpalartAllmarasBase() = default;
202 
203 
204  // Member Functions
205 
206  //- Re-read model coefficients if they have changed
207  virtual bool read();
208 
209  //- Return the effective diffusivity for nuTilda
211 
212  //- Return the (estimated) turbulent kinetic energy
213  virtual tmp<volScalarField> k() const;
214 
215  //- Return the (estimated) turbulent kinetic energy dissipation rate
216  virtual tmp<volScalarField> epsilon() const;
217 
218  //- Return the (estimated) specific dissipation rate
219  virtual tmp<volScalarField> omega() const;
220 
221  //- Return the modified kinematic viscosity
223  {
224  return nuTilda_;
225  }
226 
227  //- Correct nuTilda and related properties
228  virtual void correct();
229 };
230 
231 
232 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
233 
234 } // End namespace Foam
235 
236 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
237 
238 #ifdef NoRepository
239  #include "SpalartAllmarasBase.C"
240 #endif
241 
242 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
243 
244 #endif
245 
246 // ************************************************************************* //
tmp< volScalarField > chi() const
virtual tmp< volScalarField > omega() const
Return the (estimated) specific dissipation rate.
virtual ~SpalartAllmarasBase()=default
Destructor.
tmp< volScalarField > DnuTildaEff() const
Return the effective diffusivity for nuTilda.
tmp< volScalarField > fv1(const volScalarField &chi) const
tmp< volScalarField > Omega(const volTensorField &gradU) const
virtual tmp< volScalarField > dTilda(const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU) const =0
Length scale.
tmp< volScalarField > nuTilda() const
Return the modified kinematic viscosity.
volScalarField nuTilda_
Modified kinematic viscosity [m^2/s].
BasicTurbulenceModel::rhoField rhoField
Definition: RASModel.H:143
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 tmp< volScalarField > Stilda(const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU, const volScalarField &dTilda) const
tmp< volScalarField > r(const volScalarField &nur, const volScalarField &Stilda, const volScalarField &dTilda) const
BasicTurbulenceModel::transportModel transportModel
Definition: RASModel.H:144
virtual tmp< volScalarField > k() const
Return the (estimated) turbulent kinetic energy.
tmp< volScalarField > ft2(const volScalarField &chi) const
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
tmp< volScalarField > fv2(const volScalarField &chi, const volScalarField &fv1) const
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:63
BasicEddyViscosityModel::rhoField rhoField
virtual bool read()
Re-read model coefficients if they have changed.
tmp< volScalarField::Internal > fw(const volScalarField &Stilda, const volScalarField &dTilda) const
BasicEddyViscosityModel::transportModel transportModel
BasicEddyViscosityModel::alphaField alphaField
U
Definition: pEqn.H:72
virtual void correct()
Correct nuTilda and related properties.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const volScalarField & y_
Wall distance.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
virtual tmp< volScalarField > epsilon() const
Return the (estimated) turbulent kinetic energy dissipation rate.
Namespace for OpenFOAM.
BasicTurbulenceModel::alphaField alphaField
Definition: RASModel.H:142