CelikEtaIndex.C
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 -------------------------------------------------------------------------------
10 License
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 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "CelikEtaIndex.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace resolutionIndexModels
37 {
40 }
41 }
42 
43 
44 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 
47 Foam::resolutionIndexModels::CelikEtaIndex::eta() const
48 {
49  const auto& nu = getOrReadField<volScalarField>(nuName_);
50  tmp<volScalarField> tepsilon = epsilon();
51 
52  const dimensionedScalar epsilonMin(tepsilon.cref().dimensions(), SMALL);
53 
54  // (CKJ:Eq. 23)
55  return pow025(pow3(nu)/max(epsilonMin, tepsilon));
56 }
57 
58 
60 Foam::resolutionIndexModels::CelikEtaIndex::epsilon() const
61 {
62  const auto& kSgs = getOrReadField<volScalarField>(kName_);
63  const auto& Delta = getOrReadField<volScalarField>(deltaName_);
64  tmp<volScalarField> tnuEff = nuEff();
65 
66  // (Derived based on CKJ:Eq. 25-26, p.031102-5)
67  return tnuEff*kSgs/(Ck_*sqr(Delta));
68 }
69 
70 
72 Foam::resolutionIndexModels::CelikEtaIndex::nuEff() const
73 {
74  const auto& nu = getOrReadField<volScalarField>(nuName_);
75  const auto& nuSgs = getOrReadField<volScalarField>(nutName_);
76  tmp<volScalarField> tnuNum = nuNum();
77 
78  // (CKJ:p. 031102-3)
79  return tnuNum + nuSgs + nu;
80 }
81 
82 
84 Foam::resolutionIndexModels::CelikEtaIndex::nuNum() const
85 {
86  const auto& Delta = getOrReadField<volScalarField>(deltaName_);
87 
88  tmp<volScalarField> tkNum = kNum();
89 
90  // (CKJ:Eq. 35)
91  return sign(tkNum.cref())*Cnu_*Delta*sqrt(mag(tkNum.cref()));
92 }
93 
94 
96 Foam::resolutionIndexModels::CelikEtaIndex::kNum() const
97 {
98  const auto& kSgs = getOrReadField<volScalarField>(kName_);
99  const auto& Delta = getOrReadField<volScalarField>(deltaName_);
100 
101  tmp<volScalarField> th = cbrt(V());
102 
103  // (CKJ:Eq. 28)
104  return Cn_*sqr(th/Delta)*kSgs;
105 }
106 
107 
108 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
109 
111 (
112  const word& name,
113  const fvMesh& mesh,
114  const dictionary& dict
115 )
116 :
118  alphaEta_(),
119  m_(),
120  Cnu_(),
121  Cn_(),
122  Ck_(),
123  kName_(),
124  deltaName_(),
125  nuName_(),
126  nutName_()
127 {
128  read(dict);
129 }
130 
131 
132 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
133 
135 {
137  {
138  return false;
139  }
140 
141  // (Default values from CKJ:p. 031102-3, 031102-5)
142  alphaEta_ = dict.getOrDefault<scalar>("alphaEta", 0.05);
143  m_ = dict.getOrDefault<scalar>("m", 0.5);
144  Cnu_ = dict.getOrDefault<scalar>("Cnu", 0.1);
145  Cn_ = dict.getOrDefault<scalar>("Cn", 1.0);
146  Ck_ = dict.getOrDefault<scalar>("Ck", 0.0376);
147  kName_ = dict.getOrDefault<word>("k", "k");
148  deltaName_ = dict.getOrDefault<word>("delta", "delta");
149  nuName_ = dict.getOrDefault<word>("nu", "nu");
150  nutName_ = dict.getOrDefault<word>("nut", "nut");
151 
152  return true;
153 }
154 
155 
157 {
158  // Calculate Kolmogorov and mesh length scales
159  tmp<volScalarField> teta = eta();
160  tmp<volScalarField> th = cbrt(V());
161 
162  // Calculate index field
163  auto& index = getOrReadField<volScalarField>(resultName());
164 
165  // (CKJ:Eq. 9)
166  index = 1.0/(1.0 + alphaEta_*pow(th/teta, m_));
167  index.correctBoundaryConditions();
168 
169  return true;
170 }
171 
172 
174 {
175  const auto& index = getOrReadField<volScalarField>(resultName());
176 
177  Info<< tab << "writing field:" << index.name() << endl;
178 
179  index.write();
180 
181  return true;
182 }
183 
184 
185 // ************************************************************************* //
virtual bool write()
Write the result field.
dimensionedScalar sign(const dimensionedScalar &ds)
dictionary dict
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
virtual bool execute()
Calculate the result field.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
dimensionedScalar pow025(const dimensionedScalar &ds)
constexpr char tab
The tab &#39;\t&#39; character(0x09)
Definition: Ostream.H:49
Macros for easy insertion into run-time selection tables.
CelikEtaIndex(const word &name, const fvMesh &mesh, const dictionary &dict)
Construct from components.
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:127
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
A class for handling words, derived from Foam::string.
Definition: word.H:63
defineTypeNameAndDebug(CelikEtaIndex, 0)
dimensionedScalar cbrt(const dimensionedScalar &ds)
A base class for resolutionIndex models.
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
virtual bool read(const dictionary &dict)
Read top-level dictionary.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
addToRunTimeSelectionTable(resolutionIndexModel, CelikEtaIndex, dictionary)
messageStream Info
Information stream (stdout output on master, null elsewhere)
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Computes a single-mesh resolution index according to Celik et al.&#39;s index using Kolmogorov length sca...
volScalarField & nu
Namespace for OpenFOAM.
virtual bool read(const dictionary &dict)
Read top-level dictionary.