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-2020 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"
28 #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  velGroup_
62  (
63  refCast<const velocityGroup>
64  (
65  popBal.mesh().lookupObject<phaseModel>
66  (
67  IOobject::groupName
68  (
69  "alpha",
70  dict.lookup("velocityGroup")
71  )
72  ).dPtr()()
73  )
74  )
75 {}
76 
77 
78 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
79 
81 {
82  const volScalarField& alphat =
84  (
86  );
87 
88  const volScalarField::Boundary& alphatBf = alphat.boundaryField();
89 
90  typedef compressible::alphatWallBoilingWallFunctionFvPatchScalarField
91  alphatWallBoilingWallFunction;
92 
93  forAll(alphatBf, patchi)
94  {
95  if
96  (
97  isA<alphatWallBoilingWallFunction>(alphatBf[patchi])
98  )
99  {
100  const alphatWallBoilingWallFunction& alphatw =
101  refCast<const alphatWallBoilingWallFunction>(alphatBf[patchi]);
102 
103  const scalarField& dDep = alphatw.dDeparture();
104 
105  if (min(dDep) < velGroup_.sizeGroups().first().dSph().value())
106  {
107  Warning
108  << "Minimum departure diameter " << min(dDep)
109  << " m outside of range ["
110  << velGroup_.sizeGroups().first().dSph().value() << ", "
111  << velGroup_.sizeGroups().last().dSph().value() << "] m"
112  << " at patch " << alphatw.patch().name()
113  << endl
114  << " The nucleation rate in populationBalance "
115  << popBal_.name() << " is set to zero." << endl
116  << " Adjust discretization over property space to"
117  << " suppress this warning."
118  << endl;
119  }
120  else if (max(dDep) > velGroup_.sizeGroups().last().dSph().value())
121  {
122  Warning
123  << "Maximum departure diameter " << max(dDep)
124  << " m outside of range ["
125  << velGroup_.sizeGroups().first().dSph().value() << ", "
126  << velGroup_.sizeGroups().last().dSph().value() << "] m"
127  << " at patch " << alphatw.patch().name()
128  << endl
129  << " The nucleation rate in populationBalance "
130  << popBal_.name() << " is set to zero." << endl
131  << " Adjust discretization over property space to"
132  << " suppress this warning."
133  << endl;
134  }
135  }
136  }
137 }
138 
139 
140 void
142 (
143  volScalarField& nucleationRate,
144  const label i
145 )
146 {
147  const sizeGroup& fi = popBal_.sizeGroups()[i];
148  const phaseModel& phase = fi.phase();
149  const volScalarField& rho = phase.rho();
150 
151  const volScalarField& alphat =
153  (
155  );
156 
157  const volScalarField::Boundary& alphatBf = alphat.boundaryField();
158 
159  typedef compressible::alphatWallBoilingWallFunctionFvPatchScalarField
160  alphatWallBoilingWallFunction;
161 
162  forAll(alphatBf, patchi)
163  {
164  if
165  (
166  isA<alphatWallBoilingWallFunction>(alphatBf[patchi])
167  )
168  {
169  const alphatWallBoilingWallFunction& alphatw =
170  refCast<const alphatWallBoilingWallFunction>(alphatBf[patchi]);
171 
172  const scalarField& dmdt = alphatw.dmdtf();
173  const scalarField& dDep = alphatw.dDeparture();
174 
175  const labelList& faceCells = alphatw.patch().faceCells();
176 
177  dimensionedScalar unitLength("unitLength", dimLength, 1.0);
178 
179  forAll(alphatw, facei)
180  {
181  if (dmdt[facei] > small)
182  {
183  const label faceCelli = faceCells[facei];
184 
185  nucleationRate[faceCelli] +=
186  popBal_.eta
187  (
188  i,
189  fi.x()/pow3(fi.dSph())*pow3(dDep[facei]*unitLength)
190  ).value()
191  *dmdt[facei]/rho[faceCelli]/fi.x().value();
192  }
193  }
194  }
195  }
196 }
197 
198 
199 // ************************************************************************* //
virtual void correct()
Correct diameter independent expressions.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
const word & name() const
Return name.
Definition: IOobject.H:303
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:106
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const word & name() const
Definition: phaseModel.H:109
const populationBalanceModel & popBal_
Reference to the populationBalanceModel.
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.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
const phaseModel & continuousPhase() const
Return continuous phase.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static word groupName(Name name, const word &group)
const Type & value() const
Return const reference to value.
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 sizeGroups belonging to this populationBalance.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
messageStream Warning
dimensionedScalar pow3(const dimensionedScalar &ds)
label patchi
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
wallBoiling(const populationBalanceModel &popBal, const dictionary &dict)
const fvMesh & mesh() const
Return reference to the mesh.
Namespace for OpenFOAM.