alphatPhaseJayatillekeWallFunctionFvPatchScalarField.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) 2015-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 
27 #include "phaseSystem.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace compressible
36 {
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 scalar alphatPhaseJayatillekeWallFunctionFvPatchScalarField::maxExp_ = 50.0;
41 scalar alphatPhaseJayatillekeWallFunctionFvPatchScalarField::tolerance_ = 0.01;
42 label alphatPhaseJayatillekeWallFunctionFvPatchScalarField::maxIters_ = 10;
43 
44 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 
47 (
48  const scalarField& Prat
49 ) const
50 {
51  return 9.24*(pow(Prat, 0.75) - 1.0)*(1.0 + 0.28*exp(-0.007*Prat));
52 }
53 
54 
55 tmp<scalarField>
57 (
58  const nutWallFunctionFvPatchScalarField& nutw,
59  const scalarField& P,
60  const scalarField& Prat
61 ) const
62 {
63  tmp<scalarField> typsf(new scalarField(this->size()));
64  scalarField& ypsf = typsf.ref();
65 
66  forAll(ypsf, facei)
67  {
68  scalar ypt = 11.0;
69 
70  for (int i=0; i<maxIters_; i++)
71  {
72  const scalar f =
73  ypt - (log(nutw.E()*ypt)/nutw.kappa() + P[facei])/Prat[facei];
74  const scalar df = 1 - 1.0/(ypt*nutw.kappa()*Prat[facei]);
75  const scalar yptNew = ypt - f/df;
76 
77  if (yptNew < vSmall)
78  {
79  ypsf[facei] = 0;
80  break;
81  }
82  else if (mag(yptNew - ypt) < tolerance_)
83  {
84  ypsf[facei] = yptNew;
85  }
86  else
87  {
88  ypt = yptNew;
89  }
90  }
91 
92  ypsf[facei] = ypt;
93  }
94 
95  return typsf;
96 }
97 
98 
99 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
100 
103 (
104  const fvPatch& p,
105  const DimensionedField<scalar, volMesh>& iF
106 )
107 :
108  fixedValueFvPatchScalarField(p, iF),
109  Prt_(0.85)
110 {}
111 
112 
115 (
117  const fvPatch& p,
118  const DimensionedField<scalar, volMesh>& iF,
119  const fvPatchFieldMapper& mapper
120 )
121 :
122  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
123  Prt_(ptf.Prt_)
124 {}
125 
126 
129 (
130  const fvPatch& p,
131  const DimensionedField<scalar, volMesh>& iF,
132  const dictionary& dict
133 )
134 :
135  fixedValueFvPatchScalarField(p, iF, dict),
136  Prt_(dict.lookupOrDefault<scalar>("Prt", 0.85))
137 {}
138 
139 
140 
143 (
145 )
146 :
147  fixedValueFvPatchScalarField(awfpsf),
148  Prt_(awfpsf.Prt_)
149 {}
150 
151 
154 (
156  const DimensionedField<scalar, volMesh>& iF
157 )
158 :
159  fixedValueFvPatchScalarField(awfpsf, iF),
160  Prt_(awfpsf.Prt_)
161 {}
162 
163 
164 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
165 
166 tmp<scalarField>
168 (
169  const scalarField& prevAlphat
170 ) const
171 {
172  // Lookup the fluid model
173  const phaseSystem& fluid =
174  db().lookupObject<phaseSystem>("phaseProperties");
175 
176  const phaseModel& phase
177  (
178  fluid.phases()[internalField().group()]
179  );
180 
181  const label patchi = patch().index();
182 
183  // Retrieve turbulence properties from model
184  const phaseCompressibleMomentumTransportModel& turbModel =
185  db().lookupObject<phaseCompressibleMomentumTransportModel>
186  (
187  IOobject::groupName(momentumTransportModel::typeName, phase.name())
188  );
189 
190  const nutWallFunctionFvPatchScalarField& nutw =
191  nutWallFunctionFvPatchScalarField::nutw(turbModel, patchi);
192 
193  const scalar Cmu25 = pow025(nutw.Cmu());
194 
195  const scalarField& y = turbModel.y()[patchi];
196 
197  const tmp<scalarField> tmuw = turbModel.mu(patchi);
198  const scalarField& muw = tmuw();
199 
200  const tmp<scalarField> talphaw = phase.thermo().alpha(patchi);
201  const scalarField& alphaw = talphaw();
202 
203  const tmp<volScalarField> tk = turbModel.k();
204  const volScalarField& k = tk();
205  const fvPatchScalarField& kw = k.boundaryField()[patchi];
206 
207  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
208  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
209  const scalarField magGradUw(mag(Uw.snGrad()));
210 
211  const fvPatchScalarField& rhow = turbModel.rho().boundaryField()[patchi];
212  const fvPatchScalarField& hew =
213  phase.thermo().he().boundaryField()[patchi];
214 
215  const fvPatchScalarField& Tw =
216  phase.thermo().T().boundaryField()[patchi];
217 
218  const scalarField Tp(Tw.patchInternalField());
219 
220  // Heat flux [W/m^2] - lagging alphatw
221  const scalarField qDot
222  (
223  (prevAlphat + alphaw)*hew.snGrad()
224  );
225 
226  const scalarField uTau(Cmu25*sqrt(kw));
227 
228  const scalarField yPlus(uTau*y/(muw/rhow));
229 
230  const scalarField Pr(muw/alphaw);
231 
232  // Molecular-to-turbulent Prandtl number ratio
233  const scalarField Prat(Pr/Prt_);
234 
235  // Thermal sublayer thickness
236  const scalarField P(this->Psmooth(Prat));
237 
238  const scalarField yPlusTherm(this->yPlusTherm(nutw, P, Prat));
239 
240  tmp<scalarField> talphatConv(new scalarField(this->size()));
241  scalarField& alphatConv = talphatConv.ref();
242 
243  // Populate boundary values
244  forAll(alphatConv, facei)
245  {
246  // Evaluate new effective thermal diffusivity
247  scalar alphaEff = 0.0;
248  if (yPlus[facei] < yPlusTherm[facei])
249  {
250  const scalar A = qDot[facei]*rhow[facei]*uTau[facei]*y[facei];
251 
252  const scalar B = qDot[facei]*Pr[facei]*yPlus[facei];
253 
254  const scalar C =
255  Pr[facei]*0.5*rhow[facei]*uTau[facei]*sqr(magUp[facei]);
256 
257  alphaEff = A/(B + C + vSmall);
258  }
259  else
260  {
261  const scalar A = qDot[facei]*rhow[facei]*uTau[facei]*y[facei];
262 
263  const scalar B =
264  qDot[facei]*Prt_
265  *(1.0/nutw.kappa()*log(nutw.E()*yPlus[facei]) + P[facei]);
266 
267  const scalar magUc =
268  uTau[facei]/nutw.kappa()
269  *log(nutw.E()*yPlusTherm[facei]) - mag(Uw[facei]);
270 
271  const scalar C =
272  0.5*rhow[facei]*uTau[facei]
273  *(Prt_*sqr(magUp[facei]) + (Pr[facei] - Prt_)*sqr(magUc));
274 
275  alphaEff = A/(B + C + vSmall);
276  }
277 
278  // Update convective heat transfer turbulent thermal diffusivity
279  alphatConv[facei] = max(0.0, alphaEff - alphaw[facei]);
280  }
281 
282  return talphatConv;
283 }
284 
285 
287 {
288  if (updated())
289  {
290  return;
291  }
292 
293  operator==(calcAlphat(*this));
294 
295  fixedValueFvPatchScalarField::updateCoeffs();
296 }
297 
298 
300 (
301  Ostream& os
302 ) const
303 {
305  writeEntry(os, "Prt", Prt_);
306  writeEntry(os, "value", *this);
307 }
308 
309 
310 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
311 
313 (
316 );
317 
318 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
319 
320 } // End namespace compressible
321 } // End namespace Foam
322 
323 // ************************************************************************* //
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
dictionary dict
#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
PhaseCompressibleMomentumTransportModel< phaseModel > phaseCompressibleMomentumTransportModel
Typedef for phaseCompressibleMomentumTransportModel.
fvPatchField< vector > fvPatchVectorField
dimensionedScalar log(const dimensionedScalar &ds)
tmp< scalarField > calcAlphat(const scalarField &prevAlphat) const
Evaluate the turbulent thermal diffusivity.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensionedScalar pow025(const dimensionedScalar &ds)
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:280
Macros for easy insertion into run-time selection tables.
scalar magUp
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
scalar uTau
static const nutWallFunctionFvPatchScalarField & nutw(const momentumTransportModel &turbModel, const label patchi)
Return the nut patchField for the given wall patch.
dimensionedScalar exp(const dimensionedScalar &ds)
fvPatchField< scalar > fvPatchScalarField
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static word groupName(Name name, const word &group)
tmp< scalarField > yPlusTherm(const nutWallFunctionFvPatchScalarField &nutw, const scalarField &P, const scalarField &Prat) const
Calculate y+ at the edge of the thermal laminar sublayer.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
tmp< scalarField > Psmooth(const scalarField &Prat) const
&#39;P&#39; function
label patchi
alphatPhaseJayatillekeWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
scalar yPlus
dimensioned< scalar > mag(const dimensioned< Type > &)
volScalarField & p
makePatchTypeField(fvPatchScalarField, thermalBaffleFvPatchScalarField)
Namespace for OpenFOAM.