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-2022 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 BasicMomentumTransportModel>
62 :
63  public LESeddyViscosity<BasicMomentumTransportModel>
64 {
65 protected:
66 
67  // Protected data
68 
69  // Model constants
70 
73 
83 
84  // Fields
85 
87 
88  //- Wall distance
89  // Note: different to wall distance in parent RASModel
90  // which is for near-wall cells only
91  const volScalarField& y_;
92 
93 
94  // Protected Member Functions
95 
96  tmp<volScalarField> chi() const;
97 
99 
101  (
104  ) const;
105 
107  (
108  const volTensorField::Internal& gradU
109  ) const;
110 
112  (
117  ) const;
118 
120  (
121  const volScalarField::Internal& nur,
124  ) const;
125 
127  (
130  ) const;
131 
132  //- Length scale
134  (
137  const volTensorField::Internal& gradU
138  ) const;
139 
140  //- Cache the LES region indicator field
141  virtual void cacheLESRegion
142  (
144  ) const;
145 
146  void correctNut(const volScalarField& fv1);
147  virtual void correctNut();
148 
149 
150 public:
151 
152  typedef typename BasicMomentumTransportModel::alphaField alphaField;
153  typedef typename BasicMomentumTransportModel::rhoField rhoField;
154 
155 
156  //- Runtime type information
157  TypeName("SpalartAllmarasDES");
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 viscosity& viscosity,
171  const word& type = typeName
172  );
173 
174  //- Disallow default bitwise copy construction
175  SpalartAllmarasDES(const SpalartAllmarasDES&) = delete;
176 
177 
178  //- Destructor
179  virtual ~SpalartAllmarasDES()
180  {}
181 
182 
183  // Member Functions
184 
185  //- Read model coefficients if they have changed
186  virtual bool read();
187 
188  //- Return the effective diffusivity for nuTilda
190 
191  //- Return SGS kinetic energy
192  virtual tmp<volScalarField> k() const;
193 
195  {
196  return nuTilda_;
197  }
198 
199  //- Correct nuTilda and related properties
200  virtual void correct();
201 
202 
203  // Member Operators
204 
205  //- Disallow default bitwise assignment
206  void operator=(const SpalartAllmarasDES&) = delete;
207 };
208 
209 
210 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
211 
212 } // End namespace LESModels
213 } // End namespace Foam
214 
215 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
216 
217 #ifdef NoRepository
218  #include "SpalartAllmarasDES.C"
219 #endif
220 
221 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
222 
223 #endif
224 
225 // ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
Eddy viscosity LES SGS model base class.
BasicMomentumTransportModel::alphaField alphaField
BasicMomentumTransportModel::rhoField rhoField
SpalartAllmarasDES DES turbulence model for incompressible and compressible flows.
BasicMomentumTransportModel::alphaField alphaField
tmp< volScalarField > chi() const
tmp< volScalarField::Internal > Stilda(const volScalarField::Internal &chi, const volScalarField::Internal &fv1, const volScalarField::Internal &Omega, const volScalarField::Internal &dTilda) const
tmp< volScalarField::Internal > Omega(const volTensorField::Internal &gradU) const
void operator=(const SpalartAllmarasDES &)=delete
Disallow default bitwise assignment.
virtual tmp< volScalarField > k() const
Return SGS kinetic energy.
virtual void correct()
Correct nuTilda and related properties.
tmp< volScalarField > DnuTildaEff() const
Return the effective diffusivity for nuTilda.
tmp< volScalarField > fv1(const volScalarField &chi) const
tmp< volScalarField::Internal > fv2(const volScalarField::Internal &chi, const volScalarField::Internal &fv1) const
virtual tmp< volScalarField::Internal > dTilda(const volScalarField::Internal &chi, const volScalarField::Internal &fv1, const volTensorField::Internal &gradU) const
Length scale.
tmp< volScalarField > nuTilda() const
SpalartAllmarasDES(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity, const word &type=typeName)
Construct from components.
const volScalarField & y_
Wall distance.
virtual void cacheLESRegion(const volScalarField::Internal &dTilda) const
Cache the LES region indicator field.
TypeName("SpalartAllmarasDES")
Runtime type information.
tmp< volScalarField::Internal > fw(const volScalarField::Internal &Omega, const volScalarField::Internal &dTilda) const
tmp< volScalarField::Internal > r(const volScalarField::Internal &nur, const volScalarField::Internal &Omega, const volScalarField::Internal &dTilda) const
virtual bool read()
Read model coefficients if they have changed.
BasicMomentumTransportModel::rhoField rhoField
A class for managing temporary objects.
Definition: tmp.H:55
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
A class for handling words, derived from string.
Definition: word.H:62
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
Namespace for OpenFOAM.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488