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-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 
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace compressible
35 {
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 scalar alphatJayatillekeWallFunctionFvPatchScalarField::maxExp_ = 50.0;
40 scalar alphatJayatillekeWallFunctionFvPatchScalarField::tolerance_ = 0.01;
41 label alphatJayatillekeWallFunctionFvPatchScalarField::maxIters_ = 10;
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
48 (
49  const fvPatch& p,
51  const dictionary& dict
52 )
53 :
54  fixedValueFvPatchScalarField(p, iF, dict),
55  Prt_(dict.lookupOrDefault<scalar>("Prt", 0.85))
56 {}
57 
58 
61 (
63  const fvPatch& p,
65  const fvPatchFieldMapper& mapper
66 )
67 :
68  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
69  Prt_(ptf.Prt_)
70 {}
71 
72 
75 (
78 )
79 :
80  fixedValueFvPatchScalarField(awfpsf, iF),
81  Prt_(awfpsf.Prt_)
82 {}
83 
84 
85 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
86 
88 (
89  const scalarField& Prat
90 )
91 {
92  return 9.24*(pow(Prat, 0.75) - 1.0)*(1.0 + 0.28*exp(-0.007*Prat));
93 }
94 
95 
97 (
99  const scalarField& P,
100  const scalarField& Prat
101 )
102 {
103  tmp<scalarField> typt(new scalarField(nutw.size()));
104  scalarField& ypt = typt.ref();
105 
106  const scalar E = nutw.E();
107  const scalar kappa = nutw.kappa();
108 
109  forAll(ypt, facei)
110  {
111  ypt[facei] = 11.0;
112 
113  for (int i=0; i<maxIters_; i++)
114  {
115  const scalar f =
116  ypt[facei] - (log(E*ypt[facei])/kappa + P[facei])/Prat[facei];
117 
118  const scalar df = 1.0 - 1.0/(ypt[facei]*nutw.kappa()*Prat[facei]);
119 
120  const scalar dypt = - f/df;
121 
122  ypt[facei] += dypt;
123 
124  if (ypt[facei] < vSmall)
125  {
126  ypt[facei] = 0;
127  break;
128  }
129 
130  if (mag(dypt) < tolerance_)
131  {
132  break;
133  }
134  }
135  }
136 
137  return typt;
138 }
139 
140 
142 (
144  const scalar Prt,
145  const label patchi
146 )
147 {
148  const compressibleMomentumTransportModel& turbModel =
149  ttm.momentumTransport();
150 
153 
154  const scalar Cmu25 = pow025(nutw.Cmu());
155 
156  const scalarField& y = turbModel.y()[patchi];
157 
158  const tmp<scalarField> tnuw = turbModel.nu(patchi);
159  const scalarField& nuw = tnuw();
160 
161  const tmp<scalarField> talphaw
162  (
163  ttm.thermo().kappa().boundaryField()[patchi]
164  /ttm.thermo().Cp().boundaryField()[patchi]
165  );
166  const scalarField& alphaw = talphaw();
167 
168  const tmp<volScalarField> tk = turbModel.k();
169  const volScalarField& k = tk();
170 
171  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
172  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
173  const scalarField magGradUw(mag(Uw.snGrad()));
174 
175  const scalarField& rhow = turbModel.rho().boundaryField()[patchi];
176  const fvPatchScalarField& hew = ttm.thermo().he().boundaryField()[patchi];
177 
178  // Enthalpy gradient
179  const scalarField gradHew(hew.snGrad());
180 
181  // Molecular Prandtl number
182  const scalarField Pr(rhow*nuw/alphaw);
183 
184  // Molecular-to-turbulent Prandtl number ratio
185  const scalarField Prat(Pr/Prt);
186 
187  // Thermal sublayer thickness
188  const scalarField P
189  (
191  );
192  const scalarField yPlusTherm
193  (
195  (
196  nutw,
197  P,
198  Prat
199  )
200  );
201 
202  // Populate boundary values
203  tmp<scalarField> talphatw(new scalarField(nutw.size()));
204  scalarField& alphatw = talphatw.ref();
205  forAll(alphatw, facei)
206  {
207  const label celli = nutw.patch().faceCells()[facei];
208 
209  const scalar uTau = Cmu25*sqrt(k[celli]);
210 
211  const scalar yPlus = uTau*y[facei]/nuw[facei];
212 
213  // Evaluate new effective thermal diffusivity
214  scalar alphaEff = 0.0;
215  if (yPlus < yPlusTherm[facei])
216  {
217  const scalar A = gradHew[facei]*rhow[facei]*uTau*y[facei];
218 
219  const scalar B = gradHew[facei]*Pr[facei]*yPlus;
220 
221  const scalar C = Pr[facei]*0.5*rhow[facei]*uTau*sqr(magUp[facei]);
222 
223  alphaEff = (A - C)/(B + sign(B)*rootVSmall);
224  }
225  else
226  {
227  const scalar A = gradHew[facei]*rhow[facei]*uTau*y[facei];
228 
229  const scalar B =
230  gradHew[facei]*Prt
231  *(1.0/nutw.kappa()*log(nutw.E()*yPlus) + P[facei]);
232 
233  const scalar magUc =
234  uTau/nutw.kappa()
235  *log(nutw.E()*yPlusTherm[facei]) - mag(Uw[facei]);
236 
237  const scalar C =
238  0.5*rhow[facei]*uTau
239  *(Prt*sqr(magUp[facei]) + (Pr[facei] - Prt)*sqr(magUc));
240 
241  alphaEff = (A - C)/(B + sign(B)*rootVSmall);
242  }
243 
244  // Bounds on turbulent thermal diffusivity
245  static const scalar alphatwMin = 0;
246  const scalar alphatwMax = great*alphaw[facei]*nutw[facei]/nuw[facei];
247 
248  // Update turbulent thermal diffusivity
249  alphatw[facei] =
250  min(max(alphaEff - alphaw[facei], alphatwMin), alphatwMax);
251  }
252 
253  return talphatw;
254 }
255 
256 
258 {
259  if (updated())
260  {
261  return;
262  }
263 
265  db().lookupType<fluidThermophysicalTransportModel>
266  (
267  internalField().group()
268  );
269 
270  this->operator==(alphat(ttm, Prt_, patch().index()));
271 
273 }
274 
275 
277 {
279  writeEntry(os, "Prt", Prt_);
280  writeEntry(os, "value", *this);
281 }
282 
283 
284 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
285 
287 (
290 );
291 
292 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
293 
294 } // End namespace compressible
295 } // End namespace Foam
296 
297 // ************************************************************************* //
static const Foam::dimensionedScalar A("A", Foam::dimPressure, 611.21)
static const Foam::dimensionedScalar B("B", Foam::dimless, 18.678)
static const Foam::dimensionedScalar C("C", Foam::dimTemperature, 234.5)
scalar y
label k
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
Graphite solid properties.
Definition: C.H:51
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
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 volScalarField & he()=0
Enthalpy/Internal energy [J/kg].
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.
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:160
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.
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:231
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:172
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:164
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:204
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
const volVectorField & U() const
Access function to velocity field.
const volScalarField::Boundary & y() const
Return the near wall distance.
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.
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
const scalar uTau
const scalar magUp
label patchi
makePatchTypeField(fvPatchScalarField, alphatJayatillekeWallFunctionFvPatchScalarField)
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
Namespace for OpenFOAM.
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedScalar sign(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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar log(const dimensionedScalar &ds)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
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 > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar pow025(const dimensionedScalar &ds)
labelList f(nPoints)
dictionary dict
volScalarField & p