SpalartAllmaras.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::RASModels::SpalartAllmaras
26 
27 Description
28  Spalart-Allmaras one-eqn mixing-length model for incompressible and
29  compressible external flows.
30 
31  Reference:
32  \verbatim
33  Spalart, P.R., & Allmaras, S.R. (1994).
34  A one-equation turbulence model for aerodynamic flows.
35  La Recherche Aerospatiale, 1, 5-21.
36  \endverbatim
37 
38  The model is implemented without the trip-term and hence the ft2 term is
39  not needed.
40 
41  It is necessary to limit the Stilda generation term as the model generates
42  unphysical results if this term becomes negative which occurs for complex
43  flow. Several approaches have been proposed to limit Stilda but it is not
44  clear which is the most appropriate. Here the limiter proposed by Spalart
45  is implemented in which Stilda is clipped at Cs*Omega with the default value
46  of Cs = 0.3.
47 
48  The default model coefficients are
49  \verbatim
50  SpalartAllmarasCoeffs
51  {
52  Cb1 0.1355;
53  Cb2 0.622;
54  Cw2 0.3;
55  Cw3 2.0;
56  Cv1 7.1;
57  Cs 0.3;
58  sigmaNut 0.66666;
59  kappa 0.41;
60  }
61  \endverbatim
62 
63 SourceFiles
64  SpalartAllmaras.C
65 
66 \*---------------------------------------------------------------------------*/
67 
68 #ifndef SpalartAllmaras_H
69 #define SpalartAllmaras_H
70 
71 #include "RASModel.H"
72 #include "eddyViscosity.H"
73 
74 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
75 
76 namespace Foam
77 {
78 namespace RASModels
79 {
80 
81 /*---------------------------------------------------------------------------*\
82  Class SpalartAllmaras Declaration
83 \*---------------------------------------------------------------------------*/
84 
85 template<class BasicMomentumTransportModel>
86 class SpalartAllmaras
87 :
88  public eddyViscosity<RASModel<BasicMomentumTransportModel>>
89 {
90 
91 protected:
92 
93  // Protected data
94 
95  // Model coefficients
96 
99 
107 
108 
109  // Fields
110 
112 
113  //- Wall distance
114  // Note: different to wall distance in parent RASModel
115  // which is for near-wall cells only
117 
118 
119  // Protected Member Functions
120 
121  tmp<volScalarField> chi() const;
122 
124 
126  (
129  ) const;
130 
132  (
135  ) const;
136 
138  (
140  ) const;
141 
142  void correctNut(const volScalarField& fv1);
143  virtual void correctNut();
144 
145 
146 public:
147 
148  typedef typename BasicMomentumTransportModel::alphaField alphaField;
149  typedef typename BasicMomentumTransportModel::rhoField rhoField;
150 
151 
152  //- Runtime type information
153  TypeName("SpalartAllmaras");
154 
155 
156  // Constructors
157 
158  //- Construct from components
160  (
161  const alphaField& alpha,
162  const rhoField& rho,
163  const volVectorField& U,
164  const surfaceScalarField& alphaRhoPhi,
165  const surfaceScalarField& phi,
166  const viscosity& viscosity,
167  const word& type = typeName
168  );
169 
170  //- Disallow default bitwise copy construction
171  SpalartAllmaras(const SpalartAllmaras&) = delete;
172 
173 
174  //- Destructor
175  virtual ~SpalartAllmaras()
176  {}
177 
178 
179  // Member Functions
180 
181  //- Read RASProperties dictionary
182  virtual bool read();
183 
184  //- Return the effective diffusivity for nuTilda
186 
187  //- Return the turbulence kinetic energy
188  virtual tmp<volScalarField> k() const;
189 
190  //- Return the turbulence kinetic energy dissipation rate
191  virtual tmp<volScalarField> epsilon() const;
192 
193  //- Return the turbulence specific dissipation rate
194  virtual tmp<volScalarField> omega() const;
195 
196  //- Solve the turbulence equations and correct the turbulence viscosity
197  virtual void correct();
198 
199 
200  // Member Operators
201 
202  //- Disallow default bitwise assignment
203  void operator=(const SpalartAllmaras&) = delete;
204 };
205 
206 
207 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
208 
209 } // End namespace RASModels
210 } // End namespace Foam
211 
212 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
213 
214 #ifdef NoRepository
215  #include "SpalartAllmaras.C"
216 #endif
217 
218 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
219 
220 #endif
221 
222 // ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
Spalart-Allmaras one-eqn mixing-length model for incompressible and compressible external flows.
BasicMomentumTransportModel::alphaField alphaField
virtual ~SpalartAllmaras()
Destructor.
tmp< volScalarField > chi() const
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
tmp< volScalarField > DnuTildaEff() const
Return the effective diffusivity for nuTilda.
const volScalarField::Internal & y_
Wall distance.
TypeName("SpalartAllmaras")
Runtime type information.
void operator=(const SpalartAllmaras &)=delete
Disallow default bitwise assignment.
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 turbulence kinetic energy dissipation rate.
tmp< volScalarField::Internal > fw(const volScalarField::Internal &Stilda) const
SpalartAllmaras(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity, const word &type=typeName)
Construct from components.
tmp< volScalarField::Internal > Stilda(const volScalarField::Internal &chi, const volScalarField::Internal &fv1) const
virtual tmp< volScalarField > omega() const
Return the turbulence specific dissipation rate.
virtual bool read()
Read RASProperties dictionary.
BasicMomentumTransportModel::rhoField rhoField
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:52
BasicMomentumTransportModel::alphaField alphaField
BasicMomentumTransportModel::rhoField rhoField
A class for managing temporary objects.
Definition: tmp.H:55
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
A class for handling words, derived from string.
Definition: word.H:62
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
Namespace for OpenFOAM.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488