nutkRoughWallFunctionFvPatchScalarField.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) 2011-2016, 2019 OpenFOAM Foundation
9 Copyright (C) 2019-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
30#include "turbulenceModel.H"
31#include "fvPatchFieldMapper.H"
32#include "volFields.H"
34
35// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
36
38(
39 const scalar KsPlus,
40 const scalar Cs
41) const
42{
43 // Return fn based on non-dimensional roughness height
44
45 if (KsPlus < 90.0)
46 {
47 return pow
48 (
49 (KsPlus - 2.25)/87.75 + Cs*KsPlus,
50 sin(0.4258*(log(KsPlus) - 0.811))
51 );
52 }
53
54 return (1.0 + Cs*KsPlus);
55}
56
57
59calcNut() const
60{
61 const label patchi = patch().index();
62
63 const auto& turbModel = db().lookupObject<turbulenceModel>
64 (
66 (
68 internalField().group()
69 )
70 );
71
72 const scalarField& y = turbModel.y()[patchi];
73
74 const tmp<volScalarField> tk = turbModel.k();
75 const volScalarField& k = tk();
76
77 const tmp<scalarField> tnuw = turbModel.nu(patchi);
78 const scalarField& nuw = tnuw();
79
80 const scalar Cmu25 = pow025(wallCoeffs_.Cmu());
81 const scalar kappa = wallCoeffs_.kappa();
82 const scalar E = wallCoeffs_.E();
83
84 auto tnutw = tmp<scalarField>::New(*this);
85 auto& nutw = tnutw.ref();
86
87 forAll(nutw, facei)
88 {
89 const label celli = patch().faceCells()[facei];
90
91 const scalar uStar = Cmu25*sqrt(k[celli]);
92 const scalar yPlus = uStar*y[facei]/nuw[facei];
93 const scalar KsPlus = uStar*Ks_[facei]/nuw[facei];
94
95 scalar Edash = E;
96 if (2.25 < KsPlus)
97 {
98 Edash /= fnRough(KsPlus, Cs_[facei]);
99 }
100
101 const scalar limitingNutw = max(nutw[facei], nuw[facei]);
102
103 // To avoid oscillations limit the change in the wall viscosity
104 // which is particularly important if it temporarily becomes zero
105 nutw[facei] =
106 max
107 (
108 min
109 (
110 nuw[facei]
111 *(yPlus*kappa/log(max(Edash*yPlus, 1+1e-4)) - 1),
112 2*limitingNutw
113 ), 0.5*limitingNutw
114 );
115
116 if (debug)
117 {
118 Info<< "yPlus = " << yPlus
119 << ", KsPlus = " << KsPlus
120 << ", Edash = " << Edash
121 << ", nutw = " << nutw[facei]
122 << endl;
123 }
124 }
125
126 return tnutw;
127}
128
129
131(
132 Ostream& os
133) const
134{
135 Cs_.writeEntry("Cs", os);
136 Ks_.writeEntry("Ks", os);
137}
138
139
140// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
141
142Foam::nutkRoughWallFunctionFvPatchScalarField::
143nutkRoughWallFunctionFvPatchScalarField
144(
145 const fvPatch& p,
147)
148:
150 Ks_(p.size(), Zero),
151 Cs_(p.size(), Zero)
152{}
153
154
155Foam::nutkRoughWallFunctionFvPatchScalarField::
156nutkRoughWallFunctionFvPatchScalarField
157(
159 const fvPatch& p,
161 const fvPatchFieldMapper& mapper
162)
163:
164 nutkWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
165 Ks_(ptf.Ks_, mapper),
166 Cs_(ptf.Cs_, mapper)
167{}
168
169
170Foam::nutkRoughWallFunctionFvPatchScalarField::
171nutkRoughWallFunctionFvPatchScalarField
172(
173 const fvPatch& p,
175 const dictionary& dict
176)
177:
179 Ks_("Ks", dict, p.size()),
180 Cs_("Cs", dict, p.size())
181{}
182
183
184Foam::nutkRoughWallFunctionFvPatchScalarField::
185nutkRoughWallFunctionFvPatchScalarField
186(
188)
189:
191 Ks_(rwfpsf.Ks_),
192 Cs_(rwfpsf.Cs_)
193{}
194
195
196Foam::nutkRoughWallFunctionFvPatchScalarField::
197nutkRoughWallFunctionFvPatchScalarField
198(
201)
202:
204 Ks_(rwfpsf.Ks_),
205 Cs_(rwfpsf.Cs_)
206{}
207
208
209// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
210
212(
213 const fvPatchFieldMapper& m
214)
215{
217 Ks_.autoMap(m);
218 Cs_.autoMap(m);
219}
220
221
223(
224 const fvPatchScalarField& ptf,
225 const labelList& addr
226)
227{
229
230 const auto& nrwfpsf =
231 refCast<const nutkRoughWallFunctionFvPatchScalarField>(ptf);
232
233 Ks_.rmap(nrwfpsf.Ks_, addr);
234 Cs_.rmap(nrwfpsf.Cs_, addr);
235}
236
237
239(
240 Ostream& os
241) const
242{
244 writeLocalEntries(os);
245 writeEntry("value", os);
246}
247
248
249// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
250
251namespace Foam
252{
254 (
257 );
258}
259
260// ************************************************************************* //
scalar y
label k
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
void rmap(const atmBoundaryLayer &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
virtual bool write()
Write the output fields.
A FieldMapper for finite-volume patch fields.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
This boundary condition provides a wall-function for the turbulent viscosity (i.e....
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual scalar fnRough(const scalar KsPlus, const scalar Cs) const
Compute the roughness function.
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
void writeLocalEntries(Ostream &os) const
Write local wall function variables.
virtual tmp< scalarField > calcNut() const
Calculate the turbulent viscosity.
This boundary condition provides a wall function for the turbulent viscosity (i.e....
A class for managing temporary objects.
Definition: tmp.H:65
Abstract base class for turbulence models (RAS, LES and laminar).
static const word propertiesName
Default name of the turbulence properties dictionary.
volScalarField & p
scalar yPlus
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Definition: fvPatchField.H:676
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensionedScalar sin(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere)
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensionedScalar sqrt(const dimensionedScalar &ds)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
dimensionedScalar pow025(const dimensionedScalar &ds)
dictionary dict
const scalarField & Cs
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333