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-2021 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:
128 
129  typedef typename BasicMomentumTransportModel::alphaField alphaField;
130  typedef typename BasicMomentumTransportModel::rhoField rhoField;
131 
132 
133  //- Runtime type information
134  TypeName("dynamicKEqn");
135 
136 
137  // Constructors
138 
139  //- Construct from components
141  (
142  const alphaField& alpha,
143  const rhoField& rho,
144  const volVectorField& U,
145  const surfaceScalarField& alphaRhoPhi,
146  const surfaceScalarField& phi,
147  const viscosity& viscosity,
148  const word& type = typeName
149  );
150 
151  //- Disallow default bitwise copy construction
152  dynamicKEqn(const dynamicKEqn&) = delete;
153 
154 
155  //- Destructor
156  virtual ~dynamicKEqn()
157  {}
158 
159 
160  // Member Functions
161 
162  //- Read model coefficients if they have changed
163  virtual bool read();
164 
165  //- Return SGS kinetic energy
166  virtual tmp<volScalarField> k() const
167  {
168  return k_;
169  }
170 
171  //- Return sub-grid disipation rate
172  virtual tmp<volScalarField> epsilon() const;
173 
174  //- Return the effective diffusivity for k
175  tmp<volScalarField> DkEff() const
176  {
177  return volScalarField::New
178  (
179  "DkEff",
180  this->nut_ + this->nu()
181  );
182  }
183 
184  //- Correct Eddy-Viscosity and related properties
185  virtual void correct();
186 
187 
188  // Member Operators
189 
190  //- Disallow default bitwise assignment
191  void operator=(const dynamicKEqn&) = delete;
192 };
193 
194 
195 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
196 
197 } // End namespace LESModels
198 } // End namespace Foam
199 
200 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
201 
202 #ifdef NoRepository
203  #include "dynamicKEqn.C"
204 #endif
205 
206 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
207 
208 #endif
209 
210 // ************************************************************************* //
static const Foam::dimensionedScalar D("D", Foam::dimTemperature, 257.14)
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
Eddy viscosity LES SGS model base class.
Dynamic one equation eddy-viscosity model.
Definition: dynamicKEqn.H:76
BasicMomentumTransportModel::alphaField alphaField
Definition: dynamicKEqn.H:128
volScalarField Ce() const
Definition: dynamicKEqn.C:101
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
Definition: dynamicKEqn.H:174
volScalarField Ck(const volSymmTensorField &D, const volScalarField &KK) const
Calculate Ck by filtering the velocity field U.
Definition: dynamicKEqn.C:53
virtual tmp< volScalarField > k() const
Return SGS kinetic energy.
Definition: dynamicKEqn.H:165
virtual void correct()
Correct Eddy-Viscosity and related properties.
Definition: dynamicKEqn.C:224
TypeName("dynamicKEqn")
Runtime type information.
virtual ~dynamicKEqn()
Destructor.
Definition: dynamicKEqn.H:155
virtual tmp< volScalarField > epsilon() const
Return sub-grid disipation rate.
Definition: dynamicKEqn.C:213
virtual tmp< fvScalarMatrix > kSource() const
Definition: dynamicKEqn.C:129
tmp< volScalarField > KK() const
Return the KK field obtained by filtering the velocity field U.
Definition: dynamicKEqn.C:41
dynamicKEqn(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: dynamicKEqn.C:147
void operator=(const dynamicKEqn &)=delete
Disallow default bitwise assignment.
autoPtr< LESfilter > filterPtr_
Definition: dynamicKEqn.H:89
virtual bool read()
Read model coefficients if they have changed.
Definition: dynamicKEqn.C:197
BasicMomentumTransportModel::rhoField rhoField
Definition: dynamicKEqn.H:129
Abstract class for LES filters.
Definition: LESfilter.H:55
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Simple top-hat filter used in dynamic LES models.
Definition: simpleFilter.H:54
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))
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