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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2019-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Class
28 Foam::LESModels::dynamicKEqn
29
30Group
31 grpLESTurbulence
32
33Description
34 Dynamic one equation eddy-viscosity model
35
36 Eddy viscosity SGS model using a modeled balance equation to simulate
37 the behaviour of k in which a dynamic procedure is applied to evaluate the
38 coefficients.
39
40 Reference:
41 \verbatim
42 Kim, W and Menon, S. (1995).
43 A new dynamic one-equation subgrid-scale model for
44 large eddy simulation.
45 In 33rd Aerospace Sciences Meeting and Exhibit, Reno, NV, 1995.
46 \endverbatim
47
48 There are no default model coefficients but the filter used for KK must be
49 supplied, e.g.
50 \verbatim
51 dynamicKEqnCoeffs
52 {
53 filter simple;
54 }
55 \endverbatim
56
57SourceFiles
58 dynamicKEqn.C
59
60\*---------------------------------------------------------------------------*/
61
62#ifndef dynamicKEqn_H
63#define dynamicKEqn_H
64
65#include "LESeddyViscosity.H"
66#include "simpleFilter.H"
67
68// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
69
70namespace Foam
71{
72namespace LESModels
73{
74
75/*---------------------------------------------------------------------------*\
76 Class dynamicKEqn Declaration
77\*---------------------------------------------------------------------------*/
78
79template<class BasicTurbulenceModel>
80class dynamicKEqn
81:
82 public LESeddyViscosity<BasicTurbulenceModel>
83{
84 // Private Member Functions
85
86 //- No copy construct
87 dynamicKEqn(const dynamicKEqn&) = delete;
88
89 //- No copy assignment
90 void operator=(const dynamicKEqn&) = delete;
91
92
93protected:
94
95 // Protected data
96
97 // Fields
100
101
102 // Filters
107
108
109 // Protected Member Functions
110
111 //- Calculate Ck by filtering the velocity field U
113 (
114 const volSymmTensorField& D,
115 const volScalarField& KK
116 ) const;
117
118 //- Calculate Ce by filtering the velocity field U
120 (
121 const volSymmTensorField& D,
122 const volScalarField& KK
123 ) const;
124
125 volScalarField Ce() const;
126
127 //- Update sub-grid eddy-viscosity
128 void correctNut
129 (
130 const volSymmTensorField& D,
131 const volScalarField& KK
132 );
133
134 virtual void correctNut();
135
136 virtual tmp<fvScalarMatrix> kSource() const;
137
138
139public:
141 typedef typename BasicTurbulenceModel::alphaField alphaField;
142 typedef typename BasicTurbulenceModel::rhoField rhoField;
143 typedef typename BasicTurbulenceModel::transportModel transportModel;
144
145
146 //- Runtime type information
147 TypeName("dynamicKEqn");
148
149
150 // Constructors
151
152 //- Construct from components
154 (
155 const alphaField& alpha,
156 const rhoField& rho,
157 const volVectorField& U,
158 const surfaceScalarField& alphaRhoPhi,
159 const surfaceScalarField& phi,
160 const transportModel& transport,
161 const word& propertiesName = turbulenceModel::propertiesName,
162 const word& type = typeName
163 );
164
165
166 //- Destructor
167 virtual ~dynamicKEqn() = default;
168
169
170 // Member Functions
171
172 //- Read model coefficients if they have changed
173 virtual bool read();
174
175 //- Return SGS kinetic energy
176 virtual tmp<volScalarField> k() const
177 {
178 return k_;
179 }
180
181 //- Return the effective diffusivity for k
183 {
185 (
186 new volScalarField("DkEff", this->nut_ + this->nu())
187 );
188 }
189
190 //- Correct Eddy-Viscosity and related properties
191 virtual void correct();
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// ************************************************************************* //
surfaceScalarField & phi
Eddy viscosity LES SGS model base class.
Dynamic one equation eddy-viscosity model.
Definition: dynamicKEqn.H:82
BasicTurbulenceModel::alphaField alphaField
Definition: dynamicKEqn.H:140
volScalarField Ce() const
Definition: dynamicKEqn.C:97
BasicTurbulenceModel::rhoField rhoField
Definition: dynamicKEqn.H:141
volScalarField Ck(const volSymmTensorField &D, const volScalarField &KK) const
Calculate Ck by filtering the velocity field U.
Definition: dynamicKEqn.C:43
virtual void correct()
Correct Eddy-Viscosity and related properties.
Definition: dynamicKEqn.C:223
TypeName("dynamicKEqn")
Runtime type information.
virtual ~dynamicKEqn()=default
Destructor.
virtual tmp< fvScalarMatrix > kSource() const
Definition: dynamicKEqn.C:139
BasicTurbulenceModel::transportModel transportModel
Definition: dynamicKEqn.H:142
virtual tmp< volScalarField > k() const
Return SGS kinetic energy.
Definition: dynamicKEqn.H:175
autoPtr< LESfilter > filterPtr_
Definition: dynamicKEqn.H:104
virtual bool read()
Read model coefficients if they have changed.
Definition: dynamicKEqn.C:209
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
Definition: dynamicKEqn.H:181
Abstract class for LES filters.
Definition: LESfilter.H:58
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
Simple top-hat filter used in dynamic LES models.
Definition: simpleFilter.H:56
A class for managing temporary objects.
Definition: tmp.H:65
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
Namespace for OpenFOAM.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
volScalarField & nu
volScalarField & alpha
const dimensionedScalar & D
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73