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