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-2020 OpenCFD Ltd
10 -------------------------------------------------------------------------------
11 License
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"
32 #include "wallFvPatch.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 namespace incompressible
41 {
42 
43 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44 
47 
48 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
49 
51 {
52  if (!isA<wallFvPatch>(patch()))
53  {
55  << "Invalid wall function specification" << nl
56  << " Patch type for patch " << patch().name()
57  << " must be wall" << nl
58  << " Current patch type is " << patch().type() << nl << endl
59  << abort(FatalError);
60  }
61 }
62 
63 
65 (
66  const turbulenceModel& turbModel
67 ) const
68 {
69  const label patchi = patch().index();
70  const tmp<volScalarField> tnut = turbModel.nut();
71  const volScalarField& nut = tnut();
72 
73  if (isA<nutWallFunctionFvPatchScalarField>(nut.boundaryField()[patchi]))
74  {
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 
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 
105 (
106  const scalar P,
107  const scalar Prat
108 ) const
109 {
110  scalar ypt = 11.0;
111 
112  for (int i=0; i<maxIters_; i++)
113  {
114  scalar f = ypt - (log(E_*ypt)/kappa_ + P)/Prat;
115  scalar df = 1.0 - 1.0/(ypt*kappa_*Prat);
116  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 
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 
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  checkType();
169 }
170 
171 
174 (
175  const fvPatch& p,
177  const dictionary& dict
178 )
179 :
180  fixedValueFvPatchScalarField(p, iF, dict),
181  Prt_(dict.get<scalar>("Prt")), // force read to avoid ambiguity
182  kappa_(dict.getOrDefault<scalar>("kappa", 0.41)),
183  E_(dict.getOrDefault<scalar>("E", 9.8))
184 {
185  checkType();
186 }
187 
188 
191 (
193 )
194 :
195  fixedValueFvPatchScalarField(wfpsf),
196  Prt_(wfpsf.Prt_),
197  kappa_(wfpsf.kappa_),
198  E_(wfpsf.E_)
199 {
200  checkType();
201 }
202 
203 
206 (
209 )
210 :
211  fixedValueFvPatchScalarField(wfpsf, iF),
212  Prt_(wfpsf.Prt_),
213  kappa_(wfpsf.kappa_),
214  E_(wfpsf.E_)
215 {
216  checkType();
217 }
218 
219 
220 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
221 
223 {
224  if (updated())
225  {
226  return;
227  }
228 
229  const label patchi = patch().index();
230 
231  // Retrieve turbulence properties from model
232 
233  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
234  (
236  (
238  internalField().group()
239  )
240  );
241 
242  const scalarField yPlusp(yPlus(turbModel));
243 
244  const tmp<volScalarField> tnu = turbModel.nu();
245  const volScalarField& nu = tnu();
246  const scalarField& nuw = nu.boundaryField()[patchi];
247 
249  db().lookupObject<IOdictionary>("transportProperties");
250 
251  // Molecular Prandtl number
252  const scalar Pr
253  (
255  );
256 
257  // Populate boundary values
258  scalarField& alphatw = *this;
259  forAll(alphatw, facei)
260  {
261  scalar yPlus = yPlusp[facei];
262 
263  // Molecular-to-turbulent Prandtl number ratio
264  scalar Prat = Pr/Prt_;
265 
266  // Thermal sublayer thickness
267  scalar P = Psmooth(Prat);
268  scalar yPlusTherm = this->yPlusTherm(P, Prat);
269 
270  // Update turbulent thermal conductivity
271  if (yPlus > yPlusTherm)
272  {
273  scalar nu = nuw[facei];
274  scalar kt = nu*(yPlus/(Prt_*(log(E_*yPlus)/kappa_ + P)) - 1/Pr);
275  alphatw[facei] = max(0.0, kt);
276  }
277  else
278  {
279  alphatw[facei] = 0.0;
280  }
281  }
282 
284 }
285 
286 
288 {
290  os.writeEntry("Prt", Prt_);
291  os.writeEntry("kappa", kappa_);
292  os.writeEntry("E", E_);
293  writeEntry("value", os);
294 }
295 
296 
297 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
298 
300 (
303 );
304 
305 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
306 
307 } // End namespace incompressible
308 } // End namespace Foam
309 
310 // ************************************************************************* //
Foam::incompressible::alphatJayatillekeWallFunctionFvPatchScalarField::checkType
virtual void checkType()
Check the type of the patch.
Definition: alphatJayatillekeWallFunctionFvPatchScalarField.C:50
Foam::fvPatchField< vector >
alphatJayatillekeWallFunctionFvPatchScalarField.H
volFields.H
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
Foam::fvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:384
Foam::fvPatchField::snGrad
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:225
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::constant::atomic::group
constexpr const char *const group
Group name for atomic constants.
Definition: atomicConstants.H:52
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::turbulenceModel::nut
virtual tmp< volScalarField > nut() const =0
Return the turbulence viscosity.
transportProperties
IOdictionary transportProperties(IOobject("transportProperties", runTime.constant(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
nutWallFunctionFvPatchScalarField.H
Foam::nutWallFunctionFvPatchScalarField::yPlus
virtual tmp< scalarField > yPlus() const =0
Calculate and return the yPlus at the boundary.
wallFvPatch.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
fvPatchFieldMapper.H
Foam::incompressible::alphatJayatillekeWallFunctionFvPatchScalarField::tolerance_
static scalar tolerance_
Definition: alphatJayatillekeWallFunctionFvPatchScalarField.H:120
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261
Foam::turbulenceModel::y
const nearWallDist & y() const
Return the near wall distances.
Definition: turbulenceModel.H:159
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::incompressible::alphatJayatillekeWallFunctionFvPatchScalarField::yPlus
tmp< scalarField > yPlus(const turbulenceModel &turbModel) const
Return the patch y+.
Definition: alphatJayatillekeWallFunctionFvPatchScalarField.C:65
Foam::incompressible::alphatJayatillekeWallFunctionFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: alphatJayatillekeWallFunctionFvPatchScalarField.C:222
Foam::incompressible::alphatJayatillekeWallFunctionFvPatchScalarField
This boundary condition provides a kinematic turbulent thermal conductivity for using wall functions,...
Definition: alphatJayatillekeWallFunctionFvPatchScalarField.H:100
nu
volScalarField & nu
Definition: readMechanicalProperties.H:176
Foam::incompressible::alphatJayatillekeWallFunctionFvPatchScalarField::Psmooth
scalar Psmooth(const scalar Prat) const
`P' function
Definition: alphatJayatillekeWallFunctionFvPatchScalarField.C:96
Foam::incompressible::alphatJayatillekeWallFunctionFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: alphatJayatillekeWallFunctionFvPatchScalarField.C:287
Foam::Field< scalar >
Foam::incompressible::alphatJayatillekeWallFunctionFvPatchScalarField::alphatJayatillekeWallFunctionFvPatchScalarField
alphatJayatillekeWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: alphatJayatillekeWallFunctionFvPatchScalarField.C:140
Foam::incompressible::alphatJayatillekeWallFunctionFvPatchScalarField::Prt_
scalar Prt_
Turbulent Prandtl number.
Definition: alphatJayatillekeWallFunctionFvPatchScalarField.H:109
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Pr
dimensionedScalar Pr("Pr", dimless, laminarTransport)
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
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:105
Foam::log
dimensionedScalar log(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:262
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::incompressible::alphatJayatillekeWallFunctionFvPatchScalarField::E_
scalar E_
E coefficient.
Definition: alphatJayatillekeWallFunctionFvPatchScalarField.H:115
Foam::fvPatchField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:321
Foam::incompressible::alphatJayatillekeWallFunctionFvPatchScalarField::maxIters_
static label maxIters_
Definition: alphatJayatillekeWallFunctionFvPatchScalarField.H:121
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::incompressible::alphatJayatillekeWallFunctionFvPatchScalarField::kappa_
scalar kappa_
Von Karman constant.
Definition: alphatJayatillekeWallFunctionFvPatchScalarField.H:112
Foam::TurbulenceModel::nu
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
Definition: TurbulenceModel.H:159
Foam::nl
constexpr char nl
Definition: Ostream.H:404
f
labelList f(nPoints)
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::IncompressibleTurbulenceModel
Templated abstract base class for single-phase incompressible turbulence models.
Definition: IncompressibleTurbulenceModel.H:55
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::IOobject::groupName
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
nut
scalar nut
Definition: evaluateNearWall.H:5
Foam::incompressible::makePatchTypeField
makePatchTypeField(fvPatchScalarField, alphatJayatillekeWallFunctionFvPatchScalarField)
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::turbulenceModel::nuEff
virtual tmp< volScalarField > nuEff() const =0
Return the effective viscosity.
Foam::nutWallFunctionFvPatchScalarField
The class nutWallFunction is a base class that parents the derived boundary conditions which provide ...
Definition: nutWallFunctionFvPatchScalarField.H:251
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189
Foam::turbulenceModel::U
const volVectorField & U() const
Access function to velocity field.
Definition: turbulenceModel.H:144
y
scalar y
Definition: LISASMDCalcMethod1.H:14