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-2023 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 "fieldMapper.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().lookupType<momentumTransportModel>(internalField().group());
45 
46  const tmp<scalarField> tnuw = turbModel.nu(patchi);
47  const scalarField& nuw = tnuw();
48 
49  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
50  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
51 
52  const scalarField yPlus(this->yPlus(magUp));
53 
54  tmp<scalarField> tnutw(new scalarField(patch().size(), 0.0));
55  scalarField& nutw = tnutw.ref();
56 
57  forAll(yPlus, facei)
58  {
59  if (yPlus[facei] > yPlusLam_)
60  {
61  nutw[facei] =
62  nuw[facei]*(yPlus[facei]*kappa_/log(E_*yPlus[facei]) - 1);
63  }
64  }
65 
66  return tnutw;
67 }
68 
69 
71 (
72  const scalarField& magUp
73 ) const
74 {
75  const label patchi = patch().index();
76 
77  const momentumTransportModel& turbModel =
78  db().lookupType<momentumTransportModel>(internalField().group());
79 
80  const scalarField& y = turbModel.yb()[patchi];
81  const tmp<scalarField> tnuw = turbModel.nu(patchi);
82  const scalarField& nuw = tnuw();
83 
84  tmp<scalarField> tyPlus(new scalarField(patch().size(), 0.0));
85  scalarField& yPlus = tyPlus.ref();
86 
87  forAll(yPlus, facei)
88  {
89  const scalar Re = magUp[facei]*y[facei]/nuw[facei];
90  const scalar ryPlusLam = 1/yPlusLam_;
91 
92  int iter = 0;
93  scalar yp = yPlusLam_;
94  scalar yPlusLast = yp;
95 
96  do
97  {
98  yPlusLast = yp;
99  if (yp > yPlusLam_)
100  {
101  yp = (kappa_*Re + yp)/(1 + log(E_*yp));
102  }
103  else
104  {
105  yp = sqrt(Re);
106  }
107  } while(mag(ryPlusLam*(yp - yPlusLast)) > 0.0001 && ++iter < 20);
108 
109  yPlus[facei] = yp;
110  }
111 
112  return tyPlus;
113 }
114 
115 
116 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
117 
119 (
120  const fvPatch& p,
122  const dictionary& dict
123 )
124 :
126 {}
127 
128 
130 (
132  const fvPatch& p,
134  const fieldMapper& mapper
135 )
136 :
137  nutWallFunctionFvPatchScalarField(ptf, p, iF, mapper)
138 {}
139 
140 
142 (
143  const nutUWallFunctionFvPatchScalarField& sawfpsf,
145 )
146 :
148 {}
149 
150 
151 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
152 
154 {
155  const label patchi = patch().index();
156 
157 
158  const momentumTransportModel& turbModel =
159  db().lookupType<momentumTransportModel>(internalField().group());
160 
161  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
162  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
163 
164  return yPlus(magUp);
165 }
166 
167 
169 {
171  writeLocalEntries(os);
172  writeEntry(os, "value", *this);
173 }
174 
175 
176 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
177 
179 (
182 );
183 
184 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
185 
186 } // End namespace Foam
187 
188 // ************************************************************************* //
scalar y
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const Boundary & boundaryField() const
Return const-reference to the boundary field.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:229
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:170
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
Abstract base class for turbulence models (RAS, LES and laminar).
const volScalarField::Boundary & yb() const
Return the near wall distance.
const volVectorField & U() const
Access function to velocity field.
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions,...
virtual tmp< scalarField > yPlus() const
Calculate and return the yPlus at the boundary.
virtual tmp< scalarField > nut() const
Calculate the turbulence viscosity.
nutUWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions,...
scalar yPlusLam_
Y+ at the edge of the laminar sublayer.
virtual void writeLocalEntries(Ostream &) const
Write local wall function variables.
static const nutWallFunctionFvPatchScalarField & nutw(const momentumTransportModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
const scalar magUp
label patchi
Namespace for OpenFOAM.
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
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar log(const dimensionedScalar &ds)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
dictionary dict
volScalarField & p