alphatJayatillekeWallFunctionFvPatchScalarField.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2017-2022 OpenCFD Ltd
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
30#include "fvPatchFieldMapper.H"
31#include "volFields.H"
33#include "wallFvPatch.H"
35
36// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37
38namespace Foam
39{
40namespace compressible
41{
42
43// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44
45scalar alphatJayatillekeWallFunctionFvPatchScalarField::maxExp_ = 50.0;
46scalar alphatJayatillekeWallFunctionFvPatchScalarField::tolerance_ = 0.01;
47label alphatJayatillekeWallFunctionFvPatchScalarField::maxIters_ = 10;
48
49// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50
51void alphatJayatillekeWallFunctionFvPatchScalarField::checkType()
52{
53 if (!isA<wallFvPatch>(patch()))
54 {
56 << "Patch type for patch " << patch().name() << " must be wall\n"
57 << "Current patch type is " << patch().type() << nl
58 << exit(FatalError);
59 }
60}
61
62
63tmp<scalarField> alphatJayatillekeWallFunctionFvPatchScalarField::yPlus
64(
65 const compressible::turbulenceModel& turbModel
66) const
67{
68 const label patchi = patch().index();
69
70 const tmp<volScalarField> tnut = turbModel.nut();
71 const volScalarField& nut = tnut();
72
73 if (isA<nutWallFunctionFvPatchScalarField>(nut.boundaryField()[patchi]))
74 {
75 const auto& nutPf =
76 dynamic_cast<const nutWallFunctionFvPatchScalarField&>
77 (
78 nut.boundaryField()[patchi]
79 );
80
81 return nutPf.yPlus();
82 }
83 else
84 {
85 const scalarField& y = turbModel.y()[patchi];
86 const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
87
88 return
89 y*sqrt(turbModel.nuEff(patchi)*mag(Uw.snGrad()))
90 /turbModel.nu(patchi);
91 }
92}
93
94
95scalar alphatJayatillekeWallFunctionFvPatchScalarField::Psmooth
96(
97 const scalar Prat
98) const
99{
100 return 9.24*(pow(Prat, 0.75) - 1.0)*(1.0 + 0.28*exp(-0.007*Prat));
101}
102
103
104scalar alphatJayatillekeWallFunctionFvPatchScalarField::yPlusTherm
105(
106 const scalar P,
107 const scalar Prat
108) const
109{
110 scalar ypt = 11;
111
112 for (int iter = 0; iter < maxIters_; ++iter)
113 {
114 const scalar f = ypt - (log(E_*ypt)/kappa_ + P)/Prat;
115 const scalar df = 1.0 - 1.0/(ypt*kappa_*Prat);
116 const scalar yptNew = ypt - f/df;
117
118 if (yptNew < VSMALL)
119 {
120 return 0;
121 }
122 else if (mag(yptNew - ypt) < tolerance_)
123 {
124 return yptNew;
125 }
126 else
127 {
128 ypt = yptNew;
129 }
130 }
131
132 return ypt;
133}
134
135
136// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
137
138alphatJayatillekeWallFunctionFvPatchScalarField::
139alphatJayatillekeWallFunctionFvPatchScalarField
140(
141 const fvPatch& p,
143)
144:
145 fixedValueFvPatchScalarField(p, iF),
146 Prt_(0.85),
147 kappa_(0.41),
148 E_(9.8)
149{
150 checkType();
151}
152
153
154alphatJayatillekeWallFunctionFvPatchScalarField::
155alphatJayatillekeWallFunctionFvPatchScalarField
156(
158 const fvPatch& p,
160 const fvPatchFieldMapper& mapper
161)
162:
163 fixedValueFvPatchScalarField(ptf, p, iF, mapper),
164 Prt_(ptf.Prt_),
165 kappa_(ptf.kappa_),
166 E_(ptf.E_)
167{}
168
169
170alphatJayatillekeWallFunctionFvPatchScalarField::
171alphatJayatillekeWallFunctionFvPatchScalarField
172(
173 const fvPatch& p,
175 const dictionary& dict
176)
177:
178 fixedValueFvPatchScalarField(p, iF, dict),
179 Prt_(dict.getOrDefault<scalar>("Prt", 0.85)),
180 kappa_(dict.getOrDefault<scalar>("kappa", 0.41)),
181 E_(dict.getOrDefault<scalar>("E", 9.8))
182{
183 checkType();
184}
185
186
187alphatJayatillekeWallFunctionFvPatchScalarField::
188alphatJayatillekeWallFunctionFvPatchScalarField
189(
191)
192:
193 fixedValueFvPatchScalarField(awfpsf),
194 Prt_(awfpsf.Prt_),
195 kappa_(awfpsf.kappa_),
196 E_(awfpsf.E_)
197{
198 checkType();
199}
200
201
202alphatJayatillekeWallFunctionFvPatchScalarField::
203alphatJayatillekeWallFunctionFvPatchScalarField
204(
207)
208:
209 fixedValueFvPatchScalarField(awfpsf, iF),
210 Prt_(awfpsf.Prt_),
211 kappa_(awfpsf.kappa_),
212 E_(awfpsf.E_)
213{
214 checkType();
215}
216
217
218// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
219
221{
222 if (updated())
223 {
224 return;
225 }
226
227 const label patchi = patch().index();
228
229 // Retrieve turbulence properties from model
230 const auto& turbModel = db().lookupObject<compressible::turbulenceModel>
231 (
233 (
235 internalField().group()
236 )
237 );
238
239 const scalarField yPlusp(yPlus(turbModel));
240
241 const scalarField& y = turbModel.y()[patchi];
242
243 const tmp<scalarField> tmuw = turbModel.mu(patchi);
244 const scalarField& muw = tmuw();
245
246 const tmp<scalarField> talphaw = turbModel.alpha(patchi);
247 const scalarField& alphaw = talphaw();
248
249 const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
250 const scalarField magUp(mag(Uw.patchInternalField() - Uw));
251 const scalarField magGradUw(mag(Uw.snGrad()));
252
253 const scalarField& rhow = turbModel.rho().boundaryField()[patchi];
254 const fvPatchScalarField& hew =
255 turbModel.transport().he().boundaryField()[patchi];
256
257 scalarField& alphatw = *this;
258
259 // Heat flux [W/m2] - lagging alphatw
260 const scalarField qDot
261 (
262 turbModel.transport().alphaEff(alphatw, patchi)*hew.snGrad()
263 );
264
265 // Populate boundary values
266 forAll(alphatw, facei)
267 {
268 const scalar yPlus = yPlusp[facei];
269
270 const scalar uTau = yPlus/y[facei]*(muw[facei]/rhow[facei]);
271
272 // Molecular Prandtl number
273 const scalar Pr = muw[facei]/alphaw[facei];
274
275 // Molecular-to-turbulent Prandtl number ratio
276 const scalar Prat = Pr/Prt_;
277
278 // Thermal sublayer thickness
279 const scalar P = Psmooth(Prat);
280 const scalar yPlusTherm = this->yPlusTherm(P, Prat);
281
282 // Evaluate new effective thermal diffusivity
283 scalar alphaEff = 0;
284 if (yPlus < yPlusTherm)
285 {
286 const scalar A = qDot[facei]*rhow[facei]*uTau*y[facei];
287 const scalar B = qDot[facei]*Pr*yPlus;
288 const scalar C = Pr*0.5*rhow[facei]*uTau*sqr(magUp[facei]);
289
290 alphaEff = A/(B + C + VSMALL);
291 }
292 else
293 {
294 const scalar A = qDot[facei]*rhow[facei]*uTau*y[facei];
295 const scalar B = qDot[facei]*Prt_*(1.0/kappa_*log(E_*yPlus) + P);
296 const scalar magUc =
297 uTau/kappa_*log(E_*yPlusTherm) - mag(Uw[facei]);
298 const scalar C =
299 0.5*rhow[facei]*uTau
300 *(Prt_*sqr(magUp[facei]) + (Pr - Prt_)*sqr(magUc));
301
302 alphaEff = A/(B + C + VSMALL);
303 }
304
305 // Update turbulent thermal diffusivity
306 alphatw[facei] = max(scalar(0), alphaEff - alphaw[facei]);
307
308 if (debug)
309 {
310 Info<< " uTau = " << uTau << nl
311 << " Pr = " << Pr << nl
312 << " Prt = " << Prt_ << nl
313 << " qDot = " << qDot[facei] << nl
314 << " yPlus = " << yPlus << nl
315 << " yPlusTherm = " << yPlusTherm << nl
316 << " alphaEff = " << alphaEff << nl
317 << " alphaw = " << alphaw[facei] << nl
318 << " alphatw = " << alphatw[facei] << nl
319 << endl;
320 }
321 }
322
324}
325
326
328{
330 os.writeEntryIfDifferent<scalar>("Prt", 0.85, Prt_);
331 os.writeEntryIfDifferent<scalar>("kappa", 0.41, kappa_);
332 os.writeEntryIfDifferent<scalar>("E", 9.8, E_);
333 writeEntry("value", os);
334}
335
336
337// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
338
340(
343);
344
345// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
346
347} // End namespace compressible
348} // End namespace Foam
349
350// ************************************************************************* //
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
scalar y
Macros for easy insertion into run-time selection tables.
Graphite solid properties.
Definition: C.H:53
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:251
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
This boundary condition provides a thermal wall function for turbulent thermal diffusivity (usuallyal...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
virtual bool write()
Write the output fields.
A FieldMapper for finite-volume patch fields.
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:237
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:229
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
static const word propertiesName
Default name of the phase properties dictionary.
Definition: phaseSystem.H:290
A class for managing temporary objects.
Definition: tmp.H:65
volScalarField & p
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
scalar yPlus
scalar magUp
scalar uTau
scalar nut
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Definition: fvPatchField.H:676
volScalarField alphaEff("alphaEff", turbulence->nu()/Pr+alphat)
bool compressible
Definition: pEqn.H:2
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
messageStream Info
Information stream (stdout output on master, null elsewhere)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
fvPatchField< vector > fvPatchVectorField
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
labelList f(nPoints)
dictionary dict
dimensionedScalar Pr("Pr", dimless, laminarTransport)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333