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-2022 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 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44 
45 scalar alphatJayatillekeWallFunctionFvPatchScalarField::Psmooth
46 (
47  const scalar Prat
48 ) const
49 {
50  return 9.24*(pow(Prat, 0.75) - 1.0)*(1.0 + 0.28*exp(-0.007*Prat));
51 }
52 
53 
54 scalar alphatJayatillekeWallFunctionFvPatchScalarField::yPlusTherm
55 (
56  const nutWallFunctionFvPatchScalarField& nutw,
57  const scalar P,
58  const scalar Prat
59 ) const
60 {
61  scalar ypt = 11.0;
62 
63  for (int i=0; i<maxIters_; i++)
64  {
65  scalar f = ypt - (log(nutw.E()*ypt)/nutw.kappa() + P)/Prat;
66  scalar df = 1.0 - 1.0/(ypt*nutw.kappa()*Prat);
67  scalar yptNew = ypt - f/df;
68 
69  if (yptNew < vSmall)
70  {
71  return 0;
72  }
73  else if (mag(yptNew - ypt) < tolerance_)
74  {
75  return yptNew;
76  }
77  else
78  {
79  ypt = yptNew;
80  }
81  }
82 
83  return ypt;
84 }
85 
86 
87 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
88 
91 (
92  const fvPatch& p,
94 )
95 :
96  fixedValueFvPatchScalarField(p, iF),
97  Prt_(0.85)
98 {}
99 
100 
103 (
104  const fvPatch& p,
106  const dictionary& dict
107 )
108 :
109  fixedValueFvPatchScalarField(p, iF, dict),
110  Prt_(dict.lookupOrDefault<scalar>("Prt", 0.85))
111 {}
112 
113 
116 (
118  const fvPatch& p,
120  const fvPatchFieldMapper& mapper
121 )
122 :
123  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
124  Prt_(ptf.Prt_)
125 {}
126 
127 
130 (
133 )
134 :
135  fixedValueFvPatchScalarField(awfpsf, iF),
136  Prt_(awfpsf.Prt_)
137 {}
138 
139 
140 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
141 
143 {
144  if (updated())
145  {
146  return;
147  }
148 
149  const label patchi = patch().index();
150 
151  const thermophysicalTransportModel& ttm =
152  db().lookupObject<thermophysicalTransportModel>
153  (
155  (
156  thermophysicalTransportModel::typeName,
157  internalField().group()
158  )
159  );
160 
161  const compressibleMomentumTransportModel& turbModel =
162  ttm.momentumTransport();
163 
165  nutWallFunctionFvPatchScalarField::nutw(turbModel, patchi);
166 
167  const scalar Cmu25 = pow025(nutw.Cmu());
168 
169  const scalarField& y = turbModel.y()[patchi];
170 
171  const tmp<scalarField> tnuw = turbModel.nu(patchi);
172  const scalarField& nuw = tnuw();
173 
174  const tmp<scalarField> talphaw
175  (
176  ttm.thermo().kappa().boundaryField()[patchi]
177  /ttm.thermo().Cp().boundaryField()[patchi]
178  );
179  const scalarField& alphaw = talphaw();
180 
181  scalarField& alphatw = *this;
182 
183  const tmp<volScalarField> tk = turbModel.k();
184  const volScalarField& k = tk();
185 
186  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
187  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
188  const scalarField magGradUw(mag(Uw.snGrad()));
189 
190  const scalarField& rhow = turbModel.rho().boundaryField()[patchi];
191  const fvPatchScalarField& hew = ttm.thermo().he().boundaryField()[patchi];
192 
193  // Heat flux [W/m^2] - lagging alphatw
194  const scalarField qDot(ttm.thermo().alphaEff(alphatw, patchi)*hew.snGrad());
195 
196  // Populate boundary values
197  forAll(alphatw, facei)
198  {
199  label celli = patch().faceCells()[facei];
200 
201  scalar uTau = Cmu25*sqrt(k[celli]);
202 
203  scalar yPlus = uTau*y[facei]/nuw[facei];
204 
205  // Molecular Prandtl number
206  scalar Pr = rhow[facei]*nuw[facei]/alphaw[facei];
207 
208  // Molecular-to-turbulent Prandtl number ratio
209  scalar Prat = Pr/Prt_;
210 
211  // Thermal sublayer thickness
212  scalar P = Psmooth(Prat);
213  scalar yPlusTherm = this->yPlusTherm(nutw, P, Prat);
214 
215  // Evaluate new effective thermal diffusivity
216  scalar alphaEff = 0.0;
217  if (yPlus < yPlusTherm)
218  {
219  const scalar A = qDot[facei]*rhow[facei]*uTau*y[facei];
220 
221  const scalar B = qDot[facei]*Pr*yPlus;
222 
223  const scalar C = Pr*0.5*rhow[facei]*uTau*sqr(magUp[facei]);
224 
225  alphaEff = A/(B + C + vSmall);
226  }
227  else
228  {
229  const scalar A = qDot[facei]*rhow[facei]*uTau*y[facei];
230 
231  const scalar B =
232  qDot[facei]*Prt_*(1.0/nutw.kappa()*log(nutw.E()*yPlus) + P);
233 
234  const scalar magUc =
235  uTau/nutw.kappa()*log(nutw.E()*yPlusTherm) - mag(Uw[facei]);
236 
237  const scalar C =
238  0.5*rhow[facei]*uTau
239  *(Prt_*sqr(magUp[facei]) + (Pr - Prt_)*sqr(magUc));
240 
241  alphaEff = A/(B + C + vSmall);
242  }
243 
244  // Update turbulent thermal diffusivity
245  alphatw[facei] = max(0.0, alphaEff - alphaw[facei]);
246 
247  if (debug)
248  {
249  Info<< " uTau = " << uTau << nl
250  << " Pr = " << Pr << nl
251  << " Prt = " << Prt_ << nl
252  << " qDot = " << qDot[facei] << nl
253  << " yPlus = " << yPlus << nl
254  << " yPlusTherm = " << yPlusTherm << nl
255  << " alphaEff = " << alphaEff << nl
256  << " alphaw = " << alphaw[facei] << nl
257  << " alphatw = " << alphatw[facei] << nl
258  << endl;
259  }
260  }
261 
263 }
264 
265 
267 {
269  writeEntry(os, "Prt", Prt_);
270  writeEntry(os, "value", *this);
271 }
272 
273 
274 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
275 
277 (
280 );
281 
282 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
283 
284 } // End namespace compressible
285 } // End namespace Foam
286 
287 // ************************************************************************* //
const char *const group
Group name for atomic constants.
Graphite solid properties.
Definition: C.H:48
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:166
dimensionedScalar log(const dimensionedScalar &ds)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
virtual const fluidThermo & thermo() const =0
Access function to incompressible transport model.
const compressibleMomentumTransportModel & momentumTransport() const
const Boundary & boundaryField() const
Return const-reference to the boundary field.
virtual tmp< volScalarField > alphaEff(const volScalarField &alphat) const =0
Effective turbulent thermal diffusivity of energy.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
This boundary condition provides a thermal wall function for turbulent thermal diffusivity (usuallyal...
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
dimensionedScalar pow025(const dimensionedScalar &ds)
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:243
Macros for easy insertion into run-time selection tables.
scalar magUp
virtual volScalarField & he()=0
Enthalpy/Internal energy [J/kg].
virtual const volScalarField & Cp() const =0
Heat capacity at constant pressure [J/kg/K].
scalar uTau
alphatJayatillekeWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
static const nutWallFunctionFvPatchScalarField & nutw(const momentumTransportModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
scalar y
This boundary condition provides a turbulent kinematic viscosity condition when using wall functions...
dimensionedScalar exp(const dimensionedScalar &ds)
static word groupName(Name name, const word &group)
Foam::fvPatchFieldMapper.
virtual const volScalarField & kappa() const =0
Thermal conductivity of mixture [W/m/K].
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
static const char nl
Definition: Ostream.H:260
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:174
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
Abstract base class for thermophysical transport models (RAS, LES and laminar).
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
label patchi
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
scalar yPlus
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:216
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
A class for managing temporary objects.
Definition: PtrList.H:53
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Base class for single-phase compressible turbulence models.
makePatchTypeField(fvPatchScalarField, thermalBaffleFvPatchScalarField)
Namespace for OpenFOAM.