phaseSurfaceBoiling.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) 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::phaseSurfaceBoiling
26 
27 Description
28  Model for nucleate wall boiling on the surface of a third (solid) phase.
29 
30  This model functions very similarly to the wall boiling model (see that
31  model for references). The same sub-models are used, with exactly the same
32  specification syntax. The only difference is that the third phase must be
33  additionally specified, and that a two-resistance heat transfer model must
34  be in operation between the liquid phase and the third phase.
35 
36 Usage
37  Example usage:
38  \verbatim
39  phaseSurfaceBoiling
40  {
41  type phaseSurfaceBoiling;
42  libs ("libmultiphaseEulerFvModels.so");
43 
44  phase solid;
45 
46  // Note: Order is important. This model is one-way. It turns liquid
47  // into vapour. The phases should be specified in this order.
48  phases (water steam);
49 
50  energySemiImplicit no;
51 
52  saturationTemperature
53  {
54  type constant;
55  value 372.76;
56  }
57 
58  partitioningModel
59  {
60  type Lavieville;
61  alphaCrit 0.2;
62  }
63 
64  nucleationSiteModel
65  {
66  type LemmertChawla;
67  }
68 
69  departureDiameterModel
70  {
71  type TolubinskiKostanchuk;
72  }
73 
74  departureFrequencyModel
75  {
76  type KocamustafaogullariIshii;
77  Cf 1.18;
78  }
79  }
80  \endverbatim
81 
82 See also
83  Foam::fv::wallBoiling
84 
85 SourceFiles
86  wallBoiling.C
87 
88 \*---------------------------------------------------------------------------*/
89 
90 #ifndef phaseSurfaceBoiling_H
91 #define phaseSurfaceBoiling_H
92 
93 #include "phaseChange.H"
94 #include "nucleation.H"
95 #include "phaseSystem.H"
96 
97 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
98 
99 namespace Foam
100 {
101 
102 // Forward declaration of classes
103 namespace solvers
104 {
105  class multiphaseEuler;
106 }
107 
108 class saturationTemperatureModel;
109 
110 namespace wallBoilingModels
111 {
112  class partitioningModel;
113  class nucleationSiteModel;
116 }
117 
118 namespace fv
119 {
120 
121 /*---------------------------------------------------------------------------*\
122  Class phaseSurfaceBoiling Declaration
123 \*---------------------------------------------------------------------------*/
124 
126 :
127  public phaseChange,
128  public nucleation
129 {
130  // Private Data
131 
132  //- Solver
133  const solvers::multiphaseEuler& solver_;
134 
135  //- Reference to the phase system
136  const phaseSystem& fluid_;
137 
138  //- Reference to the liquid phase
139  const phaseModel& liquid_;
140 
141  //- Reference to the vapour phase
142  const phaseModel& vapour_;
143 
144  //- Reference to the solid phase on which the boiling occurs
145  const phaseModel& solid_;
146 
147  //- The saturation curve
148  autoPtr<saturationTemperatureModel> saturationModelPtr_;
149 
150  //- Run-time selected heat flux partitioning model
152  partitioningModel_;
153 
154  //- Run-time selected nucleation site density model
156  nucleationSiteModel_;
157 
158  //- Run-time selected bubble departure diameter model
160  departureDiameterModel_;
161 
162  //- Run-time selected bubble departure frequency model
164  departureFrequencyModel_;
165 
166  //- Counter for the evaluations of the pressure equation sources
167  mutable label pressureEquationIndex_;
168 
169  //- Surface liquid fraction
170  mutable volScalarField::Internal wetFraction_;
171 
172  //- Departure diameter
173  mutable volScalarField::Internal dDeparture_;
174 
175  //- Departure frequency
176  mutable volScalarField::Internal fDeparture_;
177 
178  //- Nucleation site density
179  mutable volScalarField::Internal nucleationSiteDensity_;
180 
181  //- Quenching heat transfer rate
182  mutable volScalarField::Internal qQuenching_;
183 
184  //- Evaporative heat transfer rate
185  mutable volScalarField::Internal qEvaporative_;
186 
187  //- Mass transfer rate
188  mutable volScalarField::Internal mDot_;
189 
190 
191  // Private Member Functions
192 
193  //- Non-virtual read
194  void readCoeffs(const dictionary& dict);
195 
196  //- Correct the phase change rate
197  void correctMDot() const;
198 
199 
200 public:
201 
202  //- Runtime type information
203  TypeName("phaseSurfaceBoiling");
204 
205 
206  // Constructors
207 
208  //- Construct from explicit source name and mesh
209  //- Construct from components
211  (
212  const word& name,
213  const word& modelType,
214  const fvMesh& mesh,
215  const dictionary& dict
216  );
217 
218 
219  // Member Functions
220 
221  // Checks
222 
223  //- Return true if the fvModel adds a source term to the given
224  // field's transport equation
225  virtual bool addsSupToField(const word& fieldName) const;
226 
227 
228  // Access
229 
230  //- Access the liquid phase
231  inline const phaseModel& liquid() const;
232 
233  //- Access the vapour phase
234  inline const phaseModel& vapour() const;
235 
236  //- Access the solid phase on which the boiling occurs
237  inline const phaseModel& solid() const;
238 
239 
240  // Evaluation
241 
242  //- Return the fraction of the latent heat that is transferred into
243  // the second phase
245 
246  //- Return the diameter of nuclei
247  virtual tmp<DimensionedField<scalar, volMesh>> d() const;
248 
249  //- Return the number rate at which nuclei are generated
251 
252  //- Return the nucleation time scale
254 
255 
256  // Sources
257 
258  //- Return the mass transfer rate
260 
261  //- Override the pressure equation to add the mass transfer rate
262  // linearised in the pressure
263  void addSup
264  (
265  const volScalarField& alpha,
266  const volScalarField& rho,
267  fvMatrix<scalar>& eqn
268  ) const;
269 
270  //- Override the phase and liquid's energy equations to add
271  // additional latent heat and quenching heat transfer terms
272  void addSup
273  (
274  const volScalarField& alpha,
275  const volScalarField& rho,
276  const volScalarField& he,
277  fvMatrix<scalar>& eqn
278  ) const;
279 
280 
281  //- Correct the fvModel
282  virtual void correct();
283 
284 
285  // IO
286 
287  //- Read source dictionary
288  virtual bool read(const dictionary& dict);
289 };
290 
291 
292 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
293 
294 } // End namespace fv
295 } // End namespace Foam
296 
297 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
298 
299 #include "phaseSurfaceBoiling.H"
300 
301 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
302 
303 #endif
304 
305 // ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
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
static const word & fieldName()
Return the name of the field associated with a source term. Special.
Definition: fvModelI.H:49
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
Model for nucleate wall boiling on the surface of a third (solid) phase.
const phaseModel & vapour() const
Access the vapour 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.
TypeName("phaseSurfaceBoiling")
Runtime type information.
virtual tmp< DimensionedField< scalar, volMesh > > mDot() const
Return the mass transfer rate.
const phaseModel & liquid() const
Access the liquid phase.
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.
const phaseModel & solid() const
Access the solid phase on which the boiling occurs.
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.
Class to represent a system of phases.
Definition: phaseSystem.H:74
Solver module for a system of any number of compressible fluid phases with a common pressure,...
A class for managing temporary objects.
Definition: tmp.H:55
Base class for bubble departure diameter models.
Base class for bubble departure frequency models.
Base class for nucleation site density models.
Base class for wall heat flux partitioning models.
A class for handling words, derived from string.
Definition: word.H:62
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
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
thermo he()
labelList fv(nPoints)
dictionary dict