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"
33 #include "wallFvPatch.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 namespace compressible
41 {
42 
43 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44 
45 scalar alphatJayatillekeWallFunctionFvPatchScalarField::maxExp_ = 50.0;
46 scalar alphatJayatillekeWallFunctionFvPatchScalarField::tolerance_ = 0.01;
47 label alphatJayatillekeWallFunctionFvPatchScalarField::maxIters_ = 10;
48 
49 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50 
51 void 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 
63 tmp<scalarField> alphatJayatillekeWallFunctionFvPatchScalarField::yPlus
64 (
65  const compressible::turbulenceModel& turbModel
66 ) const
67 {
68  const label patchi = patch().index();
69  const tmp<volScalarField> tnut = turbModel.nut();
70  const volScalarField& nut = tnut();
71 
72  if (isA<nutWallFunctionFvPatchScalarField>(nut.boundaryField()[patchi]))
73  {
74  const nutWallFunctionFvPatchScalarField& nutPf =
75  dynamic_cast<const nutWallFunctionFvPatchScalarField&>
76  (
77  nut.boundaryField()[patchi]
78  );
79 
80  return nutPf.yPlus();
81  }
82  else
83  {
84  const scalarField& y = turbModel.y()[patchi];
85  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
86 
87  return
88  y*sqrt(turbModel.nuEff(patchi)*mag(Uw.snGrad()))
89  /turbModel.nu(patchi);
90  }
91 }
92 
93 
94 scalar alphatJayatillekeWallFunctionFvPatchScalarField::Psmooth
95 (
96  const scalar Prat
97 ) const
98 {
99  return 9.24*(pow(Prat, 0.75) - 1.0)*(1.0 + 0.28*exp(-0.007*Prat));
100 }
101 
102 
103 scalar alphatJayatillekeWallFunctionFvPatchScalarField::yPlusTherm
104 (
105  const scalar P,
106  const scalar Prat
107 ) const
108 {
109  scalar ypt = 11.0;
110 
111  for (int i=0; i<maxIters_; i++)
112  {
113  scalar f = ypt - (log(E_*ypt)/kappa_ + P)/Prat;
114  scalar df = 1.0 - 1.0/(ypt*kappa_*Prat);
115  scalar yptNew = ypt - f/df;
116 
117  if (yptNew < VSMALL)
118  {
119  return 0;
120  }
121  else if (mag(yptNew - ypt) < tolerance_)
122  {
123  return yptNew;
124  }
125  else
126  {
127  ypt = yptNew;
128  }
129  }
130 
131  return ypt;
132 }
133 
134 
135 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
136 
139 (
140  const fvPatch& p,
142 )
143 :
144  fixedValueFvPatchScalarField(p, iF),
145  Prt_(0.85),
146  kappa_(0.41),
147  E_(9.8)
148 {
149  checkType();
150 }
151 
152 
155 (
157  const fvPatch& p,
159  const fvPatchFieldMapper& mapper
160 )
161 :
162  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
163  Prt_(ptf.Prt_),
164  kappa_(ptf.kappa_),
165  E_(ptf.E_)
166 {}
167 
168 
171 (
172  const fvPatch& p,
174  const dictionary& dict
175 )
176 :
177  fixedValueFvPatchScalarField(p, iF, dict),
178  Prt_(dict.getOrDefault<scalar>("Prt", 0.85)),
179  kappa_(dict.getOrDefault<scalar>("kappa", 0.41)),
180  E_(dict.getOrDefault<scalar>("E", 9.8))
181 {
182  checkType();
183 }
184 
185 
188 (
190 )
191 :
192  fixedValueFvPatchScalarField(awfpsf),
193  Prt_(awfpsf.Prt_),
194  kappa_(awfpsf.kappa_),
195  E_(awfpsf.E_)
196 {
197  checkType();
198 }
199 
200 
203 (
206 )
207 :
208  fixedValueFvPatchScalarField(awfpsf, iF),
209  Prt_(awfpsf.Prt_),
210  kappa_(awfpsf.kappa_),
211  E_(awfpsf.E_)
212 {
213  checkType();
214 }
215 
216 
217 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
218 
220 {
221  if (updated())
222  {
223  return;
224  }
225 
226  const label patchi = patch().index();
227 
228  // Retrieve turbulence properties from model
229  const compressible::turbulenceModel& turbModel =
230  db().lookupObject<compressible::turbulenceModel>
231  (
233  (
234  compressible::turbulenceModel::propertiesName,
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  scalar yPlus = yPlusp[facei];
269 
270  scalar uTau = yPlus/y[facei]*(muw[facei]/rhow[facei]);
271 
272  // Molecular Prandtl number
273  scalar Pr = muw[facei]/alphaw[facei];
274 
275  // Molecular-to-turbulent Prandtl number ratio
276  scalar Prat = Pr/Prt_;
277 
278  // Thermal sublayer thickness
279  scalar P = Psmooth(Prat);
280  scalar yPlusTherm = this->yPlusTherm(P, Prat);
281 
282  // Evaluate new effective thermal diffusivity
283  scalar alphaEff = 0.0;
284  if (yPlus < yPlusTherm)
285  {
286  scalar A = qDot[facei]*rhow[facei]*uTau*y[facei];
287  scalar B = qDot[facei]*Pr*yPlus;
288  scalar C = Pr*0.5*rhow[facei]*uTau*sqr(magUp[facei]);
289  alphaEff = A/(B + C + VSMALL);
290  }
291  else
292  {
293  scalar A = qDot[facei]*rhow[facei]*uTau*y[facei];
294  scalar B = qDot[facei]*Prt_*(1.0/kappa_*log(E_*yPlus) + P);
295  scalar magUc = uTau/kappa_*log(E_*yPlusTherm) - mag(Uw[facei]);
296  scalar C =
297  0.5*rhow[facei]*uTau
298  *(Prt_*sqr(magUp[facei]) + (Pr - Prt_)*sqr(magUc));
299  alphaEff = A/(B + C + VSMALL);
300  }
301 
302  // Update turbulent thermal diffusivity
303  alphatw[facei] = max(0.0, alphaEff - alphaw[facei]);
304 
305  if (debug)
306  {
307  Info<< " uTau = " << uTau << nl
308  << " Pr = " << Pr << nl
309  << " Prt = " << Prt_ << nl
310  << " qDot = " << qDot[facei] << nl
311  << " yPlus = " << yPlus << nl
312  << " yPlusTherm = " << yPlusTherm << nl
313  << " alphaEff = " << alphaEff << nl
314  << " alphaw = " << alphaw[facei] << nl
315  << " alphatw = " << alphatw[facei] << nl
316  << endl;
317  }
318  }
319 
321 }
322 
323 
325 {
327  os.writeEntry("Prt", Prt_);
328  os.writeEntry("kappa", kappa_);
329  os.writeEntry("E", E_);
330  writeEntry("value", os);
331 }
332 
333 
334 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
335 
337 (
340 );
341 
342 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
343 
344 } // End namespace compressible
345 } // End namespace Foam
346 
347 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::fvPatchField< vector >
volFields.H
Foam::fvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:384
Foam::compressible::alphatJayatillekeWallFunctionFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: alphatJayatillekeWallFunctionFvPatchScalarField.C:219
Foam::fvPatchField::snGrad
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:225
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
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
B
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
nutWallFunctionFvPatchScalarField.H
wallFvPatch.H
A
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
fvPatchFieldMapper.H
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261
magUp
scalar magUp
Definition: evaluateNearWall.H:10
Foam::compressible::makePatchTypeField
makePatchTypeField(fvPatchScalarField, alphatFixedDmdtWallBoilingWallFunctionFvPatchScalarField)
Foam::compressible::alphatJayatillekeWallFunctionFvPatchScalarField
This boundary condition provides a thermal wall function for turbulent thermal diffusivity (usuallyal...
Definition: alphatJayatillekeWallFunctionFvPatchScalarField.H:101
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::compressible::alphatJayatillekeWallFunctionFvPatchScalarField::write
void write(Ostream &) const
Write.
Definition: alphatJayatillekeWallFunctionFvPatchScalarField.C:324
Foam::Field< scalar >
Foam::compressible::turbulenceModel
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
Definition: turbulentFluidThermoModel.H:63
Foam::compressible::alphatJayatillekeWallFunctionFvPatchScalarField::alphatJayatillekeWallFunctionFvPatchScalarField
alphatJayatillekeWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: alphatJayatillekeWallFunctionFvPatchScalarField.C:139
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Pr
dimensionedScalar Pr("Pr", dimless, laminarTransport)
Foam::fvPatchField::patchInternalField
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:233
alphatJayatillekeWallFunctionFvPatchScalarField.H
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
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
uTau
scalar uTau
Definition: evaluateNearWall.H:14
compressible
bool compressible
Definition: pEqn.H:2
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::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::fvPatchField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:321
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::fvPatchVectorField
fvPatchField< vector > fvPatchVectorField
Definition: fvPatchFieldsFwd.H:43
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::ThermalDiffusivity::alphaEff
virtual tmp< volScalarField > alphaEff() const
Return the effective turbulent thermal diffusivity for enthalpy.
Definition: ThermalDiffusivity.H:160
f
labelList f(nPoints)
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::ThermalDiffusivity
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
Definition: phaseCompressibleTurbulenceModelFwd.H:47
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
alphaEff
volScalarField alphaEff("alphaEff", turbulence->nu()/Pr+alphat)
yPlus
scalar yPlus
Definition: evaluateNearWall.H:16
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::C
Graphite solid properties.
Definition: C.H:50
Foam::ThermalDiffusivity::alpha
virtual tmp< volScalarField > alpha() const
Return the laminar thermal diffusivity for enthalpy [kg/m/s].
Definition: ThermalDiffusivity.H:125
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
y
scalar y
Definition: LISASMDCalcMethod1.H:14