SpalartAllmarasDES.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::LESModels::SpalartAllmarasDES
26 
27 Description
28  SpalartAllmarasDES DES turbulence model for incompressible and
29  compressible flows
30 
31  Reference:
32  \verbatim
33  Spalart, P. R., Jou, W. H., Strelets, M., & Allmaras, S. R. (1997).
34  Comments on the feasibility of LES for wings, and on a hybrid
35  RANS/LES approach.
36  Advances in DNS/LES, 1, 4-8.
37  \endverbatim
38 
39 SourceFiles
40  SpalartAllmarasDES.C
41 
42 \*---------------------------------------------------------------------------*/
43 
44 #ifndef SpalartAllmarasDES_H
45 #define SpalartAllmarasDES_H
46 
47 #include "LESeddyViscosity.H"
48 
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51 namespace Foam
52 {
53 namespace LESModels
54 {
55 
56 /*---------------------------------------------------------------------------*\
57  Class SpalartAllmarasDES Declaration
58 \*---------------------------------------------------------------------------*/
59 
60 template<class BasicTurbulenceModel>
62 :
63  public LESeddyViscosity<BasicTurbulenceModel>
64 {
65  // Private Member Functions
66 
67  // Disallow default bitwise copy construct and assignment
69  void operator=(const SpalartAllmarasDES&);
70 
71 
72 protected:
73 
74  // Protected data
75 
76  // Model constants
77 
80 
90 
91  // Fields
92 
94 
95  //- Wall distance
96  // Note: different to wall distance in parent RASModel
97  // which is for near-wall cells only
98  const volScalarField& y_;
99 
100 
101  // Protected Member Functions
102 
103  tmp<volScalarField> chi() const;
104 
106 
108  (
109  const volScalarField& chi,
110  const volScalarField& fv1
111  ) const;
112 
113  tmp<volScalarField> S(const volTensorField& gradU) const;
114 
115  tmp<volScalarField> Omega(const volTensorField& gradU) const;
116 
118  (
119  const volScalarField& chi,
120  const volScalarField& fv1,
121  const volScalarField& Omega,
122  const volScalarField& dTilda
123  ) const;
124 
126  (
127  const volScalarField& nur,
128  const volScalarField& Omega,
129  const volScalarField& dTilda
130  ) const;
131 
133  (
134  const volScalarField& Omega,
135  const volScalarField& dTilda
136  ) const;
137 
138  //- Length scale
140  (
141  const volScalarField& chi,
142  const volScalarField& fv1,
143  const volTensorField& gradU
144  ) const;
145 
146  void correctNut(const volScalarField& fv1);
147  virtual void correctNut();
148 
149 
150 public:
152  typedef typename BasicTurbulenceModel::alphaField alphaField;
153  typedef typename BasicTurbulenceModel::rhoField rhoField;
154  typedef typename BasicTurbulenceModel::transportModel transportModel;
155 
156 
157  //- Runtime type information
158  TypeName("SpalartAllmarasDES");
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& propertiesName = turbulenceModel::propertiesName,
173  const word& type = typeName
174  );
175 
176 
177  //- Destructor
178  virtual ~SpalartAllmarasDES()
179  {}
180 
181 
182  // Member Functions
183 
184  //- Read model coefficients if they have changed
185  virtual bool read();
186 
187  //- Return the effective diffusivity for nuTilda
189 
190  //- Return SGS kinetic energy
191  virtual tmp<volScalarField> k() const;
194  {
195  return nuTilda_;
196  }
197 
198  //- Return the LES field indicator
199  virtual tmp<volScalarField> LESRegion() const;
200 
201  //- Correct nuTilda and related properties
202  virtual void correct();
203 };
204 
205 
206 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
207 
208 } // End namespace LESModels
209 } // End namespace Foam
210 
211 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
212 
213 #ifdef NoRepository
214  #include "SpalartAllmarasDES.C"
215 #endif
216 
217 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
218 
219 #endif
220 
221 // ************************************************************************* //
Eddy viscosity LES SGS model base class.
TypeName("SpalartAllmarasDES")
Runtime type information.
surfaceScalarField & phi
tmp< volScalarField > Stilda(const volScalarField &chi, const volScalarField &fv1, const volScalarField &Omega, const volScalarField &dTilda) const
const volScalarField & y_
Wall distance.
tmp< volScalarField > fv2(const volScalarField &chi, const volScalarField &fv1) const
tmp< volScalarField > r(const volScalarField &nur, const volScalarField &Omega, const volScalarField &dTilda) const
virtual bool read()
Read model coefficients if they have changed.
SpalartAllmarasDES DES turbulence model for incompressible and compressible flows.
BasicTurbulenceModel::transportModel transportModel
virtual tmp< volScalarField > LESRegion() const
Return the LES field indicator.
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
tmp< volScalarField > DnuTildaEff() const
Return the effective diffusivity for nuTilda.
tmp< volScalarField > S(const volTensorField &gradU) const
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
BasicTurbulenceModel::rhoField rhoField
U
Definition: pEqn.H:72
tmp< volScalarField > nuTilda() const
BasicTurbulenceModel::alphaField alphaField
A class for managing temporary objects.
Definition: PtrList.H:53
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
tmp< volScalarField > chi() const
virtual tmp< volScalarField > dTilda(const volScalarField &chi, const volScalarField &fv1, const volTensorField &gradU) const
Length scale.
virtual void correct()
Correct nuTilda and related properties.
tmp< volScalarField > fw(const volScalarField &Omega, const volScalarField &dTilda) const
tmp< volScalarField > Omega(const volTensorField &gradU) const
Namespace for OpenFOAM.
virtual tmp< volScalarField > k() const
Return SGS kinetic energy.