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 
27 #include "turbulenceModel.H"
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
30 #include "wallFvPatch.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace incompressible
38 {
39 
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 
44 
45 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
46 
48 {
49  if (!isA<wallFvPatch>(patch()))
50  {
52  << "Invalid wall function specification" << nl
53  << " Patch type for patch " << patch().name()
54  << " must be wall" << nl
55  << " Current patch type is " << patch().type() << nl << endl
56  << abort(FatalError);
57  }
58 }
59 
60 
62 (
63  const scalar Prat
64 ) const
65 {
66  return 9.24*(pow(Prat, 0.75) - 1.0)*(1.0 + 0.28*exp(-0.007*Prat));
67 }
68 
69 
71 (
72  const scalar P,
73  const scalar Prat
74 ) const
75 {
76  scalar ypt = 11.0;
77 
78  for (int i=0; i<maxIters_; i++)
79  {
80  scalar f = ypt - (log(E_*ypt)/kappa_ + P)/Prat;
81  scalar df = 1.0 - 1.0/(ypt*kappa_*Prat);
82  scalar yptNew = ypt - f/df;
83 
84  if (yptNew < VSMALL)
85  {
86  return 0;
87  }
88  else if (mag(yptNew - ypt) < tolerance_)
89  {
90  return yptNew;
91  }
92  else
93  {
94  ypt = yptNew;
95  }
96  }
97 
98  return ypt;
99 }
100 
101 
102 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
103 
106 (
107  const fvPatch& p,
109 )
110 :
111  fixedValueFvPatchScalarField(p, iF),
112  Prt_(0.85),
113  Cmu_(0.09),
114  kappa_(0.41),
115  E_(9.8)
116 {
117  checkType();
118 }
119 
120 
123 (
125  const fvPatch& p,
127  const fvPatchFieldMapper& mapper
128 )
129 :
130  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
131  Prt_(ptf.Prt_),
132  Cmu_(ptf.Cmu_),
133  kappa_(ptf.kappa_),
134  E_(ptf.E_)
135 {
136  checkType();
137 }
138 
139 
142 (
143  const fvPatch& p,
145  const dictionary& dict
146 )
147 :
148  fixedValueFvPatchScalarField(p, iF, dict),
149  Prt_(readScalar(dict.lookup("Prt"))), // force read to avoid ambiguity
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(wfpsf),
165  Prt_(wfpsf.Prt_),
166  Cmu_(wfpsf.Cmu_),
167  kappa_(wfpsf.kappa_),
168  E_(wfpsf.E_)
169 {
170  checkType();
171 }
172 
173 
176 (
179 )
180 :
181  fixedValueFvPatchScalarField(wfpsf, iF),
182  Prt_(wfpsf.Prt_),
183  Cmu_(wfpsf.Cmu_),
184  kappa_(wfpsf.kappa_),
185  E_(wfpsf.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 
204  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
205  (
207  (
210  )
211  );
212 
213  const scalar Cmu25 = pow(Cmu_, 0.25);
214  const scalarField& y = turbModel.y()[patchi];
215  const tmp<volScalarField> tnu = turbModel.nu();
216  const volScalarField& nu = tnu();
217  const scalarField& nuw = nu.boundaryField()[patchi];
218  const tmp<volScalarField> tk = turbModel.k();
219  const volScalarField& k = tk();
220 
221  const IOdictionary& transportProperties =
222  db().lookupObject<IOdictionary>("transportProperties");
223 
224  // Molecular Prandtl number
225  const scalar Pr
226  (
228  (
229  "Pr",
230  dimless,
231  transportProperties.lookup("Pr")
232  ).value()
233  );
234 
235  // Populate boundary values
236  scalarField& alphatw = *this;
237  forAll(alphatw, faceI)
238  {
239  label faceCellI = patch().faceCells()[faceI];
240 
241  // y+
242  scalar yPlus = Cmu25*sqrt(k[faceCellI])*y[faceI]/nuw[faceI];
243 
244  // Molecular-to-turbulent Prandtl number ratio
245  scalar Prat = Pr/Prt_;
246 
247  // Thermal sublayer thickness
248  scalar P = Psmooth(Prat);
249  scalar yPlusTherm = this->yPlusTherm(P, Prat);
250 
251  // Update turbulent thermal conductivity
252  if (yPlus > yPlusTherm)
253  {
254  scalar nu = nuw[faceI];
255  scalar kt = nu*(yPlus/(Prt_*(log(E_*yPlus)/kappa_ + P)) - 1/Pr);
256  alphatw[faceI] = max(0.0, kt);
257  }
258  else
259  {
260  alphatw[faceI] = 0.0;
261  }
262  }
263 
265 }
266 
267 
269 {
271  os.writeKeyword("Prt") << Prt_ << token::END_STATEMENT << nl;
272  os.writeKeyword("Cmu") << Cmu_ << token::END_STATEMENT << nl;
273  os.writeKeyword("kappa") << kappa_ << token::END_STATEMENT << nl;
274  os.writeKeyword("E") << E_ << token::END_STATEMENT << nl;
275  writeEntry("value", os);
276 }
277 
278 
279 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
280 
282 (
285 );
286 
287 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
288 
289 } // End namespace incompressible
290 } // End namespace Foam
291 
292 // ************************************************************************* //
Foam::incompressible::alphatJayatillekeWallFunctionFvPatchScalarField::checkType
virtual void checkType()
Check the type of the patch.
Definition: alphatJayatillekeWallFunctionFvPatchScalarField.C:47
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
alphatJayatillekeWallFunctionFvPatchScalarField.H
volFields.H
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
Foam::fvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:349
Foam::token::END_STATEMENT
@ END_STATEMENT
Definition: token.H:99
Foam::IOobject::groupName
static word groupName(Name name, const word &group)
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
p
p
Definition: pEqn.H:62
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
dimensionedInternalField
rDeltaT dimensionedInternalField()
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:97
wallFvPatch.H
Foam::incompressible::alphatJayatillekeWallFunctionFvPatchScalarField::Cmu_
scalar Cmu_
Cmu coefficient.
Definition: alphatJayatillekeWallFunctionFvPatchScalarField.H:115
Foam::constant::atomic::group
const char *const group
Group name for atomic constants.
Definition: atomicConstants.C:40
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
fvPatchFieldMapper.H
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:261
Foam::incompressible::alphatJayatillekeWallFunctionFvPatchScalarField::tolerance_
static scalar tolerance_
Definition: alphatJayatillekeWallFunctionFvPatchScalarField.H:126
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:252
Foam::turbulenceModel::y
const nearWallDist & y() const
Return the near wall distances.
Definition: turbulenceModel.H:157
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::incompressible::alphatJayatillekeWallFunctionFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: alphatJayatillekeWallFunctionFvPatchScalarField.C:193
Foam::incompressible::alphatJayatillekeWallFunctionFvPatchScalarField
This boundary condition provides a kinematic turbulent thermal conductivity for using wall functions,...
Definition: alphatJayatillekeWallFunctionFvPatchScalarField.H:103
nu
volScalarField & nu
Definition: readMechanicalProperties.H:179
Foam::incompressible::alphatJayatillekeWallFunctionFvPatchScalarField::Psmooth
scalar Psmooth(const scalar Prat) const
`P' function
Definition: alphatJayatillekeWallFunctionFvPatchScalarField.C:62
Foam::incompressible::alphatJayatillekeWallFunctionFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: alphatJayatillekeWallFunctionFvPatchScalarField.C:268
Foam::label
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
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::incompressible::alphatJayatillekeWallFunctionFvPatchScalarField::alphatJayatillekeWallFunctionFvPatchScalarField
alphatJayatillekeWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: alphatJayatillekeWallFunctionFvPatchScalarField.C:106
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::incompressible::alphatJayatillekeWallFunctionFvPatchScalarField::Prt_
scalar Prt_
Turbulent Prandtl number.
Definition: alphatJayatillekeWallFunctionFvPatchScalarField.H:112
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Pr
dimensionedScalar Pr("Pr", dimless, laminarTransport)
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:73
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::incompressible::alphatJayatillekeWallFunctionFvPatchScalarField::yPlusTherm
scalar yPlusTherm(const scalar P, const scalar Prat) const
Calculate y+ at the edge of the thermal laminar sublayer.
Definition: alphatJayatillekeWallFunctionFvPatchScalarField.C:71
Foam::log
dimensionedScalar log(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:253
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::turbulenceModel::k
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::incompressible::alphatJayatillekeWallFunctionFvPatchScalarField::E_
scalar E_
E coefficient.
Definition: alphatJayatillekeWallFunctionFvPatchScalarField.H:121
Foam::fvPatchField< Type >::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:296
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::incompressible::alphatJayatillekeWallFunctionFvPatchScalarField::maxIters_
static label maxIters_
Definition: alphatJayatillekeWallFunctionFvPatchScalarField.H:127
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::incompressible::alphatJayatillekeWallFunctionFvPatchScalarField::kappa_
scalar kappa_
Von Karman constant.
Definition: alphatJayatillekeWallFunctionFvPatchScalarField.H:118
Foam::readScalar
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
f
labelList f(nPoints)
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
Foam::TurbulenceModel::nu
tmp< volScalarField > nu() const
Return the laminar viscosity.
Definition: TurbulenceModel.H:160
Foam::Ostream::writeKeyword
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:59
Foam::yPlus
This function object evaluates and outputs turbulence y+ for turbulence models. The field is stored o...
Definition: yPlus.H:110
patchi
label patchi
Definition: getPatchFieldScalar.H:1
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::IncompressibleTurbulenceModel
Templated abstract base class for single-phase incompressible turbulence models.
Definition: IncompressibleTurbulenceModel.H:52
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:45
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
turbulenceModel.H
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51
y
scalar y
Definition: LISASMDCalcMethod1.H:14
Foam::incompressible::makePatchTypeField
makePatchTypeField(fvPatchScalarField, turbulentBoundaryCoupledFvPatchScalarField)