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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
29 #include "fvPatchFieldMapper.H"
31 
32 #include "phaseSystem.H"
34 #include "ThermalDiffusivity.H"
36 #include "wallFvPatch.H"
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42 namespace compressible
43 {
44 
45 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
46 
48  = 50.0;
50  = 0.01;
52  = 10;
53 
54 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
55 
57 {
58  if (!isA<wallFvPatch>(patch()))
59  {
61  << "Patch type for patch " << patch().name() << " must be wall\n"
62  << "Current patch type is " << patch().type() << nl
63  << exit(FatalError);
64  }
65 }
66 
67 
70 (
71  const scalarField& Prat
72 ) const
73 {
74  return 9.24*(pow(Prat, 0.75) - 1)*(1 + 0.28*exp(-0.007*Prat));
75 }
76 
77 
80 (
81  const scalarField& P,
82  const scalarField& Prat
83 ) const
84 {
85  tmp<scalarField> typsf(new scalarField(this->size()));
86  scalarField& ypsf = typsf.ref();
87 
88  forAll(ypsf, facei)
89  {
90  scalar ypt = 11.0;
91 
92  for (int i=0; i<maxIters_; i++)
93  {
94  scalar f = ypt - (log(E_*ypt)/kappa_ + P[facei])/Prat[facei];
95  scalar df = 1 - 1.0/(ypt*kappa_*Prat[facei]);
96  scalar yptNew = ypt - f/df;
97 
98  if (yptNew < VSMALL)
99  {
100  ypsf[facei] = 0;
101  }
102  else if (mag(yptNew - ypt) < tolerance_)
103  {
104  ypsf[facei] = yptNew;
105  }
106  else
107  {
108  ypt = yptNew;
109  }
110  }
111 
112  ypsf[facei] = ypt;
113  }
114 
115  return typsf;
116 }
117 
120 (
121  const scalarField& prevAlphat
122 ) const
123 {
124 
125  // Lookup the fluid model
126  const phaseSystem& fluid =
127  db().lookupObject<phaseSystem>("phaseProperties");
128 
129  const phaseModel& phase
130  (
131  fluid.phases()[internalField().group()]
132  );
133 
134  const label patchi = patch().index();
135 
136  // Retrieve turbulence properties from model
137  const phaseCompressibleTurbulenceModel& turbModel =
138  db().lookupObject<phaseCompressibleTurbulenceModel>
139  (
141  );
142 
143  const scalar Cmu25 = pow025(Cmu_);
144 
145  const scalarField& y = turbModel.y()[patchi];
146 
147  const tmp<scalarField> tmuw = turbModel.mu(patchi);
148  const scalarField& muw = tmuw();
149 
150  const tmp<scalarField> talphaw = phase.thermo().alpha(patchi);
151  const scalarField& alphaw = talphaw();
152 
153  const tmp<volScalarField> tk = turbModel.k();
154  const volScalarField& k = tk();
155  const fvPatchScalarField& kw = k.boundaryField()[patchi];
156 
157  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
158  const scalarField magUp(mag(Uw.patchInternalField() - Uw));
159  const scalarField magGradUw(mag(Uw.snGrad()));
160 
161  const fvPatchScalarField& rhow = turbModel.rho().boundaryField()[patchi];
162  const fvPatchScalarField& hew =
163  phase.thermo().he().boundaryField()[patchi];
164 
165  const fvPatchScalarField& Tw =
166  phase.thermo().T().boundaryField()[patchi];
167 
169 
170  // Heat flux [W/m2] - lagging alphatw
171  const scalarField qDot
172  (
173  (prevAlphat + alphaw)*hew.snGrad()
174  );
175 
176  scalarField uTau(Cmu25*sqrt(kw));
177 
178  scalarField yPlus(uTau*y/(muw/rhow));
179 
180  scalarField Pr(muw/alphaw);
181 
182  // Molecular-to-turbulent Prandtl number ratio
183  scalarField Prat(Pr/Prt_);
184 
185  // Thermal sublayer thickness
186  scalarField P(this->Psmooth(Prat));
187 
188  scalarField yPlusTherm(this->yPlusTherm(P, Prat));
189 
190  tmp<scalarField> talphatConv(new scalarField(this->size()));
191  scalarField& alphatConv = talphatConv.ref();
192 
193  // Populate boundary values
194  forAll(alphatConv, facei)
195  {
196  // Evaluate new effective thermal diffusivity
197  scalar alphaEff = 0.0;
198  if (yPlus[facei] < yPlusTherm[facei])
199  {
200  scalar A = qDot[facei]*rhow[facei]*uTau[facei]*y[facei];
201  scalar B = qDot[facei]*Pr[facei]*yPlus[facei];
202  scalar C = Pr[facei]*0.5*rhow[facei]*uTau[facei]*sqr(magUp[facei]);
203  alphaEff = A/(B + C + VSMALL);
204  }
205  else
206  {
207  scalar A = qDot[facei]*rhow[facei]*uTau[facei]*y[facei];
208  scalar B =
209  qDot[facei]*Prt_*(1.0/kappa_*log(E_*yPlus[facei]) + P[facei]);
210  scalar magUc =
211  uTau[facei]/kappa_*log(E_*yPlusTherm[facei]) - mag(Uw[facei]);
212  scalar C =
213  0.5*rhow[facei]*uTau[facei]
214  *(Prt_*sqr(magUp[facei]) + (Pr[facei] - Prt_)*sqr(magUc));
215  alphaEff = A/(B + C + VSMALL);
216  }
217 
218  // Update convective heat transfer turbulent thermal diffusivity
219  alphatConv[facei] = max(0.0, alphaEff - alphaw[facei]);
220  }
221 
222  return talphatConv;
223 }
224 
225 
226 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
227 
230 (
231  const fvPatch& p,
233 )
234 :
236  Prt_(0.85),
237  Cmu_(0.09),
238  kappa_(0.41),
239  E_(9.8)
240 {
241  checkType();
242 }
243 
244 
247 (
248  const fvPatch& p,
250  const dictionary& dict
251 )
252 :
254  Prt_(dict.lookupOrDefault<scalar>("Prt", 0.85)),
255  Cmu_(dict.lookupOrDefault<scalar>("Cmu", 0.09)),
256  kappa_(dict.lookupOrDefault<scalar>("kappa", 0.41)),
257  E_(dict.lookupOrDefault<scalar>("E", 9.8))
258 {}
259 
260 
263 (
265  const fvPatch& p,
267  const fvPatchFieldMapper& mapper
268 )
269 :
271  Prt_(ptf.Prt_),
272  Cmu_(ptf.Cmu_),
273  kappa_(ptf.kappa_),
274  E_(ptf.E_)
275 {}
276 
277 
280 (
282 )
283 :
285  Prt_(awfpsf.Prt_),
286  Cmu_(awfpsf.Cmu_),
287  kappa_(awfpsf.kappa_),
288  E_(awfpsf.E_)
289 {}
290 
291 
294 (
297 )
298 :
300  Prt_(awfpsf.Prt_),
301  Cmu_(awfpsf.Cmu_),
302  kappa_(awfpsf.kappa_),
303  E_(awfpsf.E_)
304 {}
305 
306 
307 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
308 
310 {
311  if (updated())
312  {
313  return;
314  }
315 
316  operator==(calcAlphat(*this));
317 
318  fixedValueFvPatchScalarField::updateCoeffs();
319 }
320 
321 
323 (
324  Ostream& os
325 ) const
326 {
328  os.writeEntry("Prt", Prt_);
329  os.writeEntry("Cmu", Cmu_);
330  os.writeEntry("kappa", kappa_);
331  os.writeEntry("E", E_);
332  dmdt_.writeEntry("dmdt", os);
333  writeEntry("value", os);
334 }
335 
336 
337 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
338 
340 (
343 );
344 
345 
346 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
347 
348 } // End namespace compressible
349 } // End namespace Foam
350 
351 // ************************************************************************* //
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:57
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: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:59
Foam::compressible::alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::checkType
void checkType()
Check the type of the patch.
Definition: alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.C:56
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:70
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:290
Foam::pow025
dimensionedScalar pow025(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:133
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:258
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:62
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:63
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:230
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:80
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:3
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:323
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:1015
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:355
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::compressible::alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.C:309
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::IOobject::groupName
static word groupName(StringType name, const word &group)
Create dot-delimited name.group.
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:219
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:69
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