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-2018 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 BasicTurbulenceModel>
86 class SpalartAllmaras
87 :
88  public eddyViscosity<RASModel<BasicTurbulenceModel>>
89 {
90  // Private Member Functions
91 
92  // Disallow default bitwise copy construct and assignment
94  void operator=(const SpalartAllmaras&);
95 
96 
97 protected:
98 
99  // Protected data
100 
101  // Model coefficients
113 
114 
115  // Fields
118 
119  //- Wall distance
120  // Note: different to wall distance in parent RASModel
121  // which is for near-wall cells only
122  const volScalarField& y_;
123 
124 
125  // Protected Member Functions
126 
127  tmp<volScalarField> chi() const;
128 
130 
132  (
133  const volScalarField& chi,
134  const volScalarField& fv1
135  ) const;
136 
138  (
139  const volScalarField& chi,
140  const volScalarField& fv1
141  ) const;
142 
144 
145  void correctNut(const volScalarField& fv1);
146  virtual void correctNut();
147 
148 
149 public:
151  typedef typename BasicTurbulenceModel::alphaField alphaField;
152  typedef typename BasicTurbulenceModel::rhoField rhoField;
153  typedef typename BasicTurbulenceModel::transportModel transportModel;
154 
155 
156  //- Runtime type information
157  TypeName("SpalartAllmaras");
158 
159 
160  // Constructors
161 
162  //- Construct from components
164  (
165  const alphaField& alpha,
166  const rhoField& rho,
167  const volVectorField& U,
168  const surfaceScalarField& alphaRhoPhi,
169  const surfaceScalarField& phi,
170  const transportModel& transport,
171  const word& propertiesName = turbulenceModel::propertiesName,
172  const word& type = typeName
173  );
174 
175 
176  //- Destructor
177  virtual ~SpalartAllmaras()
178  {}
179 
180 
181  // Member Functions
182 
183  //- Read RASProperties dictionary
184  virtual bool read();
185 
186  //- Return the effective diffusivity for nuTilda
188 
189  //- Return the turbulence kinetic energy
190  virtual tmp<volScalarField> k() const;
191 
192  //- Return the turbulence kinetic energy dissipation rate
193  virtual tmp<volScalarField> epsilon() const;
194 
195  //- Solve the turbulence equations and correct the turbulence viscosity
196  virtual void correct();
197 };
198 
199 
200 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
201 
202 } // End namespace RASModels
203 } // End namespace Foam
204 
205 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
206 
207 #ifdef NoRepository
208  #include "SpalartAllmaras.C"
209 #endif
210 
211 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
212 
213 #endif
214 
215 // ************************************************************************* //
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
surfaceScalarField & phi
TypeName("SpalartAllmaras")
Runtime type information.
tmp< volScalarField > Stilda(const volScalarField &chi, const volScalarField &fv1) const
BasicTurbulenceModel::alphaField alphaField
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:49
tmp< volScalarField > chi() const
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
const volScalarField & y_
Wall distance.
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from string.
Definition: word.H:59
tmp< volScalarField > DnuTildaEff() const
Return the effective diffusivity for nuTilda.
virtual bool read()
Read RASProperties dictionary.
BasicTurbulenceModel::rhoField rhoField
BasicTurbulenceModel::transportModel transportModel
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
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
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
A class for managing temporary objects.
Definition: PtrList.H:53
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
tmp< volScalarField > fv2(const volScalarField &chi, const volScalarField &fv1) const
tmp< volScalarField > fw(const volScalarField &Stilda) const
Namespace for OpenFOAM.