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-2020 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 "fvcVolumeIntegrate.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  namespace functionObjects
38  {
39  defineTypeNameAndDebug(regionSizeDistribution, 0);
41  (
42  functionObject,
43  regionSizeDistribution,
44  dictionary
45  );
46  }
47 }
48 
49 
50 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51 
52 void Foam::functionObjects::regionSizeDistribution::writeGraph
53 (
54  const coordSet& coords,
55  const word& valueName,
56  const scalarField& values
57 ) const
58 {
59  const wordList valNames(1, valueName);
60 
61  fileName outputPath = baseTimeDir();
62  mkDir(outputPath);
63 
64  OFstream str(outputPath/formatterPtr_().getFileName(coords, valNames));
65 
66  Log << " Writing distribution of " << valueName << " to " << str.name()
67  << endl;
68 
69  List<const scalarField*> valPtrs(1);
70  valPtrs[0] = &values;
71  formatterPtr_().write(coords, valNames, valPtrs, str);
72 }
73 
74 
75 void Foam::functionObjects::regionSizeDistribution::writeAlphaFields
76 (
77  const regionSplit& regions,
78  const Map<label>& patchRegions,
79  const Map<scalar>& regionVolume,
80  const volScalarField& alpha
81 ) const
82 {
83  const scalar maxDropletVol = 1.0/6.0*pow3(maxDiam_);
84 
85  // Split alpha field
86  // ~~~~~~~~~~~~~~~~~
87  // Split into
88  // - liquidCore : region connected to inlet patches
89  // - per region a volume : for all other regions
90  // - backgroundAlpha : remaining alpha
91 
92 
93  // Construct field
94  volScalarField liquidCore
95  (
96  IOobject
97  (
98  scopedName(alphaName_ + "_liquidCore"),
99  obr_.time().timeName(),
100  obr_,
102  ),
103  alpha,
105  );
106 
107  volScalarField backgroundAlpha
108  (
109  IOobject
110  (
111  scopedName(alphaName_ + "_background"),
112  obr_.time().timeName(),
113  obr_,
115  ),
116  alpha,
118  );
119 
120 
121  // Knock out any cell not in patchRegions
122  forAll(liquidCore, celli)
123  {
124  const label regioni = regions[celli];
125  if (patchRegions.found(regioni))
126  {
127  backgroundAlpha[celli] = 0;
128  }
129  else
130  {
131  liquidCore[celli] = 0;
132 
133  const scalar regionVol = regionVolume[regioni];
134  if (regionVol < maxDropletVol)
135  {
136  backgroundAlpha[celli] = 0;
137  }
138  }
139  }
140  liquidCore.correctBoundaryConditions();
141  backgroundAlpha.correctBoundaryConditions();
142 
143  if (log)
144  {
145  Info<< " Volume of liquid-core = "
146  << fvc::domainIntegrate(liquidCore).value()
147  << nl
148  << " Volume of background = "
149  << fvc::domainIntegrate(backgroundAlpha).value()
150  << endl;
151  }
152 
153  Log << " Writing liquid-core field to " << liquidCore.name() << endl;
154  liquidCore.write();
155 
156  Log<< " Writing background field to " << backgroundAlpha.name() << endl;
157  backgroundAlpha.write();
158 }
159 
160 
162 Foam::functionObjects::regionSizeDistribution::findPatchRegions
163 (
164  const regionSplit& regions
165 ) const
166 {
167  // Mark all regions starting at patches
168  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
169 
170  // Count number of patch faces (just for initial sizing)
171  const labelHashSet patchIDs(mesh_.boundaryMesh().patchSet(patchNames_));
172 
173  label nPatchFaces = 0;
174  for (const label patchi : patchIDs)
175  {
176  nPatchFaces += mesh_.boundaryMesh()[patchi].size();
177  }
178 
179 
180  Map<label> patchRegions(nPatchFaces);
181  for (const label patchi : patchIDs)
182  {
183  const polyPatch& pp = mesh_.boundaryMesh()[patchi];
184 
185  // Collect all regions on the patch
186  const labelList& faceCells = pp.faceCells();
187 
188  for (const label celli : faceCells)
189  {
190  patchRegions.insert
191  (
192  regions[celli],
193  Pstream::myProcNo() // dummy value
194  );
195  }
196  }
197 
198 
199  // Make sure all the processors have the same set of regions
200  Pstream::mapCombineGather(patchRegions, minEqOp<label>());
201  Pstream::mapCombineScatter(patchRegions);
202 
203  return patchRegions;
204 }
205 
206 
208 Foam::functionObjects::regionSizeDistribution::divide
209 (
210  const scalarField& num,
211  const scalarField& denom
212 )
213 {
214  auto tresult = tmp<scalarField>::New(num.size());
215  auto& result = tresult.ref();
216 
217  forAll(denom, i)
218  {
219  if (denom[i] != 0)
220  {
221  result[i] = num[i]/denom[i];
222  }
223  else
224  {
225  result[i] = 0.0;
226  }
227  }
228  return tresult;
229 }
230 
231 
232 void Foam::functionObjects::regionSizeDistribution::writeGraphs
233 (
234  const word& fieldName, // name of field
235  const labelList& indices, // index of bin for each region
236  const scalarField& sortedField, // per region field data
237  const scalarField& binCount, // per bin number of regions
238  const coordSet& coords // graph data for bins
239 ) const
240 {
241  if (Pstream::master())
242  {
243  // Calculate per-bin average
244  scalarField binSum(nBins_, Zero);
245  forAll(sortedField, i)
246  {
247  binSum[indices[i]] += sortedField[i];
248  }
249 
250  scalarField binAvg(divide(binSum, binCount));
251 
252  // Per bin deviation
253  scalarField binSqrSum(nBins_, Zero);
254  forAll(sortedField, i)
255  {
256  binSqrSum[indices[i]] += Foam::sqr(sortedField[i]);
257  }
258  scalarField binDev
259  (
260  sqrt(divide(binSqrSum, binCount) - Foam::sqr(binAvg))
261  );
262 
263  // Write average
264  writeGraph(coords, fieldName + "_sum", binSum);
265  // Write average
266  writeGraph(coords, fieldName + "_avg", binAvg);
267  // Write deviation
268  writeGraph(coords, fieldName + "_dev", binDev);
269  }
270 }
271 
272 
273 void Foam::functionObjects::regionSizeDistribution::writeGraphs
274 (
275  const word& fieldName, // name of field
276  const scalarField& cellField, // per cell field data
277  const regionSplit& regions, // per cell the region(=droplet)
278  const labelList& sortedRegions, // valid regions in sorted order
279  const scalarField& sortedNormalisation,
280 
281  const labelList& indices, // per region index of bin
282  const scalarField& binCount, // per bin number of regions
283  const coordSet& coords // graph data for bins
284 ) const
285 {
286  // Sum on a per-region basis. Parallel reduced.
287  Map<scalar> regionField(regionSum(regions, cellField));
288 
289  // Extract in region order
290  scalarField sortedField
291  (
292  sortedNormalisation
293  * extractData
294  (
295  sortedRegions,
296  regionField
297  )
298  );
299 
300  writeGraphs
301  (
302  fieldName, // name of field
303  indices, // index of bin for each region
304  sortedField, // per region field data
305  binCount, // per bin number of regions
306  coords // graph data for bins
307  );
308 }
309 
310 
311 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
312 
314 (
315  const word& name,
316  const Time& runTime,
317  const dictionary& dict
318 )
319 :
321  writeFile(obr_, name),
322  alphaName_(dict.get<word>("field")),
323  patchNames_(dict.get<wordRes>("patches")),
324  isoPlanes_(dict.getOrDefault("isoPlanes", false))
325 {
326  read(dict);
327 }
328 
329 
330 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
331 
333 {
336 
337  dict.readEntry("nBins", nBins_);
338  dict.readEntry("field", alphaName_);
339  dict.readEntry("threshold", threshold_);
340  dict.readEntry("maxDiameter", maxDiam_);
341  minDiam_ = 0.0;
342  dict.readIfPresent("minDiameter", minDiam_);
343  dict.readEntry("patches", patchNames_);
344  dict.readEntry("fields", fields_);
345 
346  const word format(dict.get<word>("setFormat"));
347  formatterPtr_ = writer<scalar>::New(format);
348 
349  if (dict.found(coordinateSystem::typeName_()))
350  {
351  csysPtr_.reset
352  (
353  coordinateSystem::New(obr_, dict, coordinateSystem::typeName_())
354  );
355 
356  Info<< "Transforming all vectorFields with coordinate system "
357  << csysPtr_->name() << endl;
358  }
359  else
360  {
361  csysPtr_.clear();
362  }
363 
364  if (isoPlanes_)
365  {
366  dict.readEntry("origin", origin_);
367  dict.readEntry("direction", direction_);
368  dict.readEntry("maxD", maxDiameter_);
369  dict.readEntry("nDownstreamBins", nDownstreamBins_);
370  dict.readEntry("maxDownstream", maxDownstream_);
371  direction_.normalise();
372  }
373 
374  return true;
375 }
376 
377 
379 {
380  return true;
381 }
382 
383 
385 {
386  Log << type() << " " << name() << " write:" << nl;
387 
389  if (obr_.foundObject<volScalarField>(alphaName_))
390  {
391  Log << " Looking up field " << alphaName_ << endl;
392  }
393  else
394  {
395  Info<< " Reading field " << alphaName_ << endl;
396  alphaPtr.reset
397  (
398  new volScalarField
399  (
400  IOobject
401  (
402  alphaName_,
403  mesh_.time().timeName(),
404  mesh_,
407  ),
408  mesh_
409  )
410  );
411  }
412 
413 
414  const volScalarField& alpha =
415  (
416  alphaPtr
417  ? *alphaPtr
418  : obr_.lookupObject<volScalarField>(alphaName_)
419  );
420 
421  Log << " Volume of alpha = "
422  << fvc::domainIntegrate(alpha).value()
423  << endl;
424 
425  const scalar meshVol = gSum(mesh_.V());
426  const scalar maxDropletVol = 1.0/6.0*pow3(maxDiam_);
427  const scalar delta = (maxDiam_-minDiam_)/nBins_;
428 
429  Log << " Mesh volume = " << meshVol << nl
430  << " Maximum droplet diameter = " << maxDiam_ << nl
431  << " Maximum droplet volume = " << maxDropletVol
432  << endl;
433 
434 
435  // Determine blocked faces
436  boolList blockedFace(mesh_.nFaces(), false);
437  label nBlocked = 0;
438 
439  {
440  for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
441  {
442  scalar ownVal = alpha[mesh_.faceOwner()[facei]];
443  scalar neiVal = alpha[mesh_.faceNeighbour()[facei]];
444 
445  if
446  (
447  (ownVal < threshold_ && neiVal > threshold_)
448  || (ownVal > threshold_ && neiVal < threshold_)
449  )
450  {
451  blockedFace[facei] = true;
452  nBlocked++;
453  }
454  }
455 
456  // Block coupled faces
457  forAll(alpha.boundaryField(), patchi)
458  {
459  const fvPatchScalarField& fvp = alpha.boundaryField()[patchi];
460  if (fvp.coupled())
461  {
462  tmp<scalarField> townFld(fvp.patchInternalField());
463  const scalarField& ownFld = townFld();
464 
465  tmp<scalarField> tnbrFld(fvp.patchNeighbourField());
466  const scalarField& nbrFld = tnbrFld();
467 
468  label start = fvp.patch().patch().start();
469 
470  forAll(ownFld, i)
471  {
472  scalar ownVal = ownFld[i];
473  scalar neiVal = nbrFld[i];
474 
475  if
476  (
477  (ownVal < threshold_ && neiVal > threshold_)
478  || (ownVal > threshold_ && neiVal < threshold_)
479  )
480  {
481  blockedFace[start+i] = true;
482  nBlocked++;
483  }
484  }
485  }
486  }
487  }
488 
489 
490  regionSplit regions(mesh_, blockedFace);
491 
492  Log << " Determined " << regions.nRegions()
493  << " disconnected regions" << endl;
494 
495 
496  if (debug)
497  {
498  volScalarField region
499  (
500  IOobject
501  (
502  "region",
503  mesh_.time().timeName(),
504  mesh_,
507  ),
508  mesh_,
510  );
511  Info<< " Dumping region as volScalarField to " << region.name()
512  << endl;
513 
514  forAll(regions, celli)
515  {
516  region[celli] = regions[celli];
517  }
518  region.correctBoundaryConditions();
519  region.write();
520  }
521 
522 
523  // Determine regions connected to supplied patches
524  Map<label> patchRegions(findPatchRegions(regions));
525 
526 
527  // Sum all regions
528  const scalarField alphaVol(alpha.primitiveField()*mesh_.V());
529  Map<scalar> allRegionVolume(regionSum(regions, mesh_.V()));
530  Map<scalar> allRegionAlphaVolume(regionSum(regions, alphaVol));
531  Map<label> allRegionNumCells
532  (
533  regionSum
534  (
535  regions,
536  labelField(mesh_.nCells(), 1.0)
537  )
538  );
539 
540  if (debug)
541  {
542  Info<< " " << token::TAB << "Region"
543  << token::TAB << "Volume(mesh)"
544  << token::TAB << "Volume(" << alpha.name() << "):"
545  << token::TAB << "nCells"
546  << nl;
547  scalar meshSumVol = 0.0;
548  scalar alphaSumVol = 0.0;
549  label nCells = 0;
550 
551  auto vIter = allRegionVolume.cbegin();
552  auto aIter = allRegionAlphaVolume.cbegin();
553  auto numIter = allRegionNumCells.cbegin();
554  for
555  (
556  ;
557  vIter.good() && aIter.good();
558  ++vIter, ++aIter, ++numIter
559  )
560  {
561  Info<< " " << token::TAB << vIter.key()
562  << token::TAB << vIter()
563  << token::TAB << aIter()
564  << token::TAB << numIter()
565  << nl;
566 
567  meshSumVol += vIter();
568  alphaSumVol += aIter();
569  nCells += numIter();
570  }
571  Info<< " " << token::TAB << "Total:"
572  << token::TAB << meshSumVol
573  << token::TAB << alphaSumVol
574  << token::TAB << nCells
575  << endl;
576  }
577 
578 
579  if (log)
580  {
581  Info<< " Patch connected regions (liquid core):" << nl;
582  Info<< token::TAB << " Region"
583  << token::TAB << "Volume(mesh)"
584  << token::TAB << "Volume(" << alpha.name() << "):"
585  << nl;
586 
587  forAllConstIters(patchRegions, iter)
588  {
589  const label regioni = iter.key();
590  Info<< " " << token::TAB << regioni
591  << token::TAB << allRegionVolume[regioni]
592  << token::TAB << allRegionAlphaVolume[regioni] << nl;
593 
594  }
595  Info<< endl;
596  }
597 
598  if (log)
599  {
600  Info<< " Background regions:" << nl;
601  Info<< " " << token::TAB << "Region"
602  << token::TAB << "Volume(mesh)"
603  << token::TAB << "Volume(" << alpha.name() << "):"
604  << nl;
605 
606  auto vIter = allRegionVolume.cbegin();
607  auto aIter = allRegionAlphaVolume.cbegin();
608 
609  for
610  (
611  ;
612  vIter.good() && aIter.good();
613  ++vIter, ++aIter
614  )
615  {
616  if
617  (
618  !patchRegions.found(vIter.key())
619  && vIter() >= maxDropletVol
620  )
621  {
622  Info<< " " << token::TAB << vIter.key()
623  << token::TAB << vIter()
624  << token::TAB << aIter() << nl;
625  }
626  }
627  Info<< endl;
628  }
629 
630 
631 
632  // Split alpha field
633  // ~~~~~~~~~~~~~~~~~
634  // Split into
635  // - liquidCore : region connected to inlet patches
636  // - per region a volume : for all other regions
637  // - backgroundAlpha : remaining alpha
638  writeAlphaFields(regions, patchRegions, allRegionVolume, alpha);
639 
640 
641  // Extract droplet-only allRegionVolume, i.e. delete liquid core
642  // (patchRegions) and background regions from maps.
643  // Note that we have to use mesh volume (allRegionVolume) and not
644  // allRegionAlphaVolume since background might not have alpha in it.
645  // Deleting regions where the volume-alpha-weighted is lower than
646  // threshold
647  forAllIters(allRegionVolume, vIter)
648  {
649  const label regioni = vIter.key();
650  if
651  (
652  patchRegions.found(regioni)
653  || vIter() >= maxDropletVol
654  || (allRegionAlphaVolume[regioni]/vIter() < threshold_)
655  )
656  {
657  allRegionVolume.erase(vIter);
658  allRegionAlphaVolume.erase(regioni);
659  allRegionNumCells.erase(regioni);
660  }
661  }
662 
663  if (allRegionVolume.size())
664  {
665  // Construct mids of bins for plotting
666  pointField xBin(nBins_);
667 
668  scalar x = 0.5*delta;
669  forAll(xBin, i)
670  {
671  xBin[i] = point(x, 0, 0);
672  x += delta;
673  }
674 
675  const coordSet coords("diameter", "x", xBin, mag(xBin));
676 
677 
678  // Get in region order the alpha*volume and diameter
679  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
680 
681  const labelList sortedRegions = allRegionAlphaVolume.sortedToc();
682 
683  scalarField sortedVols
684  (
685  extractData
686  (
687  sortedRegions,
688  allRegionAlphaVolume
689  )
690  );
691 
692  vectorField centroids(sortedVols.size(), Zero);
693 
694  // Check if downstream bins are calculated
695  if (isoPlanes_)
696  {
697  vectorField alphaDistance
698  (
699  (alpha.primitiveField()*mesh_.V())
700  *(mesh_.C().primitiveField() - origin_)()
701  );
702 
703  Map<vector> allRegionAlphaDistance
704  (
705  regionSum
706  (
707  regions,
708  alphaDistance
709  )
710  );
711 
712  // 2. centroid
713  vectorField sortedMoment
714  (
715  extractData
716  (
717  sortedRegions,
718  allRegionAlphaDistance
719  )
720  );
721 
722  centroids = sortedMoment/sortedVols + origin_;
723 
724  // Bin according to centroid
725  scalarField distToPlane((centroids - origin_) & direction_);
726 
727  vectorField radialDistToOrigin
728  (
729  (centroids - origin_) - (distToPlane*direction_)
730  );
731 
732  const scalar deltaX = maxDownstream_/nDownstreamBins_;
733  labelList downstreamIndices(distToPlane.size(), -1);
734  forAll(distToPlane, i)
735  {
736  if
737  (
738  (mag(radialDistToOrigin[i]) < maxDiameter_)
739  && (distToPlane[i] < maxDownstream_)
740  )
741  {
742  downstreamIndices[i] = distToPlane[i]/deltaX;
743  }
744  }
745 
746  scalarField binDownCount(nDownstreamBins_, Zero);
747  forAll(distToPlane, i)
748  {
749  if (downstreamIndices[i] != -1)
750  {
751  binDownCount[downstreamIndices[i]] += 1.0;
752  }
753  }
754 
755  // Write
756  if (Pstream::master())
757  {
758  // Construct mids of bins for plotting
759  pointField xBin(nDownstreamBins_);
760 
761  scalar x = 0.5*deltaX;
762  forAll(xBin, i)
763  {
764  xBin[i] = point(x, 0, 0);
765  x += deltaX;
766  }
767 
768  const coordSet coords("distance", "x", xBin, mag(xBin));
769  writeGraph(coords, "isoPlanes", binDownCount);
770  }
771 
772  // Write to log
773  if (log)
774  {
775  Info<< " Iso-planes Bins:" << nl
776  << " " << token::TAB << "Bin"
777  << token::TAB << "Min distance"
778  << token::TAB << "Count:"
779  << nl;
780 
781  scalar delta = 0.0;
782  forAll(binDownCount, bini)
783  {
784  Info<< " " << token::TAB << bini
785  << token::TAB << delta
786  << token::TAB << binDownCount[bini] << nl;
787  delta += deltaX;
788  }
789  Info<< endl;
790 
791  }
792  }
793 
794  // Calculate the diameters
795  scalarField sortedDiameters(sortedVols.size());
796  forAll(sortedDiameters, i)
797  {
798  sortedDiameters[i] = Foam::cbrt
799  (
800  sortedVols[i]
802  );
803  }
804 
805  // Determine the bin index for all the diameters
806  labelList indices(sortedDiameters.size());
807  forAll(sortedDiameters, i)
808  {
809  indices[i] = (sortedDiameters[i]-minDiam_)/delta;
810  }
811 
812  // Calculate the counts per diameter bin
813  scalarField binCount(nBins_, Zero);
814  forAll(sortedDiameters, i)
815  {
816  binCount[indices[i]] += 1.0;
817  }
818 
819  // Write counts
820  if (Pstream::master())
821  {
822  writeGraph(coords, "count", binCount);
823  }
824 
825  // Write to log
826  if (log)
827  {
828  Info<< " Bins:" << nl
829  << " " << token::TAB << "Bin"
830  << token::TAB << "Min diameter"
831  << token::TAB << "Count:"
832  << nl;
833 
834  scalar diam = 0.0;
835  forAll(binCount, bini)
836  {
837  Info<< " " << token::TAB << bini
838  << token::TAB << diam
839  << token::TAB << binCount[bini] << nl;
840 
841  diam += delta;
842  }
843 
844  Info<< endl;
845  }
846 
847 
848  // Write average and deviation of droplet volume.
849  writeGraphs
850  (
851  "volume", // name of field
852  indices, // per region the bin index
853  sortedVols, // per region field data
854  binCount, // per bin number of regions
855  coords // graph data for bins
856  );
857 
858  // Collect some more field
859  {
860  wordList scalarNames(obr_.names(volScalarField::typeName));
861 
862  const labelList selected(fields_.matching(scalarNames));
863 
864  for (const label fieldi : selected)
865  {
866  const word& fldName = scalarNames[fieldi];
867 
868  Log << " Scalar field " << fldName << endl;
869 
870  const scalarField& fld = obr_.lookupObject
871  <
873  >(fldName).primitiveField();
874 
875  writeGraphs
876  (
877  fldName, // name of field
878  alphaVol*fld, // per cell field data
879 
880  regions, // per cell the region(=droplet)
881  sortedRegions, // valid regions in sorted order
882  1.0/sortedVols, // per region normalisation
883 
884  indices, // index of bin for each region
885  binCount, // per bin number of regions
886  coords // graph data for bins
887  );
888  }
889  }
890  {
891  wordList vectorNames(obr_.names(volVectorField::typeName));
892 
893  const labelList selected(fields_.matching(vectorNames));
894 
895  for (const label fieldi : selected)
896  {
897  const word& fldName = vectorNames[fieldi];
898 
899  Log << " Vector field " << fldName << endl;
900 
901  vectorField fld = obr_.lookupObject
902  <
904  >(fldName).primitiveField();
905 
906  if (csysPtr_)
907  {
908  Log << "Transforming vector field " << fldName
909  << " with coordinate system "
910  << csysPtr_->name()
911  << endl;
912 
913  fld = csysPtr_->localVector(fld);
914  }
915 
916 
917  // Components
918 
919  for (direction cmp = 0; cmp < vector::nComponents; cmp++)
920  {
921  writeGraphs
922  (
923  fldName + vector::componentNames[cmp],
924  alphaVol*fld.component(cmp),// 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  // Magnitude
937  writeGraphs
938  (
939  fldName + "mag", // name of field
940  alphaVol*mag(fld), // per cell field data
941 
942  regions, // per cell the region(=droplet)
943  sortedRegions, // valid regions in sorted order
944  1.0/sortedVols, // per region normalisation
945 
946  indices, // index of bin for each region
947  binCount, // per bin number of regions
948  coords // graph data for bins
949  );
950  }
951  }
952  }
953 
954  return true;
955 }
956 
957 
958 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::Pstream::mapCombineGather
static void mapCombineGather(const List< commsStruct > &comms, Container &Values, const CombineOp &cop, const int tag, const label comm)
Definition: combineGatherScatter.C:551
Foam::functionObjects::regionSizeDistribution::write
virtual bool write()
Calculate the regionSizeDistribution and write.
Definition: regionSizeDistribution.C:384
Foam::fvPatchField< scalar >
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
regionSizeDistribution.H
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Log
#define Log
Definition: PDRblock.C:35
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::componentNames
static const char *const componentNames[]
Definition: VectorSpace.H:114
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::HashTableOps::values
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition: HashOps.H:149
Foam::fvc::domainIntegrate
dimensioned< Type > domainIntegrate(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcVolumeIntegrate.C:88
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
alphaPtr
autoPtr< volScalarField > alphaPtr
Definition: readThermalProperties.H:159
Foam::dimensioned::name
const word & name() const
Return const reference to name.
Definition: dimensionedType.C:406
Foam::Map
A HashTable to objects of type <T> with a label key.
Definition: lumpedPointController.H:69
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Am I the master process.
Definition: UPstream.H:457
Foam::functionObjects::regionSizeDistribution::regionSizeDistribution
regionSizeDistribution(const word &name, const Time &runTime, const dictionary &)
Construct for given objectRegistry and dictionary.
Definition: regionSizeDistribution.C:314
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
Foam::fvPatchField< scalar >::calculatedType
static const word & calculatedType()
Return the type of the calculated for of fvPatchField.
Definition: calculatedFvPatchField.C:35
Foam::functionObjects::fvMeshFunctionObject
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Definition: fvMeshFunctionObject.H:64
Foam::functionObjects::regionSizeDistribution::read
virtual bool read(const dictionary &)
Read the regionSizeDistribution data.
Definition: regionSizeDistribution.C:332
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::divide
void divide(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)
Foam::wordList
List< word > wordList
A List of words.
Definition: fileName.H:62
format
word format(conversionProperties.get< word >("format"))
Foam::fvPatchField::coupled
virtual bool coupled() const
Return true if this patch field is coupled.
Definition: fvPatchField.H:345
Foam::functionObjects::writeFile::read
virtual bool read(const dictionary &dict)
Read.
Definition: writeFile.C:213
Foam::Field< scalar >
Foam::regIOobject::write
virtual bool write(const bool valid=true) const
Write using setting from DB.
Definition: regIOobjectWrite.C:132
Foam::labelField
Field< label > labelField
Specialisation of Field<T> for label.
Definition: labelField.H:52
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::coordinateSystem::New
static autoPtr< coordinateSystem > New(word modelType, const objectRegistry &obr, const dictionary &dict)
Definition: coordinateSystemNew.C:84
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
Foam::fvPatchField::patchInternalField
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:233
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::regionSplit
This class separates the mesh into distinct unconnected regions, each of which is then given a label ...
Definition: regionSplit.H:140
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::regionSplit::nRegions
label nRegions() const
Return total number of regions.
Definition: regionSplit.H:294
fld
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))
Definition: gmvOutputLagrangian.H:23
forAllIters
#define forAllIters(container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:223
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::functionObjects::regionFunctionObject::read
virtual bool read(const dictionary &dict)
Read optional controls.
Definition: regionFunctionObject.C:173
Foam::log
dimensionedScalar log(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:262
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::Ostream::write
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
Foam::coordSet
Holds list of sampling positions.
Definition: coordSet.H:53
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::polyPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:361
Foam::fvPatchField::patchNeighbourField
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField on the opposite patch of a coupled patch.
Definition: fvPatchField.H:448
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=worldComm)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:463
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::nl
constexpr char nl
Definition: Ostream.H:404
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
Foam::Pstream::mapCombineScatter
static void mapCombineScatter(const List< commsStruct > &comms, Container &Values, const int tag, const label comm)
Scatter data. Reverse of combineGather.
Definition: combineGatherScatter.C:666
Foam::functionObjects::addToRunTimeSelectionTable
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
Foam::List< bool >
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::token::TAB
Tab [isspace].
Definition: token.H:123
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::wordRes
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:51
fvcVolumeIntegrate.H
Volume integrate volField creating a volField.
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Foam::direction
uint8_t direction
Definition: direction.H:52
Foam::functionObjects::defineTypeNameAndDebug
defineTypeNameAndDebug(ObukhovLength, 0)
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::fvPatch::patch
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:161
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::functionObjects::regionSizeDistribution::execute
virtual bool execute()
Do nothing.
Definition: regionSizeDistribution.C:378
Foam::functionObjects::log
Computes the natural logarithm of an input volScalarField.
Definition: log.H:227
Foam::functionObjects::writeFile
Base class for writing single files from the function objects.
Definition: writeFile.H:119
Foam::fvPatchField::patch
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:357
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::functionObjects::writeFile::baseTimeDir
fileName baseTimeDir() const
Return the base directory for the current time value.
Definition: writeFile.C:76
Foam::labelHashSet
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
Foam::cbrt
dimensionedScalar cbrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:155
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::mkDir
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition: MSwindows.C:507
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::nComponents
static constexpr direction nComponents
Number of components in this vector space.
Definition: VectorSpace.H:101
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189
Foam::IOobject::MUST_READ
Definition: IOobject.H:185
Foam::writer::New
static autoPtr< writer > New(const word &writeFormat)
Return a reference to the selected writer.
Definition: writer.C:38