Loading...
Searching...
No Matches
CelikNuIndex.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-------------------------------------------------------------------------------
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
26\*---------------------------------------------------------------------------*/
27
28#include "CelikNuIndex.H"
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34namespace Foam
35{
37{
40}
41}
42
43
44// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45
47Foam::resolutionIndexModels::CelikNuIndex::nuNum() const
48{
49 const auto& Delta = getOrReadField<volScalarField>(deltaName_);
50
51 tmp<volScalarField> tkNum = kNum();
52
53 // (CKJ:Eq. 35)
54 return sign(tkNum.cref())*Cnu_*Delta*sqrt(tkNum.cref());
55}
56
57
59Foam::resolutionIndexModels::CelikNuIndex::kNum() const
60{
61 const auto& kSgs = getOrReadField<volScalarField>(kName_);
62 const auto& Delta = getOrReadField<volScalarField>(deltaName_);
63
64 tmp<volScalarField> th = cbrt(V());
65
66 // (CKJ:Eq. 28)
67 return Cn_*sqr(th/Delta)*kSgs;
68}
69
70
71// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
72
74(
75 const word& name,
76 const fvMesh& mesh,
77 const dictionary& dict
78)
79:
81 alphaNu_(),
82 n_(),
83 Cnu_(),
84 Cn_(),
85 kName_(),
86 deltaName_(),
87 nuName_(),
88 nutName_()
90 read(dict);
91}
92
93
94// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
95
97{
99 {
100 return false;
101 }
102
103 // (Default values from CKJ:p. 031102-3, 031102-5)
104 alphaNu_ = dict.getOrDefault<scalar>("alphaNu", 0.05);
105 n_ = dict.getOrDefault<scalar>("n", 0.53);
106 Cnu_ = dict.getOrDefault<scalar>("Cnu", 0.1);
107 Cn_ = dict.getOrDefault<scalar>("Cn", 1.0);
108 kName_ = dict.getOrDefault<word>("k", "k");
109 deltaName_ = dict.getOrDefault<word>("delta", "delta");
110 nuName_ = dict.getOrDefault<word>("nu", "nu");
111 nutName_ = dict.getOrDefault<word>("nut", "nut");
112
113 return true;
114}
115
116
118{
119 // Calculate effective eddy viscosity field
120 const auto& nu = getOrReadField<volScalarField>(nuName_);
121 const auto& nuSgs = getOrReadField<volScalarField>(nutName_);
122 tmp<volScalarField> tnuNum = nuNum();
123 tmp<volScalarField> tnuEff = tnuNum + nuSgs + nu;
124
125 // Calculate index field
126 auto& index = getOrReadField<volScalarField>(resultName());
127
128 // (CKJ:Eq. 10)
129 index = 1.0/(1.0 + alphaNu_*pow(tnuEff/nu, n_));
130 index.correctBoundaryConditions();
131
132 return true;
133}
134
135
137{
138 const auto& index = getOrReadField<volScalarField>(resultName());
139
140 Info<< tab << "writing field:" << index.name() << endl;
141
142 index.write();
143
144 return true;
145}
146
147
148// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
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
A base class for resolutionIndex models.
const fvMesh & mesh() const noexcept
Return const reference to the mesh.
resolutionIndexModel(const word &name, const fvMesh &mesh, const dictionary &dict)
Construct from components.
GeoFieldType & getOrReadField(const word &fieldName) const
Return requested field from the object registry or read+register the field to the object registry.
const word & resultName() const noexcept
Return const reference to the result name.
virtual bool read(const dictionary &dict)
Read the function-object dictionary.
Computes a single-mesh resolution index according to Celik et al.'s index using effective viscosity,...
CelikNuIndex(const word &name, const fvMesh &mesh, const dictionary &dict)
Construct from components.
virtual bool read(const dictionary &dict)
Read the function-object dictionary.
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
const T & cref() const
Return const reference to the object or to the contents of a (non-null) managed pointer.
Definition tmpI.H:221
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
const expr V(m.psi().mesh().V())
A namespace for various resolutionIndex model implementations.
Namespace for OpenFOAM.
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition int32.H:127
dimensionedScalar sign(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
messageStream Info
Information stream (stdout output on master, null elsewhere).
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensionedScalar cbrt(const dimensionedScalar &ds)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
constexpr char tab
The tab '\t' character(0x09).
Definition Ostream.H:49
volScalarField & nu
dictionary dict