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 | Copyright (C) 2011-2016 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 "turbulenceModel.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 turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
44  (
46  (
48  internalField().group()
49  )
50  );
51  const scalarField& y = turbModel.y()[patchi];
52  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
53  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
54  const scalarField magGradU(mag(Uw.snGrad()));
55  const tmp<scalarField> tnuw = turbModel.nu(patchi);
56  const scalarField& nuw = tnuw();
57 
58  return
59  max
60  (
61  scalar(0),
62  sqr(magUp/(calcUPlus(magUp*y/nuw) + ROOTVSMALL))
63  /(magGradU + ROOTVSMALL)
64  - nuw
65  );
66 }
67 
68 
70 (
71  const scalarField& Rey
72 ) const
73 {
74  tmp<scalarField> tuPlus(new scalarField(patch().size(), 0.0));
75  scalarField& uPlus = tuPlus.ref();
76 
77  forAll(uPlus, facei)
78  {
79  uPlus[facei] = uPlusTable_.interpolateLog10(Rey[facei]);
80  }
81 
82  return tuPlus;
83 }
84 
85 
86 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
87 
90 (
91  const fvPatch& p,
93 )
94 :
96  uPlusTableName_("undefined-uPlusTableName"),
98  (
99  IOobject
100  (
102  patch().boundaryMesh().mesh().time().constant(),
103  patch().boundaryMesh().mesh(),
106  false
107  ),
108  false
109  )
110 {}
111 
112 
115 (
117  const fvPatch& p,
119  const fvPatchFieldMapper& mapper
120 )
121 :
122  nutWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
125 {}
126 
127 
130 (
131  const fvPatch& p,
133  const dictionary& dict
134 )
135 :
137  uPlusTableName_(dict.lookup("uPlusTable")),
139  (
140  IOobject
141  (
143  patch().boundaryMesh().mesh().time().constant(),
144  patch().boundaryMesh().mesh(),
147  false
148  ),
149  true
150  )
151 {}
152 
153 
156 (
158 )
159 :
162  uPlusTable_(wfpsf.uPlusTable_)
163 {}
164 
165 
168 (
171 )
172 :
175  uPlusTable_(wfpsf.uPlusTable_)
176 {}
177 
178 
179 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
180 
182 {
183  const label patchi = patch().index();
184 
185  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
186  (
188  (
190  internalField().group()
191  )
192  );
193  const scalarField& y = turbModel.y()[patchi];
194  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
195  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
196  const tmp<scalarField> tnuw = turbModel.nu(patchi);
197  const scalarField& nuw = tnuw();
198  const scalarField Rey(magUp*y/nuw);
199 
200  return Rey/(calcUPlus(Rey) + ROOTVSMALL);
201 }
202 
203 
205 {
207  os.writeKeyword("uPlusTable") << uPlusTableName_
208  << token::END_STATEMENT << nl;
209  writeEntry("value", os);
210 }
211 
212 
213 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
214 
216 (
219 );
220 
221 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
222 
223 } // End namespace Foam
224 
225 // ************************************************************************* //
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions...
const char *const group
Group name for atomic constants.
scalar uPlus
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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:216
scalar Rey
const volVectorField & U() const
Access function to velocity field.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
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)
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
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
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:362
Abstract base class for turbulence models (RAS, LES and laminar).
Macros for easy insertion into run-time selection tables.
scalar magUp
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
scalar y
Type interpolateLog10(scalar x) const
Interpolate - takes log10 flag into account.
dynamicFvMesh & mesh
static const word propertiesName
Default name of the turbulence properties dictionary.
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.
const nearWallDist & y() const
Return the near wall distances.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
static const char nl
Definition: Ostream.H:262
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:224
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
nutUTabulatedWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual tmp< scalarField > calcUPlus(const scalarField &Rey) const
Calculate wall u+ from table.
label patchi
Constant dispersed-phase particle diameter model.
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
virtual tmp< scalarField > calcNut() const
Calculate the turbulence viscosity.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576