Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
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"
35
36using namespace Foam::constant;
38// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39
40namespace Foam
41{
42 namespace functionObjects
43 {
46 (
50 );
51 }
52}
53
55// * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
56
57namespace Foam
58{
59
60template<class Type>
61static 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 }
73
74 return regionToSum;
75}
76
77
78static 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
95template<class Type>
96static 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
112void 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
199Foam::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
229Foam::tmp<Foam::scalarField>
230Foam::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
250void 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 (
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
324void 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))
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 (
445 (
447 (
448 alphaName_,
449 mesh_.time().timeName(),
450 mesh_,
454 ),
455 mesh_
456 )
457 );
458 }
459 const auto& alpha = talpha();
460
461 Log << " Volume of alpha = "
462 << fvc::domainIntegrate(alpha).value()
463 << endl;
464
465 const scalar meshVol = gSum(mesh_.V());
466 const scalar maxDropletVol = 1.0/6.0*mathematical::pi*pow3(maxDiam_);
467 const scalar delta = (maxDiam_-minDiam_)/nBins_;
468
469 Log << " Mesh volume = " << meshVol << nl
470 << " Maximum droplet diameter = " << maxDiam_ << nl
471 << " Maximum droplet volume = " << maxDropletVol
472 << endl;
473
474
475 // Determine blocked faces
476 boolList blockedFace(mesh_.nFaces(), false);
477 // label nBlocked = 0;
478
479 {
480 for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
481 {
482 scalar ownVal = alpha[mesh_.faceOwner()[facei]];
483 scalar neiVal = alpha[mesh_.faceNeighbour()[facei]];
484
485 if
486 (
487 (ownVal < threshold_ && neiVal > threshold_)
488 || (ownVal > threshold_ && neiVal < threshold_)
489 )
490 {
491 blockedFace[facei] = true;
492 // ++nBlocked;
493 }
494 }
495
496 // Block coupled faces
497 forAll(alpha.boundaryField(), patchi)
498 {
499 const fvPatchScalarField& fvp = alpha.boundaryField()[patchi];
500 if (fvp.coupled())
501 {
502 tmp<scalarField> townFld(fvp.patchInternalField());
503 tmp<scalarField> tnbrFld(fvp.patchNeighbourField());
504 const auto& ownFld = townFld();
505 const auto& nbrFld = tnbrFld();
506
507 label start = fvp.patch().patch().start();
508
509 forAll(ownFld, i)
510 {
511 scalar ownVal = ownFld[i];
512 scalar neiVal = nbrFld[i];
513
514 if
515 (
516 (ownVal < threshold_ && neiVal > threshold_)
517 || (ownVal > threshold_ && neiVal < threshold_)
518 )
519 {
520 blockedFace[start+i] = true;
521 // ++nBlocked;
522 }
523 }
524 }
525 }
526 }
527
528
529 regionSplit regions(mesh_, blockedFace);
530
531 Log << " Determined " << regions.nRegions()
532 << " disconnected regions" << endl;
533
534
535 if (debug)
536 {
537 volScalarField region
538 (
539 mesh_.newIOobject("region"),
540 mesh_,
542 );
543
544 Info<< " Dumping region as volScalarField to "
545 << region.name() << endl;
546
547 forAll(regions, celli)
548 {
549 region[celli] = regions[celli];
550 }
551 region.correctBoundaryConditions();
552 region.write();
553 }
554
555
556 // Determine regions connected to supplied patches
557 const labelHashSet patchRegions(findPatchRegions(regions));
558
559 // Sum all regions
560 const scalarField alphaVol(alpha.primitiveField()*mesh_.V());
561 Map<scalar> allRegionVolume(regionSum(regions, mesh_.V()));
562 Map<scalar> allRegionAlphaVolume(regionSum(regions, alphaVol));
563 Map<label> allRegionNumCells(regionSum(regions, mesh_.nCells()));
564
565 if (debug)
566 {
567 Info<< " " << token::TAB << "Region"
568 << token::TAB << "Volume(mesh)"
569 << token::TAB << "Volume(" << alpha.name() << "):"
570 << token::TAB << "nCells"
571 << nl;
572 scalar meshSumVol = 0.0;
573 scalar alphaSumVol = 0.0;
574 label nCells = 0;
575
576 auto vIter = allRegionVolume.cbegin();
577 auto aIter = allRegionAlphaVolume.cbegin();
578 auto numIter = allRegionNumCells.cbegin();
579 for
580 (
581 ;
582 vIter.good() && aIter.good();
583 ++vIter, ++aIter, ++numIter
584 )
585 {
586 Info<< " " << token::TAB << vIter.key()
587 << token::TAB << vIter()
588 << token::TAB << aIter()
589 << token::TAB << numIter()
590 << nl;
591
592 meshSumVol += vIter();
593 alphaSumVol += aIter();
594 nCells += numIter();
595 }
596 Info<< " " << token::TAB << "Total:"
597 << token::TAB << meshSumVol
598 << token::TAB << alphaSumVol
599 << token::TAB << nCells
600 << endl;
601 }
602
603
604 if (log)
605 {
606 Info<< " Patch connected regions (liquid core):" << nl;
607 Info<< token::TAB << " Region"
608 << token::TAB << "Volume(mesh)"
609 << token::TAB << "Volume(" << alpha.name() << "):"
610 << nl;
611
612 for (const label regioni : patchRegions.sortedToc())
613 {
614 Info<< " " << token::TAB << regioni
615 << token::TAB << allRegionVolume[regioni]
616 << token::TAB << allRegionAlphaVolume[regioni] << nl;
617
618 }
619 Info<< endl;
620 }
621
622 if (log)
623 {
624 Info<< " Background regions:" << nl;
625 Info<< " " << token::TAB << "Region"
626 << token::TAB << "Volume(mesh)"
627 << token::TAB << "Volume(" << alpha.name() << "):"
628 << nl;
629
630 auto vIter = allRegionVolume.cbegin();
631 auto aIter = allRegionAlphaVolume.cbegin();
632
633 for
634 (
635 ;
636 vIter.good() && aIter.good();
637 ++vIter, ++aIter
638 )
639 {
640 if
641 (
642 !patchRegions.found(vIter.key())
643 && vIter() >= maxDropletVol
644 )
645 {
646 Info<< " " << token::TAB << vIter.key()
647 << token::TAB << vIter()
648 << token::TAB << aIter() << nl;
649 }
650 }
651 Info<< endl;
652 }
653
654
655
656 // Split alpha field
657 // ~~~~~~~~~~~~~~~~~
658 // Split into
659 // - liquidCore : region connected to inlet patches
660 // - per region a volume : for all other regions
661 // - backgroundAlpha : remaining alpha
662 writeAlphaFields(regions, patchRegions, allRegionVolume, alpha);
663
664
665 // Extract droplet-only allRegionVolume, i.e. delete liquid core
666 // (patchRegions) and background regions from maps.
667 // Note that we have to use mesh volume (allRegionVolume) and not
668 // allRegionAlphaVolume since background might not have alpha in it.
669 // Deleting regions where the volume-alpha-weighted is lower than
670 // threshold
671 forAllIters(allRegionVolume, vIter)
672 {
673 const label regioni = vIter.key();
674 if
675 (
676 patchRegions.found(regioni)
677 || vIter() >= maxDropletVol
678 || (allRegionAlphaVolume[regioni]/vIter() < threshold_)
679 )
680 {
681 allRegionVolume.erase(vIter);
682 allRegionAlphaVolume.erase(regioni);
683 allRegionNumCells.erase(regioni);
684 }
685 }
686
687 if (allRegionVolume.size())
688 {
689 // Construct mids of bins for plotting
690 pointField xBin(nBins_, Zero);
691
692 {
693 scalar x = 0.5*delta;
694 for (point& p : xBin)
695 {
696 p.x() = x;
697 x += delta;
698 }
699 }
700
701 const coordSet coords("diameter", "x", xBin, mag(xBin));
702
703
704 // Get in region order the alpha*volume and diameter
705 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
706
707 const labelList sortedRegions = allRegionAlphaVolume.sortedToc();
708
709 scalarField sortedVols
710 (
711 extractData(sortedRegions, allRegionAlphaVolume)
712 );
713
714 vectorField centroids(sortedVols.size(), Zero);
715
716 // Check if downstream bins are calculated
717 if (isoPlanes_)
718 {
719 vectorField alphaDistance
720 (
721 (alpha.primitiveField()*mesh_.V())
722 *(mesh_.C().primitiveField() - origin_)()
723 );
724
725 Map<vector> allRegionAlphaDistance
726 (
728 (
729 regions,
730 alphaDistance
731 )
732 );
733
734 // 2. centroid
735 vectorField sortedMoment
736 (
737 extractData(sortedRegions, allRegionAlphaDistance)
738 );
739
740 centroids = sortedMoment/sortedVols + origin_;
741
742 // Bin according to centroid
743 scalarField distToPlane((centroids - origin_) & direction_);
744
745 vectorField radialDistToOrigin
746 (
747 (centroids - origin_) - (distToPlane*direction_)
748 );
749
750 const scalar deltaX = maxDownstream_/nDownstreamBins_;
751 labelList downstreamIndices(distToPlane.size(), -1);
752 forAll(distToPlane, i)
753 {
754 if
755 (
756 (mag(radialDistToOrigin[i]) < maxDiameter_)
757 && (distToPlane[i] < maxDownstream_)
758 )
759 {
760 downstreamIndices[i] = distToPlane[i]/deltaX;
761 }
762 }
763
764 scalarField binDownCount(nDownstreamBins_, Zero);
765 forAll(distToPlane, i)
766 {
767 if (downstreamIndices[i] != -1)
768 {
769 binDownCount[downstreamIndices[i]] += 1.0;
770 }
771 }
772
773 // Write
774 if (UPstream::master())
775 {
776 // Construct mids of bins for plotting
777 pointField xBin(nDownstreamBins_, Zero);
778
779 {
780 scalar x = 0.5*deltaX;
781 for (point& p : xBin)
782 {
783 p.x() = x;
784 x += deltaX;
785 }
786 }
787
788 const coordSet coords("distance", "x", xBin, mag(xBin));
789
790 auto& writer = formatterPtr_();
791 writer.nFields(1);
792
793 writer.open
794 (
795 coords,
796 writeFile::baseTimeDir() / (coords.name() + "_isoPlanes")
797 );
798
799 writer.write("isoPlanes", binDownCount);
800 writer.close(true);
801 }
802
803 // Write to log
804 if (log)
805 {
806 Info<< " Iso-planes Bins:" << nl
807 << " " << token::TAB << "Bin"
808 << token::TAB << "Min distance"
809 << token::TAB << "Count:"
810 << nl;
811
812 scalar delta = 0.0;
813 forAll(binDownCount, bini)
814 {
815 Info<< " " << token::TAB << bini
816 << token::TAB << delta
817 << token::TAB << binDownCount[bini] << nl;
818 delta += deltaX;
819 }
820 Info<< endl;
821
822 }
823 }
824
825 // Calculate the diameters
826 scalarField sortedDiameters(sortedVols.size());
827 forAll(sortedDiameters, i)
828 {
829 sortedDiameters[i] = Foam::cbrt
830 (
831 sortedVols[i]
833 );
834 }
835
836 // Determine the bin index for all the diameters
837 labelList indices(sortedDiameters.size());
838 forAll(sortedDiameters, i)
839 {
840 indices[i] = (sortedDiameters[i]-minDiam_)/delta;
841 }
842
843 // Calculate the counts per diameter bin
844 scalarField binCount(nBins_, Zero);
845 forAll(sortedDiameters, i)
846 {
847 binCount[indices[i]] += 1.0;
848 }
849
850 // Write counts
851 if (UPstream::master())
852 {
853 auto& writer = formatterPtr_();
854 writer.nFields(1);
855
856 writer.open
857 (
858 coords,
859 writeFile::baseTimeDir() / (coords.name() + "_count")
860 );
861
862 writer.write("count", binCount);
863 writer.close(true);
864 }
865
866 // Write to log
867 if (log)
868 {
869 Info<< " Bins:" << nl
870 << " " << token::TAB << "Bin"
871 << token::TAB << "Min diameter"
872 << token::TAB << "Count:"
873 << nl;
874
875 scalar diam = 0.0;
876 forAll(binCount, bini)
877 {
878 Info<< " " << token::TAB << bini
879 << token::TAB << diam
880 << token::TAB << binCount[bini] << nl;
881
882 diam += delta;
883 }
884
885 Info<< endl;
886 }
887
888
889 // Write average and deviation of droplet volume.
890 writeGraphs
891 (
892 "volume", // name of field
893 sortedVols, // per region field data
894
895 indices, // per region the bin index
896 binCount, // per bin number of regions
897 coords // graph data for bins
898 );
899
900 // Collect some more fields
901 {
902 for
903 (
904 const volScalarField& vfield
905 : obr_.csorted<volScalarField>(fields_)
906 )
907 {
908 const word& fldName = vfield.name();
909
910 Log << " Scalar field " << fldName << endl;
911
912 tmp<Field<scalar>> tfld(vfield.primitiveField());
913 const auto& fld = tfld();
914
915 writeGraphs
916 (
917 fldName, // name of field
918 alphaVol*fld, // per cell field data
919
920 regions, // per cell the region(=droplet)
921 sortedRegions, // valid regions in sorted order
922 1.0/sortedVols, // per region normalisation
923
924 indices, // index of bin for each region
925 binCount, // per bin number of regions
926 coords // graph data for bins
927 );
928 }
929 }
930
931 {
932 for
933 (
934 const volVectorField& vfield
935 : obr_.csorted<volVectorField>(fields_)
936 )
937 {
938 const word& fldName = vfield.name();
939
940 Log << " Vector field " << fldName << endl;
941
942 tmp<Field<vector>> tfld(vfield.primitiveField());
943
944 if (csysPtr_)
945 {
946 Log << "Transforming vector field " << fldName
947 << " with coordinate system "
948 << csysPtr_->name() << endl;
949
950 tfld = csysPtr_->localVector(tfld());
951 }
952 const auto& fld = tfld();
953
954 // Components
955
956 for (direction cmpt = 0; cmpt < vector::nComponents; ++cmpt)
957 {
958 writeGraphs
959 (
960 fldName + vector::componentNames[cmpt],
961 alphaVol*fld.component(cmpt),// per cell field data
962
963 regions, // per cell the region(=droplet)
964 sortedRegions, // valid regions in sorted order
965 1.0/sortedVols, // per region normalisation
966
967 indices, // index of bin for each region
968 binCount, // per bin number of regions
969 coords // graph data for bins
970 );
971 }
972
973 // Magnitude
974 writeGraphs
975 (
976 fldName + "mag", // name of field
977 alphaVol*mag(fld), // per cell field data
978
979 regions, // per cell the region(=droplet)
980 sortedRegions, // valid regions in sorted order
981 1.0/sortedVols, // per region normalisation
982
983 indices, // index of bin for each region
984 binCount, // per bin number of regions
985 coords // graph data for bins
986 );
987 }
988 }
989 }
990
991 return true;
992}
993
994
995// ************************************************************************* //
scalar delta
#define Log
Definition PDRblock.C:28
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
word setFormat(propsDict.getOrDefault< word >("setFormat", "vtk"))
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
vtk::lineWriter writer(edgeCentres, edgeList::null(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edgesCentres")))
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
void correctBoundaryConditions()
Correct boundary field.
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition HashTable.C:156
const_iterator cbegin() const
const_iterator set to the beginning of the HashTable
bool found(const Key &key) const
Same as contains().
Definition HashTable.H:1370
label size() const noexcept
The number of elements in table.
Definition HashTable.H:358
bool erase(const iterator &iter)
Erase an entry specified by given iterator.
Definition HashTable.C:489
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
A HashTable to objects of type <T> with a label key.
Definition Map.H:54
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
static void combineReduce(T &value, CombineOp cop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce) applying cop to inplace combine value from different processors.
static void mapCombineReduce(Container &values, CombineOp cop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Forwards to Pstream::mapReduce with an in-place cop.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
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:1697
static const char *const componentNames[]
static constexpr direction nComponents
Number of components in this vector space.
static autoPtr< coordSetWriter > New(const word &writeFormat)
Return a reference to the selected writer.
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:52
const word & name() const noexcept
The coord-set name.
Definition coordSet.H:152
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Abstract base-class for Time/database function objects.
virtual bool read(const dictionary &dict)
Read and set the function object if its data have changed.
word scopedName(const word &name) const
Return a scoped (prefixed) name.
bool log
Flag to write log into Info.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
const fvMesh & mesh_
Reference to the fvMesh.
fvMeshFunctionObject(const fvMeshFunctionObject &)=delete
No copy construct.
Computes the natural logarithm of an input volScalarField.
Definition log.H:212
const objectRegistry & obr_
Reference to the region objectRegistry.
Creates a droplet size distribution via interrogating a continuous phase fraction field.
regionSizeDistribution(const word &name, const Time &runTime, const dictionary &)
Construct for given objectRegistry and dictionary.
virtual bool read(const dictionary &dict)
Read the function-object dictionary.
virtual bool execute()
Execute the function-object operations.
virtual bool write()
Write the function-object results.
Base class for writing single files from the function objects.
Definition writeFile.H:113
writeFile(const objectRegistry &obr, const fileName &prefix, const word &name="undefined", const bool writeToFile=true, const string &ext=".dat")
Construct from objectRegistry, prefix, fileName.
Definition writeFile.C:200
virtual bool read(const dictionary &dict)
Read.
Definition writeFile.C:240
fileName baseTimeDir() const
Return the base directory for the current time value.
Definition writeFile.C:66
virtual bool coupled() const
True if the patch field is coupled.
static const word & calculatedType() noexcept
The type name for calculated patch fields.
const fvPatch & patch() const noexcept
Return the patch.
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField on the opposite patch of a coupled patch.
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch.
const polyPatch & patch() const noexcept
Return the polyPatch.
Definition fvPatch.H:202
label start() const noexcept
Return start label of this patch in the polyMesh face list.
Definition polyPatch.H:446
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
This class separates the mesh into distinct unconnected regions, each of which is then given a label ...
label nRegions() const
Return total number of regions.
A class for managing temporary objects.
Definition tmp.H:75
const T & cref() const
Return const reference to the object or to the contents of a (non-null) managed pointer.
Definition tmpI.H:221
void reset(tmp< T > &&other) noexcept
Clear existing and transfer ownership.
Definition tmpI.H:338
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
@ TAB
Tab [isspace].
Definition token.H:142
A List of wordRe with additional matching capabilities.
Definition wordRes.H:56
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
engineTime & runTime
word outputName("finiteArea-edges.obj")
Volume integrate volField creating a volField.
auto & name
constexpr scalar pi(M_PI)
Different types of constants.
Namespace for handling debugging switches.
Definition debug.C:45
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
dimensioned< Type > domainIntegrate(const GeometricField< Type, fvPatchField, volMesh > &vf)
Namespace for OpenFOAM.
List< word > wordList
List of word.
Definition fileName.H:60
Type gSum(const FieldField< Field, Type > &f)
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition int32.H:127
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
dimensionedScalar pow3(const dimensionedScalar &ds)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
void divide(DimensionedField< Type, GeoMesh > &result, const DimensionedField< Type, GeoMesh > &f1, const DimensionedField< scalar, GeoMesh > &f2)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Field< vector > vectorField
Specialisation of Field<T> for vector.
uint8_t direction
Definition direction.H:49
static Map< Type > regionSum(const regionSplit &regions, const Field< Type > &fld)
vector point
Point is a vector.
Definition point.H:37
List< bool > boolList
A List of bools.
Definition List.H:60
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
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 expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
vectorField pointField
pointField is a vectorField.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
UList< label > labelUList
A UList of labels.
Definition UList.H:75
fvPatchField< scalar > fvPatchScalarField
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
volScalarField & alpha
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define forAllIters(container, iter)
Iterate across all elements in the container object.
Definition stdFoam.H:214