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-2023 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 
114  // Protected Member Functions
115 
116  tmp<volScalarField> chi() const;
117 
119 
121  (
124  ) const;
125 
127  (
130  ) const;
131 
133  (
135  ) const;
136 
137  void correctNut(const volScalarField& fv1);
138  virtual void correctNut();
139 
140 
141 public:
142 
143  typedef typename BasicMomentumTransportModel::alphaField alphaField;
144  typedef typename BasicMomentumTransportModel::rhoField rhoField;
145 
146 
147  //- Runtime type information
148  TypeName("SpalartAllmaras");
149 
150 
151  // Constructors
152 
153  //- Construct from components
155  (
156  const alphaField& alpha,
157  const rhoField& rho,
158  const volVectorField& U,
159  const surfaceScalarField& alphaRhoPhi,
160  const surfaceScalarField& phi,
161  const viscosity& viscosity,
162  const word& type = typeName
163  );
164 
165  //- Disallow default bitwise copy construction
166  SpalartAllmaras(const SpalartAllmaras&) = delete;
167 
168 
169  //- Destructor
170  virtual ~SpalartAllmaras()
171  {}
172 
173 
174  // Member Functions
175 
176  //- Read RASProperties dictionary
177  virtual bool read();
178 
179  //- Return the effective diffusivity for nuTilda
181 
182  //- Return the turbulence kinetic energy
183  virtual tmp<volScalarField> k() const;
184 
185  //- Return the turbulence kinetic energy dissipation rate
186  virtual tmp<volScalarField> epsilon() const;
187 
188  //- Return the turbulence specific dissipation rate
189  virtual tmp<volScalarField> omega() const;
190 
191  //- Solve the turbulence equations and correct the turbulence viscosity
192  virtual void correct();
193 
194 
195  // Member Operators
196 
197  //- Disallow default bitwise assignment
198  void operator=(const SpalartAllmaras&) = delete;
199 };
200 
201 
202 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
203 
204 } // End namespace RASModels
205 } // End namespace Foam
206 
207 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
208 
209 #ifdef NoRepository
210  #include "SpalartAllmaras.C"
211 #endif
212 
213 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
214 
215 #endif
216 
217 // ************************************************************************* //
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.
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