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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2019-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Class
28  Foam::RASModels::SpalartAllmaras
29 
30 Group
31  grpRASTurbulence
32 
33 Description
34  Spalart-Allmaras one-eqn mixing-length model for incompressible and
35  compressible external flows.
36 
37  Reference:
38  \verbatim
39  Spalart, P.R., & Allmaras, S.R. (1994).
40  A one-equation turbulence model for aerodynamic flows.
41  La Recherche Aerospatiale, 1, 5-21.
42  \endverbatim
43 
44  The model is implemented without the trip-term and hence the ft2 term is
45  not needed.
46 
47  It is necessary to limit the Stilda generation term as the model generates
48  unphysical results if this term becomes negative which occurs for complex
49  flow. Several approaches have been proposed to limit Stilda but it is not
50  clear which is the most appropriate. Here the limiter proposed by Spalart
51  is implemented in which Stilda is clipped at Cs*Omega with the default value
52  of Cs = 0.3.
53 
54  The default model coefficients are
55  \verbatim
56  SpalartAllmarasCoeffs
57  {
58  Cb1 0.1355;
59  Cb2 0.622;
60  Cw2 0.3;
61  Cw3 2.0;
62  Cv1 7.1;
63  Cs 0.3;
64  sigmaNut 0.66666;
65  kappa 0.41;
66  }
67  \endverbatim
68 
69 SourceFiles
70  SpalartAllmaras.C
71 
72 \*---------------------------------------------------------------------------*/
73 
74 #ifndef SpalartAllmaras_H
75 #define SpalartAllmaras_H
76 
77 #include "RASModel.H"
78 #include "eddyViscosity.H"
79 
80 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
81 
82 namespace Foam
83 {
84 namespace RASModels
85 {
86 
87 /*---------------------------------------------------------------------------*\
88  Class SpalartAllmaras Declaration
89 \*---------------------------------------------------------------------------*/
90 
91 template<class BasicTurbulenceModel>
92 class SpalartAllmaras
93 :
94  public eddyViscosity<RASModel<BasicTurbulenceModel>>
95 {
96  // Private Member Functions
97 
98  //- No copy construct
99  SpalartAllmaras(const SpalartAllmaras&) = delete;
100 
101  //- No copy assignment
102  void operator=(const SpalartAllmaras&) = delete;
103 
104 
105 protected:
106 
107  // Protected data
108 
109  // Model coefficients
110 
113 
121 
122 
123  // Fields
124 
126 
127  //- Wall distance
128  // Note: different to wall distance in parent RASModel
129  // which is for near-wall cells only
130  const volScalarField& y_;
131 
132 
133  // Protected Member Functions
134 
135  tmp<volScalarField> chi() const;
136 
138 
140  (
141  const volScalarField& chi,
142  const volScalarField& fv1
143  ) const;
144 
146  (
147  const volScalarField& chi,
148  const volScalarField& fv1
149  ) const;
150 
152 
153  void correctNut(const volScalarField& fv1);
154  virtual void correctNut();
155 
156 
157 public:
158 
159  typedef typename BasicTurbulenceModel::alphaField alphaField;
160  typedef typename BasicTurbulenceModel::rhoField rhoField;
161  typedef typename BasicTurbulenceModel::transportModel transportModel;
162 
163 
164  //- Runtime type information
165  TypeName("SpalartAllmaras");
166 
167 
168  // Constructors
169 
170  //- Construct from components
172  (
173  const alphaField& alpha,
174  const rhoField& rho,
175  const volVectorField& U,
176  const surfaceScalarField& alphaRhoPhi,
177  const surfaceScalarField& phi,
178  const transportModel& transport,
179  const word& propertiesName = turbulenceModel::propertiesName,
180  const word& type = typeName
181  );
182 
183 
184  //- Destructor
185  virtual ~SpalartAllmaras() = default;
186 
187 
188  // Member Functions
189 
190  //- Re-read model coefficients if they have changed
191  virtual bool read();
192 
193  //- Return the effective diffusivity for nuTilda
195 
196  //- Return the turbulence kinetic energy
197  virtual tmp<volScalarField> k() const;
198 
199  //- Return the turbulence kinetic energy dissipation rate
200  virtual tmp<volScalarField> epsilon() const;
201 
202  //- Return the (estimated) specific dissipation rate
203  virtual tmp<volScalarField> omega() const;
204 
205  //- Solve the turbulence equations and correct the turbulence viscosity
206  virtual void correct();
207 };
208 
209 
210 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
211 
212 } // End namespace RASModels
213 } // End namespace Foam
214 
215 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
216 
217 #ifdef NoRepository
218  #include "SpalartAllmaras.C"
219 #endif
220 
221 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
222 
223 #endif
224 
225 // ************************************************************************* //
Foam::RASModels::SpalartAllmaras::DnuTildaEff
tmp< volScalarField > DnuTildaEff() const
Return the effective diffusivity for nuTilda.
Definition: SpalartAllmaras.C:294
Foam::RASModels::SpalartAllmaras::Cv1_
dimensionedScalar Cv1_
Definition: SpalartAllmaras.H:118
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::RASModels::SpalartAllmaras::correctNut
virtual void correctNut()
Definition: SpalartAllmaras.C:138
Foam::RASModels::SpalartAllmaras::Cw2_
dimensionedScalar Cw2_
Definition: SpalartAllmaras.H:116
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::RASModels::SpalartAllmaras::y_
const volScalarField & y_
Wall distance.
Definition: SpalartAllmaras.H:129
Foam::RASModels::SpalartAllmaras::fv1
tmp< volScalarField > fv1(const volScalarField &chi) const
Definition: SpalartAllmaras.C:52
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::RASModels::SpalartAllmaras::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: SpalartAllmaras.H:159
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
Foam::RASModels::SpalartAllmaras::k
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
Definition: SpalartAllmaras.C:304
Foam::RASModels::SpalartAllmaras::Cw1_
dimensionedScalar Cw1_
Definition: SpalartAllmaras.H:115
Foam::RASModels::SpalartAllmaras::~SpalartAllmaras
virtual ~SpalartAllmaras()=default
Destructor.
rho
rho
Definition: readInitialConditions.H:88
Foam::RASModels::SpalartAllmaras::omega
virtual tmp< volScalarField > omega() const
Return the (estimated) specific dissipation rate.
Definition: SpalartAllmaras.C:350
Foam::RASModels::SpalartAllmaras::fv2
tmp< volScalarField > fv2(const volScalarField &chi, const volScalarField &fv1) const
Definition: SpalartAllmaras.C:63
eddyViscosity.H
Foam::RASModels::SpalartAllmaras::Stilda
tmp< volScalarField > Stilda(const volScalarField &chi, const volScalarField &fv1) const
Definition: SpalartAllmaras.C:74
Foam::RASModels::SpalartAllmaras::Cb2_
dimensionedScalar Cb2_
Definition: SpalartAllmaras.H:114
Foam::RASModels::SpalartAllmaras::fw
tmp< volScalarField > fw(const volScalarField &Stilda) const
Definition: SpalartAllmaras.C:95
Foam::RASModels::SpalartAllmaras::kappa_
dimensionedScalar kappa_
Definition: SpalartAllmaras.H:111
Foam::RASModels::SpalartAllmaras::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: SpalartAllmaras.H:158
Foam::RASModels::SpalartAllmaras::Cb1_
dimensionedScalar Cb1_
Definition: SpalartAllmaras.H:113
Foam::RASModels::SpalartAllmaras::read
virtual bool read()
Re-read model coefficients if they have changed.
Definition: SpalartAllmaras.C:271
Foam::RASModels::SpalartAllmaras::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: SpalartAllmaras.H:160
Foam::RASModels::SpalartAllmaras::correct
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: SpalartAllmaras.C:372
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::RASModels::SpalartAllmaras::Cs_
dimensionedScalar Cs_
Definition: SpalartAllmaras.H:119
Foam::RASModel::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: RASModel.H:97
RASModel.H
Foam::dimensioned< scalar >
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::RASModel::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: RASModel.H:96
Foam::RASModels::SpalartAllmaras::nuTilda_
volScalarField nuTilda_
Definition: SpalartAllmaras.H:124
Foam::RASModel::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: RASModel.H:98
U
U
Definition: pEqn.H:72
Foam::RASModels::SpalartAllmaras
Spalart-Allmaras one-eqn mixing-length model for incompressible and compressible external flows.
Definition: SpalartAllmaras.H:91
Foam::RASModels::SpalartAllmaras::chi
tmp< volScalarField > chi() const
Definition: SpalartAllmaras.C:44
Foam::RASModels::SpalartAllmaras::Cw3_
dimensionedScalar Cw3_
Definition: SpalartAllmaras.H:117
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::eddyViscosity
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:55
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::RASModels::SpalartAllmaras::sigmaNut_
dimensionedScalar sigmaNut_
Definition: SpalartAllmaras.H:110
Foam::RASModels::SpalartAllmaras::epsilon
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
Definition: SpalartAllmaras.C:327
Foam::RASModels::SpalartAllmaras::TypeName
TypeName("SpalartAllmaras")
Runtime type information.