PopeIndex.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 "PopeIndex.H"
29 #include "resolutionIndexModel.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace resolutionIndexModels
37 {
38  defineTypeNameAndDebug(PopeIndex, 0);
39  addToRunTimeSelectionTable(resolutionIndexModel, PopeIndex, dictionary);
40 }
41 }
42 
43 
44 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 
47 Foam::resolutionIndexModels::PopeIndex::kNum() const
48 {
49  const auto& kSgs = getOrReadField<volScalarField>(kName_);
50  const auto& Delta = getOrReadField<volScalarField>(deltaName_);
51 
52  tmp<volScalarField> th = cbrt(V());
53 
54  // (CKJ:Eq. 28)
55  return Cn_*sqr(th/Delta)*kSgs;
56 }
57 
58 
59 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
60 
62 (
63  const word& name,
64  const fvMesh& mesh,
65  const dictionary& dict
66 )
67 :
69  includeKnum_(),
70  Cn_(),
71  UName_(),
72  UMeanName_(),
73  kName_(),
74  deltaName_()
75 {
76  read(dict);
77 }
78 
79 
80 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
81 
83 {
85  {
86  return false;
87  }
88 
89  includeKnum_ = dict.getOrDefault<bool>("includeKnum", true);
90  if (includeKnum_)
91  {
92  Cn_ = dict.getOrDefault<scalar>("Cnu", 1.0);
93  }
94  UName_ = dict.getOrDefault<word>("U", "U");
95  UMeanName_ = dict.getOrDefault<word>("UMean", "UMean");
96  kName_ = dict.getOrDefault<word>("k", "k");
97  deltaName_ = dict.getOrDefault<word>("delta", "delta");
98 
99  return true;
100 }
101 
102 
104 {
105  // Calculate resolved k field
106  const auto& U = getOrReadField<volVectorField>(UName_);
107  const auto& UMean = getOrReadField<volVectorField>(UMeanName_);
108  const volScalarField kRes(0.5*magSqr(U - UMean));
109 
110  // Calculate subgrid-scale k field
111  const auto& kSgs = getOrReadField<volScalarField>(kName_);
112 
113  // Calculate total k field
114  tmp<volScalarField> tkTot = kRes + kSgs;
115  if (includeKnum_)
116  {
117  tkTot.ref() += mag(kNum());
118  }
119 
120 
121  // Calculate index field
122  auto& index = getOrReadField<volScalarField>(resultName());
123 
124  const dimensionedScalar kMin(kSgs.dimensions(), SMALL);
125 
126  // (P:p. 560;CKJ:Eq. 11)
127  index = kRes/max(kMin, tkTot);
128  index.correctBoundaryConditions();
129 
130  return true;
131 }
132 
133 
135 {
136  const auto& index = getOrReadField<volScalarField>(resultName());
137 
138  Info<< tab << "writing field:" << index.name() << endl;
139 
140  index.write();
141 
142  return true;
143 }
144 
145 
146 // ************************************************************************* //
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
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition: tmpI.H:235
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual bool read(const dictionary &dict)
Read top-level dictionary.
Definition: PopeIndex.C:75
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
constexpr char tab
The tab &#39;\t&#39; character(0x09)
Definition: Ostream.H:49
volVectorField UMean(UMeanHeader, mesh)
Macros for easy insertion into run-time selection tables.
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:127
dynamicFvMesh & mesh
virtual bool execute()
Calculate the result field.
Definition: PopeIndex.C:96
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
PopeIndex(const word &name, const fvMesh &mesh, const dictionary &dict)
Construct from components.
Definition: PopeIndex.C:55
A class for handling words, derived from Foam::string.
Definition: word.H:63
defineTypeNameAndDebug(CelikEtaIndex, 0)
virtual bool write()
Write the result field.
Definition: PopeIndex.C:127
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.
U
Definition: pEqn.H:72
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
tmp< volScalarField > V() const
Return cell volume field.
addToRunTimeSelectionTable(resolutionIndexModel, CelikEtaIndex, dictionary)
messageStream Info
Information stream (stdout output on master, null elsewhere)
A class for managing temporary objects.
Definition: HashPtrTable.H:50
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Namespace for OpenFOAM.
virtual bool read(const dictionary &dict)
Read top-level dictionary.