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-2020 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  // Private Member Functions
91 
92  using IOobject::modelName;
93 
94 
95 protected:
96 
97  // Protected data
98 
99  // Model coefficients
111 
112 
113  // Fields
116 
117  //- Wall distance
118  // Note: different to wall distance in parent RASModel
119  // which is for near-wall cells only
121 
122 
123  // Protected Member Functions
124 
125  tmp<volScalarField> chi() const;
126 
128 
130  (
133  ) const;
134 
136  (
139  ) const;
140 
142  (
144  ) const;
145 
146  void correctNut(const volScalarField& fv1);
147  virtual void correctNut();
148 
149 
150 public:
152  typedef typename BasicMomentumTransportModel::alphaField alphaField;
153  typedef typename BasicMomentumTransportModel::rhoField rhoField;
154  typedef typename BasicMomentumTransportModel::transportModel transportModel;
155 
156 
157  //- Runtime type information
158  TypeName("SpalartAllmaras");
159 
160 
161  // Constructors
162 
163  //- Construct from components
165  (
166  const alphaField& alpha,
167  const rhoField& rho,
168  const volVectorField& U,
169  const surfaceScalarField& alphaRhoPhi,
170  const surfaceScalarField& phi,
171  const transportModel& transport,
172  const word& type = typeName
173  );
174 
175  //- Disallow default bitwise copy construction
176  SpalartAllmaras(const SpalartAllmaras&) = delete;
177 
178 
179  //- Destructor
180  virtual ~SpalartAllmaras()
181  {}
182 
183 
184  // Member Functions
185 
186  //- Read RASProperties dictionary
187  virtual bool read();
188 
189  //- Return the effective diffusivity for nuTilda
191 
192  //- Return the turbulence kinetic energy
193  virtual tmp<volScalarField> k() const;
194 
195  //- Return the turbulence kinetic energy dissipation rate
196  virtual tmp<volScalarField> epsilon() const;
197 
198  //- Solve the turbulence equations and correct the turbulence viscosity
199  virtual void correct();
200 
201 
202  // Member Operators
203 
204  //- Disallow default bitwise assignment
205  void operator=(const SpalartAllmaras&) = delete;
206 };
207 
208 
209 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
210 
211 } // End namespace RASModels
212 } // End namespace Foam
213 
214 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
215 
216 #ifdef NoRepository
217  #include "SpalartAllmaras.C"
218 #endif
219 
220 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
221 
222 #endif
223 
224 // ************************************************************************* //
BasicMomentumTransportModel::rhoField rhoField
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
TypeName("SpalartAllmaras")
Runtime type information.
tmp< volScalarField::Internal > fw(const volScalarField::Internal &Stilda) const
BasicMomentumTransportModel::alphaField alphaField
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:49
tmp< volScalarField > chi() const
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
phi
Definition: pEqn.H:104
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
tmp< volScalarField::Internal > Stilda(const volScalarField::Internal &chi, const volScalarField::Internal &fv1) const
SpalartAllmaras(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &type=typeName)
Construct from components.
A class for handling words, derived from string.
Definition: word.H:59
const volScalarField::Internal & y_
Wall distance.
tmp< volScalarField > DnuTildaEff() const
Return the effective diffusivity for nuTilda.
virtual bool read()
Read RASProperties dictionary.
tmp< volScalarField::Internal > fv2(const volScalarField::Internal &chi, const volScalarField::Internal &fv1) const
BasicMomentumTransportModel::transportModel transportModel
virtual ~SpalartAllmaras()
Destructor.
Spalart-Allmaras one-eqn mixing-length model for incompressible and compressible external flows...
U
Definition: pEqn.H:72
tmp< volScalarField > fv1(const volScalarField &chi) const
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
void operator=(const SpalartAllmaras &)=delete
Disallow default bitwise assignment.
static word modelName(Name name, const word &model)
Return the name of the object within the given model.
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.