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-2020 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 BasicMomentumTransportModel>
74 class dynamicKEqn
75 :
76  public LESeddyViscosity<BasicMomentumTransportModel>
77 {
78 protected:
79 
80  // Protected data
81 
82  // Fields
83 
85 
86 
87  // Filters
88 
92 
93 
94  // Protected Member Functions
95 
96  //- Return the KK field obtained by filtering the velocity field U
97  tmp<volScalarField> KK() const;
98 
99  //- Calculate Ck by filtering the velocity field U
101  (
102  const volSymmTensorField& D,
103  const volScalarField& KK
104  ) const;
105 
106  //- Calculate Ce by filtering the velocity field U
108  (
109  const volSymmTensorField& D,
110  const volScalarField& KK
111  ) const;
112 
113  volScalarField Ce() const;
114 
115  //- Update sub-grid eddy-viscosity
116  void correctNut
117  (
118  const volSymmTensorField& D,
119  const volScalarField& KK
120  );
121 
122  virtual void correctNut();
123 
124  virtual tmp<fvScalarMatrix> kSource() const;
125 
126 
127 public:
129  typedef typename BasicMomentumTransportModel::alphaField alphaField;
130  typedef typename BasicMomentumTransportModel::rhoField rhoField;
131  typedef typename BasicMomentumTransportModel::transportModel transportModel;
132 
133 
134  //- Runtime type information
135  TypeName("dynamicKEqn");
136 
137 
138  // Constructors
139 
140  //- Construct from components
142  (
143  const alphaField& alpha,
144  const rhoField& rho,
145  const volVectorField& U,
146  const surfaceScalarField& alphaRhoPhi,
147  const surfaceScalarField& phi,
148  const transportModel& transport,
149  const word& type = typeName
150  );
151 
152  //- Disallow default bitwise copy construction
153  dynamicKEqn(const dynamicKEqn&) = delete;
154 
155 
156  //- Destructor
157  virtual ~dynamicKEqn()
158  {}
159 
160 
161  // Member Functions
162 
163  //- Read model coefficients if they have changed
164  virtual bool read();
165 
166  //- Return SGS kinetic energy
167  virtual tmp<volScalarField> k() const
168  {
169  return k_;
170  }
171 
172  //- Return sub-grid disipation rate
173  virtual tmp<volScalarField> epsilon() const;
174 
175  //- Return the effective diffusivity for k
176  tmp<volScalarField> DkEff() const
177  {
178  return volScalarField::New
179  (
180  "DkEff",
181  this->nut_ + this->nu()
182  );
183  }
184 
185  //- Correct Eddy-Viscosity and related properties
186  virtual void correct();
187 
188 
189  // Member Operators
190 
191  //- Disallow default bitwise assignment
192  void operator=(const dynamicKEqn&) = delete;
193 };
194 
195 
196 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
197 
198 } // End namespace LESModels
199 } // End namespace Foam
200 
201 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
202 
203 #ifdef NoRepository
204  #include "dynamicKEqn.C"
205 #endif
206 
207 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
208 
209 #endif
210 
211 // ************************************************************************* //
Eddy viscosity LES SGS model base class.
virtual bool read()
Read model coefficients if they have changed.
Definition: dynamicKEqn.C:196
Abstract class for LES filters.
Definition: LESfilter.H:54
virtual tmp< volScalarField > k() const
Return SGS kinetic energy.
Definition: dynamicKEqn.H:166
virtual void correct()
Correct Eddy-Viscosity and related properties.
Definition: dynamicKEqn.C:223
void operator=(const dynamicKEqn &)=delete
Disallow default bitwise assignment.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
BasicMomentumTransportModel::transportModel transportModel
Definition: dynamicKEqn.H:130
phi
Definition: pEqn.H:104
virtual ~dynamicKEqn()
Destructor.
Definition: dynamicKEqn.H:156
volScalarField Ck(const volSymmTensorField &D, const volScalarField &KK) const
Calculate Ck by filtering the velocity field U.
Definition: dynamicKEqn.C:52
A class for handling words, derived from string.
Definition: word.H:59
autoPtr< LESfilter > filterPtr_
Definition: dynamicKEqn.H:89
volScalarField Ce() const
Definition: dynamicKEqn.C:100
TypeName("dynamicKEqn")
Runtime type information.
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:175
Simple top-hat filter used in dynamic LES models.
Definition: simpleFilter.H:50
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
virtual tmp< volScalarField > epsilon() const
Return sub-grid disipation rate.
Definition: dynamicKEqn.C:212
virtual tmp< fvScalarMatrix > kSource() const
Definition: dynamicKEqn.C:128
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
tmp< volScalarField > KK() const
Return the KK field obtained by filtering the velocity field U.
Definition: dynamicKEqn.C:40
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
A class for managing temporary objects.
Definition: PtrList.H:53
BasicMomentumTransportModel::alphaField alphaField
Definition: dynamicKEqn.H:128
dynamicKEqn(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &type=typeName)
Construct from components.
Definition: dynamicKEqn.C:146
Namespace for OpenFOAM.
BasicMomentumTransportModel::rhoField rhoField
Definition: dynamicKEqn.H:129