decompositionMethod.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) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2015-2021 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
29#include "decompositionMethod.H"
30#include "globalIndex.H"
31#include "syncTools.H"
32#include "faceSet.H"
33#include "regionSplit.H"
34#include "localPointRegion.H"
35#include "minData.H"
36#include "BitOps.H"
37#include "FaceCellWave.H"
38
39// Compatibility (MAY-2014)
44
45// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
46
47namespace Foam
48{
51
52} // End namespace Foam
53
54
55// * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
56
57namespace Foam
58{
59
60// Find named coefficents dictionary, or use default "coeffs"
61static inline const dictionary* cfindCoeffsDict
62(
63 const dictionary& dict,
64 const word& coeffsName,
65 const bool allowDefault
66)
67{
68 const dictionary* dictptr = dict.findDict(coeffsName);
69 if (!dictptr && allowDefault)
70 {
71 dictptr = dict.findDict("coeffs");
72 }
73 return dictptr;
74}
75
76} // End namespace Foam
77
78
79// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
80
82(
83 const dictionary& decompDict,
84 const word& regionName
85)
86{
87 const label nDomainsGlobal = decompDict.get<label>("numberOfSubdomains");
88
89 if (!regionName.empty())
90 {
91 const dictionary& regionDict =
92 optionalRegionDict(decompDict, regionName);
93
94 label nDomainsRegion;
95 if (regionDict.readIfPresent("numberOfSubdomains", nDomainsRegion))
96 {
97 if (nDomainsRegion >= 1 && nDomainsRegion <= nDomainsGlobal)
98 {
99 return nDomainsRegion;
100 }
101
103 << "ignoring out of range numberOfSubdomains "
104 << nDomainsRegion << " for region " << regionName
105 << nl << nl
106 << endl;
107 }
108 }
109
110 return nDomainsGlobal;
111}
112
113
115(
116 const dictionary& decompDict,
117 const word& regionName
118)
119{
120 const dictionary* dictptr = nullptr;
121 if
122 (
123 !regionName.empty()
124 && (dictptr = decompDict.findDict("regions")) != nullptr
125 )
126 {
127 dictptr = dictptr->findDict(regionName);
128 }
129 return (dictptr ? *dictptr : dictionary::null);
130}
131
132
133// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
134
135bool Foam::decompositionMethod::constraintCompat(const word& modelType) const
136{
137 bool usable = decompDict_.found(modelType);
138 if (!usable)
139 {
140 return false;
141 }
142
143 for (const auto& item : constraints_)
144 {
145 if (modelType == item.type())
146 {
147 usable = false;
148 break;
149 }
150 }
151
152 if (usable)
153 {
154 Warning
155 << nl << " Using '" << modelType
156 << "' constraint specification." << nl;
157 }
158 else
159 {
160 Warning
161 << nl << " Ignoring '" << modelType
162 << "' constraint specification - was already specified." << nl;
163 }
164
165 // The syntax changed MAY-2014
166 error::warnAboutAge("constraint keyword", 1406);
167
168 return usable;
169}
170
171
172void Foam::decompositionMethod::readConstraints()
173{
174 constraints_.clear();
175
176 const dictionary* dictptr = decompDict_.findDict("constraints");
177
178 if (dictptr)
179 {
180 for (const entry& dEntry : *dictptr)
181 {
182 if (!dEntry.isDict()) // safety
183 {
184 // Ignore or warn
185 continue;
186 }
187
188 const dictionary& dict = dEntry.dict();
189
190 if (dict.getOrDefault("enabled", true))
191 {
192 constraints_.append(decompositionConstraint::New(dict));
193 }
194 }
195 }
196
197 // Backwards compatibility (MAY-2014)
198 if (constraintCompat("preserveBaffles"))
199 {
200 constraints_.append
201 (
202 new decompositionConstraints::preserveBaffles()
203 );
204 }
205
206 if (constraintCompat("preservePatches"))
207 {
208 constraints_.append
209 (
210 new decompositionConstraints::preservePatches
211 (
212 decompDict_.get<wordRes>("preservePatches")
213 )
214 );
215 }
216
217 if (constraintCompat("preserveFaceZones"))
218 {
219 constraints_.append
220 (
221 new decompositionConstraints::preserveFaceZones
223 decompDict_.get<wordRes>("preserveFaceZones")
224 )
225 );
226 }
227
228 if (constraintCompat("singleProcessorFaceSets"))
230 constraints_.append
231 (
233 (
234 decompDict_.lookup("singleProcessorFaceSets")
235 )
236 );
237 }
238}
239
240// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
241
243(
245 const word& coeffsName,
246 int select
247)
248{
249 const bool allowDefault = !(select & selectionType::EXACT);
250
251 const dictionary* dictptr =
252 cfindCoeffsDict(dict, coeffsName, allowDefault);
253
254 if (dictptr)
255 {
256 return *dictptr;
257 }
259 // Not found
260 if (select & selectionType::MANDATORY)
261 {
263 << "'" << coeffsName << "' dictionary not found in dictionary "
264 << dict.name() << endl
266 }
268 if (select & selectionType::NULL_DICT)
269 {
270 return dictionary::null;
271 }
272
273 return dict; // Return original dictionary
274}
275
276
278(
279 const word& coeffsName,
280 int select
281) const
282{
283 const bool allowDefault = !(select & selectionType::EXACT);
284
285 const dictionary* dictptr = nullptr;
286
287 if (!decompRegionDict_.empty())
288 {
289 // Region-specific dictionary
290 dictptr = cfindCoeffsDict(decompRegionDict_, coeffsName, allowDefault);
291 }
292 if (!dictptr)
293 {
294 // General
295 dictptr = cfindCoeffsDict(decompDict_, coeffsName, allowDefault);
296 }
297
298 if (dictptr)
299 {
300 return *dictptr;
301 }
302
303 // Not found
304 if (select & selectionType::MANDATORY)
305 {
307 << "'" << coeffsName << "' dictionary not found in dictionary "
308 << decompDict_.name() << endl
310 }
311
312 if (select & selectionType::NULL_DICT)
313 {
314 return dictionary::null;
315 }
316
317 return decompDict_; // Return general dictionary
318}
319
320
321// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
322
324(
325 const dictionary& decompDict,
326 const word& regionName
327)
328:
329 decompDict_(decompDict),
330 decompRegionDict_
331 (
332 optionalRegionDict(decompDict_, regionName)
333 ),
334 nDomains_(nDomains(decompDict, regionName))
335{
336 readConstraints();
337}
338
339
340// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
341
343(
344 const dictionary& decompDict,
345 const word& regionName
346)
347{
348 word methodType(decompDict.get<word>("method"));
349
350 const dictionary& regionDict = optionalRegionDict(decompDict, regionName);
351 regionDict.readIfPresent("method", methodType);
352
353 auto* ctorPtr = dictionaryConstructorTable(methodType);
354
355 if (!ctorPtr)
356 {
358 (
359 decompDict,
360 "decompositionMethod",
361 methodType,
362 *dictionaryConstructorTablePtr_
363 ) << exit(FatalIOError);
364 }
365
366 // verbose
368 Info<< "Decomposition method " << methodType
369 << " [" << (nDomains(decompDict, regionName)) << ']';
370
371 if (!regionName.empty())
372 {
373 Info<< " (region " << regionName << ')';
374 }
375 Info<< endl;
376 }
377
378 return autoPtr<decompositionMethod>(ctorPtr(decompDict, regionName));
379}
380
381
383(
384 const polyMesh& mesh,
385 const pointField& points
386) const
388 scalarField weights(points.size(), scalar(1));
389
390 return decompose(mesh, points, weights);
391}
392
393
395(
396 const polyMesh& mesh,
397 const labelList& fineToCoarse,
398 const pointField& coarsePoints,
399 const scalarField& coarseWeights
400) const
401{
402 CompactListList<label> coarseCellCells;
403 calcCellCells
404 (
405 mesh,
406 fineToCoarse,
407 coarsePoints.size(),
408 true, // use global cell labels
409 coarseCellCells
410 );
411
412 // Decompose based on agglomerated points
413 labelList coarseDistribution
414 (
415 decompose
416 (
417 coarseCellCells.unpack(),
418 coarsePoints,
419 coarseWeights
420 )
421 );
422
423 // Rework back into decomposition for original mesh_
424 labelList fineDistribution(fineToCoarse.size());
425
426 forAll(fineDistribution, i)
427 {
428 fineDistribution[i] = coarseDistribution[fineToCoarse[i]];
429 }
430
431 return fineDistribution;
432}
433
434
436(
437 const polyMesh& mesh,
438 const labelList& fineToCoarse,
439 const pointField& coarsePoints
440) const
441{
442 scalarField weights(coarsePoints.size(), scalar(1));
443
444 return decompose
445 (
446 mesh,
447 fineToCoarse,
448 coarsePoints,
449 weights
450 );
451}
452
453
455(
456 const labelListList& globalCellCells,
457 const pointField& cc
458) const
459{
460 scalarField weights(cc.size(), scalar(1));
461
462 return decompose(globalCellCells, cc, weights);
463}
464
465
467(
468 const polyMesh& mesh,
469 const labelList& agglom,
470 const label nLocalCoarse,
471 const bool parallel,
472 CompactListList<label>& cellCells
473)
474{
475 const labelList& faceOwner = mesh.faceOwner();
476 const labelList& faceNeighbour = mesh.faceNeighbour();
478
479
480 // Create global cell numbers
481 // ~~~~~~~~~~~~~~~~~~~~~~~~~~
482
483 const globalIndex globalAgglom(nLocalCoarse, Pstream::worldComm, parallel);
484
485
486 // Get agglomerate owner on other side of coupled faces
487 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
488
489 labelList globalNeighbour(mesh.nBoundaryFaces());
490
491 for (const polyPatch& pp : patches)
492 {
493 if (pp.coupled() && (parallel || !isA<processorPolyPatch>(pp)))
494 {
495 label facei = pp.start();
496 label bFacei = pp.start() - mesh.nInternalFaces();
497
498 forAll(pp, i)
499 {
500 globalNeighbour[bFacei] = globalAgglom.toGlobal
501 (
502 agglom[faceOwner[facei]]
503 );
504
505 ++facei;
506 ++bFacei;
507 }
508 }
509 }
510
511 // Get the cell on the other side of coupled patches
512 syncTools::swapBoundaryFaceList(mesh, globalNeighbour);
513
514
515 // Count number of faces (internal + coupled)
516 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
517
518 // Number of faces per coarse cell
519 labelList nFacesPerCell(nLocalCoarse, Zero);
520
521 for (label facei = 0; facei < mesh.nInternalFaces(); ++facei)
522 {
523 const label own = agglom[faceOwner[facei]];
524 const label nei = agglom[faceNeighbour[facei]];
525
526 nFacesPerCell[own]++;
527 nFacesPerCell[nei]++;
528 }
529
530 for (const polyPatch& pp : patches)
531 {
532 if (pp.coupled() && (parallel || !isA<processorPolyPatch>(pp)))
533 {
534 label facei = pp.start();
535 label bFacei = pp.start()-mesh.nInternalFaces();
536
537 forAll(pp, i)
538 {
539 const label own = agglom[faceOwner[facei]];
540 const label globalNei = globalNeighbour[bFacei];
541
542 if
543 (
544 !globalAgglom.isLocal(globalNei)
545 || globalAgglom.toLocal(globalNei) != own
546 )
547 {
548 nFacesPerCell[own]++;
549 }
550
551 ++facei;
552 ++bFacei;
553 }
554 }
555 }
556
557
558 // Fill in offset and data
559 // ~~~~~~~~~~~~~~~~~~~~~~~
560
561 cellCells.setSize(nFacesPerCell);
562
563 nFacesPerCell = 0;
564
565 labelList& m = cellCells.m();
566 const labelList& offsets = cellCells.offsets();
567
568 // For internal faces is just offsetted owner and neighbour
569 for (label facei = 0; facei < mesh.nInternalFaces(); ++facei)
570 {
571 const label own = agglom[faceOwner[facei]];
572 const label nei = agglom[faceNeighbour[facei]];
573
574 m[offsets[own] + nFacesPerCell[own]++] = globalAgglom.toGlobal(nei);
575 m[offsets[nei] + nFacesPerCell[nei]++] = globalAgglom.toGlobal(own);
576 }
577
578 // For boundary faces is offsetted coupled neighbour
579 for (const polyPatch& pp : patches)
580 {
581 if (pp.coupled() && (parallel || !isA<processorPolyPatch>(pp)))
582 {
583 label facei = pp.start();
584 label bFacei = pp.start()-mesh.nInternalFaces();
585
586 forAll(pp, i)
587 {
588 const label own = agglom[faceOwner[facei]];
589 const label globalNei = globalNeighbour[bFacei];
590
591 if
592 (
593 !globalAgglom.isLocal(globalNei)
594 || globalAgglom.toLocal(globalNei) != own
595 )
596 {
597 m[offsets[own] + nFacesPerCell[own]++] = globalNei;
598 }
599
600 ++facei;
601 ++bFacei;
602 }
603 }
604 }
605
606
607 // Check for duplicates connections between cells
608 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
609 // Done as postprocessing step since we now have cellCells.
610
611 if (cellCells.size() == 0)
612 {
613 return;
614 }
615
616 label newIndex = 0;
617 labelHashSet nbrCells;
618
619 label startIndex = cellCells.offsets()[0];
620
621 forAll(cellCells, celli)
622 {
623 nbrCells.clear();
624 nbrCells.insert(globalAgglom.toGlobal(celli));
625
626 const label endIndex = cellCells.offsets()[celli+1];
627
628 for (label i = startIndex; i < endIndex; ++i)
629 {
630 if (nbrCells.insert(cellCells.m()[i]))
631 {
632 cellCells.m()[newIndex++] = cellCells.m()[i];
633 }
634 }
635 startIndex = endIndex;
636 cellCells.offsets()[celli+1] = newIndex;
637 }
638
639 cellCells.m().setSize(newIndex);
640
641 //forAll(cellCells, celli)
642 //{
643 // Pout<< "Original: Coarse cell " << celli << endl;
644 // forAll(mesh.cellCells()[celli], i)
645 // {
646 // Pout<< " nbr:" << mesh.cellCells()[celli][i] << endl;
647 // }
648 // Pout<< "Compacted: Coarse cell " << celli << endl;
649 // const labelUList cCells = cellCells[celli];
650 // forAll(cCells, i)
651 // {
652 // Pout<< " nbr:" << cCells[i] << endl;
653 // }
654 //}
655}
656
657
659(
660 const polyMesh& mesh,
661 const labelList& agglom,
662 const label nLocalCoarse,
663 const bool parallel,
664 CompactListList<label>& cellCells,
665 CompactListList<scalar>& cellCellWeights
666)
667{
668 const labelList& faceOwner = mesh.faceOwner();
669 const labelList& faceNeighbour = mesh.faceNeighbour();
671
672
673 // Create global cell numbers
674 // ~~~~~~~~~~~~~~~~~~~~~~~~~~
675
676 const globalIndex globalAgglom(nLocalCoarse, Pstream::worldComm, parallel);
677
678
679 // Get agglomerate owner on other side of coupled faces
680 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
681
682 labelList globalNeighbour(mesh.nBoundaryFaces());
683
684 for (const polyPatch& pp : patches)
685 {
686 if (pp.coupled() && (parallel || !isA<processorPolyPatch>(pp)))
687 {
688 label facei = pp.start();
689 label bFacei = pp.start() - mesh.nInternalFaces();
690
691 forAll(pp, i)
692 {
693 globalNeighbour[bFacei] = globalAgglom.toGlobal
694 (
695 agglom[faceOwner[facei]]
696 );
697
698 ++facei;
699 ++bFacei;
700 }
701 }
702 }
703
704 // Get the cell on the other side of coupled patches
705 syncTools::swapBoundaryFaceList(mesh, globalNeighbour);
706
707
708 // Count number of faces (internal + coupled)
709 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
710
711 // Number of faces per coarse cell
712 labelList nFacesPerCell(nLocalCoarse, Zero);
713
714 for (label facei = 0; facei < mesh.nInternalFaces(); ++facei)
715 {
716 const label own = agglom[faceOwner[facei]];
717 const label nei = agglom[faceNeighbour[facei]];
718
719 nFacesPerCell[own]++;
720 nFacesPerCell[nei]++;
721 }
722
723 for (const polyPatch& pp : patches)
724 {
725 if (pp.coupled() && (parallel || !isA<processorPolyPatch>(pp)))
726 {
727 label facei = pp.start();
728 label bFacei = pp.start() - mesh.nInternalFaces();
729
730 forAll(pp, i)
731 {
732 const label own = agglom[faceOwner[facei]];
733 const label globalNei = globalNeighbour[bFacei];
734
735 if
736 (
737 !globalAgglom.isLocal(globalNei)
738 || globalAgglom.toLocal(globalNei) != own
739 )
740 {
741 nFacesPerCell[own]++;
742 }
743
744 ++facei;
745 ++bFacei;
746 }
747 }
748 }
749
750
751 // Fill in offset and data
752 // ~~~~~~~~~~~~~~~~~~~~~~~
753
754 cellCells.setSize(nFacesPerCell);
755 cellCellWeights.setSize(nFacesPerCell);
756
757 nFacesPerCell = 0;
758
759 labelList& m = cellCells.m();
760 scalarList& w = cellCellWeights.m();
761 const labelList& offsets = cellCells.offsets();
762
763 // For internal faces is just offsetted owner and neighbour
764 for (label facei = 0; facei < mesh.nInternalFaces(); ++facei)
765 {
766 const label own = agglom[faceOwner[facei]];
767 const label nei = agglom[faceNeighbour[facei]];
768
769 const label ownIndex = offsets[own] + nFacesPerCell[own]++;
770 const label neiIndex = offsets[nei] + nFacesPerCell[nei]++;
771
772 m[ownIndex] = globalAgglom.toGlobal(nei);
773 w[ownIndex] = mag(mesh.faceAreas()[facei]);
774 m[neiIndex] = globalAgglom.toGlobal(own);
775 w[ownIndex] = mag(mesh.faceAreas()[facei]);
776 }
777
778 // For boundary faces is offsetted coupled neighbour
779 for (const polyPatch& pp : patches)
780 {
781 if (pp.coupled() && (parallel || !isA<processorPolyPatch>(pp)))
782 {
783 label facei = pp.start();
784 label bFacei = pp.start()-mesh.nInternalFaces();
785
786 forAll(pp, i)
787 {
788 const label own = agglom[faceOwner[facei]];
789 const label globalNei = globalNeighbour[bFacei];
790
791 if
792 (
793 !globalAgglom.isLocal(globalNei)
794 || globalAgglom.toLocal(globalNei) != own
795 )
796 {
797 const label ownIndex = offsets[own] + nFacesPerCell[own]++;
798 m[ownIndex] = globalNei;
799 w[ownIndex] = mag(mesh.faceAreas()[facei]);
800 }
801
802 ++facei;
803 ++bFacei;
804 }
805 }
806 }
807
808
809 // Check for duplicates connections between cells
810 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
811 // Done as postprocessing step since we now have cellCells.
812
813 if (cellCells.size() == 0)
814 {
815 return;
816 }
817
818 label newIndex = 0;
819 labelHashSet nbrCells;
820
821 label startIndex = cellCells.offsets()[0];
822
823 forAll(cellCells, celli)
824 {
825 nbrCells.clear();
826 nbrCells.insert(globalAgglom.toGlobal(celli));
827
828 const label endIndex = cellCells.offsets()[celli+1];
829
830 for (label i = startIndex; i < endIndex; ++i)
831 {
832 if (nbrCells.insert(cellCells.m()[i]))
833 {
834 cellCells.m()[newIndex] = cellCells.m()[i];
835 cellCellWeights.m()[newIndex] = cellCellWeights.m()[i];
836 newIndex++;
837 }
838 }
839 startIndex = endIndex;
840 cellCells.offsets()[celli+1] = newIndex;
841 cellCellWeights.offsets()[celli+1] = newIndex;
842 }
843
844 cellCells.m().setSize(newIndex);
845 cellCellWeights.m().setSize(newIndex);
846}
847
848
849// NOTE:
850// - alternative calcCellCells that handled explicitConnections was
851// deactivated (2014 or earlier) and finally removed APR-2018.
852
854(
855 const polyMesh& mesh,
856 const scalarField& cellWeights,
857
858 //- Whether owner and neighbour should be on same processor
859 // (takes priority over explicitConnections)
860 const boolList& blockedFace,
861
862 //- Whether whole sets of faces (and point neighbours) need to be kept
863 // on single processor
864 const PtrList<labelList>& specifiedProcessorFaces,
865 const labelList& specifiedProcessor,
866
867 //- Additional connections between boundary faces
868 const List<labelPair>& explicitConnections
869) const
870{
871 // Any weights specified?
872 const bool hasWeights = returnReduce(!cellWeights.empty(), orOp<bool>());
873
874 if (hasWeights && cellWeights.size() != mesh.nCells())
875 {
877 << "Number of weights " << cellWeights.size()
878 << " differs from number of cells " << mesh.nCells()
879 << exit(FatalError);
880 }
881
882 // Any faces not blocked?
883 const bool hasUnblocked =
885 (
886 (!blockedFace.empty() && !BitOps::all(blockedFace)),
887 orOp<bool>()
888 );
889
890
891 // Any non-mesh connections?
892 const label nConnections = returnReduce
893 (
894 explicitConnections.size(),
896 );
897
898
899 // Any processor sets?
900 label nProcSets = 0;
901 for (const labelList& procset : specifiedProcessorFaces)
902 {
903 nProcSets += procset.size();
904 }
905 reduce(nProcSets, sumOp<label>());
906
907
908 // Either do decomposition on cell centres or on agglomeration
909
910 if (!hasUnblocked && !nConnections && !nProcSets)
911 {
912 // No constraints, possibly weights
913
914 return
915 (
916 hasWeights
917 ? decompose(mesh, mesh.cellCentres(), cellWeights)
918 : decompose(mesh, mesh.cellCentres())
919 );
920 }
921
922
923 // The harder work.
924 // When we have processor sets, connections, or blocked faces.
925
926
927 // Determine local regions, separated by blockedFaces
928 regionSplit localRegion(mesh, blockedFace, explicitConnections, false);
929
930 if (debug)
931 {
932 // Only need to count unblocked faces for debugging
933 const label nUnblocked =
934 (
935 hasUnblocked
937 (
938 label(BitOps::count(blockedFace, false)),
940 )
941 : 0
942 );
943
944 Info<< "Constrained decomposition:" << nl
945 << " faces with same owner and neighbour processor : "
946 << nUnblocked << nl
947 << " baffle faces with same owner processor : "
948 << nConnections << nl
949 << " faces all on same processor : "
950 << nProcSets << nl
951 << " split into " << localRegion.nLocalRegions()
952 << " regions."
953 << endl;
954 }
955
956
957 // Gather region weights and determine region cell centres
958 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
959
960 // For the region centre, just take the first cell in the region.
961 // If we average the region centre instead, cyclics could cause
962 // the average domain centre to be outside of domain.
963
964 scalarField regionWeights(localRegion.nLocalRegions(), Zero);
965
966 pointField regionCentres(localRegion.nLocalRegions(), point::max);
967
968 if (hasWeights)
969 {
970 forAll(localRegion, celli)
971 {
972 const label regioni = localRegion[celli];
973
974 regionWeights[regioni] += cellWeights[celli];
975
976 if (regionCentres[regioni] == point::max)
977 {
978 regionCentres[regioni] = mesh.cellCentres()[celli];
979 }
980 }
981 }
982 else
983 {
984 forAll(localRegion, celli)
985 {
986 const label regioni = localRegion[celli];
987
988 regionWeights[regioni] += 1.0;
989
990 if (regionCentres[regioni] == point::max)
991 {
992 regionCentres[regioni] = mesh.cellCentres()[celli];
993 }
994 }
995 }
996
997 // Do decomposition on agglomeration
998 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
999
1000 labelList finalDecomp =
1001 decompose
1002 (
1003 mesh,
1004 localRegion,
1005 regionCentres,
1006 regionWeights
1007 );
1008
1009
1010 // Apply explicitConnections since decompose did not know about them
1011 for (const labelPair& baffle : explicitConnections)
1012 {
1013 const label f0 = baffle.first();
1014 const label f1 = baffle.second();
1015
1016 if (!blockedFace[f0] && !blockedFace[f1])
1017 {
1018 // Note: what if internal faces and owner and neighbour on
1019 // different processor?
1020 // So for now just push owner side proc
1021
1022 const label proci = finalDecomp[mesh.faceOwner()[f0]];
1023
1024 finalDecomp[mesh.faceOwner()[f1]] = proci;
1025 if (mesh.isInternalFace(f1))
1026 {
1027 finalDecomp[mesh.faceNeighbour()[f1]] = proci;
1028 }
1029 }
1030 else if (blockedFace[f0] != blockedFace[f1])
1031 {
1033 << "On explicit connection between faces " << f0
1034 << " and " << f1
1035 << " the two blockedFace status are not equal : "
1036 << blockedFace[f0] << " and " << blockedFace[f1]
1037 << exit(FatalError);
1038 }
1039 }
1040
1041
1042 // blockedFaces corresponding to processor faces need to be handled
1043 // separately since not handled by local regionSplit. We need to
1044 // walk now across coupled faces and make sure to move a whole
1045 // global region across
1046
1047 // This additionally consolidates/compacts the regions numbers globally,
1048 // since that was skipped in the previous regionSplit.
1049 if (Pstream::parRun())
1050 {
1051 // Re-do regionSplit
1052
1053 // Field on cells and faces.
1054 List<minData> cellData(mesh.nCells());
1055 List<minData> faceData(mesh.nFaces());
1056
1057 // Take over blockedFaces by seeding a negative number
1058 // (so is always less than the decomposition)
1059 label nUnblocked = 0;
1060 forAll(blockedFace, facei)
1061 {
1062 if (blockedFace[facei])
1063 {
1064 faceData[facei] = minData(-123);
1065 }
1066 else
1067 {
1068 ++nUnblocked;
1069 }
1070 }
1071
1072 // Seed unblocked faces with destination processor
1073 labelList seedFaces(nUnblocked);
1074 List<minData> seedData(nUnblocked);
1075 nUnblocked = 0;
1076
1077 forAll(blockedFace, facei)
1078 {
1079 if (!blockedFace[facei])
1080 {
1081 const label own = mesh.faceOwner()[facei];
1082 seedFaces[nUnblocked] = facei;
1083 seedData[nUnblocked] = minData(finalDecomp[own]);
1084 nUnblocked++;
1085 }
1086 }
1087
1088
1089 // Propagate information inwards
1090 FaceCellWave<minData> deltaCalc
1091 (
1092 mesh,
1093 seedFaces,
1094 seedData,
1095 faceData,
1096 cellData,
1098 );
1099
1100 // And extract
1101 forAll(finalDecomp, celli)
1102 {
1103 if (cellData[celli].valid(deltaCalc.data()))
1104 {
1105 finalDecomp[celli] = cellData[celli].data();
1106 }
1107 }
1108 }
1109
1110
1111 // For specifiedProcessorFaces rework the cellToProc to enforce
1112 // all on one processor since we can't guarantee that the input
1113 // to regionSplit was a single region.
1114 // E.g. faceSet 'a' with the cells split into two regions
1115 // by a notch formed by two walls
1116 //
1117 // \ /
1118 // \ /
1119 // ---a----+-----a-----
1120 //
1121 //
1122 // Note that reworking the cellToProc might make the decomposition
1123 // unbalanced.
1124 forAll(specifiedProcessorFaces, seti)
1125 {
1126 const labelList& set = specifiedProcessorFaces[seti];
1127
1128 label proci = specifiedProcessor[seti];
1129 if (proci == -1)
1130 {
1131 // If no processor specified - use the one from the 0th element
1132 if (set.size())
1133 {
1134 proci = finalDecomp[mesh.faceOwner()[set[0]]];
1135 }
1136 else
1137 {
1138 // Zero-sized processor (e.g. from redistributePar)
1139 proci = 0;
1140 }
1141 }
1142
1143 for (const label facei : set)
1144 {
1145 const face& f = mesh.faces()[facei];
1146 for (const label pointi : f)
1147 {
1148 const labelList& pFaces = mesh.pointFaces()[pointi];
1149 for (const label pFacei : pFaces)
1150 {
1151 finalDecomp[mesh.faceOwner()[pFacei]] = proci;
1152 if (mesh.isInternalFace(pFacei))
1153 {
1154 finalDecomp[mesh.faceNeighbour()[pFacei]] = proci;
1155 }
1156 }
1157 }
1158 }
1159 }
1160
1161
1162 if (debug && Pstream::parRun())
1163 {
1164 labelList nbrDecomp;
1165 syncTools::swapBoundaryCellList(mesh, finalDecomp, nbrDecomp);
1166
1168 for (const polyPatch& pp : patches)
1169 {
1170 if (pp.coupled())
1171 {
1172 forAll(pp, i)
1173 {
1174 const label facei = pp.start()+i;
1175 const label own = mesh.faceOwner()[facei];
1176 const label bFacei = facei-mesh.nInternalFaces();
1177
1178 if (!blockedFace[facei])
1179 {
1180 const label ownProc = finalDecomp[own];
1181 const label nbrProc = nbrDecomp[bFacei];
1182
1183 if (ownProc != nbrProc)
1184 {
1186 << "patch:" << pp.name()
1187 << " face:" << facei
1188 << " at:" << mesh.faceCentres()[facei]
1189 << " ownProc:" << ownProc
1190 << " nbrProc:" << nbrProc
1191 << exit(FatalError);
1192 }
1193 }
1194 }
1195 }
1196 }
1197 }
1198
1199 return finalDecomp;
1200}
1201
1202
1204(
1205 const polyMesh& mesh,
1206 boolList& blockedFace,
1207 PtrList<labelList>& specifiedProcessorFaces,
1208 labelList& specifiedProcessor,
1209 List<labelPair>& explicitConnections
1210) const
1211{
1212 blockedFace.setSize(mesh.nFaces());
1213 blockedFace = true;
1214
1215 specifiedProcessorFaces.clear();
1216 explicitConnections.clear();
1217
1218 for (const decompositionConstraint& decompConstraint : constraints_)
1219 {
1220 decompConstraint.add
1221 (
1222 mesh,
1223 blockedFace,
1224 specifiedProcessorFaces,
1225 specifiedProcessor,
1226 explicitConnections
1227 );
1228 }
1229}
1230
1231
1233(
1234 const polyMesh& mesh,
1235 const boolList& blockedFace,
1236 const PtrList<labelList>& specifiedProcessorFaces,
1237 const labelList& specifiedProcessor,
1238 const List<labelPair>& explicitConnections,
1239 labelList& decomposition
1240) const
1241{
1242 for (const decompositionConstraint& decompConstraint : constraints_)
1243 {
1244 decompConstraint.apply
1245 (
1246 mesh,
1247 blockedFace,
1248 specifiedProcessorFaces,
1249 specifiedProcessor,
1250 explicitConnections,
1251 decomposition
1252 );
1253 }
1254}
1255
1256
1258(
1259 const polyMesh& mesh,
1260 const scalarField& cellWeights
1261) const
1262{
1263 // Collect all constraints
1264
1265 boolList blockedFace;
1266 PtrList<labelList> specifiedProcessorFaces;
1267 labelList specifiedProcessor;
1268 List<labelPair> explicitConnections;
1269 setConstraints
1270 (
1271 mesh,
1272 blockedFace,
1273 specifiedProcessorFaces,
1274 specifiedProcessor,
1275 explicitConnections
1276 );
1277
1278
1279 // Construct decomposition method and either do decomposition on
1280 // cell centres or on agglomeration
1281
1282 labelList finalDecomp = decompose
1283 (
1284 mesh,
1285 cellWeights, // optional weights
1286 blockedFace, // any cells to be combined
1287 specifiedProcessorFaces,// any whole cluster of cells to be kept
1288 specifiedProcessor,
1289 explicitConnections // baffles
1290 );
1291
1292
1293 // Give any constraint the option of modifying the decomposition
1294
1295 applyConstraints
1296 (
1297 mesh,
1298 blockedFace,
1299 specifiedProcessorFaces,
1300 specifiedProcessor,
1301 explicitConnections,
1302 finalDecomp
1303 );
1304
1305 return finalDecomp;
1306}
1307
1308
1309// * * * * * * * * * * * * * * * Stub Functions * * * * * * * * * * * * * * //
1310
1312(
1313 const pointField& points,
1314 const scalarField& pointWeights
1315) const
1316{
1318 return labelList();
1319}
1320
1321
1323(
1324 const pointField& points
1325) const
1326{
1328 return labelList();
1329}
1330
1331
1332// ************************************************************************* //
A packed storage unstructured matrix of objects of type <T> using an offset table for access.
List< SubListType > unpack() const
Return non-compact list of lists.
const List< T > & m() const noexcept
Const access to the packed matrix of values.
label size() const noexcept
The primary size (the number of rows/sublists)
const labelList & offsets() const noexcept
Return the offset table (= size()+1)
void setSize(const label mRows)
Redimension - same as resize()
Wave propagation of information through grid. Every iteration information goes through one layer of c...
Definition: FaceCellWave.H:81
const TrackingData & data() const
Additional data to be passed into container.
Definition: FaceCellWave.H:348
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:191
void clear()
Clear all entries from table.
Definition: HashTable.C:678
void setSize(const label n)
Alias for resize()
Definition: List.H:218
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:116
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
void clear()
Clear the PtrList. Delete allocated entries and set size to zero.
Definition: PtrListI.H:97
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
bool empty() const noexcept
True if the UList is empty (ie, size() is zero)
Definition: UListI.H:427
T * data() noexcept
Return pointer to the underlying array serving as data storage.
Definition: UListI.H:237
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
static label worldComm
Default communicator (all processors)
Definition: UPstream.H:293
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
Abstract class for handling decomposition constraints.
Constraint to keep all cells connected to face or point of faceSet on a single processor.
Abstract base class for domain decomposition.
void applyConstraints(const polyMesh &mesh, const boolList &blockedFace, const PtrList< labelList > &specifiedProcessorFaces, const labelList &specifiedProcessor, const List< labelPair > &explicitConnections, labelList &finalDecomp) const
Helper: apply constraints to a decomposition.
void setConstraints(const polyMesh &mesh, boolList &blockedFace, PtrList< labelList > &specifiedProcessorFaces, labelList &specifiedProcessor, List< labelPair > &explicitConnections) const
Helper: extract constraints:
static void calcCellCells(const polyMesh &mesh, const labelList &agglom, const label nLocalCoarse, const bool global, CompactListList< label > &cellCells)
Helper: determine (local or global) cellCells from mesh.
static const dictionary & findCoeffsDict(const dictionary &dict, const word &coeffsName, int select=selectionType::DEFAULT)
static const dictionary & optionalRegionDict(const dictionary &decompDict, const word &regionName)
label nDomains() const noexcept
Number of domains.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
dictionary * findDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX)
Find and return a sub-dictionary pointer if present.
Definition: dictionaryI.H:127
const fileName & name() const noexcept
The dictionary name.
Definition: dictionaryI.H:48
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
static const dictionary null
An empty dictionary, which is also the parent for all dictionaries.
Definition: dictionary.H:394
static bool warnAboutAge(const int version) noexcept
Test if an age warning should be emitted.
Definition: error.C:55
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:68
label toGlobal(const label i) const
From local to global index.
Definition: globalIndexI.H:260
label toLocal(const label i) const
From global to local on current processor.
Definition: globalIndexI.H:325
bool isLocal(const label i) const
Is on local processor.
Definition: globalIndexI.H:244
label nTotalCells() const noexcept
Return total number of cells in decomposed mesh.
For use with FaceCellWave. Transports minimum passive data.
Definition: minData.H:61
static const complex max
complex (VGREAT,VGREAT)
Definition: complex.H:293
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1108
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1121
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1310
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1127
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
bool isInternalFace(const label faceIndex) const noexcept
Return true if given face label is internal to the mesh.
label nBoundaryFaces() const noexcept
Number of boundary faces (== nFaces - nInternalFaces)
const vectorField & faceCentres() const
label nInternalFaces() const noexcept
Number of internal faces.
const vectorField & cellCentres() const
label nCells() const noexcept
Number of mesh cells.
label nFaces() const noexcept
Number of mesh faces.
const labelListList & pointFaces() const
const vectorField & faceAreas() const
This class separates the mesh into distinct unconnected regions, each of which is then given a label ...
Definition: regionSplit.H:144
label nLocalRegions() const
Return local number of regions.
Definition: regionSplit.H:288
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &faceValues)
Swap coupled boundary face values. Uses eqOp.
Definition: syncTools.H:445
static void swapBoundaryCellList(const polyMesh &mesh, const UList< T > &cellData, List< T > &neighbourCellData)
Swap to obtain neighbour cell values for all boundary faces.
bool decompose() const noexcept
Query the decompose flag (normally off)
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
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
Foam::word regionName(Foam::polyMesh::defaultRegion)
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:478
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:517
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const pointField & points
#define WarningInFunction
Report a warning using Foam::Warning.
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of 'true' entries.
Definition: BitOps.H:78
bool all(const UList< bool > &bools)
True if all entries are 'true' or if the set is empty.
Definition: BitOps.H:85
Namespace for OpenFOAM.
static const dictionary * cfindCoeffsDict(const dictionary &dict, const word &coeffsName, const bool allowDefault)
List< label > labelList
A List of labels.
Definition: List.H:66
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
errorManip< error > abort(error &err)
Definition: errorManip.H:144
IOerror FatalIOError
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
messageStream Warning
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=cellModel::ref(cellModel::HEX);labelList hexLabels(8);label activeCells=0;labelList pointMap(nPoints);forAll(pointMap, i){ pointMap[i]=i;}for(label i=0;i< nPoints;i++){ if(f[i] > 0.0) { hexLabels[0]=i;hexLabels[1]=i1tab[i];hexLabels[2]=i3tab[i1tab[i]];hexLabels[3]=i3tab[i];hexLabels[4]=i8tab[i];hexLabels[5]=i1tab[i8tab[i]];hexLabels[6]=i3tab[i1tab[i8tab[i]]];hexLabels[7]=i3tab[i8tab[i]];cellShapes[activeCells].reset(hex, hexLabels);edgeList edges=cellShapes[activeCells].edges();forAll(edges, ei) { if(edges[ei].mag(points)< SMALL) { label start=pointMap[edges[ei].start()];while(start !=pointMap[start]) { start=pointMap[start];} label end=pointMap[edges[ei].end()];while(end !=pointMap[end]) { end=pointMap[end];} label minLabel=min(start, end);pointMap[start]=pointMap[end]=minLabel;} } cellZoning[activeCells]=idreg[i];activeCells++;}}cellShapes.setSize(activeCells);cellZoning.setSize(activeCells);forAll(cellShapes, celli){ cellShape &cs=cellShapes[celli];forAll(cs, i) { cs[i]=pointMap[cs[i]];} cs.collapse();}label bcIDs[11]={-1, 0, 2, 4, -1, 5, -1, 6, 7, 8, 9};const label nBCs=12;const word *kivaPatchTypes[nBCs]={ &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &symmetryPolyPatch::typeName, &wedgePolyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &symmetryPolyPatch::typeName, &oldCyclicPolyPatch::typeName};enum patchTypeNames{ PISTON, VALVE, LINER, CYLINDERHEAD, AXIS, WEDGE, INFLOW, OUTFLOW, PRESIN, PRESOUT, SYMMETRYPLANE, CYCLIC};const char *kivaPatchNames[nBCs]={ "piston", "valve", "liner", "cylinderHead", "axis", "wedge", "inflow", "outflow", "presin", "presout", "symmetryPlane", "cyclic"};List< SLList< face > > pFaces[nBCs]
Definition: readKivaGrid.H:235
labelList f(nPoints)
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333