populationBalanceMoments.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) 2022-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 
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace functionObjects
34 {
37  (
41  );
42 }
43 }
44 
45 const Foam::NamedEnum
46 <
48  4
49 >
51 {"integerMoment", "mean", "variance", "stdDev"};
52 
53 const Foam::NamedEnum
54 <
56  3
57 >
59 {"volume", "area", "diameter"};
60 
61 const Foam::NamedEnum
62 <
64  3
65 >
67 {
68  "numberConcentration",
69  "volumeConcentration",
70  "areaConcentration"
71 };
72 
73 const Foam::NamedEnum
74 <
76  3
77 >
79 {"arithmetic", "geometric", "notApplicable"};
80 
81 
82 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
83 
85 Foam::functionObjects::populationBalanceMoments::coordinateTypeSymbolicName()
86 {
87  word coordinateTypeSymbolicName(word::null);
88 
89  switch (coordinateType_)
90  {
92  {
93  coordinateTypeSymbolicName = "v";
94 
95  break;
96  }
98  {
99  coordinateTypeSymbolicName = "a";
100 
101  break;
102  }
104  {
105  coordinateTypeSymbolicName = "d";
106 
107  break;
108  }
109  }
110 
111  return coordinateTypeSymbolicName;
112 }
113 
114 
116 Foam::functionObjects::populationBalanceMoments::weightTypeSymbolicName()
117 {
118  word weightTypeSymbolicName(word::null);
119 
120  switch (weightType_)
121  {
122  case weightType::numberConcentration:
123  {
124  weightTypeSymbolicName = "N";
125 
126  break;
127  }
128  case weightType::volumeConcentration:
129  {
130  weightTypeSymbolicName = "V";
131 
132  break;
133  }
134  case weightType::areaConcentration:
135  {
136  weightTypeSymbolicName = "A";
137 
138  break;
139  }
140  }
141 
142  return weightTypeSymbolicName;
143 }
144 
145 
146 Foam::word Foam::functionObjects::populationBalanceMoments::defaultFldName()
147 {
148  word meanName
149  (
150  meanType_ == meanType::geometric
151  ? word(meanTypeNames_[meanType_]).capitalise()
152  : word("")
153  );
154 
155  return
156  word
157  (
159  (
160  "weighted"
161  + meanName
162  + word(momentTypeNames_[momentType_]).capitalise()
163  + "("
164  + weightTypeSymbolicName()
165  + ","
166  + coordinateTypeSymbolicName()
167  + ")",
168  popBalName_
169  )
170  );
171 }
172 
173 
175 Foam::functionObjects::populationBalanceMoments::integerMomentFldName()
176 {
177  return
178  word
179  (
181  (
182  word(momentTypeNames_[momentType_])
183  + Foam::name(order_)
184  + "("
185  + weightTypeSymbolicName()
186  + ","
187  + coordinateTypeSymbolicName()
188  + ")",
189  popBalName_
190  )
191  );
192 }
193 
194 
195 void Foam::functionObjects::populationBalanceMoments::setDimensions
196 (
198  momentType momType
199 )
200 {
201  switch (momType)
202  {
203  case momentType::integerMoment:
204  {
205  switch (coordinateType_)
206  {
207  case coordinateType::volume:
208  {
209  fld.dimensions().reset
210  (
211  pow(dimVolume, order_)/dimVolume
212  );
213 
214  break;
215  }
216  case coordinateType::area:
217  {
218  fld.dimensions().reset
219  (
220  pow(dimArea, order_)/dimVolume
221  );
222 
223  break;
224  }
225  case coordinateType::diameter:
226  {
227  fld.dimensions().reset
228  (
229  pow(dimLength, order_)/dimVolume
230  );
231 
232  break;
233  }
234  }
235 
236  switch (weightType_)
237  {
238  case weightType::volumeConcentration:
239  {
240  fld.dimensions().reset(fld.dimensions()*dimVolume);
241 
242  break;
243  }
244  case weightType::areaConcentration:
245  {
246  fld.dimensions().reset(fld.dimensions()*dimArea);
247 
248  break;
249  }
250  default:
251  {
252  break;
253  }
254  }
255 
256  break;
257  }
258  case momentType::mean:
259  {
260  switch (coordinateType_)
261  {
262  case coordinateType::volume:
263  {
264  fld.dimensions().reset(dimVolume);
265 
266  break;
267  }
268  case coordinateType::area:
269  {
270  fld.dimensions().reset(dimArea);
271 
272  break;
273  }
274  case coordinateType::diameter:
275  {
276  fld.dimensions().reset(dimLength);
277 
278  break;
279  }
280  }
281 
282  break;
283  }
284  case momentType::variance:
285  {
286  switch (coordinateType_)
287  {
288  case coordinateType::volume:
289  {
290  fld.dimensions().reset(sqr(dimVolume));
291 
292  break;
293  }
294  case coordinateType::area:
295  {
296  fld.dimensions().reset(sqr(dimArea));
297 
298  break;
299  }
300  case coordinateType::diameter:
301  {
302  fld.dimensions().reset(sqr(dimLength));
303 
304  break;
305  }
306  }
307 
308  if (meanType_ == meanType::geometric)
309  {
310  fld.dimensions().reset(dimless);
311  }
312 
313  break;
314  }
315  case momentType::stdDev:
316  {
317  switch (coordinateType_)
318  {
319  case coordinateType::volume:
320  {
321  fld.dimensions().reset(dimVolume);
322 
323  break;
324  }
325  case coordinateType::area:
326  {
327  fld.dimensions().reset(dimArea);
328 
329  break;
330  }
331  case coordinateType::diameter:
332  {
333  fld.dimensions().reset(dimLength);
334 
335  break;
336  }
337  }
338 
339  if (meanType_ == meanType::geometric)
340  {
341  fld.dimensions().reset(dimless);
342  }
343 
344  break;
345  }
346  }
347 }
348 
349 
351 Foam::functionObjects::populationBalanceMoments::totalConcentration
352 (
353  const diameterModels::populationBalanceModel& popBal
354 )
355 {
356  tmp<volScalarField> tTotalConcentration
357  (
359  (
360  "totalConcentration",
361  mesh_,
363  )
364  );
365 
366  volScalarField& totalConcentration = tTotalConcentration.ref();
367 
368  switch (weightType_)
369  {
370  case weightType::volumeConcentration:
371  {
372  totalConcentration.dimensions().reset
373  (
374  totalConcentration.dimensions()*dimVolume
375  );
376 
377  break;
378  }
379  case weightType::areaConcentration:
380  {
381  totalConcentration.dimensions().reset
382  (
383  totalConcentration.dimensions()*dimArea
384  );
385 
386  break;
387  }
388  default:
389  {
390  break;
391  }
392  }
393 
394  forAll(popBal.sizeGroups(), i)
395  {
396  const Foam::diameterModels::sizeGroup& fi = popBal.sizeGroups()[i];
397 
398  switch (weightType_)
399  {
400  case weightType::numberConcentration:
401  {
402  totalConcentration += fi*fi.phase()/fi.x();
403 
404  break;
405  }
406  case weightType::volumeConcentration:
407  {
408  totalConcentration += fi*fi.phase();
409 
410  break;
411  }
412  case weightType::areaConcentration:
413  {
414  totalConcentration += fi.a()*fi*fi.phase()/fi.x();
415 
416  break;
417  }
418  }
419  }
420 
421  return tTotalConcentration;
422 }
423 
424 
426 Foam::functionObjects::populationBalanceMoments::mean
427 (
428  const diameterModels::populationBalanceModel& popBal
429 )
430 {
431  tmp<volScalarField> tMean
432  (
434  (
435  "mean",
436  mesh_,
438  )
439  );
440 
441  volScalarField& mean = tMean.ref();
442 
443  setDimensions(mean, momentType::mean);
444 
445  volScalarField totalConcentration(this->totalConcentration(popBal));
446 
447  forAll(popBal.sizeGroups(), i)
448  {
449  const Foam::diameterModels::sizeGroup& fi = popBal.sizeGroups()[i];
450 
451  volScalarField concentration(fi*fi.phase()/fi.x());
452 
453  switch (weightType_)
454  {
455  case weightType::volumeConcentration:
456  {
457  concentration *= fi.x();
458 
459  break;
460  }
461  case weightType::areaConcentration:
462  {
463  concentration *= fi.a();
464 
465  break;
466  }
467  default:
468  {
469  break;
470  }
471  }
472 
473  switch (meanType_)
474  {
475  case meanType::geometric:
476  {
477  mean.dimensions().reset(dimless);
478 
479  switch (coordinateType_)
480  {
481  case coordinateType::volume:
482  {
483  dimensionedScalar unitVolume(dimVolume, 1);
484 
485  mean +=
486  Foam::log(fi.x()/unitVolume)
487  *concentration/totalConcentration;
488 
489  break;
490  }
491  case coordinateType::area:
492  {
493  dimensionedScalar unitArea(dimArea, 1);
494 
495  mean +=
496  Foam::log(fi.a()/unitArea)
497  *concentration/totalConcentration;
498 
499  break;
500  }
501  case coordinateType::diameter:
502  {
503  dimensionedScalar unitLength(dimLength, 1);
504 
505  mean +=
506  Foam::log(fi.d()/unitLength)
507  *concentration/totalConcentration;
508 
509  break;
510  }
511  }
512 
513  break;
514  }
515  default:
516  {
517  switch (coordinateType_)
518  {
519  case coordinateType::volume:
520  {
521  mean += fi.x()*concentration/totalConcentration;
522 
523  break;
524  }
525  case coordinateType::area:
526  {
527  mean += fi.a()*concentration/totalConcentration;
528 
529  break;
530  }
531  case coordinateType::diameter:
532  {
533  mean += fi.d()*concentration/totalConcentration;
534 
535  break;
536  }
537  }
538 
539  break;
540  }
541  }
542  }
543 
544  if (meanType_ == meanType::geometric)
545  {
546  mean = exp(mean);
547 
548  setDimensions(mean, momentType::mean);
549  }
550 
551  return tMean;
552 }
553 
554 
556 Foam::functionObjects::populationBalanceMoments::variance
557 (
558  const diameterModels::populationBalanceModel& popBal
559 )
560 {
561  tmp<volScalarField> tVariance
562  (
564  (
565  "variance",
566  mesh_,
568  )
569  );
570 
571  volScalarField& variance = tVariance.ref();
572 
573  setDimensions(variance, momentType::variance);
574 
575  volScalarField totalConcentration(this->totalConcentration(popBal));
576  volScalarField mean(this->mean(popBal));
577 
578  forAll(popBal.sizeGroups(), i)
579  {
580  const Foam::diameterModels::sizeGroup& fi = popBal.sizeGroups()[i];
581 
582  volScalarField concentration(fi*fi.phase()/fi.x());
583 
584  switch (weightType_)
585  {
586  case weightType::volumeConcentration:
587  {
588  concentration *= fi.x();
589 
590  break;
591  }
592  case weightType::areaConcentration:
593  {
594  concentration *= fi.a();
595 
596  break;
597  }
598  default:
599  {
600  break;
601  }
602  }
603 
604  switch (meanType_)
605  {
606  case meanType::geometric:
607  {
608  switch (coordinateType_)
609  {
610  case coordinateType::volume:
611  {
612  variance +=
613  sqr(Foam::log(fi.x()/mean))
614  *concentration/totalConcentration;
615 
616  break;
617  }
618  case coordinateType::area:
619  {
620  variance +=
621  sqr(Foam::log(fi.a()/mean))
622  *concentration/totalConcentration;
623 
624  break;
625  }
626  case coordinateType::diameter:
627  {
628  variance +=
629  sqr(Foam::log(fi.d()/mean))
630  *concentration/totalConcentration;
631 
632  break;
633  }
634  }
635 
636  break;
637  }
638  default:
639  {
640  switch (coordinateType_)
641  {
642  case coordinateType::volume:
643  {
644  variance +=
645  sqr(fi.x() - mean)*concentration/totalConcentration;
646 
647  break;
648  }
649  case coordinateType::area:
650  {
651  variance +=
652  sqr(fi.a() - mean)*concentration/totalConcentration;
653 
654  break;
655  }
656  case coordinateType::diameter:
657  {
658  variance +=
659  sqr(fi.d() - mean)*concentration/totalConcentration;
660 
661  break;
662  }
663  }
664 
665  break;
666  }
667  }
668  }
669 
670  return tVariance;
671 }
672 
673 
675 Foam::functionObjects::populationBalanceMoments::stdDev
676 (
677  const diameterModels::populationBalanceModel& popBal
678 )
679 {
680  switch (meanType_)
681  {
682  case meanType::geometric:
683  {
684  return exp(sqrt(this->variance(popBal)));
685  }
686  default:
687  {
688  return sqrt(this->variance(popBal));
689  }
690  }
691 }
692 
693 
694 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
695 
697 (
698  const word& name,
699  const Time& runTime,
700  const dictionary& dict
701 )
702 :
703  fvMeshFunctionObject(name, runTime, dict),
704  popBalName_(dict.lookup("populationBalance")),
705  momentType_(momentTypeNames_.read(dict.lookup("momentType"))),
706  coordinateType_(coordinateTypeNames_.read(dict.lookup("coordinateType"))),
707  weightType_
708  (
709  dict.found("weightType")
710  ? weightTypeNames_.read(dict.lookup("weightType"))
711  : weightType::numberConcentration
712  ),
713  meanType_(meanType::notApplicable),
714  order_(-1),
715  fldPtr_(nullptr)
716 {
717  read(dict);
718 }
719 
720 
721 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
722 
724 {}
725 
726 
727 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
728 
729 bool
731 {
733 
734  switch (momentType_)
735  {
736  case momentType::integerMoment:
737  {
738  order_ = dict.lookup<int>("order");
739 
740  break;
741  }
742  default:
743  {
744  meanType_ =
745  dict.found("meanType")
746  ? meanTypeNames_.read(dict.lookup("meanType"))
747  : meanType::arithmetic;
748 
749  break;
750  }
751  }
752 
753  switch (momentType_)
754  {
755  case momentType::integerMoment:
756  {
757  fldPtr_.set
758  (
759  new volScalarField
760  (
761  IOobject
762  (
763  this->integerMomentFldName(),
764  mesh_.time().name(),
765  mesh_,
768  ),
769  mesh_,
771  )
772  );
773 
774  volScalarField& integerMoment = fldPtr_();
775 
776  setDimensions(integerMoment, momentType::integerMoment);
777 
778  break;
779  }
780  case momentType::mean:
781  {
782  fldPtr_.set
783  (
784  new volScalarField
785  (
786  IOobject
787  (
788  this->defaultFldName(),
789  mesh_.time().name(),
790  mesh_,
793  ),
794  mesh_,
796  )
797  );
798 
799  setDimensions(fldPtr_(), momentType::mean);
800 
801  break;
802  }
803  case momentType::variance:
804  {
805  fldPtr_.set
806  (
807  new volScalarField
808  (
809  IOobject
810  (
811  this->defaultFldName(),
812  mesh_.time().name(),
813  mesh_,
816  ),
817  mesh_,
819  )
820  );
821 
822  setDimensions(fldPtr_(), momentType::variance);
823 
824  break;
825  }
826  case momentType::stdDev:
827  {
828  fldPtr_.set
829  (
830  new volScalarField
831  (
832  IOobject
833  (
834  this->defaultFldName(),
835  mesh_.time().name(),
836  mesh_,
839  ),
840  mesh_,
842  )
843  );
844 
845  setDimensions(fldPtr_(), momentType::stdDev);
846 
847  break;
848  }
849  }
850 
851  return true;
852 }
853 
854 
856 {
858  obr_.lookupObject<diameterModels::populationBalanceModel>
859  (
860  popBalName_
861  );
862 
863  switch (momentType_)
864  {
865  case momentType::integerMoment:
866  {
867  volScalarField& integerMoment = fldPtr_();
868 
869  integerMoment = Zero;
870 
871  forAll(popBal.sizeGroups(), i)
872  {
874  popBal.sizeGroups()[i];
875 
876  volScalarField concentration(fi*fi.phase()/fi.x());
877 
878  switch (weightType_)
879  {
880  case weightType::volumeConcentration:
881  {
882  concentration *= fi.x();
883 
884  break;
885  }
886  case weightType::areaConcentration:
887  {
888  concentration *= fi.a();
889 
890  break;
891  }
892  default:
893  {
894  break;
895  }
896  }
897 
898  switch (coordinateType_)
899  {
900  case coordinateType::volume:
901  {
902  integerMoment +=
903  pow(fi.x(), order_)*concentration;
904 
905  break;
906  }
907  case coordinateType::area:
908  {
909  integerMoment +=
910  pow(fi.a(), order_)*concentration;
911 
912  break;
913  }
914  case coordinateType::diameter:
915  {
916  integerMoment +=
917  pow(fi.d(), order_)*concentration;
918 
919  break;
920  }
921  }
922  }
923 
924  break;
925  }
926  case momentType::mean:
927  {
928  fldPtr_() = this->mean(popBal);
929 
930  break;
931  }
932  case momentType::variance:
933  {
934  fldPtr_() = this->variance(popBal);
935 
936  break;
937  }
938  case momentType::stdDev:
939  {
940  fldPtr_() = sqrt(this->variance(popBal));
941 
942  break;
943  }
944  }
945 
946  return true;
947 }
948 
949 
951 {
952  writeObject(fldPtr_->name());
953 
954  return true;
955 }
956 
957 
958 // ************************************************************************* //
bool found
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Macros for easy insertion into run-time selection tables.
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
static word groupName(Name name, const word &group)
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:55
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
Model for tracking the evolution of a dispersed phase size distribution due to coalescence (synonymou...
const UPtrList< sizeGroup > & sizeGroups() const
Return the size groups belonging to 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
const phaseModel & phase() const
Return const-reference to the phase.
Definition: sizeGroupI.H:37
const tmp< volScalarField > d() const
Return representative diameter of the sizeGroup.
Definition: sizeGroup.C:147
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Abstract base-class for Time/database functionObjects.
Specialisation of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
virtual bool read(const dictionary &)
Read optional controls.
Calculates and writes out integral (integer moments) or mean properties (mean, variance,...
static const NamedEnum< coordinateType, 3 > coordinateTypeNames_
Names of the coordinate types.
static const NamedEnum< meanType, 3 > meanTypeNames_
Names of the mean types.
static const NamedEnum< momentType, 4 > momentTypeNames_
Names of the moment types.
static const NamedEnum< weightType, 3 > weightTypeNames_
Names of the weight types.
coordinateType
Enumeration for the coordinate types.
populationBalanceMoments(const word &name, const Time &runTime, const dictionary &)
Construct from Time and dictionary.
virtual bool execute()
Calculate the moment fields.
virtual bool write()
Write the moment fields.
virtual bool read(const dictionary &)
Read the data.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
gmvFile<< "tracers "<< particles.size()<< nl;{ pointField positions(particles.size());label particlei=0;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter) { positions[particlei++]=iter().position(mesh);} for(i=0;i< pTraits< point >::nComponents;i++) { forAll(positions, particlei) { gmvFile<< component(positions[particlei], i)<< ' ';} gmvFile<< nl;}}forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.name(), lagrangian::cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
defineTypeNameAndDebug(adjustTimeStepToCombustion, 0)
addToRunTimeSelectionTable(functionObject, adjustTimeStepToCombustion, dictionary)
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
dimensionedScalar exp(const dimensionedScalar &ds)
bool read(const char *, int32_t &)
Definition: int32IO.C:85
const dimensionSet dimless
const dimensionSet dimLength
dimensionedScalar log(const dimensionedScalar &ds)
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
const dimensionSet dimVolume
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:62
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
const dimensionSet dimArea
void inv(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dictionary dict