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-2021 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 (
116  const fvPatch& p,
117  const DimensionedField<scalar, volMesh>& iF,
118  const dictionary& dict
119 )
120 :
121  fixedValueFvPatchScalarField(p, iF, dict),
122  Prt_(dict.lookupOrDefault<scalar>("Prt", 0.85))
123 {}
124 
125 
128 (
130  const fvPatch& p,
131  const DimensionedField<scalar, volMesh>& iF,
132  const fvPatchFieldMapper& mapper
133 )
134 :
135  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
136  Prt_(ptf.Prt_)
137 {}
138 
139 
142 (
144  const DimensionedField<scalar, volMesh>& iF
145 )
146 :
147  fixedValueFvPatchScalarField(awfpsf, iF),
148  Prt_(awfpsf.Prt_)
149 {}
150 
151 
152 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
153 
154 tmp<scalarField>
156 (
157  const scalarField& prevAlphat
158 ) const
159 {
160  // Lookup the fluid model
161  const phaseSystem& fluid =
162  db().lookupObject<phaseSystem>("phaseProperties");
163 
164  const phaseModel& phase
165  (
166  fluid.phases()[internalField().group()]
167  );
168 
169  const label patchi = patch().index();
170 
171  // Retrieve turbulence properties from model
174  (
175  IOobject::groupName(momentumTransportModel::typeName, phase.name())
176  );
177 
178  const nutWallFunctionFvPatchScalarField& nutw =
179  nutWallFunctionFvPatchScalarField::nutw(turbModel, patchi);
180 
181  const scalar Cmu25 = pow025(nutw.Cmu());
182 
183  const scalarField& y = turbModel.y()[patchi];
184 
185  const tmp<scalarField> tmuw = turbModel.mu(patchi);
186  const scalarField& muw = tmuw();
187 
188  const tmp<scalarField> talphaw = phase.thermo().alpha(patchi);
189  const scalarField& alphaw = talphaw();
190 
191  const tmp<volScalarField> tk = turbModel.k();
192  const volScalarField& k = tk();
193  const fvPatchScalarField& kw = k.boundaryField()[patchi];
194 
195  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
196  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
197  const scalarField magGradUw(mag(Uw.snGrad()));
198 
199  const fvPatchScalarField& rhow = turbModel.rho().boundaryField()[patchi];
200  const fvPatchScalarField& hew =
201  phase.thermo().he().boundaryField()[patchi];
202 
203  const fvPatchScalarField& Tw =
204  phase.thermo().T().boundaryField()[patchi];
205 
206  const scalarField Tp(Tw.patchInternalField());
207 
208  // Heat flux [W/m^2] - lagging alphatw
209  const scalarField qDot
210  (
211  (prevAlphat + alphaw)*hew.snGrad()
212  );
213 
214  const scalarField uTau(Cmu25*sqrt(kw));
215 
216  const scalarField yPlus(uTau*y/(muw/rhow));
217 
218  const scalarField Pr(muw/alphaw);
219 
220  // Molecular-to-turbulent Prandtl number ratio
221  const scalarField Prat(Pr/Prt_);
222 
223  // Thermal sublayer thickness
224  const scalarField P(this->Psmooth(Prat));
225 
226  const scalarField yPlusTherm(this->yPlusTherm(nutw, P, Prat));
227 
228  tmp<scalarField> talphatConv(new scalarField(this->size()));
229  scalarField& alphatConv = talphatConv.ref();
230 
231  // Populate boundary values
232  forAll(alphatConv, facei)
233  {
234  // Evaluate new effective thermal diffusivity
235  scalar alphaEff = 0.0;
236  if (yPlus[facei] < yPlusTherm[facei])
237  {
238  const scalar A = qDot[facei]*rhow[facei]*uTau[facei]*y[facei];
239 
240  const scalar B = qDot[facei]*Pr[facei]*yPlus[facei];
241 
242  const scalar C =
243  Pr[facei]*0.5*rhow[facei]*uTau[facei]*sqr(magUp[facei]);
244 
245  alphaEff = A/(B + C + vSmall);
246  }
247  else
248  {
249  const scalar A = qDot[facei]*rhow[facei]*uTau[facei]*y[facei];
250 
251  const scalar B =
252  qDot[facei]*Prt_
253  *(1.0/nutw.kappa()*log(nutw.E()*yPlus[facei]) + P[facei]);
254 
255  const scalar magUc =
256  uTau[facei]/nutw.kappa()
257  *log(nutw.E()*yPlusTherm[facei]) - mag(Uw[facei]);
258 
259  const scalar C =
260  0.5*rhow[facei]*uTau[facei]
261  *(Prt_*sqr(magUp[facei]) + (Pr[facei] - Prt_)*sqr(magUc));
262 
263  alphaEff = A/(B + C + vSmall);
264  }
265 
266  // Update convective heat transfer turbulent thermal diffusivity
267  alphatConv[facei] = max(0.0, alphaEff - alphaw[facei]);
268  }
269 
270  return talphatConv;
271 }
272 
273 
275 {
276  if (updated())
277  {
278  return;
279  }
280 
281  operator==(calcAlphat(*this));
282 
283  fixedValueFvPatchScalarField::updateCoeffs();
284 }
285 
286 
288 (
289  Ostream& os
290 ) const
291 {
293  writeEntry(os, "Prt", Prt_);
294  writeEntry(os, "value", *this);
295 }
296 
297 
298 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
299 
301 (
304 );
305 
306 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
307 
308 } // End namespace compressible
309 } // End namespace Foam
310 
311 // ************************************************************************* //
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
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:260
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
PhaseCompressibleMomentumTransportModel< dynamicTransportModel > momentumTransportModel
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.