faMeshDecomposition.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) 2016-2017 Wikki Ltd
9 Copyright (C) 2018-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
29#include "faMeshDecomposition.H"
30#include "Time.H"
31#include "dictionary.H"
32#include "labelIOList.H"
33#include "processorFaPatch.H"
34#include "faMesh.H"
35#include "OSspecific.H"
36#include "Map.H"
37#include "SLList.H"
38#include "globalMeshData.H"
39
40// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
41
42void Foam::faMeshDecomposition::distributeFaces()
43{
44 const word& polyMeshRegionName = mesh().name();
45
46 Info<< "\nCalculating distribution of finiteArea faces" << endl;
47
48 cpuTime decompositionTime;
49
50 for (label proci = 0; proci < nProcs(); ++proci)
51 {
52 Time processorDb
53 (
55 time().rootPath(),
56 time().caseName()/("processor" + Foam::name(proci))
57 );
58
59 polyMesh procFvMesh
60 (
61 IOobject
62 (
63 polyMeshRegionName,
64 processorDb.timeName(),
65 processorDb
66 )
67 );
68
69 IOobject ioAddr
70 (
71 "procAddressing",
72 "constant",
74 procFvMesh,
77 false // not registered
78 );
79
80
81 // faceProcAddressing (polyMesh)
82 ioAddr.rename("faceProcAddressing");
83 labelIOList fvFaceProcAddressing(ioAddr);
84
85 labelHashSet faceProcAddressingHash(2*fvFaceProcAddressing.size());
86
87 // If faMesh's fvPatch is a part of the global face zones, faces of that
88 // patch will be present on all processors. Because of that, looping
89 // through faceProcAddressing will decompose global faMesh faces to the
90 // very last processor regardless of where fvPatch is really decomposed.
91 // Since global faces which do not belong to specific processor are
92 // located at the end of the faceProcAddressing, cutting it at
93 // i = owner.size() will correctly decompose faMesh faces.
94 // Vanja Skuric, 2016-04-21
95 if (hasGlobalFaceZones_)
96 {
97 // owner (polyMesh)
98 ioAddr.rename("owner");
99 const label ownerSize = labelIOList(ioAddr).size();
100
101 for (int i = 0; i < ownerSize; ++i)
102 {
103 faceProcAddressingHash.insert(fvFaceProcAddressing[i]);
104 }
105 }
106 else
107 {
108 faceProcAddressingHash.insert
109 (
110 static_cast<labelList&>(fvFaceProcAddressing)
111 );
112 }
113
114 forAll(faceLabels(), facei)
115 {
116 // With +1 for lookup in faceMap with flip encoding
117 const label index = (faceLabels()[facei] + 1);
118
119 if (faceProcAddressingHash.found(index))
120 {
121 faceToProc_[facei] = proci;
122 }
123 }
124 }
125
126 Info<< "\nFinished decomposition in "
127 << decompositionTime.elapsedCpuTime()
128 << " s" << endl;
129}
130
131
132// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
133
135(
136 const polyMesh& mesh,
137 const label nProcessors,
138 const dictionary& params
139)
140:
141 faMesh(mesh),
142 nProcs_(nProcessors),
143 distributed_(false),
144 hasGlobalFaceZones_(false),
145 faceToProc_(nFaces()),
146 procFaceLabels_(nProcs_),
147 procMeshEdgesMap_(nProcs_),
148 procNInternalEdges_(nProcs_, Zero),
149 procPatchEdgeLabels_(nProcs_),
150 procPatchPointAddressing_(nProcs_),
151 procPatchEdgeAddressing_(nProcs_),
152 procEdgeAddressing_(nProcs_),
153 procFaceAddressing_(nProcs_),
154 procBoundaryAddressing_(nProcs_),
155 procPatchSize_(nProcs_),
156 procPatchStartIndex_(nProcs_),
157 procNeighbourProcessors_(nProcs_),
158 procProcessorPatchSize_(nProcs_),
159 procProcessorPatchStartIndex_(nProcs_),
160 globallySharedPoints_(),
161 cyclicParallel_(false)
162{
163 updateParameters(params);
164}
165
166
167// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
168
170(
171 const dictionary& params
172)
173{
174 params.readIfPresent("distributed", distributed_);
175 if (params.found("globalFaceZones"))
176 {
177 hasGlobalFaceZones_ = true;
178 }
179}
180
181
183{
184 // Decide which cell goes to which processor
185 distributeFaces();
186
187 const word& polyMeshRegionName = mesh().name();
188
189 Info<< "\nDistributing faces to processors" << endl;
190
191 labelList nLocalFaces(nProcs_, Zero);
192
193 // Pass 1: determine local sizes, sanity check
194
195 forAll(faceToProc_, facei)
196 {
197 const label proci = faceToProc_[facei];
198
199 if (proci < 0 || proci >= nProcs_)
200 {
202 << "Invalid processor label " << proci
203 << " for face " << facei << nl
204 << abort(FatalError);
205 }
206 else
207 {
208 ++nLocalFaces[proci];
209 }
210 }
211
212 // Adjust lengths
213 forAll(nLocalFaces, proci)
214 {
215 procFaceAddressing_[proci].resize(nLocalFaces[proci]);
216 nLocalFaces[proci] = 0; // restart list
217 }
218
219 // Pass 2: fill in local lists
220 forAll(faceToProc_, facei)
221 {
222 const label proci = faceToProc_[facei];
223 const label localFacei = nLocalFaces[proci];
224 ++nLocalFaces[proci];
225
226 procFaceAddressing_[proci][localFacei] = facei;
227 }
228
229
230 // Find processor mesh faceLabels and ...
231
232 for (label procI = 0; procI < nProcs(); procI++)
233 {
234 Time processorDb
235 (
237 time().rootPath(),
238 time().caseName()/("processor" + Foam::name(procI))
239 );
240
241 polyMesh procFvMesh
242 (
244 (
245 polyMeshRegionName,
246 processorDb.timeName(),
247 processorDb
248 )
249 );
250
251 IOobject ioAddr
252 (
253 "procAddressing",
254 "constant",
256 procFvMesh,
259 false // not registered
260 );
261
262
263 // pointProcAddressing (polyMesh)
264 ioAddr.rename("pointProcAddressing");
265 labelIOList fvPointProcAddressing(ioAddr);
266
267 Map<label> fvFaceProcAddressingHash;
268
269 {
270 // faceProcAddressing (polyMesh)
271 ioAddr.rename("faceProcAddressing");
272 labelIOList fvFaceProcAddressing(ioAddr);
273
274 fvFaceProcAddressingHash.resize(2*fvFaceProcAddressing.size());
275
276 forAll(fvFaceProcAddressing, facei)
277 {
278 fvFaceProcAddressingHash.insert
279 (
280 fvFaceProcAddressing[facei],
281 facei
282 );
283 }
284 };
285
286 const labelList& curProcFaceAddressing = procFaceAddressing_[procI];
287
288 labelList& curFaceLabels = procFaceLabels_[procI];
289
290 curFaceLabels = labelList(curProcFaceAddressing.size(), -1);
291
292 forAll(curProcFaceAddressing, faceI)
293 {
294 curFaceLabels[faceI] =
295 fvFaceProcAddressingHash.find
296 (
297 faceLabels()[curProcFaceAddressing[faceI]] + 1
298 ).val();
299 }
300
301 // Create processor finite area mesh
302 faMesh procMesh
303 (
304 procFvMesh,
305 labelList(procFaceLabels_[procI])
306 );
307
308 const uindirectPrimitivePatch& patch = this->patch();
309 const Map<label>& map = patch.meshPointMap();
310
311 EdgeMap<label> edgesHash;
312
313 const label nIntEdges = patch.nInternalEdges();
314
315 for (label edgei = 0; edgei < nIntEdges; ++edgei)
316 {
317 edgesHash.insert(patch.edges()[edgei], edgesHash.size());
318 }
319
320 forAll(boundary(), patchi)
321 {
322 // Include emptyFaPatch
323 const label size = boundary()[patchi].labelList::size();
324
325 for (label edgei=0; edgei < size; ++edgei)
326 {
327 edgesHash.insert
328 (
329 patch.edges()[boundary()[patchi][edgei]],
330 edgesHash.size()
331 );
332 }
333 }
334
335
336 const uindirectPrimitivePatch& procPatch = procMesh.patch();
337 const labelList& procMeshPoints = procPatch.meshPoints();
338 const edgeList& procEdges = procPatch.edges();
339
340 labelList& curPatchPointAddressing = procPatchPointAddressing_[procI];
341 curPatchPointAddressing.resize(procMeshPoints.size(), -1);
342
343 forAll(procMeshPoints, pointi)
344 {
345 curPatchPointAddressing[pointi] =
346 map[fvPointProcAddressing[procMeshPoints[pointi]]];
347 }
348
349 labelList& curPatchEdgeAddressing = procPatchEdgeAddressing_[procI];
350 curPatchEdgeAddressing.resize(procEdges.size(), -1);
351
352 Map<label>& curMap = procMeshEdgesMap_[procI];
353 curMap.clear();
354 curMap.resize(2*procEdges.size());
355
356 forAll(procEdges, edgeI)
357 {
358 edge curGlobalEdge(curPatchPointAddressing, procEdges[edgeI]);
359 curPatchEdgeAddressing[edgeI] = edgesHash.find(curGlobalEdge).val();
360 }
361
362 forAll(curPatchEdgeAddressing, edgeI)
363 {
364 curMap.insert(curPatchEdgeAddressing[edgeI], edgeI);
365 }
366
367 procNInternalEdges_[procI] = procPatch.nInternalEdges();
368 }
369
370
371 Info << "\nDistributing edges to processors" << endl;
372
373 // Loop through all internal edges and decide which processor they
374 // belong to. First visit all internal edges.
375
376 // set references to the original mesh
378 const edgeList& edges = this->edges();
379 const labelList& owner = edgeOwner();
380 const labelList& neighbour = edgeNeighbour();
381
382 // Memory management
383 {
384 List<SLList<label>> procEdgeList(nProcs());
385
386 forAll(procEdgeList, procI)
387 {
388 for(label i=0; i<procNInternalEdges_[procI]; i++)
389 {
390 procEdgeList[procI].append(procPatchEdgeAddressing_[procI][i]);
391 }
392 }
393
394
395 // Detect inter-processor boundaries
396
397 // Neighbour processor for each subdomain
398 List<SLList<label>> interProcBoundaries(nProcs());
399
400 // Edge labels belonging to each inter-processor boundary
401 List<SLList<SLList<label>>> interProcBEdges(nProcs());
402
403 List<SLList<label>> procPatchIndex(nProcs());
404
405 forAll(neighbour, edgeI)
406 {
407 if (faceToProc_[owner[edgeI]] != faceToProc_[neighbour[edgeI]])
408 {
409 // inter - processor patch edge found. Go through the list of
410 // inside boundaries for the owner processor and try to find
411 // this inter-processor patch.
412
413 bool interProcBouFound = false;
414
415 const label ownProc = faceToProc_[owner[edgeI]];
416 const label neiProc = faceToProc_[neighbour[edgeI]];
417
418 auto curInterProcBdrsOwnIter =
419 interProcBoundaries[ownProc].cbegin();
420
421 auto curInterProcBEdgesOwnIter =
422 interProcBEdges[ownProc].begin();
423
424 // WARNING: Synchronous SLList iterators
425
426 for
427 (
428 ;
429 curInterProcBdrsOwnIter.good()
430 && curInterProcBEdgesOwnIter.good();
431 ++curInterProcBdrsOwnIter,
432 ++curInterProcBEdgesOwnIter
433 )
434 {
435 if (curInterProcBdrsOwnIter() == neiProc)
436 {
437 // the inter - processor boundary exists. Add the face
438 interProcBouFound = true;
439
440 bool neighbourFound = false;
441
442 curInterProcBEdgesOwnIter().append(edgeI);
443
444 auto curInterProcBdrsNeiIter =
445 interProcBoundaries[neiProc].cbegin();
446
447 auto curInterProcBEdgesNeiIter =
448 interProcBEdges[neiProc].begin();
449
450 // WARNING: Synchronous SLList iterators
451
452 for
453 (
454 ;
455 curInterProcBdrsNeiIter.good()
456 && curInterProcBEdgesNeiIter.good();
457 ++curInterProcBdrsNeiIter,
458 ++curInterProcBEdgesNeiIter
459 )
460 {
461 if (curInterProcBdrsNeiIter() == ownProc)
462 {
463 // boundary found. Add the face
464 neighbourFound = true;
465
466 curInterProcBEdgesNeiIter().append(edgeI);
467 }
468
469 if (neighbourFound) break;
470 }
471
472 if (interProcBouFound && !neighbourFound)
473 {
475 ("faDomainDecomposition::decomposeMesh()")
476 << "Inconsistency in inter - "
477 << "processor boundary lists for processors "
478 << ownProc << " and " << neiProc
479 << abort(FatalError);
480 }
481 }
482
483 if (interProcBouFound) break;
484 }
485
486 if (!interProcBouFound)
487 {
488 // inter - processor boundaries do not exist and need to
489 // be created
490
491 // set the new addressing information
492
493 // owner
494 interProcBoundaries[ownProc].append(neiProc);
495 interProcBEdges[ownProc].append(SLList<label>(edgeI));
496
497 // neighbour
498 interProcBoundaries[neiProc].append(ownProc);
499 interProcBEdges[neiProc]
500 .append
501 (
502 SLList<label>(edgeI)
503 );
504 }
505 }
506 }
507
508
509 // Loop through patches. For cyclic boundaries detect inter-processor
510 // edges; for all other, add edges to the edge list and remember start
511 // and size of all patches.
512
513 // for all processors, set the size of start index and patch size
514 // lists to the number of patches in the mesh
515 forAll(procPatchSize_, procI)
516 {
517 procPatchSize_[procI].setSize(patches.size());
518 procPatchStartIndex_[procI].setSize(patches.size());
519 }
520
521 forAll(patches, patchI)
522 {
523 // Reset size and start index for all processors
524 forAll(procPatchSize_, procI)
525 {
526 procPatchSize_[procI][patchI] = 0;
527 procPatchStartIndex_[procI][patchI] =
528 procEdgeList[procI].size();
529 }
530
531 const label patchStart = patches[patchI].start();
532
533// if (!isA<cyclicFaPatch>(patches[patchI]))
534 if (true)
535 {
536 // Normal patch. Add edges to processor where the face
537 // next to the edge lives
538
539 const labelListList& eF = patch().edgeFaces();
540
541 const label size = patches[patchI].labelList::size();
542
543 labelList patchEdgeFaces(size, -1);
544
545 for(int eI=0; eI<size; eI++)
546 {
547 patchEdgeFaces[eI] = eF[patches[patchI][eI]][0];
548 }
549
550 forAll(patchEdgeFaces, edgeI)
551 {
552 const label curProc = faceToProc_[patchEdgeFaces[edgeI]];
553
554 // add the face
555 procEdgeList[curProc].append(patchStart + edgeI);
556
557 // increment the number of edges for this patch
558 procPatchSize_[curProc][patchI]++;
559 }
560 }
561 else
562 {
563 // Cyclic patch special treatment
564
565 const faPatch& cPatch = patches[patchI];
566
567 const label cycOffset = cPatch.size()/2;
568
569 // Set reference to faceCells for both patches
570 const labelList::subList firstEdgeFaces
571 (
572 cPatch.edgeFaces(),
573 cycOffset
574 );
575
576 const labelList::subList secondEdgeFaces
577 (
578 cPatch.edgeFaces(),
579 cycOffset,
580 cycOffset
581 );
582
583 forAll(firstEdgeFaces, edgeI)
584 {
585 if
586 (
587 faceToProc_[firstEdgeFaces[edgeI]]
588 != faceToProc_[secondEdgeFaces[edgeI]]
589 )
590 {
591 // This edge becomes an inter-processor boundary edge
592 // inter - processor patch edge found. Go through
593 // the list of inside boundaries for the owner
594 // processor and try to find this inter-processor
595 // patch.
596
597 cyclicParallel_ = true;
598
599 bool interProcBouFound = false;
600
601 const label ownProc =
602 faceToProc_[firstEdgeFaces[edgeI]];
603 const label neiProc =
604 faceToProc_[secondEdgeFaces[edgeI]];
605
606 auto curInterProcBdrsOwnIter =
607 interProcBoundaries[ownProc].cbegin();
608
609 auto curInterProcBEdgesOwnIter =
610 interProcBEdges[ownProc].begin();
611
612 // WARNING: Synchronous SLList iterators
613
614 for
615 (
616 ;
617 curInterProcBdrsOwnIter.good()
618 && curInterProcBEdgesOwnIter.good();
619 ++curInterProcBdrsOwnIter,
620 ++curInterProcBEdgesOwnIter
621 )
622 {
623 if (curInterProcBdrsOwnIter() == neiProc)
624 {
625 // the inter - processor boundary exists.
626 // Add the face
627 interProcBouFound = true;
628
629 bool neighbourFound = false;
630
631 curInterProcBEdgesOwnIter()
632 .append(patchStart + edgeI);
633
634 auto curInterProcBdrsNeiIter
635 = interProcBoundaries[neiProc].cbegin();
636
637 auto curInterProcBEdgesNeiIter =
638 interProcBEdges[neiProc].begin();
639
640 // WARNING: Synchronous SLList iterators
641
642 for
643 (
644 ;
645 curInterProcBdrsNeiIter.good()
646 && curInterProcBEdgesNeiIter.good();
647 ++curInterProcBdrsNeiIter,
648 ++curInterProcBEdgesNeiIter
649 )
650 {
651 if (curInterProcBdrsNeiIter() == ownProc)
652 {
653 // boundary found. Add the face
654 neighbourFound = true;
655
656 curInterProcBEdgesNeiIter()
657 .append
658 (
659 patchStart
660 + cycOffset
661 + edgeI
662 );
663 }
664
665 if (neighbourFound) break;
666 }
667
668 if (interProcBouFound && !neighbourFound)
669 {
671 (
672 "faDomainDecomposition::decomposeMesh()"
673 ) << "Inconsistency in inter-processor "
674 << "boundary lists for processors "
675 << ownProc << " and " << neiProc
676 << " in cyclic boundary matching"
677 << abort(FatalError);
678 }
679 }
680
681 if (interProcBouFound) break;
682 }
683
684 if (!interProcBouFound)
685 {
686 // inter - processor boundaries do not exist
687 // and need to be created
688
689 // set the new addressing information
690
691 // owner
692 interProcBoundaries[ownProc].append(neiProc);
693 interProcBEdges[ownProc]
694 .append(SLList<label>(patchStart + edgeI));
695
696 // neighbour
697 interProcBoundaries[neiProc].append(ownProc);
698 interProcBEdges[neiProc]
699 .append
700 (
702 (
703 patchStart
704 + cycOffset
705 + edgeI
706 )
707 );
708 }
709 }
710 else
711 {
712 // This cyclic edge remains on the processor
713 label ownProc = faceToProc_[firstEdgeFaces[edgeI]];
714
715 // add the edge
716 procEdgeList[ownProc].append(patchStart + edgeI);
717
718 // increment the number of edges for this patch
719 procPatchSize_[ownProc][patchI]++;
720
721 // Note: I cannot add the other side of the cyclic
722 // boundary here because this would violate the order.
723 // They will be added in a separate loop below
724 }
725 }
726
727 // Ordering in cyclic boundaries is important.
728 // Add the other half of cyclic edges for cyclic boundaries
729 // that remain on the processor
730 forAll(secondEdgeFaces, edgeI)
731 {
732 if
733 (
734 faceToProc_[firstEdgeFaces[edgeI]]
735 == faceToProc_[secondEdgeFaces[edgeI]]
736 )
737 {
738 // This cyclic edge remains on the processor
739 label ownProc = faceToProc_[firstEdgeFaces[edgeI]];
740
741 // add the second edge
742 procEdgeList[ownProc].append
743 (patchStart + cycOffset + edgeI);
744
745 // increment the number of edges for this patch
746 procPatchSize_[ownProc][patchI]++;
747 }
748 }
749 }
750 }
751
752 // Convert linked lists into normal lists
753 // Add inter-processor boundaries and remember start indices
754 forAll(procEdgeList, procI)
755 {
756 // Get internal and regular boundary processor faces
757 SLList<label>& curProcEdges = procEdgeList[procI];
758
759 // Get reference to processor edge addressing
760 labelList& curProcEdgeAddressing = procEdgeAddressing_[procI];
761
762 labelList& curProcNeighbourProcessors =
763 procNeighbourProcessors_[procI];
764
765 labelList& curProcProcessorPatchSize =
766 procProcessorPatchSize_[procI];
767
768 labelList& curProcProcessorPatchStartIndex =
769 procProcessorPatchStartIndex_[procI];
770
771 // calculate the size
772 label nEdgesOnProcessor = curProcEdges.size();
773
774 for (const auto& bedges : interProcBEdges[procI])
775 {
776 nEdgesOnProcessor += bedges.size();
777 }
778
779 curProcEdgeAddressing.setSize(nEdgesOnProcessor);
780
781 // Fill in the list. Calculate turning index.
782 // Turning index will be -1 only for some edges on processor
783 // boundaries, i.e. the ones where the current processor ID
784 // is in the face which is a edge neighbour.
785 // Turning index is stored as the sign of the edge addressing list
786
787 label nEdges = 0;
788
789 // Add internal and boundary edges
790 // Remember to increment the index by one such that the
791 // turning index works properly.
792 for (const label procEdgei : curProcEdges)
793 {
794 curProcEdgeAddressing[nEdges] = procEdgei;
795// curProcEdgeAddressing[nEdges] = procEdgei + 1;
796 ++nEdges;
797 }
798
799 // Add inter-processor boundary edges. At the beginning of each
800 // patch, grab the patch start index and size
801
802 curProcNeighbourProcessors.setSize
803 (
804 interProcBoundaries[procI].size()
805 );
806
807 curProcProcessorPatchSize.setSize
808 (
809 interProcBoundaries[procI].size()
810 );
811
812 curProcProcessorPatchStartIndex.setSize
813 (
814 interProcBoundaries[procI].size()
815 );
816
817 label nProcPatches = 0;
818
819 auto curInterProcBdrsIter =
820 interProcBoundaries[procI].cbegin();
821
822 auto curInterProcBEdgesIter =
823 interProcBEdges[procI].cbegin();
824
825 for
826 (
827 ;
828 curInterProcBdrsIter.good()
829 && curInterProcBEdgesIter.good();
830 ++curInterProcBdrsIter,
831 ++curInterProcBEdgesIter
832 )
833 {
834 curProcNeighbourProcessors[nProcPatches] =
835 curInterProcBdrsIter();
836
837 // Get start index for processor patch
838 curProcProcessorPatchStartIndex[nProcPatches] = nEdges;
839
840 label& curSize =
841 curProcProcessorPatchSize[nProcPatches];
842
843 curSize = 0;
844
845 // add faces for this processor boundary
846
847 for (const label edgei : *curInterProcBEdgesIter)
848 {
849 // add the edges
850
851 // Remember to increment the index by one such that the
852 // turning index works properly.
853 if (faceToProc_[owner[edgei]] == procI)
854 {
855 curProcEdgeAddressing[nEdges] = edgei;
856// curProcEdgeAddressing[nEdges] = edgei + 1;
857 }
858 else
859 {
860 // turning edge
861 curProcEdgeAddressing[nEdges] = edgei;
862// curProcEdgeAddressing[nEdges] = -(edgei + 1);
863 }
864
865 // increment the size
866 ++curSize;
867
868 ++nEdges;
869 }
870
871 ++nProcPatches;
872 }
873 }
874 }
875
876 Info << "\nCalculating processor boundary addressing" << endl;
877 // For every patch of processor boundary, find the index of the original
878 // patch. Mis-alignment is caused by the fact that patches with zero size
879 // are omitted. For processor patches, set index to -1.
880 // At the same time, filter the procPatchSize_ and procPatchStartIndex_
881 // lists to exclude zero-size patches
882 forAll(procPatchSize_, procI)
883 {
884 // Make a local copy of old lists
885 const labelList oldPatchSizes = procPatchSize_[procI];
886
887 const labelList oldPatchStarts = procPatchStartIndex_[procI];
888
889 labelList& curPatchSizes = procPatchSize_[procI];
890
891 labelList& curPatchStarts = procPatchStartIndex_[procI];
892
893 const labelList& curProcessorPatchSizes =
894 procProcessorPatchSize_[procI];
895
896 labelList& curBoundaryAddressing = procBoundaryAddressing_[procI];
897
898 curBoundaryAddressing.setSize
899 (
900 oldPatchSizes.size()
901 + curProcessorPatchSizes.size()
902 );
903
904 label nPatches = 0;
905
906 forAll(oldPatchSizes, patchI)
907 {
908 //- Do not suppress zero sized patches since make parallel
909 // actions inside patches near impossible.
910 //if (oldPatchSizes[patchI] > 0)
911 {
912 curBoundaryAddressing[nPatches] = patchI;
913
914 curPatchSizes[nPatches] = oldPatchSizes[patchI];
915
916 curPatchStarts[nPatches] = oldPatchStarts[patchI];
917
918 nPatches++;
919 }
920 }
921
922 // reset to the size of live patches
923 curPatchSizes.setSize(nPatches);
924 curPatchStarts.setSize(nPatches);
925
926 forAll(curProcessorPatchSizes, procPatchI)
927 {
928 curBoundaryAddressing[nPatches] = -1;
929
930 nPatches++;
931 }
932
933 curBoundaryAddressing.setSize(nPatches);
934 }
935
936
937 // Gather data about globally shared points
938
939 labelList globallySharedPoints_(0);
940
941 // Memory management
942 {
943 labelList pointsUsage(nPoints(), Zero);
944
945 // Globally shared points are the ones used by more than 2 processors
946 // Size the list approximately and gather the points
947 labelHashSet gSharedPoints
948 (
949 min(100, nPoints()/1000)
950 );
951
952 // Loop through all the processors and mark up points used by
953 // processor boundaries. When a point is used twice, it is a
954 // globally shared point
955
956 for (label procI = 0; procI < nProcs(); procI++)
957 {
958 // Get list of edge labels
959 const labelList& curEdgeLabels = procEdgeAddressing_[procI];
960
961 // Get start of processor faces
962 const labelList& curProcessorPatchStarts =
963 procProcessorPatchStartIndex_[procI];
964
965 const labelList& curProcessorPatchSizes =
966 procProcessorPatchSize_[procI];
967
968 // Reset the lookup list
969 pointsUsage = 0;
970
971 forAll(curProcessorPatchStarts, patchI)
972 {
973 const label curStart = curProcessorPatchStarts[patchI];
974 const label curEnd = curStart + curProcessorPatchSizes[patchI];
975
976 for
977 (
978 label edgeI = curStart;
979 edgeI < curEnd;
980 edgeI++
981 )
982 {
983 // Mark the original edge as used
984 // Remember to decrement the index by one (turning index)
985 const label curE = curEdgeLabels[edgeI];
986
987 const edge& e = edges[curE];
988
989 forAll(e, pointI)
990 {
991 if (pointsUsage[e[pointI]] == 0)
992 {
993 // Point not previously used
994 pointsUsage[e[pointI]] = patchI + 1;
995 }
996 else if (pointsUsage[e[pointI]] != patchI + 1)
997 {
998 // Point used by some other patch = global point!
999 gSharedPoints.insert(e[pointI]);
1000 }
1001 }
1002 }
1003 }
1004 }
1005
1006 // Grab the result from the hash list
1007 globallySharedPoints_ = gSharedPoints.sortedToc();
1008 }
1009
1010
1011 // Edge label for faPatches
1012
1013 for (label procI = 0; procI < nProcs(); procI++)
1014 {
1015 fileName processorCasePath
1016 (
1017 time().caseName()/("processor" + Foam::name(procI))
1018 );
1019
1020 // create a database
1021 Time processorDb
1022 (
1024 time().rootPath(),
1025 processorCasePath
1026 );
1027
1028
1029 // Read volume mesh
1030 polyMesh procFvMesh
1031 (
1032 IOobject
1033 (
1034 polyMeshRegionName,
1035 processorDb.timeName(),
1036 processorDb
1037 )
1038 );
1039
1040 // create finite area mesh
1041 faMesh procMesh
1042 (
1043 procFvMesh,
1044 labelList(procFaceLabels_[procI])
1045 );
1046
1047
1048 const labelList& curEdgeAddressing = procEdgeAddressing_[procI];
1049
1050 const labelList& curPatchStartIndex = procPatchStartIndex_[procI];
1051 const labelList& curPatchSize = procPatchSize_[procI];
1052
1053 const labelList& curProcessorPatchStartIndex =
1054 procProcessorPatchStartIndex_[procI];
1055
1056 const labelList& curProcessorPatchSize =
1057 procProcessorPatchSize_[procI];
1058
1059 labelListList& curPatchEdgeLabels = procPatchEdgeLabels_[procI];
1060 curPatchEdgeLabels.resize
1061 (
1062 curPatchSize.size()
1063 + curProcessorPatchSize.size()
1064 );
1065
1066 forAll(curPatchSize, patchI)
1067 {
1068 labelList& curEdgeLabels = curPatchEdgeLabels[patchI];
1069 curEdgeLabels.setSize(curPatchSize[patchI], -1);
1070
1071 label edgeI = 0;
1072
1073 for
1074 (
1075 int i=curPatchStartIndex[patchI];
1076 i<(curPatchStartIndex[patchI]+curPatchSize[patchI]);
1077 i++
1078 )
1079 {
1080 curEdgeLabels[edgeI] =
1081 procMeshEdgesMap_[procI][curEdgeAddressing[i]];
1082 edgeI++;
1083 }
1084 }
1085
1086 forAll(curProcessorPatchSize, patchI)
1087 {
1088 labelList& curEdgeLabels =
1089 curPatchEdgeLabels[curPatchSize.size() + patchI];
1090 curEdgeLabels.setSize(curProcessorPatchSize[patchI], -1);
1091
1092 label edgeI = 0;
1093
1094 for
1095 (
1096 int i=curProcessorPatchStartIndex[patchI];
1097 i<(curProcessorPatchStartIndex[patchI]
1098 +curProcessorPatchSize[patchI]);
1099 i++
1100 )
1101 {
1102 curEdgeLabels[edgeI] =
1103 procMeshEdgesMap_[procI][curEdgeAddressing[i]];
1104 edgeI++;
1105 }
1106 }
1107 }
1108}
1109
1110
1112{
1113 const word& polyMeshRegionName = mesh().name();
1114
1115 Info<< "\nConstructing processor FA meshes" << endl;
1116
1117 // Make a lookup map for globally shared points
1118 Map<label> sharedPointLookup(2*globallySharedPoints_.size());
1119
1120 forAll(globallySharedPoints_, pointi)
1121 {
1122 sharedPointLookup.insert(globallySharedPoints_[pointi], pointi);
1123 }
1124
1125 label maxProcFaces = 0;
1126 label totProcFaces = 0;
1127 label maxProcEdges = 0;
1128 label totProcEdges = 0;
1129 label maxProcPatches = 0;
1130 label totProcPatches = 0;
1131
1132 // Write out the meshes
1133 for (label procI = 0; procI < nProcs(); procI++)
1134 {
1135 // Create processor mesh without a boundary
1136
1137 fileName processorCasePath
1138 (
1139 time().caseName()/("processor" + Foam::name(procI))
1140 );
1141
1142 // create a database
1143 Time processorDb
1144 (
1146 time().rootPath(),
1147 processorCasePath
1148 );
1149
1150 // Read volume mesh
1151 polyMesh procFvMesh
1152 (
1153 IOobject
1154 (
1155 polyMeshRegionName,
1156 processorDb.timeName(),
1157 processorDb
1158 )
1159 );
1160
1161 labelIOList fvBoundaryProcAddressing
1162 (
1163 IOobject
1164 (
1165 "boundaryProcAddressing",
1166 "constant",
1168 procFvMesh,
1171 )
1172 );
1173
1174
1175 // Create finite area mesh
1176 faMesh procMesh
1177 (
1178 procFvMesh,
1179 labelList(procFaceLabels_[procI])
1180 );
1181
1182 // Create processor boundary patches
1183 const labelList& curBoundaryAddressing =
1184 procBoundaryAddressing_[procI];
1185
1186 const labelList& curPatchSizes = procPatchSize_[procI];
1187
1188 const labelList& curNeighbourProcessors =
1189 procNeighbourProcessors_[procI];
1190
1191 const labelList& curProcessorPatchSizes =
1192 procProcessorPatchSize_[procI];
1193
1194 const labelListList& curPatchEdgeLabels =
1195 procPatchEdgeLabels_[procI];
1196
1197 const faPatchList& meshPatches = boundary();
1198
1199 faPatchList procPatches
1200 (
1201 curPatchSizes.size() + curProcessorPatchSizes.size()
1202 );
1203
1204 label nPatches = 0;
1205
1206 forAll(curPatchSizes, patchi)
1207 {
1208 const labelList& curEdgeLabels = curPatchEdgeLabels[nPatches];
1209
1210 const label neiPolyPatchId =
1211 fvBoundaryProcAddressing.find
1212 (
1213 meshPatches[curBoundaryAddressing[patchi]]
1214 .ngbPolyPatchIndex()
1215 );
1216
1217 procPatches.set
1218 (
1219 nPatches,
1220 meshPatches[curBoundaryAddressing[patchi]].clone
1221 (
1222 procMesh.boundary(),
1223 curEdgeLabels,
1224 nPatches,
1225 neiPolyPatchId
1226 )
1227 );
1228 ++nPatches;
1229 }
1230
1231 forAll(curProcessorPatchSizes, procPatchI)
1232 {
1233 const labelList& curEdgeLabels = curPatchEdgeLabels[nPatches];
1234
1235 procPatches.set
1236 (
1237 nPatches,
1239 (
1240 curEdgeLabels,
1241 nPatches,
1242 procMesh.boundary(),
1243 -1,
1244 procI,
1245 curNeighbourProcessors[procPatchI]
1246 )
1247 );
1248
1249 ++nPatches;
1250 }
1251
1252 // Add boundary patches
1253 procMesh.addFaPatches(procPatches);
1254
1255 // Set the precision of the points data to 10
1257
1258 procMesh.write();
1259
1260 // Statistics
1261 Info<< nl << "Processor " << procI;
1262
1263 if (procMesh.nFaces())
1264 {
1265 Info<< nl << " ";
1266 }
1267 else
1268 {
1269 Info<< ": ";
1270 }
1271
1272 Info<< "Number of faces = " << procMesh.nFaces() << nl;
1273
1274 if (procMesh.nFaces())
1275 {
1276 Info<< " Number of points = " << procMesh.nPoints() << nl;
1277 }
1278
1279 totProcFaces += procMesh.nFaces();
1280 maxProcFaces = max(maxProcFaces, procMesh.nFaces());
1281
1282 label nBoundaryEdges = 0;
1283 label nProcPatches = 0;
1284 label nProcEdges = 0;
1285
1286 for (const faPatch& fap : procMesh.boundary())
1287 {
1288 const auto* ppp = isA<processorFaPatch>(fap);
1289
1290 if (ppp)
1291 {
1292 const auto& procPatch = *ppp;
1293
1294 Info<< " Number of edges shared with processor "
1295 << procPatch.neighbProcNo() << " = "
1296 << procPatch.size() << nl;
1297
1298 nProcEdges += procPatch.size();
1299 ++nProcPatches;
1300 }
1301 else
1302 {
1303 nBoundaryEdges += fap.size();
1304 }
1305 }
1306
1307 if (procMesh.nFaces() && (nBoundaryEdges || nProcEdges))
1308 {
1309 Info<< " Number of processor patches = " << nProcPatches << nl
1310 << " Number of processor edges = " << nProcEdges << nl
1311 << " Number of boundary edges = " << nBoundaryEdges << nl;
1312 }
1313
1314 totProcEdges += nProcEdges;
1315 totProcPatches += nProcPatches;
1316 maxProcEdges = max(maxProcEdges, nProcEdges);
1317 maxProcPatches = max(maxProcPatches, nProcPatches);
1318
1319 // Write the addressing information
1320
1321 IOobject ioAddr
1322 (
1323 "procAddressing",
1324 "constant",
1326 procMesh.thisDb(),
1329 false // not registered
1330 );
1331
1332 // pointProcAddressing
1333 ioAddr.rename("pointProcAddressing");
1334 IOListRef<label>(ioAddr, procPatchPointAddressing_[procI]).write();
1335
1336 // edgeProcAddressing
1337 ioAddr.rename("edgeProcAddressing");
1338 IOListRef<label>(ioAddr, procEdgeAddressing_[procI]).write();
1339
1340 // faceProcAddressing
1341 ioAddr.rename("faceProcAddressing");
1342 IOListRef<label>(ioAddr, procFaceAddressing_[procI]).write();
1343
1344 // boundaryProcAddressing
1345 ioAddr.rename("boundaryProcAddressing");
1346 IOListRef<label>(ioAddr, procBoundaryAddressing_[procI]).write();
1347 }
1348
1349
1350 // Summary stats
1351 Info<< nl
1352 << "Number of processor edges = " << (totProcEdges/2) << nl
1353 << "Max number of faces = " << maxProcFaces;
1354
1355 if (maxProcFaces != totProcFaces)
1356 {
1357 scalar avgValue = scalar(totProcFaces)/nProcs_;
1358
1359 Info<< " (" << 100.0*(maxProcFaces-avgValue)/avgValue
1360 << "% above average " << avgValue << ')';
1361 }
1362 Info<< nl;
1363
1364 Info<< "Max number of processor patches = " << maxProcPatches;
1365 if (totProcPatches)
1366 {
1367 scalar avgValue = scalar(totProcPatches)/nProcs_;
1368
1369 Info<< " (" << 100.0*(maxProcPatches-avgValue)/avgValue
1370 << "% above average " << avgValue << ')';
1371 }
1372 Info<< nl;
1373
1374 Info<< "Max number of edges between processors = " << maxProcEdges;
1375 if (totProcEdges)
1376 {
1377 scalar avgValue = scalar(totProcEdges)/nProcs_;
1378
1379 Info<< " (" << 100.0*(maxProcEdges-avgValue)/avgValue
1380 << "% above average " << avgValue << ')';
1381 }
1382 Info<< nl << endl;
1383
1384 return true;
1385}
1386
1387
1388// ************************************************************************* //
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
Non-intrusive singly-linked list.
Map from edge (expressed as its endpoints) to value. For easier forward declaration it is currently i...
Definition: EdgeMap.H:54
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:191
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition: HashTable.C:137
void resize(const label sz)
Resize the hash table for efficiency.
Definition: HashTable.C:601
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition: HashTableI.H:180
label size() const noexcept
The number of elements in table.
Definition: HashTableI.H:52
iterator find(const Key &key)
Find and return an iterator set at the hashed entry.
Definition: HashTableI.H:114
void clear()
Clear all entries from table.
Definition: HashTable.C:678
A IOList wrapper for writing external data.
Definition: IOList.H:123
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const fileName & caseName() const
Return the Time::caseName()
Definition: IOobject.C:518
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
IOobject(const IOobject &)=default
Copy construct.
virtual void rename(const word &newName)
Rename the object.
Definition: IOobject.H:497
const fileName & rootPath() const
Return the Time::rootPath()
Definition: IOobject.C:512
static unsigned int defaultPrecision() noexcept
Return the default precision.
Definition: IOstream.H:342
Template class for non-intrusive linked lists.
Definition: LList.H:79
void setSize(const label n)
Alias for resize()
Definition: List.H:218
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:175
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:139
A HashTable to objects of type <T> with a label key.
Definition: Map.H:60
A list of faces which address into the list of points.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
label nInternalEdges() const
Number of internal edges.
const labelList & meshPoints() const
Return labelList of mesh points in patch.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
const T * set(const label i) const
Definition: PtrList.H:138
A List obtained as a section of another List.
Definition: SubList.H:70
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
static word controlDictName
The default control dictionary name (normally "controlDict")
Definition: Time.H:226
iterator begin() noexcept
Return an iterator to begin traversing the UList.
Definition: UListI.H:329
const_iterator cbegin() const noexcept
Return const_iterator to begin traversing the constant UList.
Definition: UListI.H:343
label find(const T &val, label pos=0) const
Find index of the first occurrence of the value.
Definition: UList.C:212
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:66
Finite area boundary mesh.
Automatic faMesh decomposition class.
bool writeDecomposition()
Write decomposition.
void updateParameters(const dictionary &params)
Update flags based on the decomposition model settings.
label nProcs() const noexcept
Number of processor in decomposition.
void decomposeMesh()
Decompose mesh.
Finite area mesh (used for 2-D non-Euclidian finite area method) defined using a patch of faces on a ...
Definition: faMesh.H:100
const Time & time() const
Return reference to time.
Definition: faMesh.C:673
const faBoundaryMesh & boundary() const noexcept
Return constant reference to boundary mesh.
Definition: faMeshI.H:38
virtual bool write(const bool valid=true) const
Write mesh.
Definition: faMesh.C:1005
const polyMesh & mesh() const
Return access to polyMesh.
Definition: faMeshI.H:31
const uindirectPrimitivePatch & patch() const
Return constant reference to primitive patch.
Definition: faMeshI.H:128
virtual const objectRegistry & thisDb() const
Return reference to the mesh database.
Definition: faMesh.C:697
label nPoints() const noexcept
Number of local mesh points.
Definition: faMeshI.H:56
label nFaces() const noexcept
Number of patch faces.
Definition: faMeshI.H:80
void addFaPatches(faPatchList &plist, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: faMeshPatches.C:38
const labelList & faceLabels() const noexcept
Return the underlying polyMesh face labels.
Definition: faMeshI.H:122
static word meshSubDir
The mesh sub-directory name (usually "faMesh")
Definition: faMesh.H:504
Finite area patch class. Used for 2-D non-Euclidian finite area method.
Definition: faPatch.H:78
virtual label size() const
Patch size is the number of edge labels.
Definition: faPatch.H:311
const labelUList & edgeFaces() const
Return edge-face addressing.
Definition: faPatch.C:429
A class for handling file names.
Definition: fileName.H:76
const word & name() const
Return reference to name.
Definition: fvMesh.H:310
label start() const noexcept
The start label of boundary faces in the polyMesh face list.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:324
virtual bool write(const bool valid=true) const
Write using setting from DB.
A class for handling words, derived from Foam::string.
Definition: word.H:68
const label nProcPatches
const polyBoundaryMesh & patches
faceListList boundary
dynamicFvMesh & mesh
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:448
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const labelList nFaces(UPstream::listGatherValues< label >(aMesh.nFaces()))
const labelList nIntEdges(UPstream::listGatherValues< label >(aMesh.nInternalEdges()))
const labelList nProcEdges(UPstream::listGatherValues< label >(nLocalProcEdges))
const label nPatches
const labelList nEdges(UPstream::listGatherValues< label >(aMesh.nEdges()))
label nPoints
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
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
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
errorManip< error > abort(error &err)
Definition: errorManip.H:144
IOList< label > labelIOList
Label container classes.
Definition: labelIOList.H:44
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
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
cpuTimeCxx cpuTime
Selection of preferred clock mechanism for the elapsed cpu time.
Definition: cpuTimeFwd.H:43
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333