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  = 0.01;
51  = 10;
52 
53 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
54 
56 {
57  if (!isA<wallFvPatch>(patch()))
58  {
60  << "Patch type for patch " << patch().name() << " must be wall\n"
61  << "Current patch type is " << patch().type() << nl
62  << exit(FatalError);
63  }
64 }
65 
66 
69 (
70  const scalarField& Prat
71 ) const
72 {
73  return 9.24*(pow(Prat, 0.75) - 1)*(1 + 0.28*exp(-0.007*Prat));
74 }
75 
76 
79 (
80  const scalarField& P,
81  const scalarField& Prat
82 ) const
83 {
84  auto typsf = tmp<scalarField>::New(this->size());
85  auto& ypsf = typsf.ref();
86 
87  forAll(ypsf, facei)
88  {
89  scalar ypt = 11.0;
90 
91  for (int i = 0; i < maxIters_; ++i)
92  {
93  const scalar f = ypt - (log(E_*ypt)/kappa_ + P[facei])/Prat[facei];
94  const scalar df = 1.0 - 1.0/(ypt*kappa_*Prat[facei]);
95  const scalar yptNew = ypt - f/df;
96 
97  if (yptNew < VSMALL)
98  {
99  ypsf[facei] = 0;
100  }
101  else if (mag(yptNew - ypt) < tolerance_)
102  {
103  ypsf[facei] = yptNew;
104  }
105  else
106  {
107  ypt = yptNew;
108  }
109  }
110 
111  ypsf[facei] = ypt;
112  }
113 
114  return typsf;
115 }
116 
117 
120 (
121  const scalarField& prevAlphat
122 ) const
123 {
124  // Lookup the fluid model
125  const phaseSystem& fluid =
126  db().lookupObject<phaseSystem>("phaseProperties");
127 
128  const phaseModel& phase = fluid.phases()[internalField().group()];
129 
130  const label patchi = patch().index();
131 
132  // Retrieve turbulence properties from model
133  const auto& turbModel =
134  db().lookupObject<phaseCompressibleTurbulenceModel>
135  (
137  );
138 
139  const scalar Cmu25 = pow025(Cmu_);
140 
141  const scalarField& y = turbModel.y()[patchi];
142 
143  const tmp<scalarField> tmuw = turbModel.mu(patchi);
144  const scalarField& muw = tmuw();
145 
146  const tmp<scalarField> talphaw = phase.thermo().alpha(patchi);
147  const scalarField& alphaw = talphaw();
148 
149  const tmp<volScalarField> tk = turbModel.k();
150  const volScalarField& k = tk();
151  const fvPatchScalarField& kw = k.boundaryField()[patchi];
152 
153  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
154  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
155  const scalarField magGradUw(mag(Uw.snGrad()));
156 
157  const fvPatchScalarField& rhow = turbModel.rho().boundaryField()[patchi];
158  const fvPatchScalarField& hew =
159  phase.thermo().he().boundaryField()[patchi];
160 
161  // Heat flux [W/m2] - lagging alphatw
162  const scalarField qDot
163  (
164  (prevAlphat + alphaw)*hew.snGrad()
165  );
166 
167  scalarField uTau(Cmu25*sqrt(kw));
168 
169  scalarField yPlus(uTau*y/(muw/rhow));
170 
171  const scalarField Pr(muw/alphaw);
172 
173  // Molecular-to-turbulent Prandtl number ratio
174  const scalarField Prat(Pr/Prt_);
175 
176  // Thermal sublayer thickness
177  const scalarField P(this->Psmooth(Prat));
178 
179  tmp<scalarField> tyPlusTherm = this->yPlusTherm(P, Prat);
180  const scalarField& yPlusTherm = tyPlusTherm.cref();
181 
182  auto talphatConv = tmp<scalarField>::New(this->size());
183  auto& alphatConv = talphatConv.ref();
184 
185  // Populate boundary values
186  forAll(alphatConv, facei)
187  {
188  // Evaluate new effective thermal diffusivity
189  scalar alphaEff = 0.0;
190  if (yPlus[facei] < yPlusTherm[facei])
191  {
192  const scalar A = qDot[facei]*rhow[facei]*uTau[facei]*y[facei];
193  const scalar B = qDot[facei]*Pr[facei]*yPlus[facei];
194  const scalar C =
195  Pr[facei]*0.5*rhow[facei]*uTau[facei]*sqr(magUp[facei]);
196  alphaEff = A/(B + C + VSMALL);
197  }
198  else
199  {
200  const scalar A = qDot[facei]*rhow[facei]*uTau[facei]*y[facei];
201  const scalar B =
202  qDot[facei]*Prt_*(1.0/kappa_*log(E_*yPlus[facei]) + P[facei]);
203  const scalar magUc =
204  uTau[facei]/kappa_*log(E_*yPlusTherm[facei]) - mag(Uw[facei]);
205  const scalar C =
206  0.5*rhow[facei]*uTau[facei]
207  *(Prt_*sqr(magUp[facei]) + (Pr[facei] - Prt_)*sqr(magUc));
208  alphaEff = A/(B + C + VSMALL);
209  }
210 
211  // Update convective heat transfer turbulent thermal diffusivity
212  alphatConv[facei] = max(0.0, alphaEff - alphaw[facei]);
213  }
214 
215  return talphatConv;
216 }
217 
218 
219 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
220 
223 (
224  const fvPatch& p,
226 )
227 :
229  Prt_(0.85),
230  Cmu_(0.09),
231  kappa_(0.41),
232  E_(9.8)
233 {
234  checkType();
235 }
236 
237 
240 (
241  const fvPatch& p,
243  const dictionary& dict
244 )
245 :
247  Prt_(dict.getOrDefault<scalar>("Prt", 0.85)),
248  Cmu_(dict.getOrDefault<scalar>("Cmu", 0.09)),
249  kappa_(dict.getOrDefault<scalar>("kappa", 0.41)),
250  E_(dict.getOrDefault<scalar>("E", 9.8))
251 {}
252 
253 
256 (
258  const fvPatch& p,
260  const fvPatchFieldMapper& mapper
261 )
262 :
264  Prt_(ptf.Prt_),
265  Cmu_(ptf.Cmu_),
266  kappa_(ptf.kappa_),
267  E_(ptf.E_)
268 {}
269 
270 
273 (
275 )
276 :
278  Prt_(awfpsf.Prt_),
279  Cmu_(awfpsf.Cmu_),
280  kappa_(awfpsf.kappa_),
281  E_(awfpsf.E_)
282 {}
283 
284 
287 (
290 )
291 :
293  Prt_(awfpsf.Prt_),
294  Cmu_(awfpsf.Cmu_),
295  kappa_(awfpsf.kappa_),
296  E_(awfpsf.E_)
297 {}
298 
299 
300 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
301 
303 {
304  if (updated())
305  {
306  return;
307  }
308 
309  operator==(calcAlphat(*this));
310 
311  fixedValueFvPatchScalarField::updateCoeffs();
312 }
313 
314 
316 (
317  Ostream& os
318 ) const
319 {
321  os.writeEntry("Prt", Prt_);
322  os.writeEntry("Cmu", Cmu_);
323  os.writeEntry("kappa", kappa_);
324  os.writeEntry("E", E_);
325  dmdt_.writeEntry("dmdt", os);
326  writeEntry("value", os);
327 }
328 
329 
330 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
331 
333 (
336 );
337 
338 
339 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
340 
341 } // End namespace compressible
342 } // End namespace Foam
343 
344 // ************************************************************************* //
Foam::fvPatchField< scalar >
Foam::fvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:384
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:225
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:120
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:55
ThermalDiffusivity.H
Foam::compressible::alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField
This boundary condition provides a thermal wall function for turbulent thermal diffusivity (usuallyal...
Definition: alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.H:132
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:100
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:69
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::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:141
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
Foam::compressible::alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::Cmu_
scalar Cmu_
Empirical model coefficient.
Definition: alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.H:144
Foam::compressible::alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::E_
scalar E_
Wall roughness parameter.
Definition: alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.H:150
Foam::compressible::alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField
alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.C:223
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:79
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::compressible::alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.C:316
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::compressible::alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::tolerance_
static scalar tolerance_
Absolute tolerance.
Definition: alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.H:155
Foam::compressible::alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::kappa_
scalar kappa_
Von Karman constant.
Definition: alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.H:147
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::compressible::alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.C:302
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::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
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:66
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_
Maximum number of iterations.
Definition: alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.H:158
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::tmp::cref
const T & cref() const
Definition: tmpI.H:213
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