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 -------------------------------------------------------------------------------
11 License
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"
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 
59 calcNut() 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 
144 (
145  const fvPatch& p,
147 )
148 :
150  Ks_(p.size(), Zero),
151  Cs_(p.size(), Zero)
152 {}
153 
154 
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 
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 
186 (
188 )
189 :
191  Ks_(rwfpsf.Ks_),
192  Cs_(rwfpsf.Cs_)
193 {}
194 
195 
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 {
216  nutkWallFunctionFvPatchScalarField::autoMap(m);
217  Ks_.autoMap(m);
218  Cs_.autoMap(m);
219 }
220 
221 
223 (
224  const fvPatchScalarField& ptf,
225  const labelList& addr
226 )
227 {
228  nutkWallFunctionFvPatchScalarField::rmap(ptf, addr);
229 
230  const auto& nrwfpsf =
231  refCast<const nutkRoughWallFunctionFvPatchScalarField>(ptf);
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);
246 }
247 
248 
249 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
250 
251 namespace Foam
252 {
254  (
256  nutkRoughWallFunctionFvPatchScalarField
257  );
258 }
259 
260 // ************************************************************************* //
dictionary dict
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
dimensionedScalar log(const dimensionedScalar &ds)
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
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
dimensionedScalar pow025(const dimensionedScalar &ds)
This boundary condition provides a wall-function for the turbulent viscosity (i.e. nut) when using wall functions for rough walls, based on the turbulent kinetic energy (i.e. k). The condition manipulates the wall roughness parameter (i.e. E) to account for roughness effects.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
label k
Boltzmann constant.
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
constexpr const char *const group
Group name for atomic constants.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
scalar y
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
void writeLocalEntries(Ostream &os) const
Write local wall function variables.
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
This boundary condition provides a wall function for the turbulent viscosity (i.e. nut) based on the turbulent kinetic energy, (i.e. k) for for low- and high-Reynolds number applications.
fvPatchField< scalar > fvPatchScalarField
static const word propertiesName
Default name of the turbulence properties dictionary.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A FieldMapper for finite-volume patch fields.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
dimensionedScalar sin(const dimensionedScalar &ds)
nutkRoughWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
int debug
Static debugging option.
OBJstream os(runTime.globalPath()/outputName)
virtual tmp< scalarField > calcNut() const
Calculate the turbulent viscosity.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
scalar yPlus
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
messageStream Info
Information stream (stdout output on master, null elsewhere)
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
virtual scalar fnRough(const scalar KsPlus, const scalar Cs) const
Compute the roughness function.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
const scalarField & Cs
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127