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-2020 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 (
105  const fvPatch& p,
107  const fvPatchFieldMapper& mapper
108 )
109 :
110  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
111  Prt_(ptf.Prt_)
112 {}
113 
114 
117 (
118  const fvPatch& p,
120  const dictionary& dict
121 )
122 :
123  fixedValueFvPatchScalarField(p, iF, dict),
124  Prt_(dict.lookupOrDefault<scalar>("Prt", 0.85))
125 {}
126 
127 
130 (
132 )
133 :
134  fixedValueFvPatchScalarField(awfpsf),
135  Prt_(awfpsf.Prt_)
136 {}
137 
138 
141 (
144 )
145 :
146  fixedValueFvPatchScalarField(awfpsf, iF),
147  Prt_(awfpsf.Prt_)
148 {}
149 
150 
151 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
152 
154 {
155  if (updated())
156  {
157  return;
158  }
159 
160  const label patchi = patch().index();
161 
162  const thermophysicalTransportModel& ttm =
163  db().lookupObject<thermophysicalTransportModel>
164  (
166  (
167  thermophysicalTransportModel::typeName,
168  internalField().group()
169  )
170  );
171 
172  const compressibleMomentumTransportModel& turbModel =
173  ttm.momentumTransport();
174 
176  nutWallFunctionFvPatchScalarField::nutw(turbModel, patchi);
177 
178  const scalar Cmu25 = pow025(nutw.Cmu());
179 
180  const scalarField& y = turbModel.y()[patchi];
181 
182  const tmp<scalarField> tmuw = turbModel.mu(patchi);
183  const scalarField& muw = tmuw();
184 
185  const tmp<scalarField> talphaw = ttm.thermo().alpha(patchi);
186  const scalarField& alphaw = talphaw();
187 
188  scalarField& alphatw = *this;
189 
190  const tmp<volScalarField> tk = turbModel.k();
191  const volScalarField& k = tk();
192 
193  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
194  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
195  const scalarField magGradUw(mag(Uw.snGrad()));
196 
197  const scalarField& rhow = turbModel.rho().boundaryField()[patchi];
198  const fvPatchScalarField& hew = ttm.thermo().he().boundaryField()[patchi];
199 
200  // Heat flux [W/m^2] - lagging alphatw
201  const scalarField qDot(ttm.thermo().alphaEff(alphatw, patchi)*hew.snGrad());
202 
203  // Populate boundary values
204  forAll(alphatw, facei)
205  {
206  label celli = patch().faceCells()[facei];
207 
208  scalar uTau = Cmu25*sqrt(k[celli]);
209 
210  scalar yPlus = uTau*y[facei]/(muw[facei]/rhow[facei]);
211 
212  // Molecular Prandtl number
213  scalar Pr = muw[facei]/alphaw[facei];
214 
215  // Molecular-to-turbulent Prandtl number ratio
216  scalar Prat = Pr/Prt_;
217 
218  // Thermal sublayer thickness
219  scalar P = Psmooth(Prat);
220  scalar yPlusTherm = this->yPlusTherm(nutw, P, Prat);
221 
222  // Evaluate new effective thermal diffusivity
223  scalar alphaEff = 0.0;
224  if (yPlus < yPlusTherm)
225  {
226  const scalar A = qDot[facei]*rhow[facei]*uTau*y[facei];
227 
228  const scalar B = qDot[facei]*Pr*yPlus;
229 
230  const scalar C = Pr*0.5*rhow[facei]*uTau*sqr(magUp[facei]);
231 
232  alphaEff = A/(B + C + vSmall);
233  }
234  else
235  {
236  const scalar A = qDot[facei]*rhow[facei]*uTau*y[facei];
237 
238  const scalar B =
239  qDot[facei]*Prt_*(1.0/nutw.kappa()*log(nutw.E()*yPlus) + P);
240 
241  const scalar magUc =
242  uTau/nutw.kappa()*log(nutw.E()*yPlusTherm) - mag(Uw[facei]);
243 
244  const scalar C =
245  0.5*rhow[facei]*uTau
246  *(Prt_*sqr(magUp[facei]) + (Pr - Prt_)*sqr(magUc));
247 
248  alphaEff = A/(B + C + vSmall);
249  }
250 
251  // Update turbulent thermal diffusivity
252  alphatw[facei] = max(0.0, alphaEff - alphaw[facei]);
253 
254  if (debug)
255  {
256  Info<< " uTau = " << uTau << nl
257  << " Pr = " << Pr << nl
258  << " Prt = " << Prt_ << nl
259  << " qDot = " << qDot[facei] << nl
260  << " yPlus = " << yPlus << nl
261  << " yPlusTherm = " << yPlusTherm << nl
262  << " alphaEff = " << alphaEff << nl
263  << " alphaw = " << alphaw[facei] << nl
264  << " alphatw = " << alphatw[facei] << nl
265  << endl;
266  }
267  }
268 
270 }
271 
272 
274 {
276  writeEntry(os, "Prt", Prt_);
277  writeEntry(os, "value", *this);
278 }
279 
280 
281 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
282 
284 (
287 );
288 
289 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
290 
291 } // End namespace compressible
292 } // End namespace Foam
293 
294 // ************************************************************************* //
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
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
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:186
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 > &)
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 thermal turbulent diffusivity of mixture [kg/m/s].
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:61
dimensionedScalar pow025(const dimensionedScalar &ds)
virtual const volScalarField & alpha() const
Thermal diffusivity for enthalpy of mixture [kg/m/s].
Definition: basicThermo.C:494
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
Macros for easy insertion into run-time selection tables.
scalar magUp
virtual volScalarField & he()=0
Enthalpy/Internal energy [J/kg].
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.
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:194
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:229
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.
Abstract base class for turbulence models (RAS, LES and laminar).
makePatchTypeField(fvPatchScalarField, thermalBaffleFvPatchScalarField)
Namespace for OpenFOAM.