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 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  turbulence_
76  (
77  popBal_.mesh().lookupObjectRef<phaseCompressibleTurbulenceModel>
78  (
79  IOobject::groupName
80  (
81  turbulenceModel::propertiesName,
82  popBal_.continuousPhase().name()
83  )
84  )
85  )
86 {}
87 
88 
89 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
90 
92 {
93  const tmp<volScalarField> talphat(turbulence_.alphat());
94  const volScalarField::Boundary& alphatBf = talphat().boundaryField();
95 
96  typedef compressible::alphatWallBoilingWallFunctionFvPatchScalarField
97  alphatWallBoilingWallFunction;
98 
99  forAll(alphatBf, patchi)
100  {
101  if
102  (
103  isA<alphatWallBoilingWallFunction>(alphatBf[patchi])
104  )
105  {
106  const alphatWallBoilingWallFunction& alphatw =
107  refCast<const alphatWallBoilingWallFunction>(alphatBf[patchi]);
108 
109  const scalarField& dDep = alphatw.dDeparture();
110 
111  if (min(dDep) < velGroup_.sizeGroups().first().d().value())
112  {
113  Warning
114  << "Minimum departure diameter " << min(dDep)
115  << " m outside of range ["
116  << velGroup_.sizeGroups().first().d().value() << ", "
117  << velGroup_.sizeGroups().last().d().value() << "] m"
118  << " at patch " << alphatw.patch().name()
119  << endl
120  << " The nucleation rate in populationBalance "
121  << popBal_.name() << " is set to zero." << endl
122  << " Adjust discretization over property space to"
123  << " suppress this warning."
124  << endl;
125  }
126  else if (max(dDep) > velGroup_.sizeGroups().last().d().value())
127  {
128  Warning
129  << "Maximum departure diameter " << max(dDep)
130  << " m outside of range ["
131  << velGroup_.sizeGroups().first().d().value() << ", "
132  << velGroup_.sizeGroups().last().d().value() << "] m"
133  << " at patch " << alphatw.patch().name()
134  << endl
135  << " The nucleation rate in populationBalance "
136  << popBal_.name() << " is set to zero." << endl
137  << " Adjust discretization over property space to"
138  << " suppress this warning."
139  << endl;
140  }
141  }
142  }
143 }
144 
145 
146 void
148 (
149  volScalarField& nucleationRate,
150  const label i
151 )
152 {
153  const sizeGroup& fi = popBal_.sizeGroups()[i];
154  const phaseModel& phase = fi.phase();
155  const volScalarField& rho = phase.rho();
156  const tmp<volScalarField> talphat(turbulence_.alphat());
157  const volScalarField::Boundary& alphatBf = talphat().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.dmdt();
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_.gamma
187  (
188  i,
189  velGroup_.formFactor()*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
ThermalDiffusivity< PhaseCompressibleTurbulenceModel< phaseModel > > phaseCompressibleTurbulenceModel
Typedef for phaseCompressibleTurbulenceModel.
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:295
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
const Boundary & boundaryField() const
Return const-reference to the boundary field.
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:106
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
const dimensionedScalar gamma(const label i, const dimensionedScalar &v) const
Return allocation coefficient.
const populationBalanceModel & popBal_
Reference to the populationBalanceModel.
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
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 > &)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
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)
Namespace for OpenFOAM.