wallBoiling.H
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) 2024-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 Class
25  Foam::fv::wallBoiling
26 
27 Description
28  Model for nucleate wall boiling between two phases on the surface of a
29  number of wall patches.
30 
31  This model implements a version of the well-known RPI wall boiling model
32  (Kurul & Podowski, 1991). The model is based on an implementation described
33  in Peltola et al. (2019) and is similar to the model described by Peltola &
34  Pättikangas (2012).
35 
36  References:
37  \verbatim
38  Kurul, N., & Podowski, M.Z. (1991).
39  On the modeling of multidimensional effects in boiling channels.
40  ANS. Proc. National Heat Transfer Con. Minneapolis, Minnesota, USA,
41  1991.
42  ISBN: 0-89448-162-1, pp. 30-40.
43  \endverbatim
44 
45  \verbatim
46  Peltola, J., Pättikangas, T., Bainbridge, W., Lehnigk, R., Schlegel, F.
47  (2019).
48  On Development and validation of subcooled nucleate boiling models for
49  OpenFOAM Foundation Release.
50  NURETH-18 Conference Proceedings, Portland, Oregon, United States, 2019.
51  \endverbatim
52 
53  \verbatim
54  Peltola, J., & Pättikangas, T.J.H. (2012).
55  Development and validation of a boiling model for OpenFOAM multiphase
56  solver.
57  CFD4NRS-4 Conference Proceedings, Daejeon, Korea, 2012.
58  paper 59.
59  \endverbatim
60 
61 Usage
62  Example usage:
63  \verbatim
64  wallBoiling
65  {
66  type wallBoiling;
67  libs ("libmultiphaseEulerFvModels.so");
68 
69  // Note: Order is important. This model is one-way. It turns liquid
70  // into vapour. The phases should be specified in this order.
71  phases (water steam);
72 
73  energySemiImplicit no;
74 
75  saturationTemperature
76  {
77  type constant;
78  value 372.76;
79  }
80 
81  partitioningModel
82  {
83  type Lavieville;
84  alphaCrit 0.2;
85  }
86 
87  nucleationSiteModel
88  {
89  type LemmertChawla;
90  Cn 1;
91  NRef 30000000;
92  deltaTRef 10;
93  }
94 
95  departureDiameterModel
96  {
97  type TolubinskiKostanchuk;
98  dRef 0.00024;
99  dMax 0.0014;
100  dMin 1e-06;
101  }
102 
103  departureFrequencyModel
104  {
105  type KocamustafaogullariIshii;
106  Cf 1.18;
107  }
108  }
109  \endverbatim
110 
111  In addition to the above fvModel specification, wall patches on which
112  boiling is to be calculated should have an alphatBoilingWallFunction
113  boundary condition applied to the turbulent thermal diffusivity field.
114 
115 See also
116  Foam::alphatBoilingWallFunctionFvPatchScalarField
117 
118 SourceFiles
119  wallBoiling.C
120 
121 \*---------------------------------------------------------------------------*/
122 
123 #ifndef wallBoiling_H
124 #define wallBoiling_H
125 
126 #include "phaseChange.H"
127 #include "nucleation.H"
128 #include "phaseSystem.H"
129 
130 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
131 
132 namespace Foam
133 {
134 
135 // Forward declaration of classes
136 class saturationTemperatureModel;
137 
138 namespace wallBoilingModels
139 {
140  class partitioningModel;
141  class nucleationSiteModel;
142  class departureDiameterModel;
143  class departureFrequencyModel;
144 }
145 
146 class wallBoilingPhaseChangeRateFvPatchScalarField;
147 
148 namespace fv
149 {
150 
151 /*---------------------------------------------------------------------------*\
152  Class wallBoiling Declaration
153 \*---------------------------------------------------------------------------*/
154 
155 class wallBoiling
156 :
157  public phaseChange,
158  public nucleation
159 {
160 private:
161 
162  // Private classes
163 
164  //- Struct to hold properties that are held constant during the
165  // iterative solution procedure
166  struct laggedProperties;
167 
168 
169  // Private Data
170 
171  //- Reference to the phase system
172  const phaseSystem& fluid_;
173 
174  //- Reference to the liquid phase
175  const phaseModel& liquid_;
176 
177  //- Reference to the vapour phase
178  const phaseModel& vapour_;
179 
180  //- Reference to the liquid turbulent thermal diffusivity
181  const volScalarField& alphatLiquid_;
182 
183  //- Reference to the vapour turbulent thermal diffusivity
184  const volScalarField& alphatVapour_;
185 
186  //- Reference to the field associated with the pressure equation
187  const volScalarField& p_rgh_;
188 
189  //- Solution tolerance
190  scalar tolerance_;
191 
192  //- Estimate liquid temperature using logarithmic wall function?
193  bool liquidTemperatureWallFunction_;
194 
195  //- Turbulent Prandtl number
196  scalar Prt_;
197 
198  //- Bubble waiting time ratio
199  scalar bubbleWaitingTimeRatio_;
200 
201  //- The saturation curve
202  autoPtr<saturationTemperatureModel> saturationModelPtr_;
203 
204  //- Run-time selected heat flux partitioning model
206 
207  //- Run-time selected nucleation site density model
209 
210  //- Run-time selected bubble departure diameter model
212  departureDiameterModel_;
213 
214  //- Run-time selected bubble departure frequency model
216  departureFrequencyModel_;
217 
218  //- Counter for the evaluations of the pressure equation sources
219  mutable label pressureEquationIndex_;
220 
221  //- The phase change rate
222  mutable volScalarField mDot_;
223 
224 
225  // Private Member Functions
226 
227  //- Non-virtual read
228  void readCoeffs(const dictionary& dict);
229 
230  //- Return the boundary condition types for the phase change rate
231  wordList mDotBoundaryTypes() const;
232 
233  //- Get the liquid temperature patch field and the parameters
234  // associated with its boundary condition that are necessary for
235  // approximately evaluating the boundary condition's heat flux at
236  // a given wall temperature.
237  const fvPatchScalarField& getLiquidTemperaturePatchField
238  (
239  const laggedProperties& lagProps,
240  scalarField& isFixed,
241  scalarField& h,
242  scalarField& hTaPlusQa
243  ) const;
244 
245  //- Calculate the boiling for the given wall temperature. Return
246  // the total sum of all heat fluxes. Set the properties passed by
247  // non-const reference. Used by the functions below.
248  tmp<scalarField> calcBoiling
249  (
250  const laggedProperties& lagProps,
251  const scalarField& TLiquid,
252  const scalarField& wetFraction,
253  scalarField& dDeparture,
254  scalarField& fDeparture,
255  scalarField& nucleationSiteDensity,
256  scalarField& qQuenching,
257  scalarField& qEvaporative,
259  ) const;
260 
261  //- Calculate the boiling for the given wall temperature. Return
262  // the total sum of all heat fluxes. Use this to solve the balance
263  // between the heat fluxes specified by the boiling models and the
264  // temperature boundary condition without changing the stored
265  // boiling state.
266  tmp<scalarField> calcBoiling
267  (
269  const laggedProperties& lagProps,
270  const scalarField& TLiquid
271  ) const;
272 
273  //- Calculate the boiling for the given wall temperature. Return
274  // the total sum of all heat fluxes. Also set the stored boiling
275  // state. Use this after solving with the final wall temperature
276  // to set the boiling state.
277  tmp<scalarField> evaluateBoiling
278  (
280  const laggedProperties& lagProps,
281  const scalarField& TLiquid
282  ) const;
283 
284  //- Correct the phase change rate
285  void correctMDot() const;
286 
287 
288 public:
289 
290  //- Runtime type information
291  TypeName("wallBoiling");
292 
293 
294  // Constructors
295 
296  //- Construct from explicit source name and mesh
298  (
299  const word& name,
300  const word& modelType,
301  const fvMesh& mesh,
302  const dictionary& dict
303  );
304 
305 
306  // Member Functions
307 
308  // Access
309 
310  //- Access the liquid phase
311  inline const phaseModel& liquid() const;
312 
313  //- Access the vapour phase
314  inline const phaseModel& vapour() const;
315 
316  //- Access the liquid turbulent thermal diffusivity
317  inline const volScalarField& alphatLiquid() const;
318 
319  //- Access the vapour turbulent thermal diffusivity
320  inline const volScalarField& alphatVapour() const;
321 
322  //- Is the given patch boiling?
323  bool isBoiling(const label patchi) const;
324 
325 
326  // Evaluation
327 
328  //- Return the fraction of the latent heat that is transferred into
329  // the second phase
331 
332  //- Return the diameter of nuclei
333  virtual tmp<DimensionedField<scalar, volMesh>> d() const;
334 
335  //- Return the number rate at which nuclei are generated
337 
338  //- Return the nucleation time scale
340 
341 
342  // Sources
343 
344  //- Return the mass transfer rate
346 
347  //- Return the mass transfer rate for the given patch
349  (
350  const label patchi
351  ) const;
352 
353  //- Return the mass transfer rate for the given patch
355  (
356  const label patchi
357  ) const;
358 
359  //- Use phaseChange's source functions
360  using phaseChange::addSup;
361 
362  //- Override the pressure equation to add the mass transfer rate
363  // linearised in the pressure
364  void addSup
365  (
366  const volScalarField& alpha,
367  const volScalarField& rho,
368  fvMatrix<scalar>& eqn
369  ) const;
370 
371 
372  //- Correct the fvModel
373  virtual void correct();
374 
375 
376  // IO
377 
378  //- Read source dictionary
379  virtual bool read(const dictionary& dict);
380 };
381 
382 
383 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
384 
385 } // End namespace fv
386 } // End namespace Foam
387 
388 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
389 
390 #include "wallBoilingI.H"
391 
392 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
393 
394 #endif
395 
396 // ************************************************************************* //
Generic GeometricField class.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
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
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:69
const word & name() const
Return const access to the source name.
Definition: fvModelI.H:57
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:91
tmp< volScalarField::Internal > rho(const label i) const
Return the density.
Definition: massTransfer.C:92
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
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 between two phases on the surface of a number of wall patches.
Definition: wallBoiling.H:158
const phaseModel & vapour() const
Access the vapour phase.
Definition: wallBoilingI.H:36
TypeName("wallBoiling")
Runtime type information.
virtual tmp< DimensionedField< scalar, volMesh > > Lfraction() const
Return the fraction of the latent heat that is transferred into.
Definition: wallBoiling.C:810
virtual void correct()
Correct the fvModel.
Definition: wallBoiling.C:969
wallBoiling(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
Definition: wallBoiling.C:738
virtual tmp< DimensionedField< scalar, volMesh > > mDot() const
Return the mass transfer rate.
Definition: wallBoiling.C:901
bool isBoiling(const label patchi) const
Is the given patch boiling?
Definition: wallBoiling.C:799
const phaseModel & liquid() const
Access the liquid phase.
Definition: wallBoilingI.H:30
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: wallBoiling.C:979
const volScalarField & alphatVapour() const
Access the vapour turbulent thermal diffusivity.
Definition: wallBoilingI.H:48
const wallBoilingPhaseChangeRateFvPatchScalarField & mDotPf(const label patchi) const
Return the mass transfer rate for the given patch.
Definition: wallBoiling.C:908
void addSup(const volScalarField &alpha, const volScalarField &rho, const volScalarField &heOrYi, fvMatrix< scalar > &eqn) const
Use phaseChange's source functions.
Definition: phaseChange.C:551
virtual tmp< DimensionedField< scalar, volMesh > > d() const
Return the diameter of nuclei.
Definition: wallBoiling.C:824
virtual tmp< DimensionedField< scalar, volMesh > > tau() const
Return the nucleation time scale.
Definition: wallBoiling.C:893
const volScalarField & alphatLiquid() const
Access the liquid turbulent thermal diffusivity.
Definition: wallBoilingI.H:42
wallBoilingPhaseChangeRateFvPatchScalarField & mDotPfRef(const label patchi) const
Return the mass transfer rate for the given patch.
Definition: wallBoiling.C:926
virtual tmp< DimensionedField< scalar, volMesh > > nDot() const
Return the number rate at which nuclei are generated.
Definition: wallBoiling.C:861
Class to represent a system of phases.
Definition: phaseSystem.H:74
A class for managing temporary objects.
Definition: tmp.H:55
This boundary condition is used for the phase change rate field of the wall boiling fvModel....
A class for handling words, derived from string.
Definition: word.H:62
label patchi
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
const dimensionedScalar h
Planck constant.
Namespace for OpenFOAM.
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:59
labelList fv(nPoints)
dictionary dict