nutUWallFunctionFvPatchScalarField.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-2019 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 tmp<scalarField> tnuw = turbModel.nu(patchi);
52  const scalarField& nuw = tnuw();
53 
54  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
55  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
56 
57  const scalarField yPlus(this->yPlus(magUp));
58 
59  tmp<scalarField> tnutw(new scalarField(patch().size(), 0.0));
60  scalarField& nutw = tnutw.ref();
61 
62  forAll(yPlus, facei)
63  {
64  if (yPlus[facei] > yPlusLam_)
65  {
66  nutw[facei] =
67  nuw[facei]*(yPlus[facei]*kappa_/log(E_*yPlus[facei]) - 1);
68  }
69  }
70 
71  return tnutw;
72 }
73 
74 
76 (
77  const scalarField& magUp
78 ) const
79 {
80  const label patchi = patch().index();
81 
82  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
83  (
85  (
87  internalField().group()
88  )
89  );
90  const scalarField& y = turbModel.y()[patchi];
91  const tmp<scalarField> tnuw = turbModel.nu(patchi);
92  const scalarField& nuw = tnuw();
93 
94  tmp<scalarField> tyPlus(new scalarField(patch().size(), 0.0));
95  scalarField& yPlus = tyPlus.ref();
96 
97  forAll(yPlus, facei)
98  {
99  const scalar Re = magUp[facei]*y[facei]/nuw[facei];
100  const scalar ryPlusLam = 1/yPlusLam_;
101 
102  int iter = 0;
103  scalar yp = yPlusLam_;
104  scalar yPlusLast = yp;
105 
106  do
107  {
108  yPlusLast = yp;
109  if (yp > yPlusLam_)
110  {
111  yp = (kappa_*Re + yp)/(1 + log(E_*yp));
112  }
113  else
114  {
115  yp = sqrt(Re);
116  }
117  } while(mag(ryPlusLam*(yp - yPlusLast)) > 0.0001 && ++iter < 20);
118 
119  yPlus[facei] = yp;
120  }
121 
122  return tyPlus;
123 }
124 
125 
126 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
127 
129 (
130  const fvPatch& p,
132 )
133 :
135 {}
136 
137 
139 (
141  const fvPatch& p,
143  const fvPatchFieldMapper& mapper
144 )
145 :
146  nutWallFunctionFvPatchScalarField(ptf, p, iF, mapper)
147 {}
148 
149 
151 (
152  const fvPatch& p,
154  const dictionary& dict
155 )
156 :
158 {}
159 
160 
162 (
164 )
165 :
167 {}
168 
169 
171 (
172  const nutUWallFunctionFvPatchScalarField& sawfpsf,
174 )
175 :
177 {}
178 
179 
180 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
181 
183 {
184  const label patchi = patch().index();
185  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
186  (
188  (
190  internalField().group()
191  )
192  );
193  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
194  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
195 
196  return yPlus(magUp);
197 }
198 
199 
201 {
203  writeLocalEntries(os);
204  writeEntry(os, "value", *this);
205 }
206 
207 
208 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
209 
211 (
214 );
215 
216 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
217 
218 } // End namespace Foam
219 
220 // ************************************************************************* //
const char *const group
Group name for atomic constants.
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
#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
dimensionedScalar log(const dimensionedScalar &ds)
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
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:158
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions...
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
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
dimensionedScalar sqrt(const dimensionedScalar &ds)
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:280
Abstract base class for turbulence models (RAS, LES and laminar).
Macros for easy insertion into run-time selection tables.
scalar magUp
scalar y
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.
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
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:194
nutUWallFunctionFvPatchScalarField(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 > nut() const
Calculate the turbulence viscosity.
label patchi
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
scalar yPlusLam_
Y+ at the edge of the laminar sublayer.
static const nutWallFunctionFvPatchScalarField & nutw(const turbulenceModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
Namespace for OpenFOAM.