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-2021 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 tmp<scalarField> tnuw = turbModel.nu(patchi);
53  const scalarField& nuw = tnuw();
54 
55  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
56  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
57 
58  const scalarField yPlus(this->yPlus(magUp));
59 
60  tmp<scalarField> tnutw(new scalarField(patch().size(), 0.0));
61  scalarField& nutw = tnutw.ref();
62 
63  forAll(yPlus, facei)
64  {
65  if (yPlus[facei] > yPlusLam_)
66  {
67  nutw[facei] =
68  nuw[facei]*(yPlus[facei]*kappa_/log(E_*yPlus[facei]) - 1);
69  }
70  }
71 
72  return tnutw;
73 }
74 
75 
77 (
78  const scalarField& magUp
79 ) const
80 {
81  const label patchi = patch().index();
82 
83  const momentumTransportModel& turbModel =
84  db().lookupObject<momentumTransportModel>
85  (
87  (
88  momentumTransportModel::typeName,
89  internalField().group()
90  )
91  );
92  const scalarField& y = turbModel.y()[patchi];
93  const tmp<scalarField> tnuw = turbModel.nu(patchi);
94  const scalarField& nuw = tnuw();
95 
96  tmp<scalarField> tyPlus(new scalarField(patch().size(), 0.0));
97  scalarField& yPlus = tyPlus.ref();
98 
99  forAll(yPlus, facei)
100  {
101  const scalar Re = magUp[facei]*y[facei]/nuw[facei];
102  const scalar ryPlusLam = 1/yPlusLam_;
103 
104  int iter = 0;
105  scalar yp = yPlusLam_;
106  scalar yPlusLast = yp;
107 
108  do
109  {
110  yPlusLast = yp;
111  if (yp > yPlusLam_)
112  {
113  yp = (kappa_*Re + yp)/(1 + log(E_*yp));
114  }
115  else
116  {
117  yp = sqrt(Re);
118  }
119  } while(mag(ryPlusLam*(yp - yPlusLast)) > 0.0001 && ++iter < 20);
120 
121  yPlus[facei] = yp;
122  }
123 
124  return tyPlus;
125 }
126 
127 
128 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
129 
131 (
132  const fvPatch& p,
134 )
135 :
137 {}
138 
139 
141 (
142  const fvPatch& p,
144  const dictionary& dict
145 )
146 :
148 {}
149 
150 
152 (
154  const fvPatch& p,
156  const fvPatchFieldMapper& mapper
157 )
158 :
159  nutWallFunctionFvPatchScalarField(ptf, p, iF, mapper)
160 {}
161 
162 
164 (
165  const nutUWallFunctionFvPatchScalarField& sawfpsf,
167 )
168 :
170 {}
171 
172 
173 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
174 
176 {
177  const label patchi = patch().index();
178  const momentumTransportModel& turbModel =
179  db().lookupObject<momentumTransportModel>
180  (
182  (
183  momentumTransportModel::typeName,
184  internalField().group()
185  )
186  );
187  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
188  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
189 
190  return yPlus(magUp);
191 }
192 
193 
195 {
197  writeLocalEntries(os);
198  writeEntry(os, "value", *this);
199 }
200 
201 
202 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
203 
205 (
208 );
209 
210 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
211 
212 } // End namespace Foam
213 
214 // ************************************************************************* //
const char *const group
Group name for atomic constants.
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
dimensionedScalar log(const dimensionedScalar &ds)
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
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:181
dimensionedScalar sqrt(const dimensionedScalar &ds)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
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:243
const volScalarField::Boundary & y() const
Return the near wall distance.
Macros for easy insertion into run-time selection tables.
scalar magUp
static const nutWallFunctionFvPatchScalarField & nutw(const momentumTransportModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
scalar y
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
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:54
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:174
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.
Abstract base class for turbulence models (RAS, LES and laminar).
const volVectorField & U() const
Access function to velocity 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:53
scalar yPlusLam_
Y+ at the edge of the laminar sublayer.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
Namespace for OpenFOAM.