nutUWallFunctionFvPatchScalarField.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 fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
52  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
53  const tmp<scalarField> tnuw = turbModel.nu(patchi);
54  const scalarField& nuw = tnuw();
55 
56  tmp<scalarField> tyPlus = calcYPlus(magUp);
57  scalarField& yPlus = tyPlus.ref();
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.0);
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  scalar kappaRe = kappa_*magUp[facei]*y[facei]/nuw[facei];
100 
101  scalar yp = yPlusLam_;
102  scalar ryPlusLam = 1.0/yp;
103 
104  int iter = 0;
105  scalar yPlusLast = 0.0;
106 
107  do
108  {
109  yPlusLast = yp;
110  yp = (kappaRe + yp)/(1.0 + log(E_*yp));
111 
112  } while (mag(ryPlusLam*(yp - yPlusLast)) > 0.01 && ++iter < 10 );
113 
114  yPlus[facei] = max(0.0, yp);
115  }
116 
117  return tyPlus;
118 }
119 
120 
121 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
122 
124 (
125  const fvPatch& p,
127 )
128 :
130 {}
131 
132 
134 (
136  const fvPatch& p,
138  const fvPatchFieldMapper& mapper
139 )
140 :
141  nutWallFunctionFvPatchScalarField(ptf, p, iF, mapper)
142 {}
143 
144 
146 (
147  const fvPatch& p,
149  const dictionary& dict
150 )
151 :
153 {}
154 
155 
157 (
159 )
160 :
162 {}
163 
164 
166 (
167  const nutUWallFunctionFvPatchScalarField& sawfpsf,
169 )
170 :
172 {}
173 
174 
175 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
176 
178 {
179  const label patchi = patch().index();
180  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
181  (
183  (
185  internalField().group()
186  )
187  );
188  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
189  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
190 
191  return calcYPlus(magUp);
192 }
193 
194 
196 {
198  writeLocalEntries(os);
199  writeEntry("value", os);
200 }
201 
202 
203 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
204 
206 (
209 );
210 
211 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
212 
213 } // End namespace Foam
214 
215 // ************************************************************************* //
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: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
dimensionedScalar log(const dimensionedScalar &ds)
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 > &)
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions...
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
const volVectorField & U() const
Access function to velocity field.
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:65
Abstract base class for turbulence models (RAS, LES and laminar).
Macros for easy insertion into run-time selection tables.
scalar magUp
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
scalar y
const Boundary & boundaryField() const
Return const-reference to the boundary field.
static const word propertiesName
Default name of the turbulence properties dictionary.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual tmp< scalarField > calcYPlus(const scalarField &magUp) const
Calculate yPLus.
static word groupName(Name name, const word &group)
Foam::fvPatchFieldMapper.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
virtual tmp< scalarField > calcNut() const
Calculate the turbulence viscosity.
nutUWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
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:54
scalar yPlusLam_
Y+ at the edge of the laminar sublayer.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
const nearWallDist & y() const
Return the near wall distances.
Namespace for OpenFOAM.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:363