phaseSurfaceBoiling.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) 2025 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 
26 #include "phaseSurfaceBoiling.H"
27 
28 #include "multiphaseEuler.H"
29 
30 #include "diameterModel.H"
32 
34 #include "partitioningModel.H"
35 #include "nucleationSiteModel.H"
36 #include "departureDiameterModel.H"
38 
40 
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 namespace fv
46 {
49 }
50 }
51 
52 
53 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
54 
55 void Foam::fv::phaseSurfaceBoiling::readCoeffs(const dictionary& dict)
56 {
57  saturationModelPtr_.reset
58  (
60  (
61  "saturationTemperature",
62  dict
63  ).ptr()
64  );
65 
66  partitioningModel_.reset
67  (
69  (
70  dict.subDict("partitioningModel")
71  ).ptr()
72  );
73 
74  nucleationSiteModel_.reset
75  (
77  (
78  dict.subDict("nucleationSiteModel")
79  ).ptr()
80  );
81 
82  departureDiameterModel_.reset
83  (
85  (
86  dict.subDict("departureDiameterModel")
87  ).ptr()
88  );
89 
90  departureFrequencyModel_.reset
91  (
93  (
94  dict.subDict("departureFrequencyModel")
95  ).ptr()
96  );
97 }
98 
99 
100 void Foam::fv::phaseSurfaceBoiling::correctMDot() const
101 {
103 
104  static const dimensionedScalar rooVSmallT(dimTemperature, rootVSmall);
105  static const dimensionedScalar rootVSmallF(dimless/dimTime, rootVSmall);
106  static const dimensionedScalar rootVSmallH
107  (
109  rootVSmall
110  );
111 
112  const rhoThermo& solidThermo = solid_.thermo();
113  const volScalarField::Internal& solidT = solidThermo.T();
114  const rhoThermo& liquidThermo = liquid_.thermo();
115  const volScalarField::Internal& liquidT = liquidThermo.T();
116  const rhoThermo& vapourThermo = vapour_.thermo();
117 
118  // Estimate the surface temperature from the surrounding temperature and
119  // heat transfer coefficients. Note that lagged values of qEvaporative and
120  // qQuenching are used in this calculation.
121  const Pair<tmp<volScalarField>> Hs =
122  solver_.heatTransfer.Hs(liquid_, solid_, scalar(0));
123  const volScalarField::Internal& liquidH = Hs.first();
124  const volScalarField::Internal& solidH = Hs.second();
125  const volScalarField::Internal Tsurface
126  (
127  (solidH*solidT + liquidH*liquidT + qEvaporative_ + qQuenching_)
128  /max(solidH + liquidH, rootVSmallH)
129  );
130 
131  const volScalarField::Internal Tsat
132  (
133  saturationModelPtr_->Tsat(liquid_.fluidThermo().p()())
134  );
135 
137  (
138  vapourThermo.ha(liquid_.fluidThermo().p()(), Tsat)
139  - liquidThermo.ha()()
140  );
141 
142  // Wetted fraction
143  wetFraction_ =
144  partitioningModel_->wetFraction(liquid_()/max(1 - solid_(), small));
145 
146  // Bubble departure diameter
147  dDeparture_ =
148  departureDiameterModel_->dDeparture
149  (
150  liquid_,
151  vapour_,
152  solid_,
153  Tsurface,
154  Tsat,
155  L
156  );
157 
158  // Bubble departure frequency
159  fDeparture_ =
160  departureFrequencyModel_->fDeparture
161  (
162  liquid_,
163  vapour_,
164  solid_,
165  Tsurface,
166  Tsat,
167  L,
168  dDeparture_
169  );
170 
171  // Nucleation site density
172  nucleationSiteDensity_ =
173  nucleationSiteModel_->nucleationSiteDensity
174  (
175  liquid_,
176  vapour_,
177  solid_,
178  Tsurface,
179  Tsat,
180  L,
181  dDeparture_,
182  fDeparture_
183  );
184 
185  const tmp<volScalarField> tliquidRho(liquidThermo.rho());
186  const volScalarField::Internal& liquidRho = tliquidRho();
187  const tmp<volScalarField> tvapourRho(vapourThermo.rho());
188  const volScalarField::Internal& vapourRho = tvapourRho();
189  const volScalarField::Internal& liquidCp = liquidThermo.Cp();
190  const volScalarField::Internal& liquidkappa = liquidThermo.kappa();
191 
192  // Area fractions: Del Valle & Kenning (1985)
193  const volScalarField::Internal Ja
194  (
195  liquidRho*liquidCp*(Tsat - min(Tsurface, Tsat))/(liquidRho*L)
196  );
197  const volScalarField::Internal Al
198  (
199  wetFraction_*4.8*exp(min(-Ja/80, log(vGreat)))
200  );
201  const volScalarField::Internal A2
202  (
203  min(pi*sqr(dDeparture_)*nucleationSiteDensity_*Al/4, scalar(1))
204  );
205  const volScalarField::Internal A2E
206  (
207  min(pi*sqr(dDeparture_)*nucleationSiteDensity_*Al/4, scalar(5))
208  );
209 
210  // Surface area per unit volume
211  const volScalarField Av(solid_.diameter().Av());
212 
213  // Relaxation factor
214  const scalar f = mesh().solution().fieldRelaxationFactor(mDot_.member());
215 
216  // Volumetric mass source in due to the wall boiling and bulk nucleation
217  mDot_ =
218  (1 - f)*mDot_
219  + f*(1.0/6.0)*A2E*dDeparture_*vapourRho*fDeparture_*Av;
220 
221  // Evaporative heat flux
222  qEvaporative_ = mDot_*L;
223 
224  // Quenching heat transfer coefficient
225  const volScalarField::Internal hQuenching
226  (
227  2
228  *liquidkappa
229  *fDeparture_
230  *sqrt
231  (
232  0.8
233  /max(fDeparture_, rootVSmallF)
234  /(pi*(liquidkappa/liquidCp)/liquidRho)
235  )
236  );
237 
238  // Quenching heat flux
239  qQuenching_ =
240  (1 - f)*qQuenching_
241  + f*A2*hQuenching*max(Tsurface - liquidThermo.T(), rooVSmallT)*Av;
242 }
243 
244 
245 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
246 
248 (
249  const word& name,
250  const word& modelType,
251  const fvMesh& mesh,
252  const dictionary& dict
253 )
254 :
255  phaseChange(name, modelType, mesh, dict, wordList()),
256  nucleation(),
257  solver_(mesh().lookupObject<solvers::multiphaseEuler>(solver::typeName)),
258  fluid_(solver_.fluid),
259  liquid_(fluid_.phases()[phaseNames().first()]),
260  vapour_(fluid_.phases()[phaseNames().second()]),
261  solid_(fluid_.phases()[dict.lookup("phase")]),
262  saturationModelPtr_(nullptr),
263  partitioningModel_(nullptr),
264  nucleationSiteModel_(nullptr),
265  departureDiameterModel_(nullptr),
266  departureFrequencyModel_(nullptr),
267  pressureEquationIndex_(-1),
268  wetFraction_
269  (
270  IOobject
271  (
272  name + ":wetFraction",
273  mesh.time().name(),
274  mesh,
275  IOobject::READ_IF_PRESENT,
276  IOobject::AUTO_WRITE
277  ),
278  mesh,
279  dimensionedScalar(dimless, scalar(0))
280  ),
281  dDeparture_
282  (
283  IOobject
284  (
285  name + ":dDeparture",
286  mesh.time().name(),
287  mesh,
288  IOobject::READ_IF_PRESENT,
289  IOobject::AUTO_WRITE
290  ),
291  mesh,
293  ),
294  fDeparture_
295  (
296  IOobject
297  (
298  name + ":fDeparture",
299  mesh.time().name(),
300  mesh,
301  IOobject::READ_IF_PRESENT,
302  IOobject::AUTO_WRITE
303  ),
304  mesh,
305  dimensionedScalar(inv(dimTime), scalar(0))
306  ),
307  nucleationSiteDensity_
308  (
309  IOobject
310  (
311  name + ":nucleationSiteDensity",
312  mesh.time().name(),
313  mesh,
314  IOobject::READ_IF_PRESENT,
315  IOobject::AUTO_WRITE
316  ),
317  mesh,
318  dimensionedScalar(inv(dimArea), scalar(0))
319  ),
320  qQuenching_
321  (
322  IOobject
323  (
324  name + ":qQuenching",
325  mesh.time().name(),
326  mesh,
327  IOobject::READ_IF_PRESENT,
328  IOobject::AUTO_WRITE
329  ),
330  mesh,
332  ),
333  qEvaporative_
334  (
335  IOobject
336  (
337  name + ":qEvaporative",
338  mesh.time().name(),
339  mesh,
340  IOobject::NO_READ,
341  IOobject::NO_WRITE
342  ),
343  mesh,
345  ),
346  mDot_
347  (
348  IOobject
349  (
350  name + ":mDot",
351  mesh.time().name(),
352  mesh,
353  IOobject::READ_IF_PRESENT,
354  IOobject::AUTO_WRITE
355  ),
356  mesh,
358  )
359 {
360  readCoeffs(coeffs(dict));
361 }
362 
363 
364 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
365 
367 {
368  return
369  phaseChange::addsSupToField(fieldName)
370  || fieldName == solid_.thermo().he().name();
371 }
372 
373 
376 {
377  // Put all the latent heat into the liquid
378  return
380  (
381  name() + ":Lfraction",
382  mesh(),
383  dimensionedScalar(dimless, scalar(0))
384  );
385 }
386 
387 
390 {
391  return dDeparture_;
392 }
393 
394 
397 {
398  return fDeparture_;
399 }
400 
401 
404 {
406  return tmp<volScalarField::Internal>(nullptr);
407 }
408 
409 
412 {
413  return mDot_;
414 }
415 
416 
418 (
419  const volScalarField& alpha,
420  const volScalarField& rho,
421  fvMatrix<scalar>& eqn
422 ) const
423 {
424  // Pressure equation (i.e., continuity, linearised in the pressure)
425  if
426  (
427  (&alpha == &liquid_ || &alpha == &vapour_)
428  && (&rho == &liquid_.rho() || &rho == &vapour_.rho())
429  && &eqn.psi() == &solver_.p_rgh
430  )
431  {
432  // Ensure that the source is up-to date if this is the first call in
433  // the current phase loop
434  if (pressureEquationIndex_ % 2 == 0) correctMDot();
435  pressureEquationIndex_ ++;
436  }
437 
438  // Let the base class add the actual source
440 }
441 
442 
444 (
445  const volScalarField& alpha,
446  const volScalarField& rho,
447  const volScalarField& he,
448  fvMatrix<scalar>& eqn
449 ) const
450 {
451  // Liquid or solid energy equation. Apply the additional explicit latent
452  // and quenching heat transfers.
453  if (&he == &liquid_.thermo().he() || &he == &solid_.thermo().he())
454  {
455  const scalar sign = &he == &solid_.thermo().he() ? -1 : +1;
456 
457  eqn += sign*(qEvaporative_ + qQuenching_);
458  }
459 
460  // Let the base class do the other liquid-vapour transfers
461  if (&he == &liquid_.thermo().he())
462  {
464  }
465 }
466 
467 
469 {
470  // Reset the p_rgh equation solution counter
471  pressureEquationIndex_ = 0;
472 
473  // Correct the total phase change rate
474  correctMDot();
475 }
476 
477 
479 {
480  if (phaseChange::read(dict))
481  {
482  readCoeffs(coeffs(dict));
483  return true;
484  }
485  else
486  {
487  return false;
488  }
489 }
490 
491 
492 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
static tmp< DimensionedField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Mesh &mesh, const dimensionSet &, const PrimitiveField< Type > &)
Return a temporary field constructed from name, mesh,.
Generic GeometricField class.
DimensionedField< Type, GeoMesh, PrimitiveField > Internal
Type of the internal field from which this GeometricField is derived.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
void reset(T *=nullptr)
If object pointer already set, delete object and set to given.
Definition: autoPtrI.H:114
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
VolField< Type > & psi()
Definition: fvMatrix.H:289
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
const fvSolution & solution() const
Return the fvSolution.
Definition: fvMesh.C:1810
Finite volume model abstract base class.
Definition: fvModel.H:60
static const dictionary & coeffs(const word &modelType, const dictionary &)
Return the coefficients sub-dictionary for a given model type.
Definition: fvModelI.H:31
virtual void addSup(fvMatrix< scalar > &eqn) const
Add a source term to a field-less proxy equation.
Definition: massTransfer.C:225
virtual bool addsSupToField(const word &fieldName) const
Return true if the fvModel adds a source term to the given.
Definition: massTransfer.C:207
Mix-in interface for nucleation models. Provides access to properties of the nucleation process,...
Definition: nucleation.H:53
Base class for phase change models.
Definition: phaseChange.H:61
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: phaseChange.C:624
void addSup(const volScalarField &alpha, const volScalarField &rho, const volScalarField &heOrYi, fvMatrix< scalar > &eqn) const
Override the energy equation to add the phase change heat, or.
Definition: phaseChange.C:551
Model for nucleate wall boiling on the surface of a third (solid) phase.
virtual tmp< DimensionedField< scalar, volMesh > > Lfraction() const
Return the fraction of the latent heat that is transferred into.
virtual void correct()
Correct the fvModel.
virtual tmp< DimensionedField< scalar, volMesh > > mDot() const
Return the mass transfer rate.
virtual bool read(const dictionary &dict)
Read source dictionary.
virtual tmp< DimensionedField< scalar, volMesh > > d() const
Return the diameter of nuclei.
virtual tmp< DimensionedField< scalar, volMesh > > tau() const
Return the nucleation time scale.
void addSup(const volScalarField &alpha, const volScalarField &rho, fvMatrix< scalar > &eqn) const
Override the pressure equation to add the mass transfer rate.
phaseSurfaceBoiling(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
virtual bool addsSupToField(const word &fieldName) const
Return true if the fvModel adds a source term to the given.
virtual tmp< DimensionedField< scalar, volMesh > > nDot() const
Return the number rate at which nuclei are generated.
static const dimensionSet dimK
Coefficient dimensions.
static autoPtr< saturationTemperatureModel > New(const dictionary &dict)
Select with dictionary.
scalar fieldRelaxationFactor(const word &name) const
Return the relaxation factor for the given field.
Definition: solution.C:271
Abstract base class for run-time selectable region solvers.
Definition: solver.H:56
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.
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:381
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
Namespace for OpenFOAM.
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedScalar sign(const dimensionedScalar &ds)
const dimensionSet dimEnergy
const dimensionSet dimless
const dimensionSet dimLength
labelList second(const UList< labelPair > &p)
Definition: patchToPatch.C:49
const dimensionSet dimTemperature
dimensionedScalar log(const dimensionedScalar &ds)
labelList first(const UList< labelPair > &p)
Definition: patchToPatch.C:39
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
const dimensionSet dimDensity
const dimensionSet dimVolume
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:62
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
const dimensionSet dimArea
void inv(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
thermo he()
labelList f(nPoints)
labelList fv(nPoints)
dictionary dict