nutUTabulatedWallFunctionFvPatchScalarField.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 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 
39 {
40  const label patchi = patch().index();
41 
42  const auto& turbModel = db().lookupObject<turbulenceModel>
43  (
45  (
47  internalField().group()
48  )
49  );
50 
51  const scalarField& y = turbModel.y()[patchi];
52 
53  const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
54  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
55  const scalarField magGradU(mag(Uw.snGrad()));
56 
57  const tmp<scalarField> tnuw = turbModel.nu(patchi);
58  const scalarField& nuw = tnuw();
59 
60  return
61  max
62  (
63  scalar(0),
64  sqr(magUp/(calcUPlus(magUp*y/nuw) + ROOTVSMALL))
65  /(magGradU + ROOTVSMALL)
66  - nuw
67  );
68 }
69 
70 
73 (
74  const scalarField& Rey
75 ) const
76 {
77  auto tuPlus = tmp<scalarField>::New(patch().size(), Zero);
78  auto& uPlus = tuPlus.ref();
79 
80  forAll(uPlus, facei)
81  {
82  uPlus[facei] = uPlusTable_.interpolateLog10(Rey[facei]);
83  }
84 
85  return tuPlus;
86 }
87 
88 
90 (
91  Ostream& os
92 ) const
93 {
94  os.writeEntry("uPlusTable", uPlusTableName_);
95 }
96 
97 
98 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
99 
102 (
103  const fvPatch& p,
105 )
106 :
108  uPlusTableName_("undefined-uPlusTableName"),
109  uPlusTable_
110  (
111  IOobject
112  (
113  uPlusTableName_,
114  patch().boundaryMesh().mesh().time().constant(),
115  patch().boundaryMesh().mesh(),
116  IOobject::NO_READ,
117  IOobject::NO_WRITE,
118  IOobject::NO_REGISTER
119  ),
120  false
121  )
122 {}
123 
124 
127 (
129  const fvPatch& p,
131  const fvPatchFieldMapper& mapper
132 )
133 :
134  nutWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
135  uPlusTableName_(ptf.uPlusTableName_),
136  uPlusTable_(ptf.uPlusTable_)
137 {}
138 
139 
142 (
143  const fvPatch& p,
145  const dictionary& dict
146 )
147 :
149  uPlusTableName_(dict.get<word>("uPlusTable")),
150  uPlusTable_
151  (
152  IOobject
153  (
154  uPlusTableName_,
155  patch().boundaryMesh().mesh().time().constant(),
156  patch().boundaryMesh().mesh(),
157  IOobject::MUST_READ_IF_MODIFIED,
158  IOobject::NO_WRITE,
159  IOobject::NO_REGISTER
160  ),
161  true
162  )
163 {}
164 
165 
168 (
170 )
171 :
173  uPlusTableName_(wfpsf.uPlusTableName_),
174  uPlusTable_(wfpsf.uPlusTable_)
175 {}
176 
177 
180 (
183 )
184 :
186  uPlusTableName_(wfpsf.uPlusTableName_),
187  uPlusTable_(wfpsf.uPlusTable_)
188 {}
189 
190 
191 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
192 
195 {
196  const label patchi = patch().index();
197 
198  const auto& turbModel = db().lookupObject<turbulenceModel>
199  (
201  (
203  internalField().group()
204  )
205  );
206 
207  const scalarField& y = turbModel.y()[patchi];
208 
209  const fvPatchVectorField& Uw = U(turbModel).boundaryField()[patchi];
210  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
211 
212  const tmp<scalarField> tnuw = turbModel.nu(patchi);
213  const scalarField& nuw = tnuw();
214 
215  const scalarField Rey(magUp*y/nuw);
216 
217  return Rey/(calcUPlus(Rey) + ROOTVSMALL);
218 }
219 
220 
222 (
223  Ostream& os
224 ) const
225 {
227  writeLocalEntries(os);
229 }
230 
231 
232 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
233 
234 namespace Foam
235 {
237  (
239  nutUTabulatedWallFunctionFvPatchScalarField
240  );
241 }
242 
243 
244 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
This boundary condition provides a wall constraint on the turbulent viscosity (i.e. nut) based on velocity (i.e. U), for low- and high-Reynolds number applications.
scalar uPlus
dictionary dict
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
virtual const volVectorField & U(const turbulenceModel &turb) const
Helper to return the velocity field either from the turbulence model (default) or the mesh database...
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:213
scalar Rey
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface...
Definition: boundaryMesh.H:58
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
tmp< scalarField > calcUPlus(const scalarField &Rey) const
Calculate wall u+ from table.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
Abstract base class for turbulence models (RAS, LES and laminar).
Macros for easy insertion into run-time selection tables.
scalar magUp
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
constexpr const char *const group
Group name for atomic constants.
scalar y
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
The class nutWallFunction is an abstract base class that hosts calculation methods and common functi...
dynamicFvMesh & mesh
fvPatchField< scalar > fvPatchScalarField
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:63
A FieldMapper for finite-volume patch fields.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch.
Definition: fvPatchField.C:221
OBJstream os(runTime.globalPath()/outputName)
nutUTabulatedWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
U
Definition: pEqn.H:72
virtual tmp< scalarField > calcNut() const
Calculate the turbulent viscosity.
void writeLocalEntries(Ostream &) const
Write local wall function variables.
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.
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127