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-2018 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 "fvPatchFieldMapper.H"
29 
30 #include "phaseSystem.H"
32 #include "ThermalDiffusivity.H"
34 #include "wallFvPatch.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 namespace compressible
41 {
42 
43 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44 
45 scalar alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::maxExp_
46  = 50.0;
47 scalar alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::tolerance_
48  = 0.01;
49 label alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::maxIters_
50  = 10;
51 
52 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
53 
55 {
56  if (!isA<wallFvPatch>(patch()))
57  {
59  << "Patch type for patch " << patch().name() << " must be wall\n"
60  << "Current patch type is " << patch().type() << nl
61  << exit(FatalError);
62  }
63 }
64 
65 
66 tmp<scalarField>
68 (
69  const scalarField& Prat
70 ) const
71 {
72  return 9.24*(pow(Prat, 0.75) - 1)*(1 + 0.28*exp(-0.007*Prat));
73 }
74 
75 
76 tmp<scalarField>
78 (
79  const scalarField& P,
80  const scalarField& Prat
81 ) const
82 {
83  tmp<scalarField> typsf(new scalarField(this->size()));
84  scalarField& ypsf = typsf.ref();
85 
86  forAll(ypsf, facei)
87  {
88  scalar ypt = 11.0;
89 
90  for (int i=0; i<maxIters_; i++)
91  {
92  scalar f = ypt - (log(E_*ypt)/kappa_ + P[facei])/Prat[facei];
93  scalar df = 1 - 1.0/(ypt*kappa_*Prat[facei]);
94  scalar yptNew = ypt - f/df;
95 
96  if (yptNew < vSmall)
97  {
98  ypsf[facei] = 0;
99  }
100  else if (mag(yptNew - ypt) < tolerance_)
101  {
102  ypsf[facei] = yptNew;
103  }
104  else
105  {
106  ypt = yptNew;
107  }
108  }
109 
110  ypsf[facei] = ypt;
111  }
112 
113  return typsf;
114 }
115 
116 tmp<scalarField>
118 (
119  const scalarField& prevAlphat
120 ) const
121 {
122 
123  // Lookup the fluid model
124  const phaseSystem& fluid =
125  db().lookupObject<phaseSystem>("phaseProperties");
126 
127  const phaseModel& phase
128  (
129  fluid.phases()[internalField().group()]
130  );
131 
132  const label patchi = patch().index();
133 
134  // Retrieve turbulence properties from model
135  const phaseCompressibleTurbulenceModel& turbModel =
136  db().lookupObject<phaseCompressibleTurbulenceModel>
137  (
139  );
140 
141  const scalar Cmu25 = pow025(Cmu_);
142 
143  const scalarField& y = turbModel.y()[patchi];
144 
145  const tmp<scalarField> tmuw = turbModel.mu(patchi);
146  const scalarField& muw = tmuw();
147 
148  const tmp<scalarField> talphaw = phase.thermo().alpha(patchi);
149  const scalarField& alphaw = talphaw();
150 
151  const tmp<volScalarField> tk = turbModel.k();
152  const volScalarField& k = tk();
153  const fvPatchScalarField& kw = k.boundaryField()[patchi];
154 
155  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
156  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
157  const scalarField magGradUw(mag(Uw.snGrad()));
158 
159  const fvPatchScalarField& rhow = turbModel.rho().boundaryField()[patchi];
160  const fvPatchScalarField& hew =
161  phase.thermo().he().boundaryField()[patchi];
162 
163  const fvPatchScalarField& Tw =
164  phase.thermo().T().boundaryField()[patchi];
165 
166  scalarField Tp(Tw.patchInternalField());
167 
168  // Heat flux [W/m2] - lagging alphatw
169  const scalarField qDot
170  (
171  (prevAlphat + alphaw)*hew.snGrad()
172  );
173 
174  scalarField uTau(Cmu25*sqrt(kw));
175 
176  scalarField yPlus(uTau*y/(muw/rhow));
177 
178  scalarField Pr(muw/alphaw);
179 
180  // Molecular-to-turbulent Prandtl number ratio
181  scalarField Prat(Pr/Prt_);
182 
183  // Thermal sublayer thickness
184  scalarField P(this->Psmooth(Prat));
185 
186  scalarField yPlusTherm(this->yPlusTherm(P, Prat));
187 
188  tmp<scalarField> talphatConv(new scalarField(this->size()));
189  scalarField& alphatConv = talphatConv.ref();
190 
191  // Populate boundary values
192  forAll(alphatConv, facei)
193  {
194  // Evaluate new effective thermal diffusivity
195  scalar alphaEff = 0.0;
196  if (yPlus[facei] < yPlusTherm[facei])
197  {
198  scalar A = qDot[facei]*rhow[facei]*uTau[facei]*y[facei];
199  scalar B = qDot[facei]*Pr[facei]*yPlus[facei];
200  scalar C = Pr[facei]*0.5*rhow[facei]*uTau[facei]*sqr(magUp[facei]);
201  alphaEff = A/(B + C + vSmall);
202  }
203  else
204  {
205  scalar A = qDot[facei]*rhow[facei]*uTau[facei]*y[facei];
206  scalar B =
207  qDot[facei]*Prt_*(1.0/kappa_*log(E_*yPlus[facei]) + P[facei]);
208  scalar magUc =
209  uTau[facei]/kappa_*log(E_*yPlusTherm[facei]) - mag(Uw[facei]);
210  scalar C =
211  0.5*rhow[facei]*uTau[facei]
212  *(Prt_*sqr(magUp[facei]) + (Pr[facei] - Prt_)*sqr(magUc));
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  Cmu_(0.09),
236  kappa_(0.41),
237  E_(9.8)
238 {
239  checkType();
240 }
241 
242 
245 (
246  const fvPatch& p,
247  const DimensionedField<scalar, volMesh>& iF,
248  const dictionary& dict
249 )
250 :
252  Prt_(dict.lookupOrDefault<scalar>("Prt", 0.85)),
253  Cmu_(dict.lookupOrDefault<scalar>("Cmu", 0.09)),
254  kappa_(dict.lookupOrDefault<scalar>("kappa", 0.41)),
255  E_(dict.lookupOrDefault<scalar>("E", 9.8))
256 {}
257 
258 
261 (
263  const fvPatch& p,
264  const DimensionedField<scalar, volMesh>& iF,
265  const fvPatchFieldMapper& mapper
266 )
267 :
269  Prt_(ptf.Prt_),
270  Cmu_(ptf.Cmu_),
271  kappa_(ptf.kappa_),
272  E_(ptf.E_)
273 {}
274 
275 
278 (
280 )
281 :
283  Prt_(awfpsf.Prt_),
284  Cmu_(awfpsf.Cmu_),
285  kappa_(awfpsf.kappa_),
286  E_(awfpsf.E_)
287 {}
288 
289 
292 (
294  const DimensionedField<scalar, volMesh>& iF
295 )
296 :
298  Prt_(awfpsf.Prt_),
299  Cmu_(awfpsf.Cmu_),
300  kappa_(awfpsf.kappa_),
301  E_(awfpsf.E_)
302 {}
303 
304 
305 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
306 
308 {
309  if (updated())
310  {
311  return;
312  }
313 
314  operator==(calcAlphat(*this));
315 
316  fixedValueFvPatchScalarField::updateCoeffs();
317 }
318 
319 
321 (
322  Ostream& os
323 ) const
324 {
326  os.writeKeyword("Prt") << Prt_ << token::END_STATEMENT << nl;
327  os.writeKeyword("Cmu") << Cmu_ << token::END_STATEMENT << nl;
328  os.writeKeyword("kappa") << kappa_ << token::END_STATEMENT << nl;
329  os.writeKeyword("E") << E_ << token::END_STATEMENT << nl;
330  dmdt_.writeEntry("dmdt", os);
331  writeEntry("value", os);
332 }
333 
334 
335 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
336 
338 (
341 );
342 
343 
344 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
345 
346 } // End namespace compressible
347 } // End namespace Foam
348 
349 // ************************************************************************* //
dimensionedScalar Pr("Pr", dimless, laminarTransport)
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:428
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)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
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)
void writeEntry(const word &keyword, Ostream &os) const
Write the field as a dictionary entry.
Definition: Field.C:726
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:362
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)
tmp< scalarField > yPlusTherm(const scalarField &P, const scalarField &Prat) const
Calculate y+ at the edge of the thermal laminar sublayer.
static const char nl
Definition: Ostream.H:265
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
alphatPhaseChangeWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
volScalarField & C
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
label patchi
scalar yPlus
dimensioned< scalar > mag(const dimensioned< Type > &)
volScalarField & p
makePatchTypeField(fvPatchScalarField, thermalBaffleFvPatchScalarField)
Namespace for OpenFOAM.