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