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-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
30#include "regionSplit.H"
31#include "volFields.H"
32#include "fvcVolumeIntegrate.H"
34
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36
37namespace Foam
38{
39 namespace functionObjects
40 {
43 (
47 );
48 }
49}
50
51
52// * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
53
54namespace Foam
55{
56
57template<class Type>
58static Map<Type> regionSum(const regionSplit& regions, const Field<Type>& fld)
59{
60 // Per region the sum of fld
61 Map<Type> regionToSum(regions.nRegions()/Pstream::nProcs());
62
63 forAll(fld, celli)
64 {
65 const label regioni = regions[celli];
66 regionToSum(regioni, Type(Zero)) += fld[celli];
67 }
68
70
71 return regionToSum;
72}
73
74
75template<class Type>
76static List<Type> extractData(const labelUList& keys, const Map<Type>& regionData)
77{
78 List<Type> sortedData(keys.size());
79
80 forAll(keys, i)
81 {
82 sortedData[i] = regionData[keys[i]];
83 }
84 return sortedData;
85}
86
87} // End namespace Foam
88
89
90// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
91
92void Foam::functionObjects::regionSizeDistribution::writeAlphaFields
93(
94 const regionSplit& regions,
95 const Map<label>& patchRegions,
96 const Map<scalar>& regionVolume,
97 const volScalarField& alpha
98) const
99{
100 const scalar maxDropletVol = 1.0/6.0*pow3(maxDiam_);
101
102 // Split alpha field
103 // ~~~~~~~~~~~~~~~~~
104 // Split into
105 // - liquidCore : region connected to inlet patches
106 // - per region a volume : for all other regions
107 // - backgroundAlpha : remaining alpha
108
109
110 // Construct field
111 volScalarField liquidCore
112 (
113 IOobject
114 (
115 scopedName(alphaName_ + "_liquidCore"),
116 obr_.time().timeName(),
117 obr_,
119 ),
120 alpha,
121 fvPatchField<scalar>::calculatedType()
122 );
123
124 volScalarField backgroundAlpha
125 (
126 IOobject
127 (
128 scopedName(alphaName_ + "_background"),
129 obr_.time().timeName(),
130 obr_,
132 ),
133 alpha,
134 fvPatchField<scalar>::calculatedType()
135 );
136
137
138 // Knock out any cell not in patchRegions
139 forAll(liquidCore, celli)
140 {
141 const label regioni = regions[celli];
142 if (patchRegions.found(regioni))
143 {
144 backgroundAlpha[celli] = 0;
145 }
146 else
147 {
148 liquidCore[celli] = 0;
149
150 const scalar regionVol = regionVolume[regioni];
151 if (regionVol < maxDropletVol)
152 {
153 backgroundAlpha[celli] = 0;
154 }
155 }
156 }
157 liquidCore.correctBoundaryConditions();
158 backgroundAlpha.correctBoundaryConditions();
159
160 if (log)
161 {
162 Info<< " Volume of liquid-core = "
163 << fvc::domainIntegrate(liquidCore).value()
164 << nl
165 << " Volume of background = "
166 << fvc::domainIntegrate(backgroundAlpha).value()
167 << endl;
168 }
169
170 Log << " Writing liquid-core field to " << liquidCore.name() << endl;
171 liquidCore.write();
172
173 Log<< " Writing background field to " << backgroundAlpha.name() << endl;
174 backgroundAlpha.write();
175}
176
177
179Foam::functionObjects::regionSizeDistribution::findPatchRegions
180(
181 const regionSplit& regions
182) const
183{
184 // Mark all regions starting at patches
185 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
186
187 // Count number of patch faces (just for initial sizing)
188 const labelHashSet patchIDs(mesh_.boundaryMesh().patchSet(patchNames_));
189
190 label nPatchFaces = 0;
191 for (const label patchi : patchIDs)
192 {
193 nPatchFaces += mesh_.boundaryMesh()[patchi].size();
194 }
195
196
197 Map<label> patchRegions(nPatchFaces);
198 for (const label patchi : patchIDs)
199 {
200 const polyPatch& pp = mesh_.boundaryMesh()[patchi];
201
202 // Collect all regions on the patch
203 const labelList& faceCells = pp.faceCells();
204
205 for (const label celli : faceCells)
206 {
207 patchRegions.insert
208 (
209 regions[celli],
210 Pstream::myProcNo() // dummy value
211 );
212 }
213 }
214
215
216 // Make sure all the processors have the same set of regions
217 Pstream::mapCombineAllGather(patchRegions, minEqOp<label>());
218
219 return patchRegions;
220}
221
222
224Foam::functionObjects::regionSizeDistribution::divide
225(
226 const scalarField& num,
227 const scalarField& denom
228)
229{
230 auto tresult = tmp<scalarField>::New(num.size());
231 auto& result = tresult.ref();
232
233 forAll(denom, i)
234 {
235 if (denom[i] != 0)
236 {
237 result[i] = num[i]/denom[i];
238 }
239 else
240 {
241 result[i] = 0;
242 }
243 }
244 return tresult;
245}
246
247
248void Foam::functionObjects::regionSizeDistribution::writeGraphs
249(
250 const word& fieldName, // name of field
251 const scalarField& sortedField, // per region field data
252
253 const labelList& indices, // index of bin for each region
254 const scalarField& binCount, // per bin number of regions
255 const coordSet& coords // graph data for bins
256) const
257{
258 if (Pstream::master())
259 {
260 // Calculate per-bin average
261 scalarField binSum(nBins_, Zero);
262 forAll(sortedField, i)
263 {
264 binSum[indices[i]] += sortedField[i];
265 }
266
267 scalarField binAvg(divide(binSum, binCount));
268
269 // Per bin deviation
270 scalarField binSqrSum(nBins_, Zero);
271 forAll(sortedField, i)
272 {
273 binSqrSum[indices[i]] += Foam::sqr(sortedField[i]);
274 }
275 scalarField binDev
276 (
277 sqrt(divide(binSqrSum, binCount) - Foam::sqr(binAvg))
278 );
279
280
281 auto& writer = formatterPtr_();
282
283 word outputName;
284 if (writer.buffering())
285 {
286 outputName =
287 (
288 coords.name()
290 (
292 ({
293 fieldName + "_sum",
294 fieldName + "_avg",
295 fieldName + "_dev"
296 })
297 )
298 );
299 }
300 else
301 {
302 outputName = coords.name();
303 }
304
305 writer.open
306 (
307 coords,
308 (baseTimeDir() / outputName)
309 );
310
311 Log << " Writing distribution of "
312 << fieldName << " to " << writer.path() << endl;
313
314 writer.write(fieldName + "_sum", binSum);
315 writer.write(fieldName + "_avg", binAvg);
316 writer.write(fieldName + "_dev", binDev);
317 writer.close(true);
318 }
319}
320
321
322void Foam::functionObjects::regionSizeDistribution::writeGraphs
323(
324 const word& fieldName, // name of field
325 const scalarField& cellField, // per cell field data
326 const regionSplit& regions, // per cell the region(=droplet)
327 const labelList& sortedRegions, // valid regions in sorted order
328 const scalarField& sortedNormalisation,
329
330 const labelList& indices, // per region index of bin
331 const scalarField& binCount, // per bin number of regions
332 const coordSet& coords // graph data for bins
333) const
334{
335 // Sum on a per-region basis. Parallel reduced.
336 Map<scalar> regionField(regionSum(regions, cellField));
337
338 // Extract in region order
339 scalarField sortedField
340 (
341 sortedNormalisation
342 * extractData(sortedRegions, regionField)
343 );
344
345 writeGraphs
346 (
347 fieldName, // name of field
348 sortedField, // per region field data
349
350 indices, // index of bin for each region
351 binCount, // per bin number of regions
352 coords // graph data for bins
353 );
354}
355
356
357// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
358
360(
361 const word& name,
362 const Time& runTime,
363 const dictionary& dict
364)
365:
367 writeFile(obr_, name),
368 alphaName_(dict.get<word>("field")),
369 patchNames_(dict.get<wordRes>("patches")),
370 isoPlanes_(dict.getOrDefault("isoPlanes", false))
371{
372 read(dict);
373}
374
375
376// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
377
379{
382
383 dict.readEntry("nBins", nBins_);
384 dict.readEntry("field", alphaName_);
385 dict.readEntry("threshold", threshold_);
386 dict.readEntry("maxDiameter", maxDiam_);
387 minDiam_ = 0.0;
388 dict.readIfPresent("minDiameter", minDiam_);
389 dict.readEntry("patches", patchNames_);
390 dict.readEntry("fields", fields_);
391
392 const word setFormat(dict.get<word>("setFormat"));
393 formatterPtr_ = coordSetWriter::New
394 (
395 setFormat,
396 dict.subOrEmptyDict("formatOptions").optionalSubDict(setFormat)
397 );
398
399 if (dict.found(coordinateSystem::typeName_()))
400 {
401 csysPtr_.reset
402 (
403 coordinateSystem::New(obr_, dict, coordinateSystem::typeName_())
404 );
405
406 Info<< "Transforming all vectorFields with coordinate system "
407 << csysPtr_->name() << endl;
408 }
409 else
410 {
411 csysPtr_.clear();
412 }
413
414 if (isoPlanes_)
415 {
416 dict.readEntry("origin", origin_);
417 dict.readEntry("direction", direction_);
418 dict.readEntry("maxD", maxDiameter_);
419 dict.readEntry("nDownstreamBins", nDownstreamBins_);
420 dict.readEntry("maxDownstream", maxDownstream_);
421 direction_.normalise();
422 }
423
424 return true;
425}
426
427
429{
430 return true;
431}
432
433
435{
436 Log << type() << " " << name() << " write:" << nl;
437
438 tmp<volScalarField> talpha;
439 talpha.cref(obr_.cfindObject<volScalarField>(alphaName_));
440 if (talpha)
441 {
442 Log << " Looking up field " << alphaName_ << endl;
443 }
444 else
445 {
446 Info<< " Reading field " << alphaName_ << endl;
447 talpha.reset
448 (
450 (
452 (
453 alphaName_,
454 mesh_.time().timeName(),
455 mesh_,
458 ),
459 mesh_
460 )
461 );
462 }
463 const auto& alpha = talpha();
464
465 Log << " Volume of alpha = "
466 << fvc::domainIntegrate(alpha).value()
467 << endl;
468
469 const scalar meshVol = gSum(mesh_.V());
470 const scalar maxDropletVol = 1.0/6.0*pow3(maxDiam_);
471 const scalar delta = (maxDiam_-minDiam_)/nBins_;
472
473 Log << " Mesh volume = " << meshVol << nl
474 << " Maximum droplet diameter = " << maxDiam_ << nl
475 << " Maximum droplet volume = " << maxDropletVol
476 << endl;
477
478
479 // Determine blocked faces
480 boolList blockedFace(mesh_.nFaces(), false);
481 label nBlocked = 0;
482
483 {
484 for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
485 {
486 scalar ownVal = alpha[mesh_.faceOwner()[facei]];
487 scalar neiVal = alpha[mesh_.faceNeighbour()[facei]];
488
489 if
490 (
491 (ownVal < threshold_ && neiVal > threshold_)
492 || (ownVal > threshold_ && neiVal < threshold_)
493 )
494 {
495 blockedFace[facei] = true;
496 nBlocked++;
497 }
498 }
499
500 // Block coupled faces
501 forAll(alpha.boundaryField(), patchi)
502 {
503 const fvPatchScalarField& fvp = alpha.boundaryField()[patchi];
504 if (fvp.coupled())
505 {
508 const auto& ownFld = townFld();
509 const auto& nbrFld = tnbrFld();
510
511 label start = fvp.patch().patch().start();
512
513 forAll(ownFld, i)
514 {
515 scalar ownVal = ownFld[i];
516 scalar neiVal = nbrFld[i];
517
518 if
519 (
520 (ownVal < threshold_ && neiVal > threshold_)
521 || (ownVal > threshold_ && neiVal < threshold_)
522 )
523 {
524 blockedFace[start+i] = true;
525 nBlocked++;
526 }
527 }
528 }
529 }
530 }
531
532
533 regionSplit regions(mesh_, blockedFace);
534
535 Log << " Determined " << regions.nRegions()
536 << " disconnected regions" << endl;
537
538
539 if (debug)
540 {
541 volScalarField region
542 (
544 (
545 "region",
546 mesh_.time().timeName(),
547 mesh_,
550 ),
551 mesh_,
553 );
554 Info<< " Dumping region as volScalarField to " << region.name()
555 << endl;
556
557 forAll(regions, celli)
558 {
559 region[celli] = regions[celli];
560 }
562 region.write();
563 }
564
565
566 // Determine regions connected to supplied patches
567 Map<label> patchRegions(findPatchRegions(regions));
568
569
570 // Sum all regions
571 const scalarField alphaVol(alpha.primitiveField()*mesh_.V());
572 Map<scalar> allRegionVolume(regionSum(regions, mesh_.V()));
573 Map<scalar> allRegionAlphaVolume(regionSum(regions, alphaVol));
574 Map<label> allRegionNumCells
575 (
577 (
578 regions,
579 labelField(mesh_.nCells(), 1.0)
580 )
581 );
582
583 if (debug)
584 {
585 Info<< " " << token::TAB << "Region"
586 << token::TAB << "Volume(mesh)"
587 << token::TAB << "Volume(" << alpha.name() << "):"
588 << token::TAB << "nCells"
589 << nl;
590 scalar meshSumVol = 0.0;
591 scalar alphaSumVol = 0.0;
592 label nCells = 0;
593
594 auto vIter = allRegionVolume.cbegin();
595 auto aIter = allRegionAlphaVolume.cbegin();
596 auto numIter = allRegionNumCells.cbegin();
597 for
598 (
599 ;
600 vIter.good() && aIter.good();
601 ++vIter, ++aIter, ++numIter
602 )
603 {
604 Info<< " " << token::TAB << vIter.key()
605 << token::TAB << vIter()
606 << token::TAB << aIter()
607 << token::TAB << numIter()
608 << nl;
609
610 meshSumVol += vIter();
611 alphaSumVol += aIter();
612 nCells += numIter();
613 }
614 Info<< " " << token::TAB << "Total:"
615 << token::TAB << meshSumVol
616 << token::TAB << alphaSumVol
617 << token::TAB << nCells
618 << endl;
619 }
620
621
622 if (log)
623 {
624 Info<< " Patch connected regions (liquid core):" << nl;
625 Info<< token::TAB << " Region"
626 << token::TAB << "Volume(mesh)"
627 << token::TAB << "Volume(" << alpha.name() << "):"
628 << nl;
629
630 forAllConstIters(patchRegions, iter)
631 {
632 const label regioni = iter.key();
633 Info<< " " << token::TAB << regioni
634 << token::TAB << allRegionVolume[regioni]
635 << token::TAB << allRegionAlphaVolume[regioni] << nl;
636
637 }
638 Info<< endl;
639 }
640
641 if (log)
642 {
643 Info<< " Background regions:" << nl;
644 Info<< " " << token::TAB << "Region"
645 << token::TAB << "Volume(mesh)"
646 << token::TAB << "Volume(" << alpha.name() << "):"
647 << nl;
648
649 auto vIter = allRegionVolume.cbegin();
650 auto aIter = allRegionAlphaVolume.cbegin();
651
652 for
653 (
654 ;
655 vIter.good() && aIter.good();
656 ++vIter, ++aIter
657 )
658 {
659 if
660 (
661 !patchRegions.found(vIter.key())
662 && vIter() >= maxDropletVol
663 )
664 {
665 Info<< " " << token::TAB << vIter.key()
666 << token::TAB << vIter()
667 << token::TAB << aIter() << nl;
668 }
669 }
670 Info<< endl;
671 }
672
673
674
675 // Split alpha field
676 // ~~~~~~~~~~~~~~~~~
677 // Split into
678 // - liquidCore : region connected to inlet patches
679 // - per region a volume : for all other regions
680 // - backgroundAlpha : remaining alpha
681 writeAlphaFields(regions, patchRegions, allRegionVolume, alpha);
682
683
684 // Extract droplet-only allRegionVolume, i.e. delete liquid core
685 // (patchRegions) and background regions from maps.
686 // Note that we have to use mesh volume (allRegionVolume) and not
687 // allRegionAlphaVolume since background might not have alpha in it.
688 // Deleting regions where the volume-alpha-weighted is lower than
689 // threshold
690 forAllIters(allRegionVolume, vIter)
691 {
692 const label regioni = vIter.key();
693 if
694 (
695 patchRegions.found(regioni)
696 || vIter() >= maxDropletVol
697 || (allRegionAlphaVolume[regioni]/vIter() < threshold_)
698 )
699 {
700 allRegionVolume.erase(vIter);
701 allRegionAlphaVolume.erase(regioni);
702 allRegionNumCells.erase(regioni);
703 }
704 }
705
706 if (allRegionVolume.size())
707 {
708 // Construct mids of bins for plotting
709 pointField xBin(nBins_, Zero);
710
711 {
712 scalar x = 0.5*delta;
713 for (point& p : xBin)
714 {
715 p.x() = x;
716 x += delta;
717 }
718 }
719
720 const coordSet coords("diameter", "x", xBin, mag(xBin));
721
722
723 // Get in region order the alpha*volume and diameter
724 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
725
726 const labelList sortedRegions = allRegionAlphaVolume.sortedToc();
727
728 scalarField sortedVols
729 (
730 extractData(sortedRegions, allRegionAlphaVolume)
731 );
732
733 vectorField centroids(sortedVols.size(), Zero);
734
735 // Check if downstream bins are calculated
736 if (isoPlanes_)
737 {
738 vectorField alphaDistance
739 (
740 (alpha.primitiveField()*mesh_.V())
741 *(mesh_.C().primitiveField() - origin_)()
742 );
743
744 Map<vector> allRegionAlphaDistance
745 (
747 (
748 regions,
749 alphaDistance
750 )
751 );
752
753 // 2. centroid
754 vectorField sortedMoment
755 (
756 extractData(sortedRegions, allRegionAlphaDistance)
757 );
758
759 centroids = sortedMoment/sortedVols + origin_;
760
761 // Bin according to centroid
762 scalarField distToPlane((centroids - origin_) & direction_);
763
764 vectorField radialDistToOrigin
765 (
766 (centroids - origin_) - (distToPlane*direction_)
767 );
768
769 const scalar deltaX = maxDownstream_/nDownstreamBins_;
770 labelList downstreamIndices(distToPlane.size(), -1);
771 forAll(distToPlane, i)
772 {
773 if
774 (
775 (mag(radialDistToOrigin[i]) < maxDiameter_)
776 && (distToPlane[i] < maxDownstream_)
777 )
778 {
779 downstreamIndices[i] = distToPlane[i]/deltaX;
780 }
781 }
782
783 scalarField binDownCount(nDownstreamBins_, Zero);
784 forAll(distToPlane, i)
785 {
786 if (downstreamIndices[i] != -1)
787 {
788 binDownCount[downstreamIndices[i]] += 1.0;
789 }
790 }
791
792 // Write
793 if (Pstream::master())
794 {
795 // Construct mids of bins for plotting
796 pointField xBin(nDownstreamBins_, Zero);
797
798 {
799 scalar x = 0.5*deltaX;
800 for (point& p : xBin)
801 {
802 p.x() = x;
803 x += deltaX;
804 }
805 }
806
807 const coordSet coords("distance", "x", xBin, mag(xBin));
808
809 auto& writer = formatterPtr_();
810 writer.nFields(1);
811
812 writer.open
813 (
814 coords,
815 writeFile::baseTimeDir() / (coords.name() + "_isoPlanes")
816 );
817
818 writer.write("isoPlanes", binDownCount);
819 writer.close(true);
820 }
821
822 // Write to log
823 if (log)
824 {
825 Info<< " Iso-planes Bins:" << nl
826 << " " << token::TAB << "Bin"
827 << token::TAB << "Min distance"
828 << token::TAB << "Count:"
829 << nl;
830
831 scalar delta = 0.0;
832 forAll(binDownCount, bini)
833 {
834 Info<< " " << token::TAB << bini
835 << token::TAB << delta
836 << token::TAB << binDownCount[bini] << nl;
837 delta += deltaX;
838 }
839 Info<< endl;
840
841 }
842 }
843
844 // Calculate the diameters
845 scalarField sortedDiameters(sortedVols.size());
846 forAll(sortedDiameters, i)
847 {
848 sortedDiameters[i] = Foam::cbrt
849 (
850 sortedVols[i]
852 );
853 }
854
855 // Determine the bin index for all the diameters
856 labelList indices(sortedDiameters.size());
857 forAll(sortedDiameters, i)
858 {
859 indices[i] = (sortedDiameters[i]-minDiam_)/delta;
860 }
861
862 // Calculate the counts per diameter bin
863 scalarField binCount(nBins_, Zero);
864 forAll(sortedDiameters, i)
865 {
866 binCount[indices[i]] += 1.0;
867 }
868
869 // Write counts
870 if (Pstream::master())
871 {
872 auto& writer = formatterPtr_();
873 writer.nFields(1);
874
875 writer.open
876 (
877 coords,
878 writeFile::baseTimeDir() / (coords.name() + "_count")
879 );
880
881 writer.write("count", binCount);
882 writer.close(true);
883 }
884
885 // Write to log
886 if (log)
887 {
888 Info<< " Bins:" << nl
889 << " " << token::TAB << "Bin"
890 << token::TAB << "Min diameter"
891 << token::TAB << "Count:"
892 << nl;
893
894 scalar diam = 0.0;
895 forAll(binCount, bini)
896 {
897 Info<< " " << token::TAB << bini
898 << token::TAB << diam
899 << token::TAB << binCount[bini] << nl;
900
901 diam += delta;
902 }
903
904 Info<< endl;
905 }
906
907
908 // Write average and deviation of droplet volume.
909 writeGraphs
910 (
911 "volume", // name of field
912 sortedVols, // per region field data
913
914 indices, // per region the bin index
915 binCount, // per bin number of regions
916 coords // graph data for bins
917 );
918
919 // Collect some more fields
920 {
921 for
922 (
923 const word& fldName
924 : obr_.sortedNames<volScalarField>(fields_)
925 )
926 {
927 Log << " Scalar field " << fldName << endl;
928
930 (
931 obr_.lookupObject<volScalarField>(fldName).primitiveField()
932 );
933 const auto& fld = tfld();
934
935 writeGraphs
936 (
937 fldName, // name of field
938 alphaVol*fld, // per cell field data
939
940 regions, // per cell the region(=droplet)
941 sortedRegions, // valid regions in sorted order
942 1.0/sortedVols, // per region normalisation
943
944 indices, // index of bin for each region
945 binCount, // per bin number of regions
946 coords // graph data for bins
947 );
948 }
949 }
950
951 {
952 for
953 (
954 const word& fldName
955 : obr_.sortedNames<volVectorField>(fields_)
956 )
957 {
958 Log << " Vector field " << fldName << endl;
959
961 (
962 obr_.lookupObject<volVectorField>(fldName).primitiveField()
963 );
964
965 if (csysPtr_)
966 {
967 Log << "Transforming vector field " << fldName
968 << " with coordinate system "
969 << csysPtr_->name()
970 << endl;
971
972 tfld = csysPtr_->localVector(tfld());
973 }
974 const auto& fld = tfld();
975
976 // Components
977
978 for (direction cmpt = 0; cmpt < vector::nComponents; ++cmpt)
979 {
980 writeGraphs
981 (
982 fldName + vector::componentNames[cmpt],
983 alphaVol*fld.component(cmpt),// 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 // Magnitude
996 writeGraphs
997 (
998 fldName + "mag", // name of field
999 alphaVol*mag(fld), // per cell field data
1000
1001 regions, // per cell the region(=droplet)
1002 sortedRegions, // valid regions in sorted order
1003 1.0/sortedVols, // per region normalisation
1004
1005 indices, // index of bin for each region
1006 binCount, // per bin number of regions
1007 coords // graph data for bins
1008 );
1009 }
1010 }
1011 }
1012
1013 return true;
1014}
1015
1016
1017// ************************************************************************* //
scalar delta
#define Log
Definition: PDRblock.C:35
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
vtk::internalMeshWriter writer(topoMesh, topoCells, vtk::formatType::INLINE_ASCII, runTime.path()/"blockTopology")
Info<< nl<< "Wrote faMesh in vtk format: "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.mesh().time().globalPath()/"finiteArea-edges"));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Generic templated field type.
Definition: Field.H:82
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
void correctBoundaryConditions()
Correct boundary field.
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition: HashTable.C:137
const_iterator cbegin() const
const_iterator set to the beginning of the HashTable
bool found(const Key &key) const
Return true if hashed entry is found in table.
Definition: HashTableI.H:100
label size() const noexcept
The number of elements in table.
Definition: HashTableI.H:52
bool erase(const iterator &iter)
Erase an entry specified by given iterator.
Definition: HashTable.C:440
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
A HashTable to objects of type <T> with a label key.
Definition: Map.H:60
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
label nProcs() const noexcept
Number of ranks associated with PstreamBuffers.
static void mapCombineAllGather(const List< commsStruct > &comms, Container &values, const CombineOp &cop, const int tag, const label comm)
After completion all processors have the same data.
virtual bool read()
Re-read model coefficients if they have changed.
The TAB Method for Numerical Calculation of Spray Droplet Breakup.
Definition: TAB.H:69
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
static word suffix(const word &fldName, const word &fileExt=word::null)
Name suffix based on fieldName (underscore separator)
Holds list of sampling positions.
Definition: coordSet.H:56
const word & name() const noexcept
The coord-set name.
Definition: coordSet.H:134
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Abstract base-class for Time/database function objects.
word scopedName(const word &name) const
Return a scoped (prefixed) name.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Computes the natural logarithm of an input volScalarField.
Definition: log.H:230
Computes the magnitude of an input field.
Definition: mag.H:153
const objectRegistry & obr_
Reference to the region objectRegistry.
Creates a droplet size distribution via interrogating a continuous phase fraction field.
virtual bool write()
Calculate the regionSizeDistribution and write.
virtual bool read(const dictionary &)
Read the regionSizeDistribution data.
Base class for writing single files from the function objects.
Definition: writeFile.H:120
fileName baseTimeDir() const
Return the base directory for the current time value.
Definition: writeFile.C:73
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField on the opposite patch of a coupled patch.
Definition: fvPatchField.H:453
virtual bool coupled() const
Return true if this patch field is coupled.
Definition: fvPatchField.H:350
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:237
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:362
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:167
const Time & time() const noexcept
Return time registry.
static const char *const componentNames[]
Definition: bool.H:104
static constexpr direction nComponents
Number of components in bool is 1.
Definition: bool.H:98
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:364
int myProcNo() const noexcept
Return processor number.
virtual bool write(const bool valid=true) const
Write using setting from DB.
This class separates the mesh into distinct unconnected regions, each of which is then given a label ...
Definition: regionSplit.H:144
label nRegions() const
Return total number of regions.
Definition: regionSplit.H:294
splitCell * master() const
Definition: splitCell.H:113
A class for managing temporary objects.
Definition: tmp.H:65
const T & cref() const
Definition: tmpI.H:213
void reset(tmp< T > &&other) noexcept
Clear existing and transfer ownership.
Definition: tmpI.H:314
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:54
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
volScalarField & p
engineTime & runTime
word outputName("finiteArea-edges.obj")
Volume integrate volField creating a volField.
constexpr scalar pi(M_PI)
dimensioned< Type > domainIntegrate(const GeometricField< Type, fvPatchField, volMesh > &vf)
Namespace for OpenFOAM.
List< word > wordList
A List of words.
Definition: fileName.H:63
Type gSum(const FieldField< Field, Type > &f)
const dimensionSet dimless
Dimensionless.
List< label > labelList
A List of labels.
Definition: List.H:66
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
messageStream Info
Information stream (stdout output on master, null elsewhere)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensionedScalar sqrt(const dimensionedScalar &ds)
uint8_t direction
Definition: direction.H:56
static Map< Type > regionSum(const regionSplit &regions, const Field< Type > &fld)
Field< label > labelField
Specialisation of Field<T> for label.
Definition: labelField.H:52
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
static List< Type > extractData(const labelUList &keys, const Map< Type > &regionData)
dimensionedScalar cbrt(const dimensionedScalar &ds)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
void divide(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
volScalarField & alpha
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
#define forAllIters(container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:260
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition: stdFoam.H:278
word setFormat(propsDict.getOrDefault< word >("setFormat", "vtk"))