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 | Copyright (C) 2011-2016 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 Group
28  grpRASTurbulence
29 
30 Description
31  Spalart-Allmaras one-eqn mixing-length model for incompressible and
32  compressible external flows.
33 
34  Reference:
35  \verbatim
36  Spalart, P.R., & Allmaras, S.R. (1994).
37  A one-equation turbulence model for aerodynamic flows.
38  La Recherche Aerospatiale, 1, 5-21.
39  \endverbatim
40 
41  The model is implemented without the trip-term and hence the ft2 term is
42  not needed.
43 
44  It is necessary to limit the Stilda generation term as the model generates
45  unphysical results if this term becomes negative which occurs for complex
46  flow. Several approaches have been proposed to limit Stilda but it is not
47  clear which is the most appropriate. Here the limiter proposed by Spalart
48  is implemented in which Stilda is clipped at Cs*Omega with the default value
49  of Cs = 0.3.
50 
51  The default model coefficients are
52  \verbatim
53  SpalartAllmarasCoeffs
54  {
55  Cb1 0.1355;
56  Cb2 0.622;
57  Cw2 0.3;
58  Cw3 2.0;
59  Cv1 7.1;
60  Cs 0.3;
61  sigmaNut 0.66666;
62  kappa 0.41;
63  }
64  \endverbatim
65 
66 SourceFiles
67  SpalartAllmaras.C
68 
69 \*---------------------------------------------------------------------------*/
70 
71 #ifndef SpalartAllmaras_H
72 #define SpalartAllmaras_H
73 
74 #include "RASModel.H"
75 #include "eddyViscosity.H"
76 
77 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
78 
79 namespace Foam
80 {
81 namespace RASModels
82 {
83 
84 /*---------------------------------------------------------------------------*\
85  Class SpalartAllmaras Declaration
86 \*---------------------------------------------------------------------------*/
87 
88 template<class BasicTurbulenceModel>
89 class SpalartAllmaras
90 :
91  public eddyViscosity<RASModel<BasicTurbulenceModel>>
92 {
93  // Private Member Functions
94 
95  // Disallow default bitwise copy construct and assignment
97  void operator=(const SpalartAllmaras&);
98 
99 
100 protected:
101 
102  // Protected data
103 
104  // Model coefficients
116 
117 
118  // Fields
121 
122  //- Wall distance
123  // Note: different to wall distance in parent RASModel
124  // which is for near-wall cells only
125  const volScalarField& y_;
126 
127 
128  // Protected Member Functions
129 
130  tmp<volScalarField> chi() const;
131 
133 
135  (
136  const volScalarField& chi,
137  const volScalarField& fv1
138  ) const;
139 
141  (
142  const volScalarField& chi,
143  const volScalarField& fv1
144  ) const;
145 
147 
148  void correctNut(const volScalarField& fv1);
149  virtual void correctNut();
150 
151 
152 public:
154  typedef typename BasicTurbulenceModel::alphaField alphaField;
155  typedef typename BasicTurbulenceModel::rhoField rhoField;
156  typedef typename BasicTurbulenceModel::transportModel transportModel;
157 
158 
159  //- Runtime type information
160  TypeName("SpalartAllmaras");
161 
162 
163  // Constructors
164 
165  //- Construct from components
167  (
168  const alphaField& alpha,
169  const rhoField& rho,
170  const volVectorField& U,
171  const surfaceScalarField& alphaRhoPhi,
172  const surfaceScalarField& phi,
173  const transportModel& transport,
174  const word& propertiesName = turbulenceModel::propertiesName,
175  const word& type = typeName
176  );
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 
203 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
204 
205 } // End namespace RASModels
206 } // End namespace Foam
207 
208 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
209 
210 #ifdef NoRepository
211  #include "SpalartAllmaras.C"
212 #endif
213 
214 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
215 
216 #endif
217 
218 // ************************************************************************* //
surfaceScalarField & phi
U
Definition: pEqn.H:83
TypeName("SpalartAllmaras")
Runtime type information.
tmp< volScalarField > Stilda(const volScalarField &chi, const volScalarField &fv1) const
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
BasicTurbulenceModel::alphaField alphaField
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:52
tmp< volScalarField > fv2(const volScalarField &chi, const volScalarField &fv1) 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 > fv1(const volScalarField &chi) const
virtual bool read()
Read RASProperties dictionary.
BasicTurbulenceModel::rhoField rhoField
BasicTurbulenceModel::transportModel transportModel
virtual ~SpalartAllmaras()
Destructor.
Spalart-Allmaras one-eqn mixing-length model for incompressible and compressible external flows...
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:461
tmp< volScalarField > fw(const volScalarField &Stilda) const
tmp< volScalarField > DnuTildaEff() const
Return the effective diffusivity for nuTilda.
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
A class for managing temporary objects.
Definition: PtrList.H:54
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Namespace for OpenFOAM.
tmp< volScalarField > chi() const