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 | Copyright (C) 2011-2015 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 
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
31 #include "wallFvPatch.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace compressible
38 {
39 
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 
42 scalar alphatJayatillekeWallFunctionFvPatchScalarField::maxExp_ = 50.0;
43 scalar alphatJayatillekeWallFunctionFvPatchScalarField::tolerance_ = 0.01;
44 label alphatJayatillekeWallFunctionFvPatchScalarField::maxIters_ = 10;
45 
46 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47 
48 void alphatJayatillekeWallFunctionFvPatchScalarField::checkType()
49 {
50  if (!isA<wallFvPatch>(patch()))
51  {
53  (
54  "alphatJayatillekeWallFunctionFvPatchScalarField::checkType()"
55  )
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 
63 scalar alphatJayatillekeWallFunctionFvPatchScalarField::Psmooth
64 (
65  const scalar Prat
66 ) const
67 {
68  return 9.24*(pow(Prat, 0.75) - 1.0)*(1.0 + 0.28*exp(-0.007*Prat));
69 }
70 
71 
72 scalar alphatJayatillekeWallFunctionFvPatchScalarField::yPlusTherm
73 (
74  const scalar P,
75  const scalar Prat
76 ) const
77 {
78  scalar ypt = 11.0;
79 
80  for (int i=0; i<maxIters_; i++)
81  {
82  scalar f = ypt - (log(E_*ypt)/kappa_ + P)/Prat;
83  scalar df = 1.0 - 1.0/(ypt*kappa_*Prat);
84  scalar yptNew = ypt - f/df;
85 
86  if (yptNew < VSMALL)
87  {
88  return 0;
89  }
90  else if (mag(yptNew - ypt) < tolerance_)
91  {
92  return yptNew;
93  }
94  else
95  {
96  ypt = yptNew;
97  }
98  }
99 
100  return ypt;
101 }
102 
103 
104 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
105 
108 (
109  const fvPatch& p,
111 )
112 :
113  fixedValueFvPatchScalarField(p, iF),
114  Prt_(0.85),
115  Cmu_(0.09),
116  kappa_(0.41),
117  E_(9.8)
118 {
119  checkType();
120 }
121 
122 
125 (
127  const fvPatch& p,
129  const fvPatchFieldMapper& mapper
130 )
131 :
132  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
133  Prt_(ptf.Prt_),
134  Cmu_(ptf.Cmu_),
135  kappa_(ptf.kappa_),
136  E_(ptf.E_)
137 {}
138 
139 
142 (
143  const fvPatch& p,
145  const dictionary& dict
146 )
147 :
148  fixedValueFvPatchScalarField(p, iF, dict),
149  Prt_(dict.lookupOrDefault<scalar>("Prt", 0.85)),
150  Cmu_(dict.lookupOrDefault<scalar>("Cmu", 0.09)),
151  kappa_(dict.lookupOrDefault<scalar>("kappa", 0.41)),
152  E_(dict.lookupOrDefault<scalar>("E", 9.8))
153 {
154  checkType();
155 }
156 
157 
160 (
162 )
163 :
164  fixedValueFvPatchScalarField(awfpsf),
165  Prt_(awfpsf.Prt_),
166  Cmu_(awfpsf.Cmu_),
167  kappa_(awfpsf.kappa_),
168  E_(awfpsf.E_)
169 {
170  checkType();
171 }
172 
173 
176 (
179 )
180 :
181  fixedValueFvPatchScalarField(awfpsf, iF),
182  Prt_(awfpsf.Prt_),
183  Cmu_(awfpsf.Cmu_),
184  kappa_(awfpsf.kappa_),
185  E_(awfpsf.E_)
186 {
187  checkType();
188 }
189 
190 
191 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
192 
194 {
195  if (updated())
196  {
197  return;
198  }
199 
200  const label patchi = patch().index();
201 
202  // Retrieve turbulence properties from model
203  const compressible::turbulenceModel& turbModel =
204  db().lookupObject<compressible::turbulenceModel>
205  (
207  (
208  compressible::turbulenceModel::propertiesName,
210  )
211  );
212 
213  const scalar Cmu25 = pow025(Cmu_);
214 
215  const scalarField& y = turbModel.y()[patchi];
216 
217  const tmp<scalarField> tmuw = turbModel.mu(patchi);
218  const scalarField& muw = tmuw();
219 
220  const tmp<scalarField> talphaw = turbModel.alpha(patchi);
221  const scalarField& alphaw = talphaw();
222 
223  scalarField& alphatw = *this;
224 
225  const tmp<volScalarField> tk = turbModel.k();
226  const volScalarField& k = tk();
227 
228  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
229  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
230  const scalarField magGradUw(mag(Uw.snGrad()));
231 
232  const scalarField& rhow = turbModel.rho().boundaryField()[patchi];
233  const fvPatchScalarField& hew =
234  turbModel.transport().he().boundaryField()[patchi];
235 
236  // Heat flux [W/m2] - lagging alphatw
237  const scalarField qDot
238  (
239  turbModel.transport().alphaEff(alphatw, patchi)*hew.snGrad()
240  );
241 
242  // Populate boundary values
243  forAll(alphatw, faceI)
244  {
245  label faceCellI = patch().faceCells()[faceI];
246 
247  scalar uTau = Cmu25*sqrt(k[faceCellI]);
248 
249  scalar yPlus = uTau*y[faceI]/(muw[faceI]/rhow[faceI]);
250 
251  // Molecular Prandtl number
252  scalar Pr = muw[faceI]/alphaw[faceI];
253 
254  // Molecular-to-turbulent Prandtl number ratio
255  scalar Prat = Pr/Prt_;
256 
257  // Thermal sublayer thickness
258  scalar P = Psmooth(Prat);
259  scalar yPlusTherm = this->yPlusTherm(P, Prat);
260 
261  // Evaluate new effective thermal diffusivity
262  scalar alphaEff = 0.0;
263  if (yPlus < yPlusTherm)
264  {
265  scalar A = qDot[faceI]*rhow[faceI]*uTau*y[faceI];
266  scalar B = qDot[faceI]*Pr*yPlus;
267  scalar C = Pr*0.5*rhow[faceI]*uTau*sqr(magUp[faceI]);
268  alphaEff = A/(B + C + VSMALL);
269  }
270  else
271  {
272  scalar A = qDot[faceI]*rhow[faceI]*uTau*y[faceI];
273  scalar B = qDot[faceI]*Prt_*(1.0/kappa_*log(E_*yPlus) + P);
274  scalar magUc = uTau/kappa_*log(E_*yPlusTherm) - mag(Uw[faceI]);
275  scalar C =
276  0.5*rhow[faceI]*uTau
277  *(Prt_*sqr(magUp[faceI]) + (Pr - Prt_)*sqr(magUc));
278  alphaEff = A/(B + C + VSMALL);
279  }
280 
281  // Update turbulent thermal diffusivity
282  alphatw[faceI] = max(0.0, alphaEff - alphaw[faceI]);
283 
284  if (debug)
285  {
286  Info<< " uTau = " << uTau << nl
287  << " Pr = " << Pr << nl
288  << " Prt = " << Prt_ << nl
289  << " qDot = " << qDot[faceI] << nl
290  << " yPlus = " << yPlus << nl
291  << " yPlusTherm = " << yPlusTherm << nl
292  << " alphaEff = " << alphaEff << nl
293  << " alphaw = " << alphaw[faceI] << nl
294  << " alphatw = " << alphatw[faceI] << nl
295  << endl;
296  }
297  }
298 
300 }
301 
302 
304 {
306  os.writeKeyword("Prt") << Prt_ << token::END_STATEMENT << nl;
307  os.writeKeyword("Cmu") << Cmu_ << token::END_STATEMENT << nl;
308  os.writeKeyword("kappa") << kappa_ << token::END_STATEMENT << nl;
309  os.writeKeyword("E") << E_ << token::END_STATEMENT << nl;
310  writeEntry("value", os);
311 }
312 
313 
314 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
315 
317 (
320 );
321 
322 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
323 
324 } // End namespace compressible
325 } // End namespace Foam
326 
327 // ************************************************************************* //
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensionedScalar Pr("Pr", dimless, laminarTransport)
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:207
scalar magUp
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
This boundary condition provides a thermal wall function for turbulent thermal diffusivity (usuallyal...
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:301
const char *const group
Group name for atomic constants.
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::fvPatchFieldMapper.
dimensionedScalar exp(const dimensionedScalar &ds)
Graphite solid properties.
Definition: C.H:57
messageStream Info
makePatchTypeField(fvPatchScalarField, thermalBaffleFvPatchScalarField)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:354
Namespace for OpenFOAM.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
static word groupName(Name name, const word &group)
dimensionedScalar log(const dimensionedScalar &ds)
static const char nl
Definition: Ostream.H:260
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:215
virtual tmp< volScalarField > alpha() const
Return the laminar thermal diffusivity for enthalpy [kg/m/s].
virtual tmp< volScalarField > alphaEff() const
Return the effective turbulent thermal diffusivity for enthalpy.
scalar uTau
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define forAll(list, i)
Definition: UList.H:421
label k
Boltzmann constant.
label patchi
dimensionedScalar pow025(const dimensionedScalar &ds)
alphatJayatillekeWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Macros for easy insertion into run-time selection tables.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:59
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
rDeltaT dimensionedInternalField()
scalar yPlus
bool compressible
Definition: pEqn.H:40
error FatalError
scalar y
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Evaluates and outputs turbulence y+ for models. Values written to time directories as field &#39;yPlus&#39;...
Definition: yPlus.H:63
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
A class for managing temporary objects.
Definition: PtrList.H:118
const volScalarField & alphaEff
Definition: setAlphaEff.H:49