Loading...
Searching...
No Matches
CelikNuIndex.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::CelikNuIndex
28
29Description
30 Computes a single-mesh resolution index according to Celik et al.'s index
31 using effective viscosity, which is used as a LES/DES quality/post
32 verification metric that does not require any experimental or DNS data.
33
34 \f[
35 \Gamma_{Celik,\nu}(\mathbf{x}, t) =
36 \frac{1}{1 + \alpha_\nu \left(\frac{\nu_{eff}}{\nu}\right)^n}
37 \f]
38
39 with
40
41 \f[
42 \nu_{eff} = \nu_{num} + \nu_{sgs} + \nu
43 \f]
44
45 \f[
46 \nu_{num} = {sgn}(k_{num}) C_\nu \Delta \sqrt{k_{num}}
47 \f]
48
49 \f[
50 k_{num} = C_n \left(\frac{h}{\Delta}\right)^2 k_{sgs}
51 \f]
52
53 where
54 \vartable
55 \Gamma_{Celik,\nu}(\mathbf{x}, t) | Celik et al.'s index [-]
56 \alpha_\nu | Empirical constant [-]
57 \nu_{eff} | Effective eddy viscosity [m^2/s]
58 \nu_{num} | Numerical eddy viscosity [m^2/s]
59 \nu_{sgs} | Subgrid-scale eddy viscosity [m^2/s]
60 \nu | Kinematic viscosity of fluid [m^2/s]
61 n | Empirical exponent [-]
62 k_{num} | Numerical turbulent kinetic energy [m^2/s^2]
63 C_\nu | Empirical constant [-]
64 \Delta | Filter length scale [m]
65 C_n | Empirical constant [-]
66 h | Characteristic length scale with \f$h = V^{1/3} \f$ [m]
67 V | Cell volume [m^3]
68 k_{sgs} | Subgrid-scale turbulent kinetic energy [m^2/s^2]
69 \endvartable
70
71 Criterion for acceptable-quality resolution:
72
73 \f[
74 \Gamma_{Celik,\nu}(\mathbf{x}) \geq 0.8
75 \f]
76
77 Required fields:
78 \verbatim
79 k | Subgrid scale turbulent kinetic energy [m^2/s^2]
80 delta | Filter length [m]
81 nu | Kinematic viscosity [m^2/s]
82 nut | Subgrid-scale eddy viscosity [m^2/s]
83 \endverbatim
84
85 References:
86 \verbatim
87 Governing equations (tag:CCY):
88 Celik, I. B., Cehreli Z. N., Yavuz I. (2005).
89 Index of resolution quality for large eddy simulations.
90 Journal of Fluids Engineering. 127:949–958.
91 DOI:10.1115/1.1990201
92
93 Governing equations (tag:CKJ):
94 Celik, I., Klein, M., & Janicka, J. (2009).
95 Assessment measures for engineering LES applications.
96 Journal of fluids engineering, 131(3).
97 DOI:10.1115/1.3059703
98 \endverbatim
99
100Usage
101 Minimal example by using \c system/controlDict.functions:
102 \verbatim
103 resolutionIndexFO
104 {
105 // Inherited entries
106 ...
107 model CelikNuIndex;
108
109 // Optional entries
110 alphaNu <scalar>;
111 n <scalar>;
112 Cnu <scalar>;
113 Cn <scalar>;
114 k <word>;
115 delta <word>;
116 nu <word>;
117 nut <word>;
118 }
119 \endverbatim
120
121 where the entries mean:
122 \table
123 Property | Description | Type | Reqd | Deflt
124 model | Model name: CelikNuIndex | word | yes | -
125 alphaNu | Empirical constant | scalar | no | 0.05
126 n | Empirical exponent | scalar | no | 0.53
127 Cnu | Empirical constant | scalar | no | 0.1
128 Cn | Empirical constant | scalar | no | 1.0
129 k | Name of subgrid-scale turbulent kinetic energy field <!--
130 --> | word | no | k
131 delta | Name of filter field | word | no | delta
132 nu | Name of kinematic viscosity field | word | no | nu
133 nut | Name of subgrid-scale eddy viscosity field | word | no | nut
134 \endtable
135
136SourceFiles
137 CelikNuIndex.C
138
139\*---------------------------------------------------------------------------*/
140
141#ifndef Foam_resolutionIndexModels_CelikNuIndex_H
142#define Foam_resolutionIndexModels_CelikNuIndex_H
143
144#include "resolutionIndexModel.H"
145
146// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
147
148namespace Foam
149{
150namespace resolutionIndexModels
151{
152
153/*---------------------------------------------------------------------------*\
154 Class CelikNuIndex Declaration
155\*---------------------------------------------------------------------------*/
156
157class CelikNuIndex
158:
159 public resolutionIndexModel
160{
161 // Private Data
162
163 //- Empirical constant
164 scalar alphaNu_;
165
166 //- Empirical exponent
167 scalar n_;
168
169 //- Empirical constant
170 scalar Cnu_;
171
172 //- Empirical constant
173 scalar Cn_;
174
175 //- Name of subgrid-scale turbulent kinetic energy field
176 word kName_;
177
178 //- Name of filter field
179 word deltaName_;
180
181 //- Name of kinematic viscosity field
182 word nuName_;
183
184 //- Name of subgrid-scale eddy viscosity field
185 word nutName_;
186
187
188 // Private Member Functions
189
190 //- Return numerical eddy viscosity field
191 tmp<volScalarField> nuNum() const;
192
193 //- Return numerical turbulent kinetic energy field
194 tmp<volScalarField> kNum() const;
195
196
197public:
198
199 //- Runtime type information
200 TypeName("CelikNuIndex");
201
202
203 // Constructors
204
205 //- Construct from components
207 (
208 const word& name,
209 const fvMesh& mesh,
210 const dictionary& dict
211 );
212
213 //- No copy construct
214 CelikNuIndex(const CelikNuIndex&) = delete;
215
216 //- No copy assignment
217 void operator=(const CelikNuIndex&) = delete;
218
219
220 // Destructor
221 virtual ~CelikNuIndex() = default;
222
223
224 // Member Functions
225
226 //- Read the function-object dictionary
227 virtual bool read(const dictionary& dict);
228
229 //- Execute the function-object operations
230 virtual bool execute();
231
232 //- Write the function-object results
233 virtual bool write();
234};
235
236
237// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
238
239} // End namespace resolutionIndexModels
240} // End namespace Foam
241
242// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
243
244#endif
245
246// ************************************************************************* //
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 Celik et al.'s index using effective viscosity,...
CelikNuIndex(const CelikNuIndex &)=delete
No copy construct.
CelikNuIndex(const word &name, const fvMesh &mesh, const dictionary &dict)
Construct from components.
virtual bool read(const dictionary &dict)
Read the function-object dictionary.
TypeName("CelikNuIndex")
Runtime type information.
void operator=(const CelikNuIndex &)=delete
No copy assignment.
virtual bool execute()
Execute the function-object operations.
virtual bool write()
Write the function-object results.
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