alphatWallBoilingWallFunctionFvPatchScalarField.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2015-2018 OpenFOAM Foundation
9  Copyright (C) 2018-2021 OpenCFD Ltd
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Class
28  Foam::compressible::alphatWallBoilingWallFunctionFvPatchScalarField
29 
30 Description
31  A thermal wall function for simulation of boiling wall.
32 
33  This alpha wall function can handle the following regimes:
34  - single phase
35  - subcooled nucleate wall boiling
36  - transitional boiling
37  - film boiling
38 
39  The wall function uses a partition method to transfer heat either
40  to the liquid or vapor phase. At the moment, this function works
41  in a wall temperature fixed mode, i.e. there is no consideration
42  for the sudden change of heat transfer coefficient (htc) after
43  reaching TDBN (deviation from nucleate boiling temperature).
44 
45  References:
46  \verbatim
47  Srinivasan, V., Moon, K. M., Greif, D.,
48  Wang, D. M., & Kim, M. H. (2010).
49  Numerical simulation of immersion quenching
50  process of an engine cylinder head.
51  Applied Mathematical Modelling, 34(8), 2111-2128.
52  DOI:10.1016/j.apm.2009.10.023
53  \endverbatim
54 
55 
56  For the single phase non-boiling regime the standard
57  \c JayatillekeWallFunction is used.
58 
59  For the sub-cool nucleate boiling regime the following runtime
60  selectable submodels are used:
61  - nucleation site density
62  - bubble departure frequency
63  - bubble departure diameter
64 
65  Implements a version of the well-known RPI wall boiling model
66  (Kurul & Podowski, 1991). The model implementation is similar to the model
67  described by Peltola & Pättikangas (2012) but has been extended with the
68  wall heat flux partitioning models.
69 
70  References:
71  \verbatim
72  Kurul, N., & Podowski, M. Z. (1991).
73  On the modeling of multidimensional effects in boiling channels.
74  Proceedings of the 27th National Heat Transfer Conference.
75  Minneapolis, Minn, USA, July 28-31, 1991.
76  ISBN: 0-89448-162-1, pp. 30-40
77 
78  Peltola, J., & Pättikangas, T. (2012).
79  Development and validation of a boiling model
80  for OpenFOAM multiphase solver.
81  Proceedings of the CFD4NRS-4. p. 59.
82  Daejeon, Democratic People's Republic of Korea, September 10-12, 2012.
83  \endverbatim
84 
85 
86  Alternatively a correlation can be used instead of the RPI wall boiling model.
87  If the keyword nucleatingModel a model is provided the BC uses it
88  instead of the RPI model.
89 
90  The transition boiling regime flux (TBF) is modelled following
91  a temperature based linear interpolation between the critical heat flux
92  (CHF) and the minimum heat flux (MHF) in such a way that when the wall
93  temperature is between the range of TDBN and the Leidenfrost temperature
94  (TLeiden) a linear interpolation is used between CHF and MHF.
95 
96  Thus, the following models are required:
97  - LeidenfrostModel
98  - CHFModel
99  - CHFSubCoolModel
100  - MHFModel
101  - TDNBModel
102  - filmBoilingModel
103 
104  The linear interpolation is as follows:
105 
106  \f[
107  TBF = CHF*\phi + (1 - \phi)*MHF
108  \f]
109 
110  with
111  \f[
112  \phi = w_p*(T_w - T_{DNB})/(T_{Leiden} - T_{DNB})
113  \f]
114 
115  where:
116  \vartable
117  w_p | Model constant
118  T_w | Wall temperature [K]
119  \endvartable
120 
121 
122  The film boiling regime is applied when \f$T_w\f$ is larger than
123  \f$T_{Leiden}\f$. In this regime the correlation from the
124  \c filmBoilingModel is used for calculating the cht from the wall.
125 
126  The \c filmBoilingModel is needed in the vapor field in order to calculate
127  the heat transfer to the vapor phase in film boiling regime.
128 
129 
130 Usage
131  Example of the boundary condition specification:
132  \verbatim
133  <patchName>
134  {
135  // Mandatory entries
136  type compressible::alphatWallBoilingWallFunction;
137  phaseType <word>;
138  otherPhase <word>;
139  relax <Function1<scalar>>;
140 
141  partitioningModel
142  {
143  type Lavieville;
144  alphaCrit 0.2;
145  }
146 
147  // Conditional entries
148 
149  // Option-1: phaseType=vapor
150 
151  // Optional entries
152  LeidenfrostModel
153  {
154  type Spiegler;
155  Tcrit 647;
156  }
157 
158  filmBoilingModel
159  {
160  type Bromley;
161  }
162 
163 
164  // Option-2: phaseType=liquid
165  nucleationSiteModel
166  {
167  type LemmertChawla;
168  }
169 
170  departureDiamModel
171  {
172  type TolubinskiKostanchuk;
173  }
174 
175  departureFreqModel
176  {
177  type Cole;
178  }
179 
180  // Optional entries
181  LeidenfrostModel
182  {
183  type Spiegler;
184  Tcrit 647;
185  }
186 
187  CHFModel
188  {
189  type Zuber;
190  }
191 
192  CHFSubCoolModel
193  {
194  type HuaXu;
195  Kburn 0.5;
196  }
197 
198  MHFModel
199  {
200  type Jeschar;
201  Kmhf 1;
202  }
203 
204  TDNBModel
205  {
206  type Schroeder;
207  }
208 
209  filmBoilingModel
210  {
211  type Bromley;
212  }
213 
214  dDep <scalarField>;
215  K <scalar>;
216  wp <scalar>;
217  qQuenching <scalarField>;
218 
219 
220  // Optional entries
221  alphatConv <scalarField>;
222 
223  //Inherited entries
224  ...
225  \endverbatim
226 
227  where the entries mean:
228  \table
229  Property | Description | Type | Reqd | Deflt
230  type | compressible::alphatWallBoilingWallFunction | word | yes | -
231  phaseType | Name of phase type | word | yes | -
232  otherPhase | Name of other phase | word | yes | -
233  relax | Relaxation factor for dmdt | Function1<scalar> <!--
234  --> | yes | -
235  alphatConv | Convective turbulent thermal diffusivity <!--
236  --> | scalarField | no | 0
237  partitioningModel | Run-time selected heat flux partitioning model <!--
238  --> | dict | yes | -
239  \endtable
240 
241  Options for the \c phaseType and \c otherPhase entries:
242  \verbatim
243  vapor | Vapor phase
244  liquid | Liquid phase
245  \endverbatim
246 
247  when \c phaseType=liquid:
248  \table
249  Property | Description | Type | Reqd | Deflt
250  nucleationSiteModel | Nucleation site density model | dict | yes | -
251  departureDiamModel | Bubble departure diameter model <!--
252  --> | dict | yes | -
253  departureFreqModel | Bubble departure frequency model | dict | yes | -
254  LeidenfrostModel | Leidenfrost temperature model | dict | no | -
255  CHFModel | Critical heat flux model | dict | no | -
256  CHFSubCoolModel | CHF sub-cool model | dict | no | -
257  MHFModel | Minium heat flux model | dict | no | -
258  TDNBModel | Departure from nulceate boiling model | dict | no | -
259  filmBoilingModel | Film boiling model | dict | no | -
260  K | Model constant for area of bubbles | scalar | no | 4.0
261  wp | Wetting parameter for transient boiling | scalar | no | 1.0
262  \endtable
263 
264  The inherited entries are elaborated in:
265  -\link alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.H\endlink
266 
267 Notes
268  - Runtime selectabale submodels may require model specific entries
269  - \c phaseType and \c otherPhase entries should be the opposite of each other.
270 
271 See also
272  Foam::alphatPhaseChangeJayatillekeWallFunctionFvPatchField
273 
274 SourceFiles
275  alphatWallBoilingWallFunctionFvPatchScalarField.C
276 
277 \*---------------------------------------------------------------------------*/
278 
279 #ifndef compressible_alphatWallBoilingWallFunctionFvPatchScalarField_H
280 #define compressible_alphatWallBoilingWallFunctionFvPatchScalarField_H
281 
282 #include "Function1.H"
283 
285 #include "partitioningModel.H"
286 #include "nucleationSiteModel.H"
287 #include "departureDiameterModel.H"
288 #include "departureFrequencyModel.H"
289 #include "nucleateFluxModel.H"
290 
291 #include "LeidenfrostModel.H"
292 #include "filmBoilingModel.H"
293 #include "CHFModel.H"
294 #include "CHFSubCoolModel.H"
295 #include "MHFModel.H"
296 #include "TDNBModel.H"
297 
298 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
299 
300 namespace Foam
301 {
302 namespace compressible
303 {
304 
305 /*---------------------------------------------------------------------------*\
306  Class alphatWallBoilingWallFunctionFvPatchScalarField Declaration
307 \*---------------------------------------------------------------------------*/
308 
309 class alphatWallBoilingWallFunctionFvPatchScalarField
310 :
311  public alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField
312 {
313 public:
314 
315  // Public Enumerations
316 
317  //- Enumeration listing the possible operational modes
318  enum phaseType
319  {
320  vaporPhase,
322  };
323 
324 
325 private:
326 
327  // Private Data
328 
329  //- Enumeration of regimes per face
330  enum regimeType
331  {
332  subcool,
333  transient,
334  film,
335  nonBoiling
336  };
337 
338  //- Name of the other phase (vapor/liquid phase)
339  word otherPhaseName_;
340 
341  //- Names of heat source types
342  static const Enum<phaseType> phaseTypeNames_;
343 
344  //- Heat source type
345  phaseType phaseType_;
346 
347  //- Relaxation factor for dmdt
348  autoPtr<Function1<scalar>> relax_;
349 
350  //- Patch face area by cell volume
351  scalarField AbyV_;
352 
353  // Sub-cooling nucleating boiling
354 
355  //- Convective turbulent thermal diffusivity
356  scalarField alphatConv_;
357 
358  //- Departure diameter field
359  scalarField dDep_;
360 
361  //- Quenching surface heat flux
362  scalarField qq_;
363 
364  //- Model constant for area of bubbles
365  scalar K_;
366 
367  //- Run-time selected heat flux partitioning model
368  autoPtr<wallBoilingModels::partitioningModel>
369  partitioningModel_;
370 
371  //- Run-time selected nucleation site density model
372  autoPtr<wallBoilingModels::nucleationSiteModel>
373  nucleationSiteModel_;
374 
375  //- Run-time selected bubble departure diameter model
376  autoPtr<wallBoilingModels::departureDiameterModel>
377  departureDiamModel_;
378 
379  //- Run-time selected bubble departure frequency model
380  autoPtr<wallBoilingModels::departureFrequencyModel>
381  departureFreqModel_;
382 
383  //- Run-time sub-cooling heat flux correlatiom
384  autoPtr<wallBoilingModels::nucleateFluxModel>
385  nucleatingModel_;
386 
387 
388  // Film boiling model
389 
390  //- Run-time selected for filmBoiling model
391  autoPtr<wallBoilingModels::filmBoilingModel>
392  filmBoilingModel_;
393 
394 
395  // Transition boiling model
396 
397  //- Run-time selected for Leidenfrost temperature
398  autoPtr<wallBoilingModels::LeidenfrostModel>
399  LeidenfrostModel_;
400 
401  //- Run-time selected for CHF
402  autoPtr<wallBoilingModels::CHFModel> CHFModel_;
403 
404  //- Run-time selected for CHF sub-cool
405  autoPtr<wallBoilingModels::CHFSubCoolModel> CHFSoobModel_;
406 
407  //- Run-time selected for MHF
408  autoPtr<wallBoilingModels::MHFModel> MHFModel_;
409 
410  //- Run-time selected for MHF
411  autoPtr<wallBoilingModels::TDNBModel> TDNBModel_;
412 
413  //- Wetting parameter for transient boiling
414  scalar wp_;
415 
416  //- Use Liquid temperature at y+=250
417  bool liquidTatYplus_;
418 
419  //- Face regime
420  labelField regimeTypes_;
421 
422 
423 public:
424 
425  //- Runtime type information
426  TypeName("compressible::alphatWallBoilingWallFunction");
427 
428 
429  // Constructors
430 
431  //- Construct from patch and internal field
433  (
434  const fvPatch&,
435  const DimensionedField<scalar, volMesh>&
436  );
437 
438  //- Construct from patch, internal field and dictionary
440  (
441  const fvPatch&,
443  const dictionary&
444  );
445 
446  //- Construct by mapping given
447  //- alphatWallBoilingWallFunctionFvPatchScalarField
448  //- onto a new patch
450  (
452  const fvPatch&,
454  const fvPatchFieldMapper&
455  );
456 
457  //- Construct as copy
459  (
461  );
462 
463  //- Construct and return a clone
464  virtual tmp<fvPatchScalarField> clone() const
465  {
467  (
469  );
470  }
471 
472  //- Construct as copy setting internal field reference
474  (
477  );
478 
479  //- Construct and return a clone setting internal field reference
481  (
483  ) const
484  {
486  (
488  );
489  }
490 
491 
492  // Member Functions
493 
495 
496  //- Is there phase change mass transfer for this phasePair
497  virtual bool activePhasePair(const phasePairKey&) const;
498 
499  //- Return the rate of phase-change for specific phase pair
500  virtual const scalarField& dmdt(const phasePairKey&) const;
501 
502  //- Return the rate of phase-change for specific phase pair
503  virtual const scalarField& mDotL(const phasePairKey&) const;
504 
505  //- Return the departure diameter field
506  const scalarField& dDeparture() const
507  {
508  return dDep_;
509  }
510 
511  //- Return the quenching surface heat flux [W/m2]
512  const scalarField& qq() const
513  {
514  return qq_;
515  }
516 
517  //- Return the evaporation surface heat flux [W/m2]
518  tmp<scalarField> qe() const
519  {
520  return mDotL_/AbyV_;
521  }
522 
523  //- Return const reference to the face regime
524  const labelField& regimeTypes() const noexcept
525  {
526  return regimeTypes_;
527  }
528 
529  // Evaluation
530 
531  //- Update the coefficients associated with the patch field
532  virtual void updateCoeffs();
533 
534 
535  // I-O
536 
537  //- Write
538  virtual void write(Ostream&) const;
539 };
540 
541 
542 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
543 
544 } // End namespace compressible
545 } // End namespace Foam
546 
547 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
548 
549 #endif
550 
551 // ************************************************************************* //
nucleationSiteModel.H
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::Enum< phaseType >
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::compressible::alphatPhaseChangeWallFunctionFvPatchScalarField::mDotL
virtual const scalarField & mDotL() const
Return the enthalpy source due to phase-change.
Definition: alphatPhaseChangeWallFunctionFvPatchScalarField.H:174
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
LeidenfrostModel.H
Foam::compressible::alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField
This boundary condition provides a thermal wall function for turbulent thermal diffusivity (usuallyal...
Definition: alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.H:132
Function1.H
alphatPhaseChangeJayatillekeWallFunctionFvPatchScalarField.H
partitioningModel.H
departureDiameterModel.H
Foam::compressible::alphatWallBoilingWallFunctionFvPatchScalarField::activePhasePair
virtual bool activePhasePair(const phasePairKey &) const
Is there phase change mass transfer for this phasePair.
Definition: alphatWallBoilingWallFunctionFvPatchScalarField.C:392
Foam::compressible::alphatWallBoilingWallFunctionFvPatchScalarField
A thermal wall function for simulation of boiling wall.
Definition: alphatWallBoilingWallFunctionFvPatchScalarField.H:438
TDNBModel.H
Foam::compressible::alphatWallBoilingWallFunctionFvPatchScalarField::vaporPhase
Definition: alphatWallBoilingWallFunctionFvPatchScalarField.H:449
Foam::Field< scalar >
Foam::labelField
Field< label > labelField
Specialisation of Field<T> for label.
Definition: labelField.H:52
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::phasePairKey
An ordered or unorder pair of phase names. Typically specified as follows.
Definition: phasePairKey.H:65
Foam::compressible::alphatWallBoilingWallFunctionFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: alphatWallBoilingWallFunctionFvPatchScalarField.C:1224
nucleateFluxModel.H
Foam::compressible::alphatWallBoilingWallFunctionFvPatchScalarField::dDeparture
const scalarField & dDeparture() const
Return the departure diameter field.
Definition: alphatWallBoilingWallFunctionFvPatchScalarField.H:635
CHFSubCoolModel.H
compressible
bool compressible
Definition: pEqn.H:2
Foam::compressible::alphatWallBoilingWallFunctionFvPatchScalarField::alphatWallBoilingWallFunctionFvPatchScalarField
alphatWallBoilingWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: alphatWallBoilingWallFunctionFvPatchScalarField.C:70
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::compressible::alphatWallBoilingWallFunctionFvPatchScalarField::qq
const scalarField & qq() const
Return the quenching surface heat flux [W/m2].
Definition: alphatWallBoilingWallFunctionFvPatchScalarField.H:641
CHFModel.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::compressible::alphatPhaseChangeWallFunctionFvPatchScalarField::mDotL_
scalarField mDotL_
Latent heat of the phase-change.
Definition: alphatPhaseChangeWallFunctionFvPatchScalarField.H:112
Foam::compressible::alphatWallBoilingWallFunctionFvPatchScalarField::clone
virtual tmp< fvPatchScalarField > clone() const
Construct and return a clone.
Definition: alphatWallBoilingWallFunctionFvPatchScalarField.H:593
Foam::compressible::alphatWallBoilingWallFunctionFvPatchScalarField::liquidPhase
Definition: alphatWallBoilingWallFunctionFvPatchScalarField.H:450
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::compressible::alphatWallBoilingWallFunctionFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: alphatWallBoilingWallFunctionFvPatchScalarField.C:433
Foam::compressible::alphatWallBoilingWallFunctionFvPatchScalarField::TypeName
TypeName("compressible::alphatWallBoilingWallFunction")
Runtime type information.
departureFrequencyModel.H
Foam::compressible::alphatWallBoilingWallFunctionFvPatchScalarField::phaseType
phaseType
Enumeration listing the possible operational modes.
Definition: alphatWallBoilingWallFunctionFvPatchScalarField.H:447
Foam::compressible::alphatWallBoilingWallFunctionFvPatchScalarField::qe
tmp< scalarField > qe() const
Return the evaporation surface heat flux [W/m2].
Definition: alphatWallBoilingWallFunctionFvPatchScalarField.H:647
Foam::compressible::alphatWallBoilingWallFunctionFvPatchScalarField::regimeTypes
const labelField & regimeTypes() const noexcept
Return const reference to the face regime.
Definition: alphatWallBoilingWallFunctionFvPatchScalarField.H:653
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::compressible::alphatPhaseChangeWallFunctionFvPatchScalarField::dmdt
virtual const scalarField & dmdt() const
Return the rate of phase-change.
Definition: alphatPhaseChangeWallFunctionFvPatchScalarField.H:168
MHFModel.H
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
filmBoilingModel.H
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54