wallBoilingHeatTransfer.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2019-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 #include "phaseSystem.H"
30 #include "diameterModel.H"
31 
32 using namespace Foam::constant::mathematical;
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace heatTransferModels
39 {
42  (
46  );
47 }
48 }
49 
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
54 (
55  const dictionary& dict,
56  const phaseInterface& interface,
57  const bool registerObject
58 )
59 :
60  heatTransferModel(dict, interface, registerObject),
61  interface_
62  (
63  interface.modelCast<wallBoilingHeatTransfer, dispersedPhaseInterface>()
64  ),
65  otherInterface_
66  (
67  interface
69  .otherInterface()
70  ),
71  vapourPhaseName_(dict.lookup("vapourPhase")),
72  liquidPhaseName_(dict.lookup("liquidPhase")),
73  heatTransferModel_
74  (
76  (
77  dict.subDict("heatTransferModel"),
78  interface,
79  false,
80  false
81  )
82  ),
83  relax_(dict.lookupOrDefault<scalar>("relax", scalar(1))),
84  seedFraction_(dict.lookupOrDefault<scalar>("nucleationSeedFraction", 1e-4)),
85  partitioningModel_(nullptr),
86  nucleationSiteModel_(nullptr),
87  departureDiamModel_(nullptr),
88  departureFreqModel_(nullptr),
89  wetFraction_
90  (
91  IOobject
92  (
93  IOobject::groupName("fWallBoiling", interface_.name()),
94  interface_.mesh().time().timeName(),
95  interface_.mesh(),
96  IOobject::READ_IF_PRESENT,
97  IOobject::AUTO_WRITE
98  ),
99  interface_.mesh(),
101  ),
102  dDep_
103  (
104  IOobject
105  (
106  IOobject::groupName("departureDiameter", interface_.name()),
107  interface_.mesh().time().timeName(),
108  interface_.mesh(),
109  IOobject::READ_IF_PRESENT,
110  IOobject::AUTO_WRITE
111  ),
112  interface_.mesh(),
114  ),
115  fDep_
116  (
117  IOobject
118  (
119  IOobject::groupName("departureFrequency", interface_.name()),
120  interface_.mesh().time().timeName(),
121  interface_.mesh(),
122  IOobject::READ_IF_PRESENT,
123  IOobject::AUTO_WRITE
124  ),
125  interface_.mesh(),
127  ),
128  nucleationSiteDensity_
129  (
130  IOobject
131  (
132  IOobject::groupName("nucleationSites", interface_.name()),
133  interface_.mesh().time().timeName(),
134  interface_.mesh(),
135  IOobject::READ_IF_PRESENT,
136  IOobject::AUTO_WRITE
137  ),
138  interface_.mesh(),
140  ),
141  dmdtf_
142  (
143  IOobject
144  (
145  IOobject::groupName(typedName("dmdtf"), interface_.name()),
146  interface_.mesh().time().timeName(),
147  interface_.mesh(),
148  IOobject::READ_IF_PRESENT,
149  IOobject::AUTO_WRITE
150  ),
151  interface_.mesh(),
153  ),
154  qq_
155  (
156  IOobject
157  (
158  IOobject::groupName("qq", interface_.name()),
159  interface_.mesh().time().timeName(),
160  interface_.mesh(),
161  IOobject::READ_IF_PRESENT,
162  IOobject::AUTO_WRITE
163  ),
164  interface_.mesh(),
166  ),
167  Tsurface_
168  (
169  IOobject
170  (
171  IOobject::groupName("Tsurface", interface_.name()),
172  interface_.mesh().time().timeName(),
173  interface_.mesh(),
174  IOobject::READ_IF_PRESENT,
175  IOobject::AUTO_WRITE
176  ),
177  interface_.mesh(),
179  ),
180  K_
181  (
182  IOobject
183  (
184  IOobject::groupName(typedName("K"), interface_.name()),
185  interface_.mesh().time().timeName(),
186  interface_.mesh(),
187  IOobject::READ_IF_PRESENT,
188  IOobject::AUTO_WRITE
189  ),
190  heatTransferModel_->K()
191  )
192 {
193  partitioningModel_ =
195  (
196  dict.subDict("partitioningModel")
197  );
198 
199  nucleationSiteModel_ =
201  (
202  dict.subDict("nucleationSiteModel")
203  );
204 
205  departureDiamModel_ =
207  (
208  dict.subDict("departureDiameterModel")
209  );
210 
211  departureFreqModel_ =
213  (
214  dict.subDict("departureFrequencyModel")
215  );
216 }
217 
218 
219 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
220 
223 {}
224 
225 
226 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
227 
230 (
231  const scalar residualAlpha
232 ) const
233 {
234  const phaseSystem& fluid = interface_.fluid();
235 
236  const phaseModel& liquid = interface_.continuous();
237  const phaseModel& vapour = fluid.phases()[vapourPhaseName_];
238  const phaseModel& solid = interface_.dispersed();
239 
240  const rhoThermo& lThermo = liquid.thermo();
241  const rhoThermo& vThermo = vapour.thermo();
242  const rhoThermo& sThermo = solid.thermo();
243 
244  // Estimate the surface temperature from the surrounding temperature and
245  // heat transfer coefficients. Note that a lagged value of K is used in
246  // this calculation.
247  {
248  static const dimensionedScalar smallK(dimK, rootVSmall);
249 
250  const heatTransferModel& otherModel =
251  fluid.lookupInterfacialModel<heatTransferModel>(otherInterface_());
252  const volScalarField otherK(otherModel.K());
253 
254  Tsurface_ =
255  (K_*lThermo.T() + otherK*sThermo.T())
256  /max(K_ + otherK, smallK);
257  }
258 
259  const volScalarField Tsat
260  (
262  <
264  >
265  (phaseInterface(liquid, vapour))
266  .Tsat(liquid.thermo().p())()
267  );
268 
269  const volScalarField L
270  (
271  vThermo.ha(lThermo.p(), Tsat) - lThermo.ha()
272  );
273 
274  // Wetted fraction
275  wetFraction_ =
276  partitioningModel_->wetFraction(liquid/max(1 - solid, small));
277 
278  // Bubble departure diameter
279  dDep_ =
280  departureDiamModel_->dDeparture
281  (
282  liquid,
283  vapour,
284  solid,
285  Tsurface_,
286  Tsat,
287  L
288  );
289 
290  // Bubble departure frequency
291  fDep_ =
292  departureFreqModel_->fDeparture
293  (
294  liquid,
295  vapour,
296  solid,
297  Tsurface_,
298  Tsat,
299  L,
300  dDep_
301  );
302 
303  // Nucleation site density
304  nucleationSiteDensity_ =
305  nucleationSiteModel_->nucleationSiteDensity
306  (
307  liquid,
308  vapour,
309  solid,
310  Tsurface_,
311  Tsat,
312  L,
313  dDep_,
314  fDep_
315  );
316 
317  const tmp<volScalarField> tlrho(lThermo.rho());
318  const volScalarField& lrho = tlrho();
319  const tmp<volScalarField> tvrho(vThermo.rho());
320  const volScalarField& vrho = tvrho();
321  const volScalarField& lCp = lThermo.Cp();
322  const volScalarField& lkappa = lThermo.kappa();
323  //const volScalarField lalpha(lkappa/lCp);
324 
325  // Area fractions
326 
327  // Del Valle & Kenning (1985)
328  const volScalarField Ja
329  (
330  lrho*lCp*(Tsat - min(Tsurface_, Tsat))/(lrho*L)
331  );
332 
333  const volScalarField Al
334  (
335  wetFraction_*4.8*exp(min(-Ja/80, log(vGreat)))
336  );
337 
338  volScalarField A2
339  (
340  min(pi*sqr(dDep_)*nucleationSiteDensity_*Al/4, scalar(1))
341  );
342  const scalarField A1(max(1 - A2, scalar(1e-4)));
343  volScalarField A2E
344  (
345  min(pi*sqr(dDep_)*nucleationSiteDensity_*Al/4, scalar(5))
346  );
347 
348  const volScalarField Av(solid.dPtr()->Av());
349 
350  // Volumetric mass source in due to the wall boiling and bulk nucleation
351  dmdtf_ =
352  relax_*(1.0/6.0)*A2E*dDep_*vrho*fDep_*Av
353  + (1 - relax_)*dmdtf_;
354 
355  const dimensionedScalar zeroT(dimTemperature, small);
356  const dimensionedScalar smallT(dimTemperature, small);
357  const dimensionedScalar smallF(dimless/dimTime, small);
358 
359  // Quenching heat transfer coefficient
360  const volScalarField hQ
361  (
362  2*lkappa*fDep_
363  *sqrt((0.8/max(fDep_, smallF))/(pi*(lkappa/lCp)/lrho))
364  );
365 
366  // Volumetric quenching heat flux
367  qq_ =
368  (1 - relax_)*qq_
369  + relax_*(A2*hQ*max(Tsurface_ - lThermo.T(), zeroT))*Av;
370 
371  // Heat transfer coefficient
372  K_ =
373  heatTransferModel_->K(residualAlpha)
374  + (dmdtf_*L + qq_)/max(Tsurface_ - liquid.thermo().T(), smallT);
375 
376  return K_;
377 }
378 
379 
382 {
383  const phaseModel& liquid = interface_.continuous();
384  const phaseSystem& fluid = liquid.fluid();
385  const phaseModel& vapour(fluid.phases()[vapourPhaseName_]);
386 
387  if (phaseInterface == phaseInterfaceKey(vapour, liquid))
388  {
389  return true;
390  }
391  else
392  {
393  return false;
394  }
395 }
396 
397 
399 activePhaseInterface() const
400 {
401  const phaseModel& liquid = interface_.continuous();
402  const phaseSystem& fluid = liquid.fluid();
403  const phaseModel& vapour(fluid.phases()[vapourPhaseName_]);
404 
405  return phaseInterfaceKey(vapour, liquid);
406 }
407 
408 
411 {
412  return dmdtf_;
413 }
414 
415 
417 writeData(Ostream& os) const
418 {
419  return os.good();
420 }
421 
422 
423 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:330
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
virtual const volScalarField & T() const =0
Temperature [K].
virtual const volScalarField & Cp() const =0
Heat capacity at constant pressure [J/kg/K].
virtual tmp< volScalarField > ha() const =0
Absolute enthalpy [J/kg].
virtual const volScalarField & kappa() const =0
Thermal conductivity of mixture [W/m/K].
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Class to represent a interface between phases where one phase is considered dispersed within the othe...
virtual volScalarField & p()=0
Pressure [Pa].
Model for heat transfer between phases.
tmp< volScalarField > K() const
The heat transfer function K used in the enthalpy equation.
bool writeData(Ostream &os) const
Dummy write for regIOobject.
const volScalarField & dmdtf() const
Return the rate of phase-change.
phaseInterfaceKey activePhaseInterface() const
Key for the phase change phaseInterface.
wallBoilingHeatTransfer(const dictionary &dict, const phaseInterface &interface, const bool registerObject)
Construct from components.
Wrapper around saturationTemperatureModel to facilitate convenient construction on interfaces.
Generic thermophysical properties class for a liquid in which the functions and coefficients for each...
Definition: liquid.H:53
Word-pair based class used for keying interface models in hash tables.
Class to represent an interface between phases. Derivations can further specify the configuration of ...
virtual const rhoThermo & thermo() const =0
Return the thermophysical model.
const autoPtr< diameterModel > & dPtr() const
Return const-reference to diameterModel of the phase.
Definition: phaseModel.C:151
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:73
const phaseModelList & phases() const
Return the phase models.
Definition: phaseSystemI.H:41
const ModelType & lookupInterfacialModel(const phaseInterface &interface) const
Return a sub model for an interface.
Base-class for fluid thermodynamic properties based on density.
Definition: rhoThermo.H:55
virtual tmp< volScalarField > rho() const =0
Density [kg/m^3].
Class to represent a certain side of an interface between phases.
A class for managing temporary objects.
Definition: tmp.H:55
static autoPtr< departureDiameterModel > New(const dictionary &dict)
Select null constructed.
static autoPtr< departureFrequencyModel > New(const dictionary &dict)
Select null constructed.
static autoPtr< nucleationSiteModel > New(const dictionary &dict)
Select null constructed.
static autoPtr< partitioningModel > New(const dictionary &dict)
Select null constructed.
K
Definition: pEqn.H:75
word timeName
Definition: getTimeIndex.H:3
mathematical constants.
Namespace for OpenFOAM.
const doubleScalar e
Definition: doubleScalar.H:105
dimensionedScalar exp(const dimensionedScalar &ds)
const dimensionSet dimEnergy
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimless
const dimensionSet dimLength
const dimensionSet dimTemperature
dimensionedScalar log(const dimensionedScalar &ds)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
dimensionedScalar sqrt(const dimensionedScalar &ds)
const dimensionSet dimDensity
const dimensionSet dimVolume
defineTypeNameAndDebug(combustionModel, 0)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
word typedName(Name name)
Return the name of the object within the given type.
Definition: typeInfo.H:149
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
const dimensionSet dimArea
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
dictionary dict