SpalartAllmaras.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) 2019-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
27Class
28 Foam::RASModels::SpalartAllmaras
29
30Group
31 grpRASTurbulence
32
33Description
34 Spalart-Allmaras one-transport-equation linear-eddy-viscosity turbulence
35 closure model for incompressible and compressible external flows.
36
37 Required fields
38 \verbatim
39 nuTilda | Modified kinematic viscosity [m2/s]
40 \endverbatim
41
42 References:
43 \verbatim
44 Standard model:
45 Spalart, P.R., & Allmaras, S.R. (1994).
46 A one-equation turbulence model for aerodynamic flows.
47 La Recherche Aerospatiale, 1, 5-21.
48
49 Standard model without trip and ft2 terms (tag:R):
50 Rumsey, C. (2020).
51 The Spalart-Allmaras Turbulence Model.
52 Spalart-Allmaras One-Equation Model without ft2 Term (SA-noft2).
53 https://turbmodels.larc.nasa.gov/spalart.html#sanoft2
54 (Retrieved:12-01-2021).
55
56 Estimation expression for k and epsilon (tag:B), Eq. 4.50:
57 Bourgoin, A. (2019).
58 Bathymetry induced turbulence modelling the
59 Alderney Race site: regional approach with TELEMAC-LES.
60 Normandie Université.
61
62 Estimation expressions for omega (tag:P):
63 Pope, S. B. (2000).
64 Turbulent flows.
65 Cambridge, UK: Cambridge Univ. Press
66 DOI:10.1017/CBO9780511840531
67 \endverbatim
68
69Usage
70 Example by using \c constant/turbulenceProperties:
71 \verbatim
72 RAS
73 {
74 // Mandatory entries (unmodifiable)
75 RASModel SpalartAllmaras;
76
77 // Optional entries (runtime modifiable)
78 turbulence on;
79 printCoeffs on;
80
81 SpalartAllmarasCoeffs
82 {
83 sigmaNut 0.66666;
84 kappa 0.41;
85 Cb1 0.1355;
86 Cb2 0.622;
87 Cw2 0.3;
88 Cw3 2.0;
89 Cv1 7.1;
90 Cs 0.3;
91 }
92 }
93 \endverbatim
94
95Note
96 - The model is implemented without the trip-term since the model has almost
97 always been used in fully turbulent applications rather than those where
98 laminar-turbulent transition occurs.
99 - It has been argued that the \c ft2 term is not needed in the absence of the
100 trip-term, hence \c ft2 term is also not implementated.
101 - The \c Stilda generation term should never be allowed to be zero or negative
102 to avoid potential numerical issues and unphysical results for complex
103 flows. To this end, a limiter proposed by Spalart (R:Note-1(b)) is applied
104 onto \c Stilda where \c Stilda is clipped at \c Cs*Omega with the default
105 value of \c Cs=0.3.
106 - The model does not produce \c k, \c epsilon or \c omega. Nevertheless,
107 these quantities can be estimated by using an approximate expressions for
108 turbulent kinetic energy and dissipation rate reported in (B:Eq. 4.50).
109
110SourceFiles
111 SpalartAllmaras.C
112
113\*---------------------------------------------------------------------------*/
114
115#ifndef SpalartAllmaras_H
116#define SpalartAllmaras_H
117
118#include "RASModel.H"
119#include "eddyViscosity.H"
120
121// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
122
123namespace Foam
124{
125namespace RASModels
126{
127
128/*---------------------------------------------------------------------------*\
129 Class SpalartAllmaras Declaration
130\*---------------------------------------------------------------------------*/
131
132template<class BasicTurbulenceModel>
133class SpalartAllmaras
134:
135 public eddyViscosity<RASModel<BasicTurbulenceModel>>
136{
137 // Private Member Functions
138
139 //- No copy construct
140 SpalartAllmaras(const SpalartAllmaras&) = delete;
141
142 //- No copy assignment
143 void operator=(const SpalartAllmaras&) = delete;
144
145
146protected:
147
148 // Protected Data
149
150 // Model coefficients
161
162
163 // Fields
164
165 //- Modified kinematic viscosity [m2/s]
167
168 //- Wall distance
169 // Note: different to wall distance in parent RASModel
170 // which is for near-wall cells only
172
173
174 // Protected Member Functions
175
176 tmp<volScalarField> chi() const;
177
179
181 (
184 ) const;
185
187
189 (
191 ) const;
192
193
194 //- Update nut with the latest available nuTilda
195 virtual void correctNut();
196
197
198public:
200 typedef typename BasicTurbulenceModel::alphaField alphaField;
201 typedef typename BasicTurbulenceModel::rhoField rhoField;
202 typedef typename BasicTurbulenceModel::transportModel transportModel;
203
204
205 //- Runtime type information
206 TypeName("SpalartAllmaras");
207
208
209 // Constructors
210
211 //- Construct from components
213 (
214 const alphaField& alpha,
215 const rhoField& rho,
216 const volVectorField& U,
217 const surfaceScalarField& alphaRhoPhi,
218 const surfaceScalarField& phi,
219 const transportModel& transport,
220 const word& propertiesName = turbulenceModel::propertiesName,
221 const word& type = typeName
222 );
223
224
225 //- Destructor
226 virtual ~SpalartAllmaras() = default;
227
228
229 // Member Functions
230
231 //- Re-read model coefficients if they have changed
232 virtual bool read();
233
234 //- Return the effective diffusivity for nuTilda
236
237 //- Return the (estimated) turbulent kinetic energy
238 virtual tmp<volScalarField> k() const;
239
240 //- Return the (estimated) turbulent kinetic energy dissipation rate
241 virtual tmp<volScalarField> epsilon() const;
242
243 //- Return the (estimated) specific dissipation rate
244 virtual tmp<volScalarField> omega() const;
245
246 //- Solve the turbulence equations and correct the turbulent viscosity
247 virtual void correct();
248};
249
250
251// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
252
253} // End namespace RASModels
254} // End namespace Foam
255
256// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
257
258#ifdef NoRepository
259 #include "SpalartAllmaras.C"
260#endif
261
262// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
263
264#endif
265
266// ************************************************************************* //
surfaceScalarField & phi
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Spalart-Allmaras one-transport-equation linear-eddy-viscosity turbulence closure model for incompress...
BasicTurbulenceModel::alphaField alphaField
tmp< volScalarField > chi() const
BasicTurbulenceModel::rhoField rhoField
virtual tmp< volScalarField > k() const
Return the (estimated) turbulent kinetic energy.
virtual void correct()
Solve the turbulence equations and correct the turbulent viscosity.
tmp< volScalarField > DnuTildaEff() const
Return the effective diffusivity for nuTilda.
const volScalarField::Internal & y_
Wall distance.
TypeName("SpalartAllmaras")
Runtime type information.
volScalarField nuTilda_
Modified kinematic viscosity [m2/s].
tmp< volScalarField > fv1(const volScalarField &chi) const
tmp< volScalarField::Internal > fv2(const volScalarField::Internal &chi, const volScalarField::Internal &fv1) const
virtual tmp< volScalarField > epsilon() const
Return the (estimated) turbulent kinetic energy dissipation rate.
virtual void correctNut()
Update nut with the latest available nuTilda.
BasicTurbulenceModel::transportModel transportModel
tmp< volScalarField::Internal > fw(const volScalarField::Internal &Stilda) const
virtual tmp< volScalarField > omega() const
Return the (estimated) specific dissipation rate.
tmp< volScalarField::Internal > Stilda() const
virtual bool read()
Re-read model coefficients if they have changed.
virtual ~SpalartAllmaras()=default
Destructor.
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:58
A class for managing temporary objects.
Definition: tmp.H:65
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
Namespace for OpenFOAM.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
volScalarField & alpha
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73