FreeFOAM The Cross-Platform CFD Toolkit
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) 1991-2010 OpenCFD Ltd.
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::SpalartAllmaras
26 
27 Description
28  SpalartAllmaras DES (SA + LES) turbulence model for incompressible flows
29 
30 SourceFiles
31  SpalartAllmaras.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef SpalartAllmaras_H
36 #define SpalartAllmaras_H
37 
39 #include <finiteVolume/volFields.H>
40 #include <finiteVolume/wallDist.H>
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 namespace incompressible
47 {
48 namespace LESModels
49 {
50 
51 /*---------------------------------------------------------------------------*\
52  Class SpalartAllmaras Declaration
53 \*---------------------------------------------------------------------------*/
54 
56 :
57  public LESModel
58 {
59  // Private member functions
60 
61  //- Update sub-grid scale fields
62  void updateSubGridScaleFields();
63 
64  // Disallow default bitwise copy construct and assignment
66  SpalartAllmaras& operator=(const SpalartAllmaras&);
67 
68 
69 protected:
70 
71  // Protected data
72 
75 
76 
77  // Model constants
78 
88 
89 
90  // Fields
91 
95 
96 
97  // Protected member functions
98 
99  virtual tmp<volScalarField> fv1() const;
100  virtual tmp<volScalarField> fv2() const;
101  virtual tmp<volScalarField> fv3() const;
102  virtual tmp<volScalarField> S(const volTensorField& gradU) const;
103 
105  (
106  const volScalarField& S,
107  const volScalarField& dTilda
108  ) const;
109 
110  virtual tmp<volScalarField> r
111  (
112  const volScalarField& visc,
113  const volScalarField& S,
114  const volScalarField& dTilda
115  ) const;
116 
117  virtual tmp<volScalarField> fw
118  (
119  const volScalarField& S,
120  const volScalarField& dTilda
121  ) const;
122 
123  //- Length scale
124  virtual tmp<volScalarField> dTilda(const volScalarField& S) const;
125 
126 
127 public:
128 
129  //- Runtime type information
130  TypeName("SpalartAllmaras");
131 
132 
133  // Constructors
134 
135  //- Construct from components
137  (
138  const volVectorField& U,
139  const surfaceScalarField& phi,
141  const word& modelName = typeName
142  );
143 
144 
145  //- Destructor
147  {}
148 
149 
150  // Member Functions
151 
152  //- Return SGS kinetic energy
153  virtual tmp<volScalarField> k() const;
154 
155  //- Return sub-grid disipation rate
156  virtual tmp<volScalarField> epsilon() const;
157 
159  {
160  return nuTilda_;
161  }
162 
163  //- Return SGS viscosity
164  virtual tmp<volScalarField> nuSgs() const
165  {
166  return nuSgs_;
167  }
168 
169  //- Return the sub-grid stress tensor.
170  virtual tmp<volSymmTensorField> B() const;
171 
172  //- Return the effective sub-grid turbulence stress tensor
173  // including the laminar stress
174  virtual tmp<volSymmTensorField> devBeff() const;
175 
176  //- Return the deviatoric part of the divergence of Beff
177  // i.e. the additional term in the filtered NSE.
179 
180  //- Correct nuTilda and related properties
181  virtual void correct(const tmp<volTensorField>& gradU);
182 
183  //- Read LESProperties dictionary
184  virtual bool read();
185 };
186 
187 
188 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
189 
190 } // End namespace LESModels
191 } // End namespace incompressible
192 } // End namespace Foam
193 
194 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
195 
196 #endif
197 
198 // ************************ vim: set sw=4 sts=4 et: ************************ //