dynamicKEqn.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::dynamicKEqn
26 
27 Description
28  Dynamic one equation eddy-viscosity model
29 
30  Eddy viscosity SGS model using a modeled balance equation to simulate
31  the behaviour of k in which a dynamic procedure is applied to evaluate the
32  coefficients.
33 
34  Reference:
35  \verbatim
36  Kim, W and Menon, S. (1995).
37  A new dynamic one-equation subgrid-scale model for
38  large eddy simulation.
39  In 33rd Aerospace Sciences Meeting and Exhibit, Reno, NV, 1995.
40  \endverbatim
41 
42  There are no default model coefficients but the filter used for KK must be
43  supplied, e.g.
44  \verbatim
45  dynamicKEqnCoeffs
46  {
47  filter simple;
48  }
49  \endverbatim
50 
51 SourceFiles
52  dynamicKEqn.C
53 
54 \*---------------------------------------------------------------------------*/
55 
56 #ifndef dynamicKEqn_H
57 #define dynamicKEqn_H
58 
59 #include "LESeddyViscosity.H"
60 #include "simpleFilter.H"
61 
62 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
63 
64 namespace Foam
65 {
66 namespace LESModels
67 {
68 
69 /*---------------------------------------------------------------------------*\
70  Class dynamicKEqn Declaration
71 \*---------------------------------------------------------------------------*/
72 
73 template<class BasicTurbulenceModel>
74 class dynamicKEqn
75 :
76  public LESeddyViscosity<BasicTurbulenceModel>
77 {
78  // Private Member Functions
79 
80  // Disallow default bitwise copy construct and assignment
81  dynamicKEqn(const dynamicKEqn&);
82  void operator=(const dynamicKEqn&);
83 
84 
85 protected:
86 
87  // Protected data
88 
89  // Fields
90 
92 
93 
94  // Filters
95 
99 
100 
101  // Protected Member Functions
102 
103  //- Calculate Ck by filtering the velocity field U
105  (
106  const volSymmTensorField& D,
107  const volScalarField& KK
108  ) const;
109 
110  //- Calculate Ce by filtering the velocity field U
112  (
113  const volSymmTensorField& D,
114  const volScalarField& KK
115  ) const;
116 
117  volScalarField Ce() const;
118 
119  //- Update sub-grid eddy-viscosity
120  void correctNut
121  (
122  const volSymmTensorField& D,
123  const volScalarField& KK
124  );
125 
126  virtual void correctNut();
127 
128  virtual tmp<fvScalarMatrix> kSource() const;
129 
130 
131 public:
133  typedef typename BasicTurbulenceModel::alphaField alphaField;
134  typedef typename BasicTurbulenceModel::rhoField rhoField;
135  typedef typename BasicTurbulenceModel::transportModel transportModel;
136 
137 
138  //- Runtime type information
139  TypeName("dynamicKEqn");
140 
141 
142  // Constructors
143 
144  //- Construct from components
146  (
147  const alphaField& alpha,
148  const rhoField& rho,
149  const volVectorField& U,
150  const surfaceScalarField& alphaRhoPhi,
151  const surfaceScalarField& phi,
152  const transportModel& transport,
153  const word& propertiesName = turbulenceModel::propertiesName,
154  const word& type = typeName
155  );
156 
157 
158  //- Destructor
159  virtual ~dynamicKEqn()
160  {}
161 
162 
163  // Member Functions
164 
165  //- Read model coefficients if they have changed
166  virtual bool read();
167 
168  //- Return SGS kinetic energy
169  virtual tmp<volScalarField> k() const
170  {
171  return k_;
172  }
173 
174  //- Return sub-grid disipation rate
175  virtual tmp<volScalarField> epsilon() const;
176 
177  //- Return the effective diffusivity for k
178  tmp<volScalarField> DkEff() const
179  {
180  return tmp<volScalarField>
181  (
182  new volScalarField("DkEff", this->nut_ + this->nu())
183  );
184  }
185 
186  //- Correct Eddy-Viscosity and related properties
187  virtual void correct();
188 };
189 
190 
191 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
192 
193 } // End namespace LESModels
194 } // End namespace Foam
195 
196 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
197 
198 #ifdef NoRepository
199  #include "dynamicKEqn.C"
200 #endif
201 
202 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
203 
204 #endif
205 
206 // ************************************************************************* //
Eddy viscosity LES SGS model base class.
virtual bool read()
Read model coefficients if they have changed.
Definition: dynamicKEqn.C:200
Abstract class for LES filters.
Definition: LESfilter.H:54
virtual tmp< volScalarField > k() const
Return SGS kinetic energy.
Definition: dynamicKEqn.H:168
surfaceScalarField & phi
virtual void correct()
Correct Eddy-Viscosity and related properties.
Definition: dynamicKEqn.C:237
BasicTurbulenceModel::rhoField rhoField
Definition: dynamicKEqn.H:133
virtual ~dynamicKEqn()
Destructor.
Definition: dynamicKEqn.H:158
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
volScalarField Ck(const volSymmTensorField &D, const volScalarField &KK) const
Calculate Ck by filtering the velocity field U.
Definition: dynamicKEqn.C:40
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from string.
Definition: word.H:59
autoPtr< LESfilter > filterPtr_
Definition: dynamicKEqn.H:96
BasicTurbulenceModel::transportModel transportModel
Definition: dynamicKEqn.H:134
BasicTurbulenceModel::alphaField alphaField
Definition: dynamicKEqn.H:132
volScalarField Ce() const
Definition: dynamicKEqn.C:88
TypeName("dynamicKEqn")
Runtime type information.
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
Dynamic one equation eddy-viscosity model.
Definition: dynamicKEqn.H:73
U
Definition: pEqn.H:72
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
Definition: dynamicKEqn.H:177
Simple top-hat filter used in dynamic LES models.
Definition: simpleFilter.H:50
virtual tmp< volScalarField > epsilon() const
Return sub-grid disipation rate.
Definition: dynamicKEqn.C:216
virtual tmp< fvScalarMatrix > kSource() const
Definition: dynamicKEqn.C:130
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
A class for managing temporary objects.
Definition: PtrList.H:53
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
volScalarField & nu
Namespace for OpenFOAM.