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-2024 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 
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace compressible
35 {
36 
37 // * * * * * * * * * * * * Private Static Data Members * * * * * * * * * * * //
38 
39 const scalar alphatJayatillekeWallFunctionFvPatchScalarField::tolerance_ = 0.01;
40 const label alphatJayatillekeWallFunctionFvPatchScalarField::maxIters_ = 10;
41 
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
47 (
48  const fvPatch& p,
50  const dictionary& dict
51 )
52 :
53  fixedValueFvPatchScalarField(p, iF, dict),
54  Prt_(dict.lookupOrDefault<scalar>("Prt", dimless, 0.85))
55 {}
56 
57 
60 (
62  const fvPatch& p,
64  const fieldMapper& mapper
65 )
66 :
67  fixedValueFvPatchScalarField(ptf, p, iF, mapper, false),
68  Prt_(ptf.Prt_)
69 {
70  mapper(*this, ptf, [&](){ return this->patchInternalField(); });
71 }
72 
73 
76 (
79 )
80 :
81  fixedValueFvPatchScalarField(awfpsf, iF),
82  Prt_(awfpsf.Prt_)
83 {}
84 
85 
86 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
87 
89 (
90  const fvPatchScalarField& ptf,
91  const fieldMapper& mapper
92 )
93 {
94  mapper(*this, ptf, [&](){ return this->patchInternalField(); });
95 }
96 
97 
99 (
100  const scalarField& Prat
101 )
102 {
103  return 9.24*(pow(Prat, 0.75) - 1)*(1 + 0.28*exp(-0.007*Prat));
104 }
105 
106 
108 (
110  const scalarField& P,
111  const scalarField& Prat
112 )
113 {
114  tmp<scalarField> typt(new scalarField(nutw.size()));
115  scalarField& ypt = typt.ref();
116 
117  const scalar E = nutw.E();
118  const scalar kappa = nutw.kappa();
119 
120  forAll(ypt, facei)
121  {
122  ypt[facei] = 11;
123 
124  for (int i=0; i<maxIters_; i++)
125  {
126  const scalar f =
127  ypt[facei] - (log(E*ypt[facei])/kappa + P[facei])/Prat[facei];
128 
129  const scalar df = 1 - 1/(ypt[facei]*nutw.kappa()*Prat[facei]);
130 
131  const scalar dypt = - f/df;
132 
133  ypt[facei] += dypt;
134 
135  if (ypt[facei] < vSmall)
136  {
137  ypt[facei] = 0;
138  break;
139  }
140 
141  if (mag(dypt) < tolerance_)
142  {
143  break;
144  }
145  }
146  }
147 
148  return typt;
149 }
150 
151 
153 (
155  const scalar Prt,
156  const label patchi
157 )
158 {
159  const compressibleMomentumTransportModel& turbModel =
160  ttm.momentumTransport();
161 
164 
165  const scalar E = nutw.E();
166  const scalar kappa = nutw.kappa();
167 
168  const tmp<scalarField> tnuw = turbModel.nu(patchi);
169  const scalarField& nuw = tnuw();
170 
171  const scalarField alphaw
172  (
173  ttm.thermo().kappa().boundaryField()[patchi]
174  /ttm.thermo().Cp().boundaryField()[patchi]
175  );
176 
177  const scalarField& rhow = turbModel.rho().boundaryField()[patchi];
178 
179  // Molecular Prandtl number
180  const scalarField Pr(rhow*nuw/alphaw);
181 
182  // Molecular-to-turbulent Prandtl number ratio
183  const scalarField Prat(Pr/Prt);
184 
185  // Momentum sublayer thickness
186  const scalarField yPlus(nutw.yPlus());
187 
188  // Thermal sublayer thickness
189  const scalarField P
190  (
192  );
193  const scalarField yPlusTherm
194  (
196  (
197  nutw,
198  P,
199  Prat
200  )
201  );
202 
203  // Populate boundary values
204  tmp<scalarField> talphatw(new scalarField(nutw.size(), scalar(0)));
205  scalarField& alphatw = talphatw.ref();
206  forAll(alphatw, facei)
207  {
208  if (yPlus[facei] > yPlusTherm[facei])
209  {
210  const scalar Tplus = Prt*(log(E*yPlus[facei])/kappa + P[facei]);
211 
212  const scalar alphaByAlphaEff = Tplus/Pr[facei]/yPlus[facei];
213 
214  alphatw[facei] =
215  alphaw[facei]*max(1/alphaByAlphaEff - 1, scalar(0));
216  }
217  }
218 
219  return talphatw;
220 }
221 
222 
224 {
225  if (updated())
226  {
227  return;
228  }
229 
231  db().lookupType<fluidThermophysicalTransportModel>
232  (
233  internalField().group()
234  );
235 
236  this->operator==(alphat(ttm, Prt_, patch().index()));
237 
239 }
240 
241 
243 {
245  writeEntry(os, "Prt", Prt_);
246  writeEntry(os, "value", *this);
247 }
248 
249 
250 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
251 
253 (
256 );
257 
258 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
259 
260 } // End namespace compressible
261 } // End namespace Foam
262 
263 // ************************************************************************* //
#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
virtual const volScalarField & Cp() const =0
Heat capacity at constant pressure [J/kg/K].
virtual const volScalarField & kappa() const =0
Thermal conductivity of mixture [W/m/K].
Base class for single-phase compressible turbulence models.
const rhoField & rho() const
Return the density field.
This boundary condition provides a thermal wall function for turbulent thermal diffusivity (usuallyal...
static tmp< scalarField > alphat(const fluidThermophysicalTransportModel &ttm, const scalar Prt, const label patchi)
Calculate the turbulent thermal diffusivity.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void map(const fvPatchScalarField &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
static tmp< scalarField > yPlusTherm(const nutWallFunctionFvPatchScalarField &nutw, const scalarField &P, const scalarField &Prat)
Calculate y+ at the edge of the thermal laminar sublayer.
alphatJayatillekeWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
static tmp< scalarField > P(const scalarField &Prat)
Calculate the smoothing function.
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 for fluid thermophysical transport models RAS, LES and laminar.
const compressibleMomentumTransportModel & momentumTransport() const
Access function to momentum transport model.
virtual const fluidThermo & thermo() const =0
Access function to fluid thermophysical properties.
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 void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:202
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions,...
static const nutWallFunctionFvPatchScalarField & nutw(const momentumTransportModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
virtual tmp< scalarField > yPlus() const =0
Calculate and return the yPlus at the boundary.
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 yPlus
label patchi
makePatchTypeField(fvPatchScalarField, alphatJayatillekeWallFunctionFvPatchScalarField)
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
Namespace for OpenFOAM.
dimensionedScalar exp(const dimensionedScalar &ds)
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
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
const dimensionSet dimless
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
labelList f(nPoints)
dictionary dict
volScalarField & p