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-2024 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  heatTransferModel_
73  (
75  (
76  dict.subDict("heatTransferModel"),
77  interface,
78  false,
79  false
80  )
81  ),
82  relax_(dict.lookupOrDefault<scalar>("relax", scalar(1))),
83  partitioningModel_(nullptr),
84  nucleationSiteModel_(nullptr),
85  departureDiamModel_(nullptr),
86  departureFreqModel_(nullptr),
87  wetFraction_
88  (
89  IOobject
90  (
91  IOobject::groupName("fWallBoiling", interface_.name()),
92  interface_.mesh().time().name(),
93  interface_.mesh(),
94  IOobject::READ_IF_PRESENT,
95  IOobject::AUTO_WRITE
96  ),
97  interface_.mesh(),
99  ),
100  dDep_
101  (
102  IOobject
103  (
104  IOobject::groupName("departureDiameter", interface_.name()),
105  interface_.mesh().time().name(),
106  interface_.mesh(),
107  IOobject::READ_IF_PRESENT,
108  IOobject::AUTO_WRITE
109  ),
110  interface_.mesh(),
112  ),
113  fDep_
114  (
115  IOobject
116  (
117  IOobject::groupName("departureFrequency", interface_.name()),
118  interface_.mesh().time().name(),
119  interface_.mesh(),
120  IOobject::READ_IF_PRESENT,
121  IOobject::AUTO_WRITE
122  ),
123  interface_.mesh(),
125  ),
126  nucleationSiteDensity_
127  (
128  IOobject
129  (
130  IOobject::groupName("nucleationSites", interface_.name()),
131  interface_.mesh().time().name(),
132  interface_.mesh(),
133  IOobject::READ_IF_PRESENT,
134  IOobject::AUTO_WRITE
135  ),
136  interface_.mesh(),
138  ),
139  dmdtf_
140  (
141  IOobject
142  (
143  IOobject::groupName(typedName("dmdtf"), interface_.name()),
144  interface_.mesh().time().name(),
145  interface_.mesh(),
146  IOobject::READ_IF_PRESENT,
147  IOobject::AUTO_WRITE
148  ),
149  interface_.mesh(),
151  ),
152  qq_
153  (
154  IOobject
155  (
156  IOobject::groupName("qq", interface_.name()),
157  interface_.mesh().time().name(),
158  interface_.mesh(),
159  IOobject::READ_IF_PRESENT,
160  IOobject::AUTO_WRITE
161  ),
162  interface_.mesh(),
164  ),
165  Tsurface_
166  (
167  IOobject
168  (
169  IOobject::groupName("Tsurface", interface_.name()),
170  interface_.mesh().time().name(),
171  interface_.mesh(),
172  IOobject::READ_IF_PRESENT,
173  IOobject::AUTO_WRITE
174  ),
175  interface_.mesh(),
177  ),
178  K_
179  (
180  IOobject
181  (
182  IOobject::groupName(typedName("K"), interface_.name()),
183  interface_.mesh().time().name(),
184  interface_.mesh(),
185  IOobject::READ_IF_PRESENT,
186  IOobject::AUTO_WRITE
187  ),
188  heatTransferModel_->K()
189  )
190 {
191  partitioningModel_ =
193  (
194  dict.subDict("partitioningModel")
195  );
196 
197  nucleationSiteModel_ =
199  (
200  dict.subDict("nucleationSiteModel")
201  );
202 
203  departureDiamModel_ =
205  (
206  dict.subDict("departureDiameterModel")
207  );
208 
209  departureFreqModel_ =
211  (
212  dict.subDict("departureFrequencyModel")
213  );
214 }
215 
216 
217 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
218 
221 {}
222 
223 
224 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
225 
228 (
229  const scalar residualAlpha
230 ) const
231 {
232  const phaseSystem& fluid = interface_.fluid();
233 
234  const phaseModel& liquid = interface_.continuous();
235  const phaseModel& vapour = fluid.phases()[vapourPhaseName_];
236  const phaseModel& solid = interface_.dispersed();
237 
238  const rhoThermo& lThermo = liquid.thermo();
239  const rhoThermo& vThermo = vapour.thermo();
240  const rhoThermo& sThermo = solid.thermo();
241 
242  // Estimate the surface temperature from the surrounding temperature and
243  // heat transfer coefficients. Note that a lagged value of K is used in
244  // this calculation.
245  {
246  static const dimensionedScalar smallK(dimK, rootVSmall);
247 
248  const heatTransferModel& otherModel =
249  fluid.lookupInterfacialModel<heatTransferModel>(otherInterface_());
250  const volScalarField otherK(otherModel.K());
251 
252  Tsurface_ =
253  (K_*lThermo.T() + otherK*sThermo.T())
254  /max(K_ + otherK, smallK);
255  }
256 
257  const volScalarField Tsat
258  (
260  <
262  >
263  (phaseInterface(liquid, vapour))
264  .Tsat(liquid.fluidThermo().p())()
265  );
266 
267  const volScalarField L
268  (
269  vThermo.ha(liquid.fluidThermo().p(), Tsat) - lThermo.ha()
270  );
271 
272  // Wetted fraction
273  wetFraction_ =
274  partitioningModel_->wetFraction(liquid/max(1 - solid, small));
275 
276  // Bubble departure diameter
277  dDep_ =
278  departureDiamModel_->dDeparture
279  (
280  liquid,
281  vapour,
282  solid,
283  Tsurface_,
284  Tsat,
285  L
286  );
287 
288  // Bubble departure frequency
289  fDep_ =
290  departureFreqModel_->fDeparture
291  (
292  liquid,
293  vapour,
294  solid,
295  Tsurface_,
296  Tsat,
297  L,
298  dDep_
299  );
300 
301  // Nucleation site density
302  nucleationSiteDensity_ =
303  nucleationSiteModel_->nucleationSiteDensity
304  (
305  liquid,
306  vapour,
307  solid,
308  Tsurface_,
309  Tsat,
310  L,
311  dDep_,
312  fDep_
313  );
314 
315  const tmp<volScalarField> tlrho(lThermo.rho());
316  const volScalarField& lrho = tlrho();
317  const tmp<volScalarField> tvrho(vThermo.rho());
318  const volScalarField& vrho = tvrho();
319  const volScalarField& lCp = lThermo.Cp();
320  const volScalarField& lkappa = lThermo.kappa();
321  //const volScalarField lalpha(lkappa/lCp);
322 
323  // Area fractions
324 
325  // Del Valle & Kenning (1985)
326  const volScalarField Ja
327  (
328  lrho*lCp*(Tsat - min(Tsurface_, Tsat))/(lrho*L)
329  );
330 
331  const volScalarField Al
332  (
333  wetFraction_*4.8*exp(min(-Ja/80, log(vGreat)))
334  );
335 
336  volScalarField A2
337  (
338  min(pi*sqr(dDep_)*nucleationSiteDensity_*Al/4, scalar(1))
339  );
340  const scalarField A1(max(1 - A2, scalar(1e-4)));
341  volScalarField A2E
342  (
343  min(pi*sqr(dDep_)*nucleationSiteDensity_*Al/4, scalar(5))
344  );
345 
346  const volScalarField Av(solid.diameter().Av());
347 
348  // Volumetric mass source in due to the wall boiling and bulk nucleation
349  dmdtf_ =
350  relax_*(1.0/6.0)*A2E*dDep_*vrho*fDep_*Av
351  + (1 - relax_)*dmdtf_;
352 
353  const dimensionedScalar zeroT(dimTemperature, small);
354  const dimensionedScalar smallT(dimTemperature, small);
355  const dimensionedScalar smallF(dimless/dimTime, small);
356 
357  // Quenching heat transfer coefficient
358  const volScalarField hQ
359  (
360  2*lkappa*fDep_
361  *sqrt((0.8/max(fDep_, smallF))/(pi*(lkappa/lCp)/lrho))
362  );
363 
364  // Volumetric quenching heat flux
365  qq_ =
366  (1 - relax_)*qq_
367  + relax_*(A2*hQ*max(Tsurface_ - lThermo.T(), zeroT))*Av;
368 
369  // Heat transfer coefficient
370  K_ =
371  heatTransferModel_->K(residualAlpha)
372  + (dmdtf_*L + qq_)/max(Tsurface_ - liquid.thermo().T(), smallT);
373 
374  return K_;
375 }
376 
377 
380 {
381  const phaseModel& liquid = interface_.continuous();
382  const phaseSystem& fluid = liquid.fluid();
383  const phaseModel& vapour(fluid.phases()[vapourPhaseName_]);
384 
385  if (phaseInterface == phaseInterfaceKey(vapour, liquid))
386  {
387  return true;
388  }
389  else
390  {
391  return false;
392  }
393 }
394 
395 
397 flipSign() const
398 {
399  const phaseModel& liquid = interface_.continuous();
400  const phaseSystem& fluid = liquid.fluid();
401  const phaseModel& vapour(fluid.phases()[vapourPhaseName_]);
402 
403  return vapour.name() != phaseInterfaceKey(vapour, liquid).first();
404 }
405 
406 
409 {
410  return dmdtf_;
411 }
412 
413 
415 writeData(Ostream& os) const
416 {
417  return os.good();
418 }
419 
420 
421 // ************************************************************************* //
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
const Type & first() const
Return first.
Definition: Pair.H:98
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].
virtual tmp< volScalarField > Av() const =0
Return the surface area per unit volume.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Class to represent a interface between phases where one phase is considered dispersed within the othe...
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.
bool activePhaseInterface(const phaseInterfaceKey &) const
Is there phase change mass transfer for this phaseInterface.
bool flipSign() const
True if the sign of dmdtf should be changed.
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 diameterModel & diameter() const
Return a reference to the diameterModel of the phase.
Definition: phaseModel.C:187
const word & name() const
Return the name of this phase.
Definition: phaseModel.C:145
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:102
const ModelType & lookupInterfacialModel(const phaseInterface &interface) const
Return a sub model for an interface.
Base-class for thermodynamic properties based on density.
Definition: rhoThermo.H:54
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
mathematical constants.
Namespace for OpenFOAM.
const doubleScalar e
Definition: doubleScalar.H:106
dimensionedScalar exp(const dimensionedScalar &ds)
addToRunTimeSelectionTable(polyPatch, mergedCyclicPolyPatch, word)
const dimensionSet dimEnergy
dimensionedSymmTensor sqr(const dimensionedVector &dv)
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
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:176
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
const dimensionSet dimArea
dictionary dict