velocityGroup.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) 2017-2023 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 "velocityGroup.H"
28 #include "populationBalanceModel.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace diameterModels
35 {
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
43 
44 Foam::tmp<Foam::volScalarField> Foam::diameterModels::velocityGroup::dsm() const
45 {
46  tmp<volScalarField> tInvDsm
47  (
49  (
50  "invDsm",
51  phase().mesh(),
53  )
54  );
55 
56  volScalarField& invDsm = tInvDsm.ref();
57 
58  forAll(sizeGroups_, i)
59  {
60  const sizeGroup& fi = sizeGroups_[i];
61 
62  invDsm += fi.a()*fi/fi.x();
63  }
64 
65  return 6/tInvDsm;
66 }
67 
68 
69 Foam::tmp<Foam::volScalarField> Foam::diameterModels::velocityGroup::N() const
70 {
71  tmp<volScalarField> tN
72  (
74  (
75  "N",
76  phase().mesh(),
78  )
79  );
80 
81  volScalarField& N = tN.ref();
82 
83  forAll(sizeGroups_, i)
84  {
85  N += phase()*sizeGroups_[i]/sizeGroups_[i].x();
86  }
87 
88  return tN;
89 }
90 
91 
93 Foam::diameterModels::velocityGroup::fSum() const
94 {
95  tmp<volScalarField> tsumSizeGroups
96  (
98  (
99  "sumSizeGroups",
100  phase().mesh(),
102  )
103  );
104 
105  volScalarField& sumSizeGroups = tsumSizeGroups.ref();
106 
107  forAll(sizeGroups_, i)
108  {
109  sumSizeGroups += sizeGroups_[i];
110  }
111 
112  return tsumSizeGroups;
113 }
114 
115 
116 void Foam::diameterModels::velocityGroup::scale()
117 {
118  Info<< "Scaling sizeGroups for velocityGroup " << phase().name() << endl;
119 
120  forAll(sizeGroups_, i)
121  {
122  sizeGroups_[i].max(0);
123  };
124 
125  const volScalarField fSum(this->fSum());
126 
127  forAll(sizeGroups_, i)
128  {
129  sizeGroups_[i] /= fSum;
130 
131  sizeGroups_[i].correctBoundaryConditions();
132  };
133 }
134 
135 
136 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
137 
139 (
140  const dictionary& diameterProperties,
141  const phaseModel& phase
142 )
143 :
144  diameterModel(diameterProperties, phase),
145  popBalName_(diameterProperties.lookup("populationBalance")),
146  popBalPtr_(nullptr),
147  sizeGroups_
148  (
149  diameterProperties.lookup("sizeGroups"),
150  sizeGroup::iNew
151  (
152  *this,
153  populationBalanceModel::groups::New
154  (
155  popBalName_,
156  phase.mesh()
157  ).nSizeGroups()
158  )
159  ),
160  d_(IOobject::groupName("d", phase.name()), dsm())
161 {
163  (
164  popBalName_,
165  phase.mesh()
166  ).insert(*this);
167 }
168 
169 
170 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
171 
173 {}
174 
175 
176 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
177 
180 {
181  if (popBalPtr_ == nullptr)
182  {
183  popBalPtr_ =
184  &phase().mesh().lookupObject<populationBalanceModel>(popBalName_);
185  }
186 
187  return *popBalPtr_;
188 }
189 
190 
192 {
193  return d_;
194 }
195 
196 
198 {
200  (
202  (
203  "a",
204  phase().mesh(),
206  )
207  );
208 
209  volScalarField& a = tA.ref();
210 
211  forAll(sizeGroups_, i)
212  {
213  const sizeGroup& fi = sizeGroups_[i];
214 
215  a += fi.a()*fi/fi.x();
216  }
217 
218  return phase()*a;
219 }
220 
221 
223 {
224  const populationBalanceModel& popBal = this->popBal();
225 
226  if (!popBal.solveOnFinalIterOnly() || popBal.fluid().pimple().finalIter())
227  {
228  forAll(sizeGroups_, i)
229  {
230  sizeGroups_[i].correct();
231  }
232 
233  if
234  (
235  phase().mesh().solution().solverDict(popBalName_)
236  .lookupOrDefault<Switch>
237  (
238  "scale",
239  true
240  )
241  )
242  {
243  scale();
244  }
245 
246  volScalarField::Internal fSum(this->fSum());
247 
248  Info<< phase().name() << " sizeGroups-sum volume fraction, min, max = "
249  << fSum.weightedAverage(phase().mesh().V()).value()
250  << ' ' << min(fSum).value()
251  << ' ' << max(fSum).value()
252  << endl;
253 
254  d_ = dsm();
255 
256  Info<< this->phase().name() << " Sauter mean diameter, min, max = "
257  << d_.weightedAverage(d_.mesh().V()).value()
258  << ' ' << min(d_).value()
259  << ' ' << max(d_).value()
260  << endl;
261  }
262 }
263 
264 
266 (
268 )
269 {
271 
272  return true;
273 }
274 
275 
276 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const Mesh & mesh() const
Return mesh.
dimensioned< Type > weightedAverage(const DimensionedField< scalar, GeoMesh > &) const
Calculate and return weighted average.
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
Abstract base-class for dispersed-phase particle diameter models.
Definition: diameterModel.H:52
const phaseModel & phase() const
Return the phase.
virtual bool read(const dictionary &phaseProperties)
Read phaseProperties dictionary.
Definition: diameterModel.C:62
static groups & New(const word &popBalName, const objectRegistry &db)
Lookup in the registry or construct new.
void insert(velocityGroup &group)
Insert a velocity group into the table.
Model for tracking the evolution of a dispersed phase size distribution due to coalescence (synonymou...
const phaseSystem & fluid() const
Return reference to the phaseSystem.
Switch solveOnFinalIterOnly() const
Solve on final pimple iteration only.
const fvMesh & mesh() const
Return reference to the mesh.
Single size class fraction field representing a fixed particle volume as defined by the user through ...
Definition: sizeGroup.H:102
const dimensionedScalar & x() const
Return representative volume of the sizeGroup.
Definition: sizeGroupI.H:58
const tmp< volScalarField > a() const
Return representative surface area of the sizeGroup.
Definition: sizeGroup.C:135
Computes the Sauter mean diameter based on a user specified size distribution, defined in terms of si...
Definition: velocityGroup.H:87
virtual void correct()
Correct the model.
const populationBalanceModel & popBal() const
Return the populationBalance this velocityGroup belongs to.
virtual bool read(const dictionary &diameterProperties)
Read diameterProperties dictionary.
virtual tmp< volScalarField > d() const
Get the diameter field.
velocityGroup(const dictionary &diameterProperties, const phaseModel &phase)
Construct from dictionary and phase.
virtual ~velocityGroup()
Destructor.
virtual tmp< volScalarField > Av() const
Get the surface area per unit volume field.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
Helper class to manage multi-specie phase properties.
const pimpleNoLoopControl & pimple() const
Return pimpleNoLoopControl.
Definition: phaseSystemI.H:95
bool finalIter() const
Flag to indicate the final iteration.
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:51
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
addToRunTimeSelectionTable(diameterModel, constant, dictionary)
defineTypeNameAndDebug(constant, 0)
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimless
messageStream Info
const dimensionSet dimLength
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimVolume
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.