nutUTabulatedWallFunctionFvPatchScalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2020 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 #include "momentumTransportModel.H"
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
38 
40 {
41  const label patchi = patch().index();
42 
43  const momentumTransportModel& turbModel =
44  db().lookupObject<momentumTransportModel>
45  (
47  (
48  momentumTransportModel::typeName,
49  internalField().group()
50  )
51  );
52  const scalarField& y = turbModel.y()[patchi];
53  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
54  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
55  const scalarField magGradU(mag(Uw.snGrad()));
56  const tmp<scalarField> tnuw = turbModel.nu(patchi);
57  const scalarField& nuw = tnuw();
58 
59  return
60  max
61  (
62  scalar(0),
63  sqr(magUp/(calcUPlus(magUp*y/nuw) + rootVSmall))
64  /(magGradU + rootVSmall)
65  - nuw
66  );
67 }
68 
69 
71 (
72  const scalarField& Rey
73 ) const
74 {
75  tmp<scalarField> tuPlus(new scalarField(patch().size(), 0.0));
76  scalarField& uPlus = tuPlus.ref();
77 
78  forAll(uPlus, facei)
79  {
80  uPlus[facei] = uPlusTable_.interpolateLog10(Rey[facei]);
81  }
82 
83  return tuPlus;
84 }
85 
86 
87 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
88 
91 (
92  const fvPatch& p,
94 )
95 :
97  uPlusTableName_("undefined-uPlusTableName"),
99  (
100  IOobject
101  (
103  patch().boundaryMesh().mesh().time().constant(),
104  patch().boundaryMesh().mesh(),
107  false
108  ),
109  false
110  )
111 {}
112 
113 
116 (
118  const fvPatch& p,
120  const fvPatchFieldMapper& mapper
121 )
122 :
123  nutWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
126 {}
127 
128 
131 (
132  const fvPatch& p,
134  const dictionary& dict
135 )
136 :
138  uPlusTableName_(dict.lookup("uPlusTable")),
140  (
141  IOobject
142  (
144  patch().boundaryMesh().mesh().time().constant(),
145  patch().boundaryMesh().mesh(),
148  false
149  ),
150  true
151  )
152 {}
153 
154 
157 (
159 )
160 :
163  uPlusTable_(wfpsf.uPlusTable_)
164 {}
165 
166 
169 (
172 )
173 :
176  uPlusTable_(wfpsf.uPlusTable_)
177 {}
178 
179 
180 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
181 
183 {
184  const label patchi = patch().index();
185 
186  const momentumTransportModel& turbModel =
187  db().lookupObject<momentumTransportModel>
188  (
190  (
191  momentumTransportModel::typeName,
192  internalField().group()
193  )
194  );
195  const scalarField& y = turbModel.y()[patchi];
196  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
197  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
198  const tmp<scalarField> tnuw = turbModel.nu(patchi);
199  const scalarField& nuw = tnuw();
200  const scalarField Rey(magUp*y/nuw);
201 
202  return Rey/(calcUPlus(Rey) + rootVSmall);
203 }
204 
205 
207 {
209  writeEntry(os, "uPlusTable", uPlusTableName_);
210  writeEntry(os, "value", *this);
211 }
212 
213 
214 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
215 
217 (
220 );
221 
222 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
223 
224 } // End namespace Foam
225 
226 // ************************************************************************* //
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions...
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
const char *const group
Group name for atomic constants.
scalar uPlus
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:186
scalar Rey
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
const Boundary & boundaryField() const
Return const-reference to the boundary field.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
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:59
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
const nearWallDist & y() const
Return the near wall distances.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:280
Macros for easy insertion into run-time selection tables.
scalar magUp
scalar y
Type interpolateLog10(scalar x) const
Interpolate - takes log10 flag into account.
dynamicFvMesh & mesh
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static word groupName(Name name, const word &group)
Foam::fvPatchFieldMapper.
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:194
nutUTabulatedWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
virtual tmp< scalarField > calcUPlus(const scalarField &Rey) const
Calculate wall u+ from table.
Abstract base class for turbulence models (RAS, LES and laminar).
const volVectorField & U() const
Access function to velocity field.
label patchi
virtual tmp< scalarField > nut() const
Calculate the turbulence viscosity.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
dimensioned< scalar > mag(const dimensioned< Type > &)
nutWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
A class for managing temporary objects.
Definition: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:812