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-2025 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 
45 Foam::diameterModels::velocityGroup::dsm() const
46 {
47  tmp<volScalarField> tsumFi
48  (
50  (
51  "sumFi",
52  phase().mesh(),
54  )
55  );
56  tmp<volScalarField> tsumFiAbyV
57  (
59  (
60  "sumFiAbyV",
61  phase().mesh(),
63  )
64  );
65 
66  volScalarField& sumFi = tsumFi.ref();
67  volScalarField& sumFiAbyV = tsumFiAbyV.ref();
68 
69  forAll(sizeGroups_, i)
70  {
71  const sizeGroup& fi = sizeGroups_[i];
72 
73  sumFi += max(fi, rootVSmall);
74  sumFiAbyV += max(fi, rootVSmall)*fi.a()/fi.x();
75  }
76 
77  return 6*sumFi/tsumFiAbyV;
78 }
79 
80 
82 Foam::diameterModels::velocityGroup::fSum() const
83 {
84  tmp<volScalarField> tsumSizeGroups
85  (
87  (
88  "sumSizeGroups",
89  phase().mesh(),
91  )
92  );
93 
94  volScalarField& sumSizeGroups = tsumSizeGroups.ref();
95 
96  forAll(sizeGroups_, i)
97  {
98  sumSizeGroups += sizeGroups_[i];
99  }
100 
101  return tsumSizeGroups;
102 }
103 
104 
105 void Foam::diameterModels::velocityGroup::scale()
106 {
107  Info<< "Scaling sizeGroups for velocityGroup " << phase().name() << endl;
108 
109  forAll(sizeGroups_, i)
110  {
111  sizeGroups_[i].max(0);
112  };
113 
114  const volScalarField fSum(this->fSum());
115 
116  forAll(sizeGroups_, i)
117  {
118  sizeGroups_[i] /= fSum;
119 
120  sizeGroups_[i].correctBoundaryConditions();
121  };
122 }
123 
124 
125 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
126 
128 (
129  const dictionary& diameterProperties,
130  const phaseModel& phase
131 )
132 :
133  diameterModel(diameterProperties, phase),
134  popBalName_(diameterProperties.lookup("populationBalance")),
135  popBalPtr_(nullptr),
136  sizeGroups_
137  (
138  diameterProperties.lookup("sizeGroups"),
139  sizeGroup::iNew
140  (
141  *this,
142  populationBalanceModel::groups::New
143  (
144  popBalName_,
145  phase.mesh()
146  ).nSizeGroups()
147  )
148  ),
149  d_
150  (
151  IOobject
152  (
153  IOobject::groupName("d", phase.name()),
154  phase.time().name(),
155  phase.mesh(),
156  IOobject::NO_READ,
157  IOobject::AUTO_WRITE
158  ),
159  dsm()
160  )
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 {
199  tmp<volScalarField> tsumFiAbyV
200  (
202  (
203  "sumFiAbyV",
204  phase().mesh(),
206  )
207  );
208  volScalarField& sumFiAbyV = tsumFiAbyV.ref();
209 
210  forAll(sizeGroups_, i)
211  {
212  const sizeGroup& fi = sizeGroups_[i];
213 
214  sumFiAbyV += fi*fi.a()/fi.x();
215  }
216 
217  return phase()*tsumFiAbyV;
218 }
219 
220 
222 {
223  const populationBalanceModel& popBal = this->popBal();
224 
225  if (!popBal.solveOnFinalIterOnly() || popBal.fluid().pimple().finalIter())
226  {
227  forAll(sizeGroups_, i)
228  {
229  sizeGroups_[i].correct();
230  }
231 
232  if (popBal.solverDict().lookupOrDefault<Switch>("scale", true))
233  {
234  scale();
235  }
236 
237  volScalarField::Internal fSum(this->fSum());
238 
239  Info<< phase().name() << " sizeGroups-sum volume fraction, min, max = "
240  << fSum.weightedAverage(phase().mesh().V()).value()
241  << ' ' << min(fSum).value()
242  << ' ' << max(fSum).value()
243  << endl;
244 
245  d_ = dsm();
246 
247  Info<< this->phase().name() << " Sauter mean diameter, min, max = "
248  << d_.weightedAverage(d_.mesh().V()).value()
249  << ' ' << min(d_).value()
250  << ' ' << max(d_).value()
251  << endl;
252  }
253 }
254 
255 
257 (
259 )
260 {
262 
263  return true;
264 }
265 
266 
267 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
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, PrimitiveField2 > &) const
Calculate and return weighted average.
Generic GeometricField class.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, 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
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
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.
bool solveOnFinalIterOnly() const
Solve on final pimple iteration only.
const fvMesh & mesh() const
Return reference to the mesh.
const dictionary & solverDict() const
Return solution settings dictionary for this populationBalance.
Single size class fraction field representing a fixed particle volume as defined by the user through ...
Definition: sizeGroup.H:96
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:140
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 keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
T lookupOrDefault(const word &, const T &, const bool writeDefault=writeOptionalEntries > 0) const
Find and return a T, if not found return the given default.
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:97
bool finalIter() const
Flag to indicate the final iteration.
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:197
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
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:258
const dimensionSet dimless
messageStream Info
const dimensionSet dimLength
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:62
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
void inv(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
tmp< DimensionedField< TypeR, GeoMesh, Field > > New(const tmp< DimensionedField< TypeR, GeoMesh, Field >> &tdf1, const word &name, const dimensionSet &dimensions)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.