regionSizeDistribution.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2013-2016 OpenFOAM Foundation
9  Copyright (C) 2016-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "regionSizeDistribution.H"
30 #include "regionSplit.H"
31 #include "volFields.H"
32 #include "fvcVolumeIntegrate.H"
33 #include "mathematicalConstants.H"
35 
36 using namespace Foam::constant;
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42  namespace functionObjects
43  {
44  defineTypeNameAndDebug(regionSizeDistribution, 0);
46  (
47  functionObject,
48  regionSizeDistribution,
49  dictionary
50  );
51  }
52 }
53 
54 
55 // * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
56 
57 namespace Foam
58 {
59 
60 template<class Type>
61 static Map<Type> regionSum(const regionSplit& regions, const Field<Type>& fld)
62 {
63  // Per region the sum of fld
64  Map<Type> regionToSum(regions.nRegions()/UPstream::nProcs());
65 
66  forAll(fld, celli)
67  {
68  const label regioni = regions[celli];
69  regionToSum(regioni, Type(Zero)) += fld[celli];
70  }
71 
73 
74  return regionToSum;
75 }
76 
77 
78 static Map<label> regionSum(const regionSplit& regions, const label nCells)
79 {
80  // Per region the sum of fld
81  Map<label> regionToSum(regions.nRegions()/UPstream::nProcs());
82 
83  for (label celli = 0; celli < nCells; ++celli)
84  {
85  const label regioni = regions[celli];
86  ++regionToSum(regioni);
87  }
88 
90 
91  return regionToSum;
92 }
93 
94 
95 template<class Type>
96 static List<Type> extractData(const labelUList& keys, const Map<Type>& regionData)
97 {
98  List<Type> sortedData(keys.size());
99 
100  forAll(keys, i)
101  {
102  sortedData[i] = regionData[keys[i]];
103  }
104  return sortedData;
105 }
106 
107 } // End namespace Foam
108 
109 
110 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
111 
112 void Foam::functionObjects::regionSizeDistribution::writeAlphaFields
113 (
114  const regionSplit& regions,
115  const labelHashSet& keepRegions,
116  const Map<scalar>& regionVolume,
117  const volScalarField& alpha
118 ) const
119 {
120  const scalar maxDropletVol = 1.0/6.0*mathematical::pi*pow3(maxDiam_);
121 
122  // Split alpha field
123  // ~~~~~~~~~~~~~~~~~
124  // Split into
125  // - liquidCore : region connected to inlet patches
126  // - per region a volume : for all other regions
127  // - backgroundAlpha : remaining alpha
128 
129 
130  // Construct field
131  volScalarField liquidCore
132  (
133  IOobject
134  (
135  scopedName(alphaName_ + "_liquidCore"),
136  obr_.time().timeName(),
137  obr_,
139  ),
140  alpha,
142  );
143 
144  volScalarField backgroundAlpha
145  (
146  IOobject
147  (
148  scopedName(alphaName_ + "_background"),
149  obr_.time().timeName(),
150  obr_,
152  ),
153  alpha,
155  );
156 
157 
158  // Knock out any cell not in keepRegions (patch regions)
159  forAll(liquidCore, celli)
160  {
161  const label regioni = regions[celli];
162  if (keepRegions.found(regioni))
163  {
164  backgroundAlpha[celli] = 0;
165  }
166  else
167  {
168  liquidCore[celli] = 0;
169 
170  const scalar regionVol = regionVolume[regioni];
171  if (regionVol < maxDropletVol)
172  {
173  backgroundAlpha[celli] = 0;
174  }
175  }
176  }
177  liquidCore.correctBoundaryConditions();
178  backgroundAlpha.correctBoundaryConditions();
179 
180  if (log)
181  {
182  Info<< " Volume of liquid-core = "
183  << fvc::domainIntegrate(liquidCore).value()
184  << nl
185  << " Volume of background = "
186  << fvc::domainIntegrate(backgroundAlpha).value()
187  << endl;
188  }
189 
190  Log << " Writing liquid-core field to " << liquidCore.name() << endl;
191  liquidCore.write();
192 
193  Log<< " Writing background field to " << backgroundAlpha.name() << endl;
194  backgroundAlpha.write();
195 }
196 
197 
199 Foam::functionObjects::regionSizeDistribution::findPatchRegions
200 (
201  const regionSplit& regions
202 ) const
203 {
204  // Mark all regions starting at patches
205  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
206 
207  labelHashSet patchRegions(2*regions.nRegions());
208 
209  labelHashSet patchSet(mesh_.boundaryMesh().patchSet(patchNames_));
210 
211  for (const label patchi : patchSet)
212  {
213  const polyPatch& pp = mesh_.boundaryMesh()[patchi];
214 
215  // All regions connected to the patch
216  for (const label celli : pp.faceCells())
217  {
218  patchRegions.insert(regions[celli]);
219  }
220  }
221 
222  // Ensure all processors have the same set of regions
223  Pstream::combineReduce(patchRegions, plusEqOp<labelHashSet>());
224 
225  return patchRegions;
226 }
227 
228 
230 Foam::functionObjects::regionSizeDistribution::divide
231 (
232  const scalarField& num,
233  const scalarField& denom
234 )
235 {
236  auto tresult = tmp<scalarField>::New(num.size(), Zero);
237  auto& result = tresult.ref();
238 
239  forAll(denom, i)
240  {
241  if (ROOTVSMALL < Foam::mag(denom[i]))
242  {
243  result[i] = num[i]/denom[i];
244  }
245  }
246  return tresult;
247 }
248 
249 
250 void Foam::functionObjects::regionSizeDistribution::writeGraphs
251 (
252  const word& fieldName, // name of field
253  const scalarField& sortedField, // per region field data
254 
255  const labelList& indices, // index of bin for each region
256  const scalarField& binCount, // per bin number of regions
257  const coordSet& coords // graph data for bins
258 ) const
259 {
260  if (UPstream::master())
261  {
262  // Calculate per-bin average
263  scalarField binSum(nBins_, Zero);
264  forAll(sortedField, i)
265  {
266  binSum[indices[i]] += sortedField[i];
267  }
268 
269  scalarField binAvg(divide(binSum, binCount));
270 
271  // Per bin deviation
272  scalarField binSqrSum(nBins_, Zero);
273  forAll(sortedField, i)
274  {
275  binSqrSum[indices[i]] += Foam::sqr(sortedField[i]);
276  }
277  scalarField binDev
278  (
279  sqrt(divide(binSqrSum, binCount) - Foam::sqr(binAvg))
280  );
281 
282 
283  auto& writer = formatterPtr_();
284 
285  word outputName;
286  if (writer.buffering())
287  {
288  outputName =
289  (
290  coords.name()
292  (
293  wordList
294  ({
295  fieldName + "_sum",
296  fieldName + "_avg",
297  fieldName + "_dev"
298  })
299  )
300  );
301  }
302  else
303  {
304  outputName = coords.name();
305  }
306 
307  writer.open
308  (
309  coords,
310  (baseTimeDir() / outputName)
311  );
312 
313  Log << " Writing distribution of "
314  << fieldName << " to " << writer.path() << endl;
315 
316  writer.write(fieldName + "_sum", binSum);
317  writer.write(fieldName + "_avg", binAvg);
318  writer.write(fieldName + "_dev", binDev);
319  writer.close(true);
320  }
321 }
322 
323 
324 void Foam::functionObjects::regionSizeDistribution::writeGraphs
325 (
326  const word& fieldName, // name of field
327  const scalarField& cellField, // per cell field data
328  const regionSplit& regions, // per cell the region(=droplet)
329  const labelList& sortedRegions, // valid regions in sorted order
330  const scalarField& sortedNormalisation,
331 
332  const labelList& indices, // per region index of bin
333  const scalarField& binCount, // per bin number of regions
334  const coordSet& coords // graph data for bins
335 ) const
336 {
337  // Sum on a per-region basis. Parallel reduced.
338  Map<scalar> regionField(regionSum(regions, cellField));
339 
340  // Extract in region order
341  scalarField sortedField
342  (
343  sortedNormalisation
344  * extractData(sortedRegions, regionField)
345  );
346 
347  writeGraphs
348  (
349  fieldName, // name of field
350  sortedField, // per region field data
351 
352  indices, // index of bin for each region
353  binCount, // per bin number of regions
354  coords // graph data for bins
355  );
356 }
357 
358 
359 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
360 
362 (
363  const word& name,
364  const Time& runTime,
365  const dictionary& dict
366 )
367 :
369  writeFile(obr_, name),
370  alphaName_(dict.get<word>("field")),
371  patchNames_(dict.get<wordRes>("patches")),
372  isoPlanes_(dict.getOrDefault("isoPlanes", false))
373 {
374  read(dict);
375 }
376 
377 
378 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
379 
381 {
384 
385  dict.readEntry("nBins", nBins_);
386  dict.readEntry("field", alphaName_);
387  dict.readEntry("threshold", threshold_);
388  dict.readEntry("maxDiameter", maxDiam_);
389  minDiam_ = 0.0;
390  dict.readIfPresent("minDiameter", minDiam_);
391  dict.readEntry("patches", patchNames_);
392  dict.readEntry("fields", fields_);
393 
394  const word setFormat(dict.get<word>("setFormat"));
395  formatterPtr_ = coordSetWriter::New
396  (
397  setFormat,
398  dict.subOrEmptyDict("formatOptions").optionalSubDict(setFormat)
399  );
400 
401  csysPtr_ = coordinateSystem::NewIfPresent(obr_, dict);
402 
403  if (csysPtr_)
404  {
405  Info<< "Transforming all vectorFields with coordinate system "
406  << csysPtr_->name() << endl;
407  }
408 
409  if (isoPlanes_)
410  {
411  dict.readEntry("origin", origin_);
412  dict.readEntry("direction", direction_);
413  dict.readEntry("maxD", maxDiameter_);
414  dict.readEntry("nDownstreamBins", nDownstreamBins_);
415  dict.readEntry("maxDownstream", maxDownstream_);
416  direction_.normalise();
417  }
418 
419  return true;
420 }
421 
424 {
425  return true;
426 }
427 
428 
430 {
431  Log << type() << " " << name() << " write:" << nl;
432 
433  tmp<volScalarField> talpha;
434  talpha.cref(obr_.cfindObject<volScalarField>(alphaName_));
435  if (talpha)
436  {
437  Log << " Looking up field " << alphaName_ << endl;
438  }
439  else
440  {
441  Info<< " Reading field " << alphaName_ << endl;
442  talpha.reset
443  (
444  new volScalarField
445  (
446  IOobject
447  (
448  alphaName_,
449  mesh_.time().timeName(),
450  mesh_,
453  ),
454  mesh_
455  )
456  );
457  }
458  const auto& alpha = talpha();
459 
460  Log << " Volume of alpha = "
461  << fvc::domainIntegrate(alpha).value()
462  << endl;
463 
464  const scalar meshVol = gSum(mesh_.V());
465  const scalar maxDropletVol = 1.0/6.0*mathematical::pi*pow3(maxDiam_);
466  const scalar delta = (maxDiam_-minDiam_)/nBins_;
467 
468  Log << " Mesh volume = " << meshVol << nl
469  << " Maximum droplet diameter = " << maxDiam_ << nl
470  << " Maximum droplet volume = " << maxDropletVol
471  << endl;
472 
473 
474  // Determine blocked faces
475  boolList blockedFace(mesh_.nFaces(), false);
476  // label nBlocked = 0;
477 
478  {
479  for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
480  {
481  scalar ownVal = alpha[mesh_.faceOwner()[facei]];
482  scalar neiVal = alpha[mesh_.faceNeighbour()[facei]];
483 
484  if
485  (
486  (ownVal < threshold_ && neiVal > threshold_)
487  || (ownVal > threshold_ && neiVal < threshold_)
488  )
489  {
490  blockedFace[facei] = true;
491  // ++nBlocked;
492  }
493  }
494 
495  // Block coupled faces
496  forAll(alpha.boundaryField(), patchi)
497  {
498  const fvPatchScalarField& fvp = alpha.boundaryField()[patchi];
499  if (fvp.coupled())
500  {
501  tmp<scalarField> townFld(fvp.patchInternalField());
502  tmp<scalarField> tnbrFld(fvp.patchNeighbourField());
503  const auto& ownFld = townFld();
504  const auto& nbrFld = tnbrFld();
505 
506  label start = fvp.patch().patch().start();
507 
508  forAll(ownFld, i)
509  {
510  scalar ownVal = ownFld[i];
511  scalar neiVal = nbrFld[i];
512 
513  if
514  (
515  (ownVal < threshold_ && neiVal > threshold_)
516  || (ownVal > threshold_ && neiVal < threshold_)
517  )
518  {
519  blockedFace[start+i] = true;
520  // ++nBlocked;
521  }
522  }
523  }
524  }
525  }
526 
527 
528  regionSplit regions(mesh_, blockedFace);
529 
530  Log << " Determined " << regions.nRegions()
531  << " disconnected regions" << endl;
532 
533 
534  if (debug)
535  {
536  volScalarField region
537  (
538  IOobject
539  (
540  "region",
541  mesh_.time().timeName(),
542  mesh_,
546  ),
547  mesh_,
549  );
550  Info<< " Dumping region as volScalarField to "
551  << region.name() << endl;
552 
553  forAll(regions, celli)
554  {
555  region[celli] = regions[celli];
556  }
557  region.correctBoundaryConditions();
558  region.write();
559  }
560 
561 
562  // Determine regions connected to supplied patches
563  const labelHashSet patchRegions(findPatchRegions(regions));
564 
565  // Sum all regions
566  const scalarField alphaVol(alpha.primitiveField()*mesh_.V());
567  Map<scalar> allRegionVolume(regionSum(regions, mesh_.V()));
568  Map<scalar> allRegionAlphaVolume(regionSum(regions, alphaVol));
569  Map<label> allRegionNumCells(regionSum(regions, mesh_.nCells()));
570 
571  if (debug)
572  {
573  Info<< " " << token::TAB << "Region"
574  << token::TAB << "Volume(mesh)"
575  << token::TAB << "Volume(" << alpha.name() << "):"
576  << token::TAB << "nCells"
577  << nl;
578  scalar meshSumVol = 0.0;
579  scalar alphaSumVol = 0.0;
580  label nCells = 0;
581 
582  auto vIter = allRegionVolume.cbegin();
583  auto aIter = allRegionAlphaVolume.cbegin();
584  auto numIter = allRegionNumCells.cbegin();
585  for
586  (
587  ;
588  vIter.good() && aIter.good();
589  ++vIter, ++aIter, ++numIter
590  )
591  {
592  Info<< " " << token::TAB << vIter.key()
593  << token::TAB << vIter()
594  << token::TAB << aIter()
595  << token::TAB << numIter()
596  << nl;
597 
598  meshSumVol += vIter();
599  alphaSumVol += aIter();
600  nCells += numIter();
601  }
602  Info<< " " << token::TAB << "Total:"
603  << token::TAB << meshSumVol
604  << token::TAB << alphaSumVol
605  << token::TAB << nCells
606  << endl;
607  }
608 
609 
610  if (log)
611  {
612  Info<< " Patch connected regions (liquid core):" << nl;
613  Info<< token::TAB << " Region"
614  << token::TAB << "Volume(mesh)"
615  << token::TAB << "Volume(" << alpha.name() << "):"
616  << nl;
617 
618  for (const label regioni : patchRegions.sortedToc())
619  {
620  Info<< " " << token::TAB << regioni
621  << token::TAB << allRegionVolume[regioni]
622  << token::TAB << allRegionAlphaVolume[regioni] << nl;
623 
624  }
625  Info<< endl;
626  }
627 
628  if (log)
629  {
630  Info<< " Background regions:" << nl;
631  Info<< " " << token::TAB << "Region"
632  << token::TAB << "Volume(mesh)"
633  << token::TAB << "Volume(" << alpha.name() << "):"
634  << nl;
635 
636  auto vIter = allRegionVolume.cbegin();
637  auto aIter = allRegionAlphaVolume.cbegin();
638 
639  for
640  (
641  ;
642  vIter.good() && aIter.good();
643  ++vIter, ++aIter
644  )
645  {
646  if
647  (
648  !patchRegions.found(vIter.key())
649  && vIter() >= maxDropletVol
650  )
651  {
652  Info<< " " << token::TAB << vIter.key()
653  << token::TAB << vIter()
654  << token::TAB << aIter() << nl;
655  }
656  }
657  Info<< endl;
658  }
659 
660 
661 
662  // Split alpha field
663  // ~~~~~~~~~~~~~~~~~
664  // Split into
665  // - liquidCore : region connected to inlet patches
666  // - per region a volume : for all other regions
667  // - backgroundAlpha : remaining alpha
668  writeAlphaFields(regions, patchRegions, allRegionVolume, alpha);
669 
670 
671  // Extract droplet-only allRegionVolume, i.e. delete liquid core
672  // (patchRegions) and background regions from maps.
673  // Note that we have to use mesh volume (allRegionVolume) and not
674  // allRegionAlphaVolume since background might not have alpha in it.
675  // Deleting regions where the volume-alpha-weighted is lower than
676  // threshold
677  forAllIters(allRegionVolume, vIter)
678  {
679  const label regioni = vIter.key();
680  if
681  (
682  patchRegions.found(regioni)
683  || vIter() >= maxDropletVol
684  || (allRegionAlphaVolume[regioni]/vIter() < threshold_)
685  )
686  {
687  allRegionVolume.erase(vIter);
688  allRegionAlphaVolume.erase(regioni);
689  allRegionNumCells.erase(regioni);
690  }
691  }
692 
693  if (allRegionVolume.size())
694  {
695  // Construct mids of bins for plotting
696  pointField xBin(nBins_, Zero);
697 
698  {
699  scalar x = 0.5*delta;
700  for (point& p : xBin)
701  {
702  p.x() = x;
703  x += delta;
704  }
705  }
706 
707  const coordSet coords("diameter", "x", xBin, mag(xBin));
708 
709 
710  // Get in region order the alpha*volume and diameter
711  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
712 
713  const labelList sortedRegions = allRegionAlphaVolume.sortedToc();
714 
715  scalarField sortedVols
716  (
717  extractData(sortedRegions, allRegionAlphaVolume)
718  );
719 
720  vectorField centroids(sortedVols.size(), Zero);
721 
722  // Check if downstream bins are calculated
723  if (isoPlanes_)
724  {
725  vectorField alphaDistance
726  (
727  (alpha.primitiveField()*mesh_.V())
728  *(mesh_.C().primitiveField() - origin_)()
729  );
730 
731  Map<vector> allRegionAlphaDistance
732  (
733  regionSum
734  (
735  regions,
736  alphaDistance
737  )
738  );
739 
740  // 2. centroid
741  vectorField sortedMoment
742  (
743  extractData(sortedRegions, allRegionAlphaDistance)
744  );
745 
746  centroids = sortedMoment/sortedVols + origin_;
747 
748  // Bin according to centroid
749  scalarField distToPlane((centroids - origin_) & direction_);
750 
751  vectorField radialDistToOrigin
752  (
753  (centroids - origin_) - (distToPlane*direction_)
754  );
755 
756  const scalar deltaX = maxDownstream_/nDownstreamBins_;
757  labelList downstreamIndices(distToPlane.size(), -1);
758  forAll(distToPlane, i)
759  {
760  if
761  (
762  (mag(radialDistToOrigin[i]) < maxDiameter_)
763  && (distToPlane[i] < maxDownstream_)
764  )
765  {
766  downstreamIndices[i] = distToPlane[i]/deltaX;
767  }
768  }
769 
770  scalarField binDownCount(nDownstreamBins_, Zero);
771  forAll(distToPlane, i)
772  {
773  if (downstreamIndices[i] != -1)
774  {
775  binDownCount[downstreamIndices[i]] += 1.0;
776  }
777  }
778 
779  // Write
780  if (UPstream::master())
781  {
782  // Construct mids of bins for plotting
783  pointField xBin(nDownstreamBins_, Zero);
784 
785  {
786  scalar x = 0.5*deltaX;
787  for (point& p : xBin)
788  {
789  p.x() = x;
790  x += deltaX;
791  }
792  }
793 
794  const coordSet coords("distance", "x", xBin, mag(xBin));
795 
796  auto& writer = formatterPtr_();
797  writer.nFields(1);
798 
799  writer.open
800  (
801  coords,
802  writeFile::baseTimeDir() / (coords.name() + "_isoPlanes")
803  );
804 
805  writer.write("isoPlanes", binDownCount);
806  writer.close(true);
807  }
808 
809  // Write to log
810  if (log)
811  {
812  Info<< " Iso-planes Bins:" << nl
813  << " " << token::TAB << "Bin"
814  << token::TAB << "Min distance"
815  << token::TAB << "Count:"
816  << nl;
817 
818  scalar delta = 0.0;
819  forAll(binDownCount, bini)
820  {
821  Info<< " " << token::TAB << bini
822  << token::TAB << delta
823  << token::TAB << binDownCount[bini] << nl;
824  delta += deltaX;
825  }
826  Info<< endl;
827 
828  }
829  }
830 
831  // Calculate the diameters
832  scalarField sortedDiameters(sortedVols.size());
833  forAll(sortedDiameters, i)
834  {
835  sortedDiameters[i] = Foam::cbrt
836  (
837  sortedVols[i]
839  );
840  }
841 
842  // Determine the bin index for all the diameters
843  labelList indices(sortedDiameters.size());
844  forAll(sortedDiameters, i)
845  {
846  indices[i] = (sortedDiameters[i]-minDiam_)/delta;
847  }
848 
849  // Calculate the counts per diameter bin
850  scalarField binCount(nBins_, Zero);
851  forAll(sortedDiameters, i)
852  {
853  binCount[indices[i]] += 1.0;
854  }
855 
856  // Write counts
857  if (UPstream::master())
858  {
859  auto& writer = formatterPtr_();
860  writer.nFields(1);
861 
862  writer.open
863  (
864  coords,
865  writeFile::baseTimeDir() / (coords.name() + "_count")
866  );
867 
868  writer.write("count", binCount);
869  writer.close(true);
870  }
871 
872  // Write to log
873  if (log)
874  {
875  Info<< " Bins:" << nl
876  << " " << token::TAB << "Bin"
877  << token::TAB << "Min diameter"
878  << token::TAB << "Count:"
879  << nl;
880 
881  scalar diam = 0.0;
882  forAll(binCount, bini)
883  {
884  Info<< " " << token::TAB << bini
885  << token::TAB << diam
886  << token::TAB << binCount[bini] << nl;
887 
888  diam += delta;
889  }
890 
891  Info<< endl;
892  }
893 
894 
895  // Write average and deviation of droplet volume.
896  writeGraphs
897  (
898  "volume", // name of field
899  sortedVols, // per region field data
900 
901  indices, // per region the bin index
902  binCount, // per bin number of regions
903  coords // graph data for bins
904  );
905 
906  // Collect some more fields
907  {
908  for
909  (
910  const volScalarField& vfield
911  : obr_.csorted<volScalarField>(fields_)
912  )
913  {
914  const word& fldName = vfield.name();
915 
916  Log << " Scalar field " << fldName << endl;
917 
918  tmp<Field<scalar>> tfld(vfield.primitiveField());
919  const auto& fld = tfld();
920 
921  writeGraphs
922  (
923  fldName, // name of field
924  alphaVol*fld, // per cell field data
925 
926  regions, // per cell the region(=droplet)
927  sortedRegions, // valid regions in sorted order
928  1.0/sortedVols, // per region normalisation
929 
930  indices, // index of bin for each region
931  binCount, // per bin number of regions
932  coords // graph data for bins
933  );
934  }
935  }
936 
937  {
938  for
939  (
940  const volVectorField& vfield
941  : obr_.csorted<volVectorField>(fields_)
942  )
943  {
944  const word& fldName = vfield.name();
945 
946  Log << " Vector field " << fldName << endl;
947 
948  tmp<Field<vector>> tfld(vfield.primitiveField());
949 
950  if (csysPtr_)
951  {
952  Log << "Transforming vector field " << fldName
953  << " with coordinate system "
954  << csysPtr_->name() << endl;
955 
956  tfld = csysPtr_->localVector(tfld());
957  }
958  const auto& fld = tfld();
959 
960  // Components
961 
962  for (direction cmpt = 0; cmpt < vector::nComponents; ++cmpt)
963  {
964  writeGraphs
965  (
966  fldName + vector::componentNames[cmpt],
967  alphaVol*fld.component(cmpt),// per cell field data
968 
969  regions, // per cell the region(=droplet)
970  sortedRegions, // valid regions in sorted order
971  1.0/sortedVols, // per region normalisation
972 
973  indices, // index of bin for each region
974  binCount, // per bin number of regions
975  coords // graph data for bins
976  );
977  }
978 
979  // Magnitude
980  writeGraphs
981  (
982  fldName + "mag", // name of field
983  alphaVol*mag(fld), // per cell field data
984 
985  regions, // per cell the region(=droplet)
986  sortedRegions, // valid regions in sorted order
987  1.0/sortedVols, // per region normalisation
988 
989  indices, // index of bin for each region
990  binCount, // per bin number of regions
991  coords // graph data for bins
992  );
993  }
994  }
995  }
996 
997  return true;
998 }
999 
1000 
1001 // ************************************************************************* //
Different types of constants.
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
scalar delta
This class separates the mesh into distinct unconnected regions, each of which is then given a label ...
Definition: regionSplit.H:136
dictionary dict
void divide(DimensionedField< Type, GeoMesh > &result, const DimensionedField< Type, GeoMesh > &f1, const DimensionedField< scalar, GeoMesh > &f2)
virtual bool write()
Calculate the regionSizeDistribution and write.
uint8_t direction
Definition: direction.H:46
dimensionedScalar log(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
vtk::lineWriter writer(edgeCentres, edgeList::null(), fileName(aMesh.time().globalPath()/"finiteArea-edgesCentres"))
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
Tab [isspace].
Definition: token.H:129
dimensionedSymmTensor sqr(const dimensionedVector &dv)
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
engineTime & runTime
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
const T & cref() const
Return const reference to the object or to the contents of a (non-null) managed pointer.
Definition: tmpI.H:221
dimensioned< Type > domainIntegrate(const GeometricField< Type, fvPatchField, volMesh > &vf)
::Foam::direction nComponents(const expressions::valueTypeCode) noexcept
The number of components associated with given valueTypeCode.
Definition: exprTraits.C:40
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
regionSizeDistribution(const word &name, const Time &runTime, const dictionary &)
Construct for given objectRegistry and dictionary.
Macros for easy insertion into run-time selection tables.
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:127
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
static word suffix(const word &fldName, const word &fileExt=word::null)
Name suffix based on fieldName (underscore separator)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
static autoPtr< coordSetWriter > New(const word &writeFormat)
Return a reference to the selected writer.
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run. ...
Definition: UPstream.H:1065
word outputName("finiteArea-edges.obj")
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
Type gSum(const FieldField< Field, Type > &f)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
Generic templated field type.
Definition: Field.H:62
fvPatchField< scalar > fvPatchScalarField
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static Map< Type > regionSum(const regionSplit &regions, const Field< Type > &fld)
dimensionedScalar cbrt(const dimensionedScalar &ds)
static List< Type > extractData(const labelUList &keys, const Map< Type > &regionData)
#define forAllIters(container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:336
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:53
constexpr scalar pi(M_PI)
Volume integrate volField creating a volField.
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
virtual bool read(const dictionary &)
Read the regionSizeDistribution data.
static void combineReduce(const List< commsStruct > &comms, T &value, const CombineOp &cop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce) applying cop to inplace combine value from different processors...
void read(Istream &, label &val, const dictionary &)
In-place read with dictionary lookup.
int debug
Static debugging option.
defineTypeNameAndDebug(combustionModel, 0)
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const word & name() const noexcept
Return const reference to name.
List< word > wordList
List of word.
Definition: fileName.H:59
dimensionedScalar pow3(const dimensionedScalar &ds)
vector point
Point is a vector.
Definition: point.H:37
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1082
Nothing to be read.
static const word & calculatedType() noexcept
The type name for calculated patch fields.
Definition: fvPatchField.H:204
#define Log
Definition: PDRblock.C:28
label nRegions() const
Return total number of regions.
Definition: regionSplit.H:328
messageStream Info
Information stream (stdout output on master, null elsewhere)
Field< vector > vectorField
Specialisation of Field<T> for vector.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
word setFormat(propsDict.getOrDefault< word >("setFormat", "vtk"))
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
static void mapCombineReduce(Container &values, const CombineOp &cop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce) applying cop to inplace combine map values from different processo...
Base class for writing single files from the function objects.
Definition: writeFile.H:112
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
List< bool > boolList
A List of bools.
Definition: List.H:60
Do not request registration (bool: false)
void reset(tmp< T > &&other) noexcept
Clear existing and transfer ownership.
Definition: tmpI.H:338
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
A HashTable to objects of type <T> with a label key.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127