wallBoiling.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) 2018-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 
26 #include "wallBoiling.H"
29 #include "phaseSystem.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace diameterModels
36 {
37 namespace nucleationModels
38 {
41  (
45  );
46 }
47 }
48 }
49 
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
55 (
56  const populationBalanceModel& popBal,
57  const dictionary& dict
58 )
59 :
60  nucleationModel(popBal, dict)
61 {}
62 
63 
64 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
65 
67 {
68  const volScalarField& alphat =
69  popBal_.mesh().lookupObject<volScalarField>
70  (
71  IOobject::groupName("alphat", popBal_.continuousPhase().name())
72  );
73 
74  const volScalarField::Boundary& alphatBf = alphat.boundaryField();
75 
77  alphatWallBoilingWallFunction;
78 
79  forAll(alphatBf, patchi)
80  {
81  if (isA<alphatWallBoilingWallFunction>(alphatBf[patchi]))
82  {
83  const alphatWallBoilingWallFunction& alphatw =
84  refCast<const alphatWallBoilingWallFunction>(alphatBf[patchi]);
85 
86  const scalarField& dDep = alphatw.dDeparture();
87 
88  if (min(dDep) < velGroup_.sizeGroups().first().dSph().value())
89  {
90  Warning
91  << "Minimum departure diameter " << min(dDep)
92  << " m outside of range ["
93  << velGroup_.sizeGroups().first().dSph().value() << ", "
94  << velGroup_.sizeGroups().last().dSph().value() << "] m"
95  << " at patch " << alphatw.patch().name()
96  << endl
97  << " The nucleation rate in populationBalance "
98  << popBal_.name() << " is set to zero." << endl
99  << " Adjust discretisation over property space to"
100  << " suppress this warning."
101  << endl;
102  }
103  else if (max(dDep) > velGroup_.sizeGroups().last().dSph().value())
104  {
105  Warning
106  << "Maximum departure diameter " << max(dDep)
107  << " m outside of range ["
108  << velGroup_.sizeGroups().first().dSph().value() << ", "
109  << velGroup_.sizeGroups().last().dSph().value() << "] m"
110  << " at patch " << alphatw.patch().name()
111  << endl
112  << " The nucleation rate in populationBalance "
113  << popBal_.name() << " is set to zero." << endl
114  << " Adjust discretisation over property space to"
115  << " suppress this warning."
116  << endl;
117  }
118  }
119  }
120 }
121 
122 
124 (
125  volScalarField& nucleationRate,
126  const label i
127 )
128 {
129  const sizeGroup& fi = popBal_.sizeGroups()[i];
130  const phaseModel& phase = fi.phase();
131  const volScalarField& rho = phase.rho();
132 
133  const volScalarField& alphat =
134  popBal_.mesh().lookupObject<volScalarField>
135  (
136  IOobject::groupName("alphat", popBal_.continuousPhase().name())
137  );
138 
139  const volScalarField::Boundary& alphatBf = alphat.boundaryField();
140 
142  alphatWallBoilingWallFunction;
143 
144  forAll(alphatBf, patchi)
145  {
146  if (!isA<alphatWallBoilingWallFunction>(alphatBf[patchi])) continue;
147 
148  const alphatWallBoilingWallFunction& alphatw =
149  refCast<const alphatWallBoilingWallFunction>(alphatBf[patchi]);
150 
151  const scalarField& dmdt = alphatw.dmdtf();
152  const scalarField& dDep = alphatw.dDeparture();
153 
154  const labelList& faceCells = alphatw.patch().faceCells();
155 
156  dimensionedScalar unitLength("unitLength", dimLength, 1);
157 
158  forAll(alphatw, facei)
159  {
160  if (dmdt[facei] > small)
161  {
162  const label faceCelli = faceCells[facei];
163 
164  nucleationRate[faceCelli] +=
165  popBal_.eta
166  (
167  i,
168  fi.x()/pow3(fi.dSph())*pow3(dDep[facei]*unitLength)
169  ).value()
170  *dmdt[facei]/rho[faceCelli]/fi.x().value();
171  }
172  }
173  }
174 }
175 
176 
177 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
const Mesh & mesh() const
Return mesh.
Generic GeometricBoundaryField class.
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
static word groupName(Name name, const word &group)
virtual const fileName & name() const
Return the name of the stream.
Definition: IOstream.H:294
A thermal wall function for simulation of subcooled nucleate wall boiling with runtime selectable sub...
Base class for nucleation models.
Wall-boiling model which requires a velocityGroup (i.e. phase) to be specified in which the nucleatio...
Definition: wallBoiling.H:60
virtual void precompute()
Precompute diameter independent expressions.
Definition: wallBoiling.C:66
virtual void addToNucleationRate(volScalarField &nucleationRate, const label i)
Add to nucleationRate.
Definition: wallBoiling.C:124
wallBoiling(const populationBalanceModel &popBal, const dictionary &dict)
Definition: wallBoiling.C:55
Model for tracking the evolution of a dispersed phase size distribution due to coalescence (synonymou...
Single size class fraction field representing a fixed particle volume as defined by the user through ...
Definition: sizeGroup.H:102
const dimensionedScalar & dSph() const
Return representative spherical diameter of the sizeGroup.
Definition: sizeGroupI.H:51
const dimensionedScalar & x() const
Return representative volume of the sizeGroup.
Definition: sizeGroupI.H:58
const phaseModel & phase() const
Return const-reference to the phase.
Definition: sizeGroupI.H:37
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
const Type & value() const
Return const reference to value.
virtual const volScalarField & rho() const =0
Return the density field.
label patchi
addToRunTimeSelectionTable(nucleationModel, reactionDriven, dictionary)
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimLength
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
messageStream Warning
dictionary dict