alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.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) 2015-2019 OpenFOAM Foundation
9  Copyright (C) 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"
32 
33 #include "phaseSystem.H"
35 #include "ThermalDiffusivity.H"
37 #include "wallFvPatch.H"
38 
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43 namespace compressible
44 {
45 
46 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
47 
49  = 50.0;
51  = 0.01;
53  = 10;
54 
55 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
56 
58 {
59  if (!isA<wallFvPatch>(patch()))
60  {
62  << "Patch type for patch " << patch().name() << " must be wall\n"
63  << "Current patch type is " << patch().type() << nl
64  << exit(FatalError);
65  }
66 }
67 
68 
71 (
72  const scalarField& Prat
73 ) const
74 {
75  return 9.24*(pow(Prat, 0.75) - 1)*(1 + 0.28*exp(-0.007*Prat));
76 }
77 
78 
81 (
82  const scalarField& P,
83  const scalarField& Prat
84 ) const
85 {
86  tmp<scalarField> typsf(new scalarField(this->size()));
87  scalarField& ypsf = typsf.ref();
88 
89  forAll(ypsf, facei)
90  {
91  scalar ypt = 11.0;
92 
93  for (int i=0; i<maxIters_; i++)
94  {
95  scalar f = ypt - (log(E_*ypt)/kappa_ + P[facei])/Prat[facei];
96  scalar df = 1 - 1.0/(ypt*kappa_*Prat[facei]);
97  scalar yptNew = ypt - f/df;
98 
99  if (yptNew < VSMALL)
100  {
101  ypsf[facei] = 0;
102  }
103  else if (mag(yptNew - ypt) < tolerance_)
104  {
105  ypsf[facei] = yptNew;
106  }
107  else
108  {
109  ypt = yptNew;
110  }
111  }
112 
113  ypsf[facei] = ypt;
114  }
115 
116  return typsf;
117 }
118 
121 (
122  const scalarField& prevAlphat
123 ) const
124 {
125 
126  // Lookup the fluid model
127  const phaseSystem& fluid =
128  db().lookupObject<phaseSystem>("phaseProperties");
129 
130  const phaseModel& phase
131  (
132  fluid.phases()[internalField().group()]
133  );
134 
135  const label patchi = patch().index();
136 
137  // Retrieve turbulence properties from model
138  const phaseCompressibleTurbulenceModel& turbModel =
139  db().lookupObject<phaseCompressibleTurbulenceModel>
140  (
142  );
143 
144  const scalar Cmu25 = pow025(Cmu_);
145 
146  const scalarField& y = turbModel.y()[patchi];
147 
148  const tmp<scalarField> tmuw = turbModel.mu(patchi);
149  const scalarField& muw = tmuw();
150 
151  const tmp<scalarField> talphaw = phase.thermo().alpha(patchi);
152  const scalarField& alphaw = talphaw();
153 
154  const tmp<volScalarField> tk = turbModel.k();
155  const volScalarField& k = tk();
156  const fvPatchScalarField& kw = k.boundaryField()[patchi];
157 
158  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
159  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
160  const scalarField magGradUw(mag(Uw.snGrad()));
161 
162  const fvPatchScalarField& rhow = turbModel.rho().boundaryField()[patchi];
163  const fvPatchScalarField& hew =
164  phase.thermo().he().boundaryField()[patchi];
165 
166  const fvPatchScalarField& Tw =
167  phase.thermo().T().boundaryField()[patchi];
168 
170 
171  // Heat flux [W/m2] - lagging alphatw
172  const scalarField qDot
173  (
174  (prevAlphat + alphaw)*hew.snGrad()
175  );
176 
177  scalarField uTau(Cmu25*sqrt(kw));
178 
179  scalarField yPlus(uTau*y/(muw/rhow));
180 
181  scalarField Pr(muw/alphaw);
182 
183  // Molecular-to-turbulent Prandtl number ratio
184  scalarField Prat(Pr/Prt_);
185 
186  // Thermal sublayer thickness
187  scalarField P(this->Psmooth(Prat));
188 
189  scalarField yPlusTherm(this->yPlusTherm(P, Prat));
190 
191  tmp<scalarField> talphatConv(new scalarField(this->size()));
192  scalarField& alphatConv = talphatConv.ref();
193 
194  // Populate boundary values
195  forAll(alphatConv, facei)
196  {
197  // Evaluate new effective thermal diffusivity
198  scalar alphaEff = 0.0;
199  if (yPlus[facei] < yPlusTherm[facei])
200  {
201  scalar A = qDot[facei]*rhow[facei]*uTau[facei]*y[facei];
202  scalar B = qDot[facei]*Pr[facei]*yPlus[facei];
203  scalar C = Pr[facei]*0.5*rhow[facei]*uTau[facei]*sqr(magUp[facei]);
204  alphaEff = A/(B + C + VSMALL);
205  }
206  else
207  {
208  scalar A = qDot[facei]*rhow[facei]*uTau[facei]*y[facei];
209  scalar B =
210  qDot[facei]*Prt_*(1.0/kappa_*log(E_*yPlus[facei]) + P[facei]);
211  scalar magUc =
212  uTau[facei]/kappa_*log(E_*yPlusTherm[facei]) - mag(Uw[facei]);
213  scalar C =
214  0.5*rhow[facei]*uTau[facei]
215  *(Prt_*sqr(magUp[facei]) + (Pr[facei] - Prt_)*sqr(magUc));
216  alphaEff = A/(B + C + VSMALL);
217  }
218 
219  // Update convective heat transfer turbulent thermal diffusivity
220  alphatConv[facei] = max(0.0, alphaEff - alphaw[facei]);
221  }
222 
223  return talphatConv;
224 }
225 
226 
227 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
228 
231 (
232  const fvPatch& p,
234 )
235 :
237  Prt_(0.85),
238  Cmu_(0.09),
239  kappa_(0.41),
240  E_(9.8)
241 {
242  checkType();
243 }
244 
245 
248 (
249  const fvPatch& p,
251  const dictionary& dict
252 )
253 :
255  Prt_(dict.getOrDefault<scalar>("Prt", 0.85)),
256  Cmu_(dict.getOrDefault<scalar>("Cmu", 0.09)),
257  kappa_(dict.getOrDefault<scalar>("kappa", 0.41)),
258  E_(dict.getOrDefault<scalar>("E", 9.8))
259 {}
260 
261 
264 (
266  const fvPatch& p,
268  const fvPatchFieldMapper& mapper
269 )
270 :
272  Prt_(ptf.Prt_),
273  Cmu_(ptf.Cmu_),
274  kappa_(ptf.kappa_),
275  E_(ptf.E_)
276 {}
277 
278 
281 (
283 )
284 :
286  Prt_(awfpsf.Prt_),
287  Cmu_(awfpsf.Cmu_),
288  kappa_(awfpsf.kappa_),
289  E_(awfpsf.E_)
290 {}
291 
292 
295 (
298 )
299 :
301  Prt_(awfpsf.Prt_),
302  Cmu_(awfpsf.Cmu_),
303  kappa_(awfpsf.kappa_),
304  E_(awfpsf.E_)
305 {}
306 
307 
308 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
309 
311 {
312  if (updated())
313  {
314  return;
315  }
316 
317  operator==(calcAlphat(*this));
318 
319  fixedValueFvPatchScalarField::updateCoeffs();
320 }
321 
322 
324 (
325  Ostream& os
326 ) const
327 {
329  os.writeEntry("Prt", Prt_);
330  os.writeEntry("Cmu", Cmu_);
331  os.writeEntry("kappa", kappa_);
332  os.writeEntry("E", E_);
333  dmdt_.writeEntry("dmdt", os);
334  writeEntry("value", os);
335 }
336 
337 
338 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
339 
341 (
344 );
345 
346 
347 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
348 
349 } // End namespace compressible
350 } // End namespace Foam
351 
352 // ************************************************************************* //
Foam::fvPatchField< scalar >
Foam::fvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:364
Foam::phaseModel
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:54
Foam::fvPatchField::snGrad
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:217
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::compressible::alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::calcAlphat
tmp< scalarField > calcAlphat(const scalarField &prevAlphat) const
Update turbulent thermal diffusivity.
Definition: alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.C:121
compressibleTurbulenceModel.H
Foam::phase
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:54
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::compressible::alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::checkType
void checkType()
Check the type of the patch.
Definition: alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.C:57
ThermalDiffusivity.H
Foam::compressible::alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField
This boundary condition provides a thermal wall function for turbulent thermal diffusivity (usuallyal...
Definition: alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.H:103
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
B
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.H
wallFvPatch.H
Foam::compressible::alphatPhaseChangeWallFunctionFvPatchScalarField
Abstract base-class for all alphatWallFunctions supporting phase-change.
Definition: alphatPhaseChangeWallFunctionFvPatchScalarField.H:57
A
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
fvPatchFieldMapper.H
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261
fluid
twoPhaseSystem & fluid
Definition: setRegionFluidFields.H:3
Foam::compressible::alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::Psmooth
tmp< scalarField > Psmooth(const scalarField &Prat) const
'P' function
Definition: alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.C:71
magUp
scalar magUp
Definition: evaluateNearWall.H:10
PhaseCompressibleTurbulenceModel.H
Foam::compressible::makePatchTypeField
makePatchTypeField(fvPatchScalarField, alphatFixedDmdtWallBoilingWallFunctionFvPatchScalarField)
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::pow025
dimensionedScalar pow025(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:133
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:228
Foam::Field< scalar >
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Foam::compressible::alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::Prt_
scalar Prt_
Turbulent Prandtl number.
Definition: alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.H:113
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:225
Foam::compressible::alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::Cmu_
scalar Cmu_
Cmu coefficient.
Definition: alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.H:116
Foam::compressible::alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::E_
scalar E_
E coefficient.
Definition: alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.H:122
Foam::compressible::alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField
alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.C:231
Foam::compressible::alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::yPlusTherm
tmp< scalarField > yPlusTherm(const scalarField &P, const scalarField &Prat) const
Calculate y+ at the edge of the thermal laminar sublayer.
Definition: alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.C:81
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::compressible::alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::maxExp_
static scalar maxExp_
Definition: alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.H:126
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:121
Foam::log
dimensionedScalar log(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:262
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::compressible::alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.C:324
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::phase::name
const word & name() const
Definition: phase.H:111
Foam::GeometricField::T
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field)
Definition: GeometricField.C:1046
Foam::compressible::alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::tolerance_
static scalar tolerance_
Definition: alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.H:127
Foam::compressible::alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::kappa_
scalar kappa_
Von Karman constant.
Definition: alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.H:119
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::compressible::alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.C:310
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)
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
alphaEff
volScalarField alphaEff("alphaEff", turbulence->nu()/Pr+alphat)
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:232
yPlus
scalar yPlus
Definition: evaluateNearWall.H:16
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::phaseSystem
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:63
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
Foam::C
Graphite solid properties.
Definition: C.H:50
Foam::compressible::alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::maxIters_
static label maxIters_
Definition: alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.H:128
Foam::GeometricField< scalar, fvPatchField, volMesh >
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
y
scalar y
Definition: LISASMDCalcMethod1.H:14