All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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-2022 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 {
39  defineTypeNameAndDebug(wallBoiling, 0);
41  (
42  nucleationModel,
43  wallBoiling,
44  dictionary
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 =
70  (
72  );
73 
74  const volScalarField::Boundary& alphatBf = alphat.boundaryField();
75 
76  typedef compressible::alphatWallBoilingWallFunctionFvPatchScalarField
77  alphatWallBoilingWallFunction;
78 
79  forAll(alphatBf, patchi)
80  {
81  if
82  (
83  isA<alphatWallBoilingWallFunction>(alphatBf[patchi])
84  )
85  {
86  const alphatWallBoilingWallFunction& alphatw =
87  refCast<const alphatWallBoilingWallFunction>(alphatBf[patchi]);
88 
89  const scalarField& dDep = alphatw.dDeparture();
90 
91  if (min(dDep) < velGroup_.sizeGroups().first().dSph().value())
92  {
93  Warning
94  << "Minimum departure diameter " << min(dDep)
95  << " m outside of range ["
96  << velGroup_.sizeGroups().first().dSph().value() << ", "
97  << velGroup_.sizeGroups().last().dSph().value() << "] m"
98  << " at patch " << alphatw.patch().name()
99  << endl
100  << " The nucleation rate in populationBalance "
101  << popBal_.name() << " is set to zero." << endl
102  << " Adjust discretisation over property space to"
103  << " suppress this warning."
104  << endl;
105  }
106  else if (max(dDep) > velGroup_.sizeGroups().last().dSph().value())
107  {
108  Warning
109  << "Maximum departure diameter " << max(dDep)
110  << " m outside of range ["
111  << velGroup_.sizeGroups().first().dSph().value() << ", "
112  << velGroup_.sizeGroups().last().dSph().value() << "] m"
113  << " at patch " << alphatw.patch().name()
114  << endl
115  << " The nucleation rate in populationBalance "
116  << popBal_.name() << " is set to zero." << endl
117  << " Adjust discretisation over property space to"
118  << " suppress this warning."
119  << endl;
120  }
121  }
122  }
123 }
124 
125 
126 void
128 (
129  volScalarField& nucleationRate,
130  const label i
131 )
132 {
133  const sizeGroup& fi = popBal_.sizeGroups()[i];
134  const phaseModel& phase = fi.phase();
135  const volScalarField& rho = phase.rho();
136 
137  const volScalarField& alphat =
139  (
141  );
142 
143  const volScalarField::Boundary& alphatBf = alphat.boundaryField();
144 
145  typedef compressible::alphatWallBoilingWallFunctionFvPatchScalarField
146  alphatWallBoilingWallFunction;
147 
148  forAll(alphatBf, patchi)
149  {
150  if
151  (
152  isA<alphatWallBoilingWallFunction>(alphatBf[patchi])
153  )
154  {
155  const alphatWallBoilingWallFunction& alphatw =
156  refCast<const alphatWallBoilingWallFunction>(alphatBf[patchi]);
157 
158  const scalarField& dmdt = alphatw.dmdtf();
159  const scalarField& dDep = alphatw.dDeparture();
160 
161  const labelList& faceCells = alphatw.patch().faceCells();
162 
163  dimensionedScalar unitLength("unitLength", dimLength, 1);
164 
165  forAll(alphatw, facei)
166  {
167  if (dmdt[facei] > small)
168  {
169  const label faceCelli = faceCells[facei];
170 
171  nucleationRate[faceCelli] +=
172  popBal_.eta
173  (
174  i,
175  fi.x()/pow3(fi.dSph())*pow3(dDep[facei]*unitLength)
176  ).value()
177  *dmdt[facei]/rho[faceCelli]/fi.x().value();
178  }
179  }
180  }
181  }
182 }
183 
184 
185 // ************************************************************************* //
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
const word & name() const
Return name.
Definition: IOobject.H:315
const velocityGroup & velGroup_
Velocity group in which the nucleation occurs.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const word & name() const
Definition: phaseModel.H:110
const populationBalanceModel & popBal_
Reference to the populationBalanceModel.
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
Type of the boundary field.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Macros for easy insertion into run-time selection tables.
const dimensionedScalar eta(const label i, const dimensionedScalar &v) const
Return allocation coefficient.
const dimensionSet dimLength
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
const phaseModel & continuousPhase() const
Return continuous phase.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static word groupName(Name name, const word &group)
virtual void precompute()
Precompute diameter independent expressions.
const Type & value() const
Return const reference to value.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
List< label > labelList
A List of labels.
Definition: labelList.H:56
virtual void addToNucleationRate(volScalarField &nucleationRate, const label i)
Add to nucleationRate.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
const UPtrList< sizeGroup > & sizeGroups() const
Return the size groups belonging to this populationBalance.
messageStream Warning
dimensionedScalar pow3(const dimensionedScalar &ds)
label patchi
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
wallBoiling(const populationBalanceModel &popBal, const dictionary &dict)
const PtrList< sizeGroup > & sizeGroups() const
Return sizeGroups belonging to this velocityGroup.
const fvMesh & mesh() const
Return reference to the mesh.
Namespace for OpenFOAM.