alphatJayatillekeWallFunctionFvPatchScalarField.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"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace incompressible
35 {
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 scalar alphatJayatillekeWallFunctionFvPatchScalarField::tolerance_ = 0.01;
40 label alphatJayatillekeWallFunctionFvPatchScalarField::maxIters_ = 10;
41 
42 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
43 
45 (
46  const scalar Prat
47 ) const
48 {
49  return 9.24*(pow(Prat, 0.75) - 1.0)*(1.0 + 0.28*exp(-0.007*Prat));
50 }
51 
52 
54 (
56  const scalar P,
57  const scalar Prat
58 ) const
59 {
60  scalar ypt = 11.0;
61 
62  for (int i=0; i<maxIters_; i++)
63  {
64  scalar f = ypt - (log(nutw.E()*ypt)/nutw.kappa() + P)/Prat;
65  scalar df = 1.0 - 1.0/(ypt*nutw.kappa()*Prat);
66  scalar yptNew = ypt - f/df;
67 
68  if (yptNew < vSmall)
69  {
70  return 0;
71  }
72  else if (mag(yptNew - ypt) < tolerance_)
73  {
74  return yptNew;
75  }
76  else
77  {
78  ypt = yptNew;
79  }
80  }
81 
82  return ypt;
83 }
84 
85 
86 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
87 
90 (
91  const fvPatch& p,
93 )
94 :
95  fixedValueFvPatchScalarField(p, iF),
96  Prt_(0.85)
97 {}
98 
99 
102 (
104  const fvPatch& p,
106  const fvPatchFieldMapper& mapper
107 )
108 :
109  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
110  Prt_(ptf.Prt_)
111 {}
112 
113 
116 (
117  const fvPatch& p,
119  const dictionary& dict
120 )
121 :
122  fixedValueFvPatchScalarField(p, iF, dict),
123  Prt_(readScalar(dict.lookup("Prt"))) // force read to avoid ambiguity
124 {}
125 
126 
129 (
131 )
132 :
133  fixedValueFvPatchScalarField(wfpsf),
134  Prt_(wfpsf.Prt_)
135 {}
136 
137 
140 (
143 )
144 :
145  fixedValueFvPatchScalarField(wfpsf, iF),
146  Prt_(wfpsf.Prt_)
147 {}
148 
149 
150 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
151 
153 {
154  if (updated())
155  {
156  return;
157  }
158 
159  const label patchi = patch().index();
160 
161  // Retrieve turbulence properties from model
162 
163  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
164  (
166  (
168  internalField().group()
169  )
170  );
171 
173  nutWallFunctionFvPatchScalarField::nutw(turbModel, patchi);
174 
175  const scalar Cmu25 = pow(nutw.Cmu(), 0.25);
176  const scalarField& y = turbModel.y()[patchi];
177  const tmp<volScalarField> tnu = turbModel.nu();
178  const volScalarField& nu = tnu();
179  const scalarField& nuw = nu.boundaryField()[patchi];
180  const tmp<volScalarField> tk = turbModel.k();
181  const volScalarField& k = tk();
182 
183  const IOdictionary& transportProperties =
184  db().lookupObject<IOdictionary>("transportProperties");
185 
186  // Molecular Prandtl number
187  const scalar Pr
188  (
190  (
191  "Pr",
192  dimless,
193  transportProperties.lookup("Pr")
194  ).value()
195  );
196 
197  // Populate boundary values
198  scalarField& alphatw = *this;
199  forAll(alphatw, facei)
200  {
201  const label celli = patch().faceCells()[facei];
202 
203  const scalar yPlus = Cmu25*sqrt(k[celli])*y[facei]/nuw[facei];
204 
205  // Molecular-to-turbulent Prandtl number ratio
206  const scalar Prat = Pr/Prt_;
207 
208  // Thermal sublayer thickness
209  const scalar P = Psmooth(Prat);
210  const scalar yPlusTherm = this->yPlusTherm(nutw, P, Prat);
211 
212  // Update turbulent thermal conductivity
213  if (yPlus > yPlusTherm)
214  {
215  const scalar nu = nuw[facei];
216  const scalar kt =
217  nu*(yPlus/(Prt_*(log(nutw.E()*yPlus)/nutw.kappa() + P)) - 1/Pr);
218 
219  alphatw[facei] = max(0.0, kt);
220  }
221  else
222  {
223  alphatw[facei] = 0.0;
224  }
225  }
226 
228 }
229 
230 
232 {
234  writeEntry(os, "Prt", Prt_);
235  writeEntry(os, "value", *this);
236 }
237 
238 
239 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
240 
242 (
245 );
246 
247 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
248 
249 } // End namespace incompressible
250 } // End namespace Foam
251 
252 // ************************************************************************* //
const char *const group
Group name for atomic constants.
#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)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
const Boundary & boundaryField() const
Return const-reference to the boundary field.
alphatJayatillekeWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
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
label k
Boltzmann constant.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:280
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Macros for easy insertion into run-time selection tables.
scalar yPlusTherm(const nutWallFunctionFvPatchScalarField &nutw, const scalar P, const scalar Prat) const
Calculate y+ at the edge of the thermal laminar sublayer.
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:52
Templated abstract base class for single-phase incompressible turbulence models.
scalar y
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions...
dimensionedScalar exp(const dimensionedScalar &ds)
static const word propertiesName
Default name of the turbulence properties dictionary.
static word groupName(Name name, const word &group)
This boundary condition provides a kinematic turbulent thermal conductivity for using wall functions...
Foam::fvPatchFieldMapper.
makePatchTypeField(fvPatchScalarField, alphatJayatillekeWallFunctionFvPatchScalarField)
const nearWallDist & y() const
Return the near wall distances.
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:68
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
labelList f(nPoints)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
label patchi
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
scalar yPlus
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:229
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
dimensioned< scalar > mag(const dimensioned< Type > &)
A class for managing temporary objects.
Definition: PtrList.H:53
static const nutWallFunctionFvPatchScalarField & nutw(const turbulenceModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
volScalarField & nu
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:583
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.