Smagorinsky.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-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
27Class
28 Foam::LESModels::Smagorinsky
29
30Group
31 grpLESTurbulence
32
33Description
34 The Smagorinsky SGS model.
35
36 Reference:
37 \verbatim
38 Smagorinsky, J. (1963).
39 General circulation experiments with the primitive equations: I.
40 The basic experiment*.
41 Monthly weather review, 91(3), 99-164.
42 \endverbatim
43
44 The form of the Smagorinsky model implemented is obtained from the
45 k-equation model assuming local equilibrium which provides estimates of both
46 k and epsilon separate from the sub-grid scale viscosity:
47
48 \verbatim
49 B = 2/3*k*I - 2*nuSgs*dev(D)
50
51 where
52
53 D = symm(grad(U));
54 k from D:B + Ce*k^3/2/delta = 0
55 nuSgs = Ck*sqrt(k)*delta
56 \endverbatim
57
58 The default model coefficients are
59 \verbatim
60 SmagorinskyCoeffs
61 {
62 Ck 0.094;
63 Ce 1.048;
64 }
65 \endverbatim
66
67See also
68 Foam::LESModels::kEqn
69
70SourceFiles
71 Smagorinsky.C
72
73\*---------------------------------------------------------------------------*/
74
75#ifndef Smagorinsky_H
76#define Smagorinsky_H
77
78#include "LESModel.H"
79#include "LESeddyViscosity.H"
80
81// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
82
83namespace Foam
84{
85namespace LESModels
86{
87
88/*---------------------------------------------------------------------------*\
89 Class Smagorinsky Declaration
90\*---------------------------------------------------------------------------*/
91
92template<class BasicTurbulenceModel>
93class Smagorinsky
94:
95 public LESeddyViscosity<BasicTurbulenceModel>
96{
97 // Private Member Functions
98
99 //- No copy construct
100 Smagorinsky(const Smagorinsky&) = delete;
101
102 //- No copy assignment
103 void operator=(const Smagorinsky&) = delete;
104
105
106protected:
107
108 // Protected data
111
112
113 // Protected Member Functions
114
115 //- Return SGS kinetic energy
116 // calculated from the given velocity gradient
117 tmp<volScalarField> k(const tmp<volTensorField>& gradU) const;
118
119 //- Update the SGS eddy viscosity
120 virtual void correctNut();
121
122
123public:
125 typedef typename BasicTurbulenceModel::alphaField alphaField;
126 typedef typename BasicTurbulenceModel::rhoField rhoField;
127 typedef typename BasicTurbulenceModel::transportModel transportModel;
128
129
130 //- Runtime type information
131 TypeName("Smagorinsky");
132
133
134 // Constructors
135
136 //- Construct from components
138 (
139 const alphaField& alpha,
140 const rhoField& rho,
141 const volVectorField& U,
142 const surfaceScalarField& alphaRhoPhi,
143 const surfaceScalarField& phi,
144 const transportModel& transport,
145 const word& propertiesName = turbulenceModel::propertiesName,
146 const word& type = typeName
147 );
148
149
150 //- Destructor
151 virtual ~Smagorinsky() = default;
152
153
154 // Member Functions
155
156 //- Read model coefficients if they have changed
157 virtual bool read();
158
159 //- Return SGS kinetic energy
160 virtual tmp<volScalarField> k() const
161 {
162 return k(fvc::grad(this->U_));
163 }
164
165 //- Correct Eddy-Viscosity and related properties
166 virtual void correct();
167};
168
169
170// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
171
172} // End namespace LESModels
173} // End namespace Foam
174
175// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
176
177#ifdef NoRepository
178 #include "Smagorinsky.C"
179#endif
180
181// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
182
183#endif
184
185// ************************************************************************* //
surfaceScalarField & phi
Eddy viscosity LES SGS model base class.
The Smagorinsky SGS model.
Definition: Smagorinsky.H:95
BasicTurbulenceModel::alphaField alphaField
Definition: Smagorinsky.H:124
BasicTurbulenceModel::rhoField rhoField
Definition: Smagorinsky.H:125
virtual ~Smagorinsky()=default
Destructor.
virtual void correct()
Correct Eddy-Viscosity and related properties.
Definition: Smagorinsky.C:143
TypeName("Smagorinsky")
Runtime type information.
virtual void correctNut()
Update the SGS eddy viscosity.
Definition: Smagorinsky.C:70
BasicTurbulenceModel::transportModel transportModel
Definition: Smagorinsky.H:126
virtual tmp< volScalarField > k() const
Return SGS kinetic energy.
Definition: Smagorinsky.H:159
virtual bool read()
Read model coefficients if they have changed.
Definition: Smagorinsky.C:129
A class for managing temporary objects.
Definition: tmp.H:65
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
Namespace for OpenFOAM.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
volScalarField & alpha
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73