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