Loading...
Searching...
No Matches
PopeIndex.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) 2022 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26Class
27 Foam::resolutionIndexModels::PopeIndex
28
29Description
30 Computes a single-mesh resolution index according to Pope's index,
31 which is used as a LES/DES quality/post-verification metric that does
32 not require any experimental or DNS data.
33
34 \f[
35 \Gamma_{Pope}(\mathbf{x}, t) = \frac{k_{res}}{k_{tot}}
36 \f]
37
38 with
39
40 \f[
41 k_{tot} = k_{res} + k_{sgs} + |k_{num}|
42 \f]
43
44 where
45 \vartable
46 \Gamma_{Pope}(\mathbf{x}, t) | Pope's index [-]
47 k_{tot} | Total turbulent kinetic energy [m^2/s^2]
48 k_{res} | Resolved turbulent kinetic energy [m^2/s^2]
49 k_{sgs} | Subgrid-scale turbulent kinetic energy [m^2/s^2]
50 k_{num} | Numerical turbulent kinetic energy [m^2/s^2]
51 \endvartable
52
53 Inclusion of \f$k_{num}\f$ is optional, and set as \c true by default:
54
55 \f[
56 k_{num} = C_n \left(\frac{h}{\Delta}\right)^2 k_{sgs}
57 \f]
58
59 where
60 \vartable
61 C_n | Empirical constant [-]
62 h | Characteristic length scale with \f$h = V^{1/3} \f$ [m]
63 V | Cell volume [m^3]
64 \Delta | Filter length scale [m]
65 \endvartable
66
67 Typical criterion for acceptable-quality resolution:
68
69 \f[
70 \Gamma_{Pope}(\mathbf{x}) \geq 0.8
71 \f]
72
73 Required fields:
74 \verbatim
75 U | Velocity [m/s]
76 UMean | Mean velocity [m/s]
77 k | Subgrid-scale turbulent kinetic energy [m^2/s^2]
78 delta | Filter length [m]
79 \endverbatim
80
81 References:
82 \verbatim
83 Governing equation (tag:P):
84 Pope, S. B. (2000).
85 Turbulent flows.
86 Cambridge, UK: Cambridge Univ. Press
87 DOI:10.1017/CBO9780511840531
88
89 Governing equations for the denominator kNum term (tag:CKJ):
90 Celik, I., Klein, M., & Janicka, J. (2009).
91 Assessment measures for engineering LES applications.
92 Journal of fluids engineering, 131(3).
93 DOI:10.1115/1.3059703
94 \endverbatim
95
96Usage
97 Minimal example by using \c system/controlDict.functions:
98 \verbatim
99 resolutionIndexFO
100 {
101 // Inherited entries
102 ...
103 model PopeIndex;
104
105 // Optional entries
106 U <word>;
107 UMean <word>;
108 k <word>;
109 delta <word>;
110 includeKnum <bool>;
111
112 // Conditional entries
113 // when includeKnum = true
114 Cn <scalar>;
115 }
116 \endverbatim
117
118 where the entries mean:
119 \table
120 Property | Description | Type | Reqd | Deflt
121 model | Model name: PopeIndex | word | yes | -
122 U | Name of velocity field | word | no | U
123 UMean | Name of mean velocity field | word | no | UMean
124 k | Name of subgrid-scale turbulent kinetic energy field <!--
125 --> | word | no | k
126 delta | Name of filter field | word | no | delta
127 includeKnum | Flag to include numerical k field | bool | no | true
128 Cn | Empirical constant | choice | no | 1.0
129 \endtable
130
131Note
132 - Some studies measured \f$\Gamma_{Pope} > 1\f$ with \f$k_{res}\f$ comparisons
133 between a LES and a corresponding filtered DNS. Nonphysical results may
134 occur.
135
136SourceFiles
137 PopeIndex.C
138
139\*---------------------------------------------------------------------------*/
140
141#ifndef Foam_resolutionIndexModels_PopeIndex_H
142#define Foam_resolutionIndexModels_PopeIndex_H
143
144#include "resolutionIndexModel.H"
145
146// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
147
148namespace Foam
149{
150namespace resolutionIndexModels
151{
152
153/*---------------------------------------------------------------------------*\
154 Class PopeIndex Declaration
155\*---------------------------------------------------------------------------*/
156
157class PopeIndex
158:
159 public resolutionIndexModel
160{
161 // Private Data
162
163 //- Flag to include numerical turbulent kinetic energy field
164 bool includeKnum_;
165
166 //- Empirical constant in numerical turbulent kinetic energy estimation
167 scalar Cn_;
168
169 //- Name of velocity field
170 word UName_;
171
172 //- Name of mean velocity field
173 word UMeanName_;
174
175 //- Name of subgrid-scale turbulent kinetic energy field
176 word kName_;
177
178 //- Name of filter field
179 word deltaName_;
180
181
182 // Private Member Functions
183
184 //- Return numerical turbulent kinetic energy field
185 tmp<volScalarField> kNum() const;
186
187
188public:
189
190 //- Runtime type information
191 TypeName("PopeIndex");
192
193
194 // Constructors
195
196 //- Construct from components
198 (
199 const word& name,
200 const fvMesh& mesh,
201 const dictionary& dict
202 );
203
204 //- No copy construct
205 PopeIndex(const PopeIndex&) = delete;
206
207 //- No copy assignment
208 void operator=(const PopeIndex&) = delete;
209
210
211 // Destructor
212 virtual ~PopeIndex() = default;
213
214
215 // Member Functions
216
217 //- Read the function-object dictionary
218 virtual bool read(const dictionary& dict);
219
220 //- Execute the function-object operations
221 virtual bool execute();
222
223 //- Write the function-object results
224 virtual bool write();
225};
226
227
228// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
229
230} // End namespace resolutionIndexModels
231} // End namespace Foam
232
233// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
234
235#endif
236
237// ************************************************************************* //
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const fvMesh & mesh() const noexcept
Return const reference to the mesh.
resolutionIndexModel(const word &name, const fvMesh &mesh, const dictionary &dict)
Construct from components.
Computes a single-mesh resolution index according to Pope's index, which is used as a LES/DES quality...
Definition PopeIndex.H:241
PopeIndex(const word &name, const fvMesh &mesh, const dictionary &dict)
Construct from components.
Definition PopeIndex.C:55
void operator=(const PopeIndex &)=delete
No copy assignment.
virtual bool read(const dictionary &dict)
Read the function-object dictionary.
Definition PopeIndex.C:75
TypeName("PopeIndex")
Runtime type information.
virtual bool execute()
Execute the function-object operations.
Definition PopeIndex.C:96
PopeIndex(const PopeIndex &)=delete
No copy construct.
virtual bool write()
Write the function-object results.
Definition PopeIndex.C:127
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
auto & name
A namespace for various resolutionIndex model implementations.
Namespace for OpenFOAM.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68