sizeDistribution.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-2019 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 "sizeDistribution.H"
27 #include "sizeGroup.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace functionObjects
35 {
36  defineTypeNameAndDebug(sizeDistribution, 0);
37  addToRunTimeSelectionTable(functionObject, sizeDistribution, dictionary);
38 }
39 }
40 
41 template<>
42 const char*
44 <
46  2
47 >::names[] = {"cellZone", "all"};
48 
49 template<>
50 const char*
52 <
54  4
55 >::names[] =
56  {
57  "numberDensity",
58  "volumeDensity",
59  "numberConcentration",
60  "moments"
61  };
62 
63 template<>
64 const char*
66 <
68  2
69 >::names[] = {"diameter", "volume"};
70 
71 const Foam::NamedEnum
72 <
74  2
76 
77 const Foam::NamedEnum
78 <
80  4
82 
83 const Foam::NamedEnum
84 <
86  2
88 
89 
90 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
91 
93 (
94  const dictionary& dict
95 )
96 {
97  switch (functionType_)
98  {
99  case ftNdf:
100  {
101  break;
102  }
103 
104  case ftVdf:
105  {
106  break;
107  }
108 
109  case ftNc:
110  {
111  break;
112  }
113 
114  case ftMom:
115  {
116  break;
117  }
118 
119  default:
120  {
122  << "Unknown functionType. Valid types are:"
123  << functionTypeNames_ << nl << exit(FatalIOError);
124  }
125  }
126 
127  switch (abszissaType_)
128  {
129  case atDiameter:
130  {
131  break;
132  }
133 
134  case atVolume:
135  {
136  break;
137  }
138 
139  default:
140  {
142  << "Unknown abszissaType. Valid types are:"
143  << abszissaTypeNames_ << nl << exit(FatalIOError);
144  }
145  }
146 
147  setCellZoneCells();
148 
149  if (nCells_ == 0)
150  {
152  << type() << " " << name() << ": "
153  << selectionModeTypeNames_[selectionModeType_]
154  << "(" << selectionModeTypeName_ << "):" << nl
155  << " Selection has no cells" << exit(FatalIOError);
156  }
157 
158  volume_ = volume();
159 
160  Info<< type() << " " << name() << ":"
161  << selectionModeTypeNames_[selectionModeType_]
162  << "(" << selectionModeTypeName_ << "):" << nl
163  << " total cells = " << nCells_ << nl
164  << " total volume = " << volume_
165  << nl << endl;
166 }
167 
168 
170 {
171  switch (selectionModeType_)
172  {
173  case rtCellZone:
174  {
175  dict().lookup("cellZone") >> selectionModeTypeName_;
176 
177  label zoneId =
178  mesh().cellZones().findZoneID(selectionModeTypeName_);
179 
180  if (zoneId < 0)
181  {
183  << "Unknown cellZone name: " << selectionModeTypeName_
184  << ". Valid cellZone names are: "
185  << mesh().cellZones().names()
186  << nl << exit(FatalIOError);
187  }
188 
189  cellId_ = mesh().cellZones()[zoneId];
190  nCells_ = returnReduce(cellId_.size(), sumOp<label>());
191  break;
192  }
193 
194  case rtAll:
195  {
196  cellId_ = identity(mesh().nCells());
197  nCells_ = returnReduce(cellId_.size(), sumOp<label>());
198  break;
199  }
200 
201  default:
202  {
204  << "Unknown selectionMode type. Valid selectionMode types are:"
205  << selectionModeTypeNames_ << nl << exit(FatalIOError);
206  }
207  }
208 }
209 
210 
212 {
213  return gSum(filterField(mesh().V()));
214 }
215 
216 
218 {
219  List<scalarField> allValues(Pstream::nProcs());
220 
221  allValues[Pstream::myProcNo()] = field;
222 
223  Pstream::gatherList(allValues);
224 
225  if (Pstream::master())
226  {
227  field =
228  ListListOps::combine<scalarField>
229  (
230  allValues,
231  accessOp<scalarField>()
232  );
233  }
234 }
235 
236 
239 (
240  const scalarField& field
241 ) const
242 {
243  return tmp<scalarField>(new scalarField(field, cellId_));
244 }
245 
246 
248 (
249  const label i
250 )
251 {
252  OFstream& file = this->file();
253 
254  switch (functionType_)
255  {
256  case ftNdf:
257  {
258  writeHeader(file, "Number density function");
259  break;
260  }
261 
262  case ftVdf:
263  {
264  writeHeader(file, "Volume density function");
265  break;
266  }
267 
268  case ftNc:
269  {
270  writeHeader(file, "Number concentration");
271  break;
272  }
273 
274  case ftMom:
275  {
276  writeHeader(file, "Moments");
277  break;
278  }
279  }
280 
281  switch (abszissaType_)
282  {
283  case atVolume:
284  {
285  writeCommented(file, "Time/volume");
286  break;
287  }
288 
289  case atDiameter:
290  {
291  writeCommented(file, "Time/diameter");
292  break;
293  }
294  }
295 
296  switch (functionType_)
297  {
298  case ftMom:
299  {
300  for (label i = 0; i <= momentOrder_; i++)
301  {
302  file() << tab << i;
303  }
304 
305  break;
306  }
307  default:
308  {
309  forAll(popBal_.sizeGroups(), sizeGroupi)
310  {
311  const diameterModels::sizeGroup& fi =
312  popBal_.sizeGroups()[sizeGroupi];
313 
314  switch (abszissaType_)
315  {
316  case atDiameter:
317  {
318  file() << tab << fi.d().value();
319 
320  break;
321  }
322 
323  case atVolume:
324  {
325  file() << tab << fi.x().value();
326 
327  break;
328  }
329  }
330  }
331 
332  break;
333  }
334  }
335 
336  file << endl;
337 }
338 
339 
340 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
341 
343 (
344  const word& name,
345  const Time& runTime,
346  const dictionary& dict
347 )
348 :
349  fvMeshFunctionObject(name, runTime, dict),
350  logFiles(obr_, name),
351  dict_(dict),
352  selectionModeType_
353  (
354  selectionModeTypeNames_.read(dict.lookup("selectionMode"))
355  ),
356  selectionModeTypeName_(word::null),
357  functionType_(functionTypeNames_.read(dict.lookup("functionType"))),
358  abszissaType_(abszissaTypeNames_.read(dict.lookup("abszissaType"))),
359  nCells_(0),
360  cellId_(),
361  volume_(0.0),
362  writeVolume_(dict.lookupOrDefault("writeVolume", false)),
363  popBal_
364  (
365  obr_.lookupObject<Foam::diameterModels::populationBalanceModel>
366  (
367  dict.lookup("populationBalance")
368  )
369  ),
370  N_(popBal_.sizeGroups().size()),
371  momentOrder_(dict.lookupOrDefault<label>("momentOrder", 0)),
372  normalize_(dict.lookupOrDefault("normalize", false)),
373  sumN_(0.0),
374  sumV_(0.0)
375 {
376  read(dict);
377  resetName(name);
378 }
379 
380 
381 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
382 
384 {}
385 
386 
387 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
388 
389 bool Foam::functionObjects::sizeDistribution::read(const dictionary& dict)
390 {
391  if (dict != dict_)
392  {
393  dict_ = dict;
394  }
395 
397 
398  initialise(dict);
399 
400  return true;
401 }
402 
403 
405 {
406  return true;
407 }
408 
409 
411 {
412  logFiles::write();
413 
414  if (Pstream::master())
415  {
416  writeTime(file());
417  }
418 
419  Log << type() << " " << name() << " write" << nl;
420 
421  scalarField V(filterField(mesh().V()));
422  combineFields(V);
423 
424  sumN_ = 0;
425  sumV_ = 0;
426 
427  forAll(N_, i)
428  {
429  const Foam::diameterModels::sizeGroup& fi = popBal_.sizeGroups()[i];
430 
431  const volScalarField& alpha = fi.VelocityGroup().phase();
432 
433  scalarField Ni(fi*alpha/fi.x());
434  scalarField values(filterField(Ni));
435  scalarField V(filterField(mesh().V()));
436 
437  // Combine onto master
438  combineFields(values);
439  combineFields(V);
440 
441  if (Pstream::master())
442  {
443  // Calculate volume-averaged number concentration
444  N_[i] = sum(V*values)/sum(V);
445  }
446 
447  sumN_ += N_[i];
448 
449  sumV_ += N_[i]*fi.x().value();
450  }
451 
452  if (Pstream::master())
453  {
454  switch (functionType_)
455  {
456  case ftMom:
457  {
458  for (label m = 0; m <= momentOrder_; m++)
459  {
460  scalar result(0.0);
461 
462  forAll(N_, i)
463  {
465  popBal_.sizeGroups()[i];
466 
467  switch (abszissaType_)
468  {
469  case atVolume:
470  {
471  result += pow(fi.x().value(), m)*N_[i];
472 
473  break;
474  }
475 
476  case atDiameter:
477  {
478  result += pow(fi.d().value(), m)*N_[i];
479 
480  break;
481  }
482  }
483  }
484 
485  file() << tab << result;
486  }
487 
488  break;
489  }
490 
491  default:
492  {
493  forAll(popBal_.sizeGroups(), i)
494  {
496  popBal_.sizeGroups()[i];
497 
498  scalar result(0.0);
499  scalar delta(0.0);
500 
501  switch (abszissaType_)
502  {
503  case atVolume:
504  {
505  delta = popBal_.v()[i+1].value()
506  - popBal_.v()[i].value();
507 
508  break;
509  }
510 
511  case atDiameter:
512  {
513  const scalar& formFactor =
514  fi.VelocityGroup().formFactor().value();
515 
516  delta =
517  pow
518  (
519  popBal_.v()[i+1].value()
520  /formFactor,
521  1.0/3.0
522  )
523  - pow
524  (
525  popBal_.v()[i].value()
526  /formFactor,
527  1.0/3.0
528  );
529 
530  break;
531  }
532  }
533 
534  switch (functionType_)
535  {
536  case ftNdf:
537  {
538  if (normalize_ == true)
539  {
540  result = N_[i]/delta/sumN_;
541  }
542  else
543  {
544  result = N_[i]/delta;
545  }
546 
547  break;
548  }
549 
550  case ftVdf:
551  {
552  if (normalize_ == true)
553  {
554  result = N_[i]*fi.x().value()/delta/sumV_;
555  }
556  else
557  {
558  result = N_[i]*fi.x().value()/delta;
559  }
560 
561  break;
562  }
563 
564  case ftNc:
565  {
566  if (normalize_ == true)
567  {
568  result = N_[i]/sumN_;
569  }
570  else
571  {
572  result = N_[i];
573  }
574 
575  break;
576  }
577 
578  default:
579  {
580  break;
581  }
582  }
583 
584  file()<< tab << result;
585  }
586  }
587  }
588  }
589 
590  if (Pstream::master())
591  {
592  file()<< endl;
593  }
594 
595  Log << endl;
596 
597  return true;
598 }
599 
600 
601 // ************************************************************************* //
scalar delta
tmp< scalarField > filterField(const scalarField &field) const
Filter field according to cellIds.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual bool write()
Write function.
Definition: logFiles.C:180
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
virtual bool read(const dictionary &dict)
Read the input data.
const dimensionedScalar & x() const
Return representative volume of the sizeGroup.
Definition: sizeGroupI.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const word & name() const
Return the name of this functionObject.
static const char tab
Definition: Ostream.H:264
virtual ~sizeDistribution()
Destructor.
addToRunTimeSelectionTable(functionObject, Qdot, dictionary)
static const NamedEnum< functionTypes, 4 > functionTypeNames_
Function type names.
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
void combineFields(scalarField &field)
Combine fields from all processor domains into single field.
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
static const NamedEnum< abszissaTypes, 2 > abszissaTypeNames_
Abszissa type names.
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
virtual bool read(const dictionary &dict)
Read from dictionary.
void initialise(const dictionary &dict)
Initialise, e.g. cell addressing.
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:51
Macros for easy insertion into run-time selection tables.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
virtual bool read(const dictionary &)
Read optional controls.
bool read(const char *, int32_t &)
Definition: int32IO.C:85
stressControl lookup("compactNormalStress") >> compactNormalStress
const phaseModel & phase() const
Return the phase.
dynamicFvMesh & mesh
Type gSum(const FieldField< Field, Type > &f)
This class represents a single sizeGroup belonging to a velocityGroup. The main property of a sizeGro...
Definition: sizeGroup.H:94
virtual void writeFileHeader(const label i)
Output file header information.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
scalar volume() const
Calculate and return volume of the evaluated cell zone.
virtual const word & type() const =0
Runtime type information.
const Type & value() const
Return const reference to value.
static const char nl
Definition: Ostream.H:265
const dimensionedScalar & formFactor() const
Return the form factor.
selectionModeTypes
Selection mode type enumeration.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
defineTypeNameAndDebug(Qdot, 0)
sizeDistribution(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
void setCellZoneCells()
Set cells to evaluate based on a cell zone.
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
static const NamedEnum< selectionModeTypes, 2 > selectionModeTypeNames_
Selection mode type names.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:331
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
#define Log
Report write to Foam::Info if the local log switch is true.
abszissaTypes
abszissa type enumeration
messageStream Info
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
const dimensionedScalar & d() const
Return representative diameter of the sizeGroup.
Definition: sizeGroupI.H:50
A class for managing temporary objects.
Definition: PtrList.H:53
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
functionTypes
Function type enumeration.
const velocityGroup & VelocityGroup() const
Return const-reference to the velocityGroup.
Definition: sizeGroupI.H:43
Namespace for OpenFOAM.
IOerror FatalIOError