Smagorinsky.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::Smagorinsky
26 
27 Description
28  The Smagorinsky SGS model.
29 
30  Reference:
31  \verbatim
32  Smagorinsky, J. (1963).
33  General circulation experiments with the primitive equations: I.
34  The basic experiment*.
35  Monthly weather review, 91(3), 99-164.
36  \endverbatim
37 
38  The form of the Smagorinsky model implemented is obtained from the
39  k-equation model assuming local equilibrium which provides estimates of both
40  k and epsilon separate from the sub-grid scale viscosity:
41 
42  \verbatim
43  B = 2/3*k*I - 2*nuSgs*dev(D)
44 
45  where
46 
47  D = symm(grad(U));
48  k from D:B + Ce*k^3/2/delta = 0
49  nuSgs = Ck*sqrt(k)*delta
50  \endverbatim
51 
52  The default model coefficients are
53  \verbatim
54  SmagorinskyCoeffs
55  {
56  Ck 0.094;
57  Ce 1.048;
58  }
59  \endverbatim
60 
61 See also
62  Foam::LESModels::kEqn
63 
64 SourceFiles
65  Smagorinsky.C
66 
67 \*---------------------------------------------------------------------------*/
68 
69 #ifndef Smagorinsky_H
70 #define Smagorinsky_H
71 
72 #include "LESModel.H"
73 #include "LESeddyViscosity.H"
74 
75 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
76 
77 namespace Foam
78 {
79 namespace LESModels
80 {
81 
82 /*---------------------------------------------------------------------------*\
83  Class Smagorinsky Declaration
84 \*---------------------------------------------------------------------------*/
85 
86 template<class BasicMomentumTransportModel>
87 class Smagorinsky
88 :
89  public LESeddyViscosity<BasicMomentumTransportModel>
90 {
91 protected:
92 
93  // Protected Member Functions
94 
95  //- Return SGS kinetic energy
96  // calculated from the given velocity gradient
97  tmp<volScalarField> k(const tmp<volTensorField>& gradU) const;
98 
99  //- Update the SGS eddy viscosity
100  virtual void correctNut();
101 
102 
103 public:
104 
105  typedef typename BasicMomentumTransportModel::alphaField alphaField;
106  typedef typename BasicMomentumTransportModel::rhoField rhoField;
107 
108 
109  //- Runtime type information
110  TypeName("Smagorinsky");
111 
112 
113  // Constructors
114 
115  //- Construct from components
117  (
118  const alphaField& alpha,
119  const rhoField& rho,
120  const volVectorField& U,
121  const surfaceScalarField& alphaRhoPhi,
122  const surfaceScalarField& phi,
123  const viscosity& viscosity,
124  const word& type = typeName
125  );
126 
127  //- Disallow default bitwise copy construction
128  Smagorinsky(const Smagorinsky&) = delete;
129 
130 
131  //- Destructor
132  virtual ~Smagorinsky()
133  {}
134 
135 
136  // Member Functions
137 
138  //- Read model coefficients if they have changed
139  virtual bool read();
140 
141  //- Return SGS kinetic energy
142  virtual tmp<volScalarField> k() const
143  {
144  return k(fvc::grad(this->U_));
145  }
146 
147  //- Correct Eddy-Viscosity and related properties
148  virtual void correct();
149 
150 
151  // Member Operators
152 
153  //- Disallow default bitwise assignment
154  void operator=(const Smagorinsky&) = delete;
155 };
156 
157 
158 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
159 
160 } // End namespace LESModels
161 } // End namespace Foam
162 
163 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
164 
165 #ifdef NoRepository
166  #include "Smagorinsky.C"
167 #endif
168 
169 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
170 
171 #endif
172 
173 // ************************************************************************* //
Generic GeometricField class.
Eddy viscosity LES SGS model base class.
BasicMomentumTransportModel::alphaField alphaField
BasicMomentumTransportModel::rhoField rhoField
The Smagorinsky SGS model.
Definition: Smagorinsky.H:89
BasicMomentumTransportModel::alphaField alphaField
Definition: Smagorinsky.H:104
Smagorinsky(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.
Definition: Smagorinsky.C:74
virtual tmp< volScalarField > k() const
Return SGS kinetic energy.
Definition: Smagorinsky.H:141
virtual void correct()
Correct Eddy-Viscosity and related properties.
Definition: Smagorinsky.C:112
virtual ~Smagorinsky()
Destructor.
Definition: Smagorinsky.H:131
TypeName("Smagorinsky")
Runtime type information.
virtual void correctNut()
Update the SGS eddy viscosity.
Definition: Smagorinsky.C:60
void operator=(const Smagorinsky &)=delete
Disallow default bitwise assignment.
virtual bool read()
Read model coefficients if they have changed.
Definition: Smagorinsky.C:105
BasicMomentumTransportModel::rhoField rhoField
Definition: Smagorinsky.H:105
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))
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
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