domainDecompositionDecompose.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "domainDecomposition.H"
27 #include "decompositionMethod.H"
28 #include "IOobjectList.H"
29 #include "cyclicFvPatch.H"
30 #include "processorCyclicFvPatch.H"
33 
34 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
35 
36 void Foam::domainDecomposition::mark
37 (
38  const labelList& zoneElems,
39  const label zoneI,
40  labelList& elementToZone
41 )
42 {
43  forAll(zoneElems, i)
44  {
45  label pointi = zoneElems[i];
46 
47  if (elementToZone[pointi] == -1)
48  {
49  // First occurrence
50  elementToZone[pointi] = zoneI;
51  }
52  else if (elementToZone[pointi] >= 0)
53  {
54  // Multiple zones
55  elementToZone[pointi] = -2;
56  }
57  }
58 }
59 
60 
61 void Foam::domainDecomposition::addInterProcFace
62 (
63  const label facei,
64  const label ownerProc,
65  const label nbrProc,
66  const label subPatchID,
67  List<Map<label>>& nbrToInterPatch,
68  List<DynamicList<DynamicList<label>>>& interPatchFaces,
69  List<labelListList>& subPatchIDs,
70  List<labelListList>& subPatchStarts
71 ) const
72 {
73  // Create new interproc patches, if they do not already exist
74  label toNbrProcPatchi = -1, toOwnerProcPatchi = -1;
75  if (!nbrToInterPatch[ownerProc].found(nbrProc))
76  {
77  toNbrProcPatchi = nbrToInterPatch[ownerProc].size();
78  nbrToInterPatch[ownerProc].insert(nbrProc, toNbrProcPatchi);
79  interPatchFaces[ownerProc].append(DynamicList<label>());
80  subPatchIDs[ownerProc].append(labelList(1, subPatchID));
81  subPatchStarts[ownerProc].append(labelList(1, label(0)));
82 
83  if (facei != -1 && completeMesh().isInternalFace(facei))
84  {
85  toOwnerProcPatchi = nbrToInterPatch[nbrProc].size();
86  nbrToInterPatch[nbrProc].insert(ownerProc, toOwnerProcPatchi);
87  interPatchFaces[nbrProc].append(DynamicList<label>());
88  subPatchIDs[nbrProc].append(labelList(1, subPatchID));
89  subPatchStarts[nbrProc].append(labelList(1, label(0)));
90  }
91  }
92  else
93  {
94  toNbrProcPatchi = nbrToInterPatch[ownerProc][nbrProc];
95 
96  if (facei != -1 && completeMesh().isInternalFace(facei))
97  {
98  toOwnerProcPatchi = nbrToInterPatch[nbrProc][ownerProc];
99  }
100  }
101 
102  // If the sub patch has changed then add new sub-patch entries
103  if (subPatchIDs[ownerProc][toNbrProcPatchi].last() != subPatchID)
104  {
105  subPatchIDs[ownerProc][toNbrProcPatchi].append(subPatchID);
106  subPatchStarts[ownerProc][toNbrProcPatchi].append
107  (
108  interPatchFaces[ownerProc][toNbrProcPatchi].size()
109  );
110 
111  if (facei != -1 && completeMesh().isInternalFace(facei))
112  {
113  subPatchIDs[nbrProc][toOwnerProcPatchi].append(subPatchID);
114  subPatchStarts[nbrProc][toOwnerProcPatchi].append
115  (
116  interPatchFaces[nbrProc][toOwnerProcPatchi].size()
117  );
118  }
119  }
120 
121  // Add face to the inter-proc patches. Note use of turning index.
122  if (facei != -1)
123  {
124  interPatchFaces[ownerProc][toNbrProcPatchi].append(facei + 1);
125 
126  if (completeMesh().isInternalFace(facei))
127  {
128  interPatchFaces[nbrProc][toOwnerProcPatchi].append(- facei - 1);
129  }
130  }
131 }
132 
133 
134 Foam::labelList Foam::domainDecomposition::distributeCells()
135 {
136  Info<< "\nCalculating distribution of cells" << endl;
137 
138  cpuTime decompositionTime;
139 
140  const dictionary decomposeParDict =
142 
143  scalarField cellWeights;
144  if (decomposeParDict.found("weightField"))
145  {
146  const word weightName = decomposeParDict.lookup("weightField");
147 
148  volScalarField weights
149  (
150  IOobject
151  (
152  weightName,
153  completeMesh().time().timeName(),
154  completeMesh(),
157  ),
158  completeMesh()
159  );
160  cellWeights = weights.primitiveField();
161  }
162 
163  const labelList result =
164  decompositionMethod::NewDecomposer(decomposeParDict)->decompose
165  (
166  completeMesh(),
167  cellWeights
168  );
169 
170  Info<< "\nFinished decomposition in "
171  << decompositionTime.elapsedCpuTime()
172  << " s" << endl;
173 
174  return result;
175 }
176 
177 
178 template<class BinaryOp>
179 inline void Foam::domainDecomposition::processInterCyclics
180 (
181  const labelList& cellToProc,
182  const polyBoundaryMesh& patches,
183  List<DynamicList<DynamicList<label>>>& interPatchFaces,
184  List<Map<label>>& procNbrToInterPatch,
185  List<labelListList>& subPatchIDs,
186  List<labelListList>& subPatchStarts,
187  bool owner,
188  BinaryOp bop
189 ) const
190 {
191  // Processor boundaries from split cyclics
192  forAll(patches, patchi)
193  {
194  if (isA<nonConformalCyclicPolyPatch>(patches[patchi]))
195  {
196  const nonConformalCyclicPolyPatch& nccPp =
197  refCast<const nonConformalCyclicPolyPatch>(patches[patchi]);
198 
199  if (nccPp.owner() != owner) continue;
200 
201  const polyPatch& origPp = nccPp.origPatch();
202  const polyPatch& nbrOrigPp = nccPp.nbrPatch().origPatch();
203 
204  const labelUList& origPatchFaceCells = origPp.faceCells();
205  const labelUList& nbrOrigPatchFaceCells = nbrOrigPp.faceCells();
206 
207  // Add all possible interfaces between processors that contain the
208  // original patches
209  forAll(origPatchFaceCells, origFacei)
210  {
211  const label ownerProc =
212  cellToProc[origPatchFaceCells[origFacei]];
213 
214  forAll(nbrOrigPatchFaceCells, nbrOrigFacei)
215  {
216  const label nbrProc =
217  cellToProc[nbrOrigPatchFaceCells[nbrOrigFacei]];
218 
219  if (bop(ownerProc, nbrProc))
220  {
221  addInterProcFace
222  (
223  -1,
224  ownerProc,
225  nbrProc,
226  patchi,
227  procNbrToInterPatch,
228  interPatchFaces,
229  subPatchIDs,
230  subPatchStarts
231  );
232  }
233  }
234  }
235  }
236  else if (isA<cyclicPolyPatch>(patches[patchi]))
237  {
238  const cyclicPolyPatch& cPp =
239  refCast<const cyclicPolyPatch>(patches[patchi]);
240 
241  if (cPp.owner() != owner) continue;
242 
243  const labelUList& patchFaceCells = cPp.faceCells();
244  const labelUList& nbrPatchFaceCells = cPp.nbrPatch().faceCells();
245 
246  // Add faces with different owner and neighbour processors
247  forAll(patchFaceCells, facei)
248  {
249  const label ownerProc = cellToProc[patchFaceCells[facei]];
250  const label nbrProc = cellToProc[nbrPatchFaceCells[facei]];
251 
252  if (bop(ownerProc, nbrProc))
253  {
254  addInterProcFace
255  (
256  cPp.start() + facei,
257  ownerProc,
258  nbrProc,
259  patchi,
260  procNbrToInterPatch,
261  interPatchFaces,
262  subPatchIDs,
263  subPatchStarts
264  );
265  }
266  }
267  }
268  }
269 }
270 
271 
272 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
273 
275 {
276  // Decide which cell goes to which processor
277  const labelList cellToProc = distributeCells();
278 
279  // Distribute the cells according to the given processor label
280 
281  // calculate the addressing information for the original mesh
282  Info<< "\nCalculating original mesh data" << endl;
283 
284  // set references to the original mesh
285  const polyBoundaryMesh& patches = completeMesh().boundaryMesh();
286  const faceList& fcs = completeMesh().faces();
287  const labelList& owner = completeMesh().faceOwner();
288  const labelList& neighbour = completeMesh().faceNeighbour();
289 
290  // loop through the list of processor labels for the cell and add the
291  // cell shape to the list of cells for the appropriate processor
292 
293  Info<< "\nDistributing cells to processors" << endl;
294 
295  // Cells per processor
296  procCellAddressing_ = invertOneToMany(nProcs(), cellToProc);
297 
298  Info<< "\nDistributing faces to processors" << endl;
299 
300  // Loop through all internal faces and decide which processor they belong to
301  // First visit all internal faces. If cells at both sides belong to the
302  // same processor, the face is an internal face. If they are different,
303  // it belongs to both processors.
304  List<DynamicList<label>> dynProcFaceAddressing(nProcs());
305 
306  // Internal faces
307  forAll(neighbour, facei)
308  {
309  if (cellToProc[owner[facei]] == cellToProc[neighbour[facei]])
310  {
311  // Face internal to processor. Notice no turning index.
312  dynProcFaceAddressing[cellToProc[owner[facei]]].append(facei+1);
313  }
314  }
315 
316  // for all processors, set the size of start index and patch size
317  // lists to the number of patches in the mesh
318  labelListList procPatchSize(nProcs());
319  labelListList procPatchStartIndex(nProcs());
320  forAll(procPatchSize, proci)
321  {
322  procPatchSize[proci].setSize(patches.size());
323  procPatchStartIndex[proci].setSize(patches.size());
324  }
325 
326  // Patch faces
327  forAll(patches, patchi)
328  {
329  // Reset size and start index for all processors
330  forAll(procPatchSize, proci)
331  {
332  procPatchSize[proci][patchi] = 0;
333  procPatchStartIndex[proci][patchi] =
334  dynProcFaceAddressing[proci].size();
335  }
336 
337  const label patchStart = patches[patchi].start();
338 
339  if (!isA<cyclicPolyPatch>(patches[patchi]))
340  {
341  // Normal patch. Add faces to processor where the cell
342  // next to the face lives.
343 
344  const labelUList& patchFaceCells = patches[patchi].faceCells();
345 
346  forAll(patchFaceCells, facei)
347  {
348  const label curProc = cellToProc[patchFaceCells[facei]];
349 
350  // add the face without turning index
351  dynProcFaceAddressing[curProc].append(patchStart+facei+1);
352 
353  // increment the number of faces for this patch
354  procPatchSize[curProc][patchi]++;
355  }
356  }
357  else
358  {
359  // Cyclic patch. Add if the opposite sides are on the same
360  // processor.
361 
362  const cyclicPolyPatch& pp =
363  refCast<const cyclicPolyPatch>(patches[patchi]);
364 
365  const labelUList& patchFaceCells = pp.faceCells();
366  const labelUList& nbrPatchFaceCells = pp.nbrPatch().faceCells();
367 
368  forAll(patchFaceCells, facei)
369  {
370  const label curProc = cellToProc[patchFaceCells[facei]];
371  const label nbrProc = cellToProc[nbrPatchFaceCells[facei]];
372 
373  if (curProc == nbrProc)
374  {
375  // add the face without turning index
376  dynProcFaceAddressing[curProc].append(patchStart+facei+1);
377 
378  // increment the number of faces for this patch
379  procPatchSize[curProc][patchi]++;
380  }
381  }
382  }
383  }
384 
385  // Done internal bits of the new mesh and the ordinary patches.
386 
387  // Per processor, from neighbour processor to the inter-processor patch
388  // that communicates with that neighbour
389  List<Map<label>> procNbrToInterPatch(nProcs());
390 
391  // Per processor the faces per inter-processor patch
392  List<DynamicList<DynamicList<label>>> interPatchFaces(nProcs());
393 
394  // Sub-patch ID-s and starts for inter-processor patches
395  List<labelListList> subPatchIDs(nProcs());
396  List<labelListList> subPatchStarts(nProcs());
397 
398  // Processor boundaries from internal faces
399  forAll(neighbour, facei)
400  {
401  label ownerProc = cellToProc[owner[facei]];
402  label nbrProc = cellToProc[neighbour[facei]];
403 
404  if (ownerProc != nbrProc)
405  {
406  // inter - processor patch face found.
407  addInterProcFace
408  (
409  facei,
410  ownerProc,
411  nbrProc,
412  -1,
413  procNbrToInterPatch,
414  interPatchFaces,
415  subPatchIDs,
416  subPatchStarts
417  );
418  }
419  }
420 
421  // Special handling needed for the case that multiple processor cyclic
422  // patches are created on each local processor domain, e.g. if a 3x3 case
423  // is decomposed using the decomposition:
424  //
425  // | 1 | 0 | 2 |
426  // cyclic left | 2 | 0 | 1 | cyclic right
427  // | 2 | 0 | 1 |
428  //
429  // - processors 1 and 2 will both have pieces of both cyclic left- and
430  // right sub-patches present
431  // - the interface patch faces are stored in a single list, where each
432  // sub-patch is referenced into the list using a patch start index and
433  // size
434  // - if the patches are in order (in the boundary file) of left, right
435  // - processor 1 will send: left, right
436  // - processor 1 will need to receive in reverse order: right, left
437  // - similarly for processor 2
438  // - the sub-patches are therefore generated in 4 passes of the patch lists
439  // 1. add faces from owner patch where local proc i < nbr proc i
440  // 2. add faces from nbr patch where local proc i < nbr proc i
441  // 3. add faces from nbr patch where local proc i > nbr proc i
442  // 4. add faces from owner patch where local proc i > nbr proc i
443 
444  processInterCyclics
445  (
446  cellToProc,
447  patches,
448  interPatchFaces,
449  procNbrToInterPatch,
450  subPatchIDs,
451  subPatchStarts,
452  true,
453  lessOp<label>()
454  );
455 
456  processInterCyclics
457  (
458  cellToProc,
459  patches,
460  interPatchFaces,
461  procNbrToInterPatch,
462  subPatchIDs,
463  subPatchStarts,
464  false,
465  lessOp<label>()
466  );
467 
468  processInterCyclics
469  (
470  cellToProc,
471  patches,
472  interPatchFaces,
473  procNbrToInterPatch,
474  subPatchIDs,
475  subPatchStarts,
476  false,
478  );
479 
480  processInterCyclics
481  (
482  cellToProc,
483  patches,
484  interPatchFaces,
485  procNbrToInterPatch,
486  subPatchIDs,
487  subPatchStarts,
488  true,
490  );
491 
492  // Sort inter-proc patch by neighbour
493  labelListList procNeighbourProcessors(nProcs());
494  labelListList procProcessorPatchSize(nProcs());
495  labelListList procProcessorPatchStartIndex(nProcs());
496  List<labelListList> procProcessorPatchSubPatchIDs(nProcs());
497  List<labelListList> procProcessorPatchSubPatchStarts(nProcs());
498  labelList order;
499  forAll(procNbrToInterPatch, proci)
500  {
501  label nInterfaces = procNbrToInterPatch[proci].size();
502 
503  procNeighbourProcessors[proci].setSize(nInterfaces);
504  procProcessorPatchSize[proci].setSize(nInterfaces);
505  procProcessorPatchStartIndex[proci].setSize(nInterfaces);
506  procProcessorPatchSubPatchIDs[proci].setSize(nInterfaces);
507  procProcessorPatchSubPatchStarts[proci].setSize(nInterfaces);
508 
509  // Get sorted neighbour processors
510  const Map<label>& curNbrToInterPatch = procNbrToInterPatch[proci];
511  labelList nbrs = curNbrToInterPatch.toc();
512 
513  sortedOrder(nbrs, order);
514 
515  DynamicList<DynamicList<label>>& curInterPatchFaces =
516  interPatchFaces[proci];
517 
518  forAll(nbrs, i)
519  {
520  const label nbrProc = nbrs[i];
521  const label interPatch = curNbrToInterPatch[nbrProc];
522 
523  procNeighbourProcessors[proci][i] = nbrProc;
524  procProcessorPatchSize[proci][i] =
525  curInterPatchFaces[interPatch].size();
526  procProcessorPatchStartIndex[proci][i] =
527  dynProcFaceAddressing[proci].size();
528 
529  // Add size as last element to substarts and transfer
530  subPatchStarts[proci][interPatch].append
531  (
532  curInterPatchFaces[interPatch].size()
533  );
534  procProcessorPatchSubPatchIDs[proci][i].transfer
535  (
536  subPatchIDs[proci][interPatch]
537  );
538  procProcessorPatchSubPatchStarts[proci][i].transfer
539  (
540  subPatchStarts[proci][interPatch]
541  );
542 
543  // And add all the face labels for interPatch
544  DynamicList<label>& interPatchFaces =
545  curInterPatchFaces[interPatch];
546 
547  forAll(interPatchFaces, j)
548  {
549  dynProcFaceAddressing[proci].append(interPatchFaces[j]);
550  }
551 
552  interPatchFaces.clearStorage();
553  }
554 
555  curInterPatchFaces.clearStorage();
556  dynProcFaceAddressing[proci].shrink();
557  }
558 
559  // Transfer face addressing to non-dynamic member data
560  forAll(dynProcFaceAddressing, proci)
561  {
562  procFaceAddressing_[proci].transfer(dynProcFaceAddressing[proci]);
563  }
564 
565  if (debug)
566  {
567  forAll(procPatchStartIndex, proci)
568  {
569  Info<< "Processor:" << proci << endl;
570  Info<< " total faces:" << procFaceAddressing_[proci].size()
571  << endl;
572 
573  const labelList& curProcPatchStartIndex =
574  procPatchStartIndex[proci];
575 
576  forAll(curProcPatchStartIndex, patchi)
577  {
578  Info<< " patch:" << patchi
579  << "\tstart:" << curProcPatchStartIndex[patchi]
580  << "\tsize:" << procPatchSize[proci][patchi]
581  << endl;
582  }
583  }
584  Info<< endl;
585 
586  forAll(procNeighbourProcessors, proci)
587  {
588  Info<< "Processor " << proci << endl;
589  forAll(procNeighbourProcessors[proci], i)
590  {
591  Info<< " nbr:" << procNeighbourProcessors[proci][i] << endl;
592  Info<< " size:" << procProcessorPatchSize[proci][i] << endl;
593  Info<< " start:" << procProcessorPatchStartIndex[proci][i]
594  << endl;
595  }
596  }
597  Info<< endl;
598  }
599 
600  if (debug > 1)
601  {
602  forAll(procFaceAddressing_, proci)
603  {
604  Info<< "Processor:" << proci << endl;
605  Info<< " faces:" << procFaceAddressing_[proci] << endl;
606  }
607  }
608 
609  Info<< "\nDistributing points to processors" << endl;
610 
611  // For every processor, loop through the list of faces for the processor.
612  // For every face, loop through the list of points and mark the point as
613  // used for the processor. Collect the list of used points for the
614  // processor.
615 
616  forAll(procPointAddressing_, proci)
617  {
619 
620  // Get reference to list of used faces
621  const labelList& procFaceLabels = procFaceAddressing_[proci];
622 
623  forAll(procFaceLabels, facei)
624  {
625  // Because of the turning index, some labels may be negative
626  const labelList& facePoints = fcs[mag(procFaceLabels[facei]) - 1];
627 
628  forAll(facePoints, pointi)
629  {
630  // Mark the point as used
631  pointLabels[facePoints[pointi]] = true;
632  }
633  }
634 
635  // Collect the used points
636  labelList& procPointLabels = procPointAddressing_[proci];
637 
638  procPointLabels.setSize(pointLabels.size());
639 
640  label nUsedPoints = 0;
641 
642  forAll(pointLabels, pointi)
643  {
644  if (pointLabels[pointi])
645  {
646  procPointLabels[nUsedPoints] = pointi;
647 
648  nUsedPoints++;
649  }
650  }
651 
652  // Reset the size of used points
653  procPointLabels.setSize(nUsedPoints);
654  }
655 
656  Info<< "\nConstructing processor meshes" << endl;
657 
658  // Mark point/faces/cells that are in zones.
659  // -1 : not in zone
660  // -2 : in multiple zones
661  // >= 0 : in single given zone
662  // This will give direct lookup of elements that are in a single zone
663  // and we'll only have to revert back to searching through all zones
664  // for the duplicate elements
665 
666  // Point zones
667  labelList pointToZone(completeMesh().points().size(), -1);
668  forAll(completeMesh().pointZones(), zoneI)
669  {
670  mark(completeMesh().pointZones()[zoneI], zoneI, pointToZone);
671  }
672 
673  // Face zones
674  labelList faceToZone(completeMesh().faces().size(), -1);
675  forAll(completeMesh().faceZones(), zoneI)
676  {
677  mark(completeMesh().faceZones()[zoneI], zoneI, faceToZone);
678  }
679 
680  // Cell zones
681  labelList cellToZone(completeMesh().nCells(), -1);
682  forAll(completeMesh().cellZones(), zoneI)
683  {
684  mark(completeMesh().cellZones()[zoneI], zoneI, cellToZone);
685  }
686 
687  // Initialise information for reporting
688  label maxProcCells = 0;
689  label totProcFaces = 0;
690  label maxProcPatches = 0;
691  label totProcPatches = 0;
692  label maxProcFaces = 0;
693 
694  // Generate the meshes
695  for (label proci = 0; proci < nProcs(); proci++)
696  {
697  // Create processor points
698  const labelList& curPointLabels = procPointAddressing_[proci];
699  const pointField& meshPoints = completeMesh().points();
700  labelList pointLookup(completeMesh().nPoints(), -1);
701  pointField procPoints(curPointLabels.size());
702  forAll(curPointLabels, pointi)
703  {
704  procPoints[pointi] = meshPoints[curPointLabels[pointi]];
705  pointLookup[curPointLabels[pointi]] = pointi;
706  }
707 
708  // Create processor faces
709  const labelList& curFaceLabels = procFaceAddressing_[proci];
710  const faceList& meshFaces = completeMesh().faces();
711  labelList faceLookup(completeMesh().nFaces(), -1);
712  faceList procFaces(curFaceLabels.size());
713  forAll(curFaceLabels, facei)
714  {
715  // Mark the original face as used
716  // Remember to decrement the index by one (turning index)
717  label curF = mag(curFaceLabels[facei]) - 1;
718 
719  faceLookup[curF] = facei;
720 
721  // get the original face
722  labelList origFaceLabels;
723 
724  if (curFaceLabels[facei] >= 0)
725  {
726  // face not turned
727  origFaceLabels = meshFaces[curF];
728  }
729  else
730  {
731  origFaceLabels = meshFaces[curF].reverseFace();
732  }
733 
734  // translate face labels into local point list
735  face& procFaceLabels = procFaces[facei];
736 
737  procFaceLabels.setSize(origFaceLabels.size());
738 
739  forAll(origFaceLabels, pointi)
740  {
741  procFaceLabels[pointi] = pointLookup[origFaceLabels[pointi]];
742  }
743  }
744 
745  // Create processor cells
746  const labelList& curCellLabels = procCellAddressing_[proci];
747  const cellList& meshCells = completeMesh().cells();
748  cellList procCells(curCellLabels.size());
749  forAll(curCellLabels, celli)
750  {
751  const labelList& origCellLabels = meshCells[curCellLabels[celli]];
752  cell& curCell = procCells[celli];
753  curCell.setSize(origCellLabels.size());
754  forAll(origCellLabels, cellFacei)
755  {
756  curCell[cellFacei] = faceLookup[origCellLabels[cellFacei]];
757  }
758  }
759 
760  // Create a processor mesh without a boundary
761  procMeshes_.set
762  (
763  proci,
764  new fvMesh
765  (
766  IOobject
767  (
768  regionName_,
769  completeMesh().facesInstance(),
770  runTimes_.procTimes()[proci]
771  ),
772  move(procPoints),
773  move(procFaces),
774  move(procCells)
775  )
776  );
777  fvMesh& procMesh = procMeshes_[proci];
778  procMesh.setPointsInstance(completeMesh().pointsInstance());
779 
780  // Create processor boundary patch information
781  const labelList& curPatchSizes = procPatchSize[proci];
782  const labelList& curPatchStarts = procPatchStartIndex[proci];
783  const labelList& curNeighbourProcessors =
784  procNeighbourProcessors[proci];
785  const labelList& curProcessorPatchSizes =
786  procProcessorPatchSize[proci];
787  const labelList& curProcessorPatchStarts =
788  procProcessorPatchStartIndex[proci];
789  const labelListList& curSubPatchIDs =
790  procProcessorPatchSubPatchIDs[proci];
791  const labelListList& curSubStarts =
792  procProcessorPatchSubPatchStarts[proci];
793  const polyPatchList& meshPatches = completeMesh().boundaryMesh();
794 
795  // Count the number of inter-proc patches
796  label nInterProcPatches = 0;
797  forAll(curSubPatchIDs, procPatchi)
798  {
799  nInterProcPatches += curSubPatchIDs[procPatchi].size();
800  }
801  List<polyPatch*> procPatches
802  (
803  curPatchSizes.size() + nInterProcPatches,
804  nullptr
805  );
806 
807  label nPatches = 0;
808 
809  // Map existing non-proc patches
810  forAll(curPatchSizes, patchi)
811  {
812  const polyPatch& meshPatch = meshPatches[patchi];
813 
814  procPatches[nPatches] =
815  meshPatch.clone
816  (
817  procMesh.boundaryMesh(),
818  nPatches,
820  (
821  curFaceLabels,
822  curPatchSizes[patchi],
823  curPatchStarts[patchi]
824  ) - 1 - meshPatch.start(),
825  curPatchStarts[patchi]
826  ).ptr();
827 
828  nPatches++;
829  }
830 
831  // Create new inter-proc patches
832  forAll(curProcessorPatchSizes, procPatchi)
833  {
834  const labelList& subPatchID = curSubPatchIDs[procPatchi];
835  const labelList& subStarts = curSubStarts[procPatchi];
836 
837  label curStart = curProcessorPatchStarts[procPatchi];
838 
839  forAll(subPatchID, i)
840  {
841  const label size =
842  i < subPatchID.size() - 1
843  ? subStarts[i+1] - subStarts[i]
844  : curProcessorPatchSizes[procPatchi] - subStarts[i];
845 
846  if (subPatchID[i] == -1)
847  {
848  // Processor patch from internal faces
849  procPatches[nPatches] =
851  (
852  size,
853  curStart,
854  nPatches,
855  procMesh.boundaryMesh(),
856  proci,
857  curNeighbourProcessors[procPatchi]
858  );
859  }
860  else if
861  (
862  isA<nonConformalCyclicPolyPatch>
863  (
864  completeMesh().boundaryMesh()[subPatchID[i]]
865  )
866  )
867  {
868  // Non-conformal processor cyclic patch from split
869  // non-conformal cyclic patch
870  const nonConformalCyclicPolyPatch& nccPp =
871  refCast<const nonConformalCyclicPolyPatch>
872  (
873  completeMesh().boundaryMesh()[subPatchID[i]]
874  );
875 
876  procPatches[nPatches] =
878  (
879  size,
880  curStart,
881  nPatches,
882  procMesh.boundaryMesh(),
883  proci,
884  curNeighbourProcessors[procPatchi],
885  nccPp.name(),
886  nccPp.origPatch().name()
887  );
888  }
889  else if
890  (
891  isA<cyclicPolyPatch>
892  (
893  completeMesh().boundaryMesh()[subPatchID[i]]
894  )
895  )
896  {
897  // Processor cyclic patch from split cyclic faces
898  const cyclicPolyPatch& cPp =
899  refCast<const cyclicPolyPatch>
900  (
901  completeMesh().boundaryMesh()[subPatchID[i]]
902  );
903 
904  procPatches[nPatches] =
906  (
907  size,
908  curStart,
909  nPatches,
910  procMesh.boundaryMesh(),
911  proci,
912  curNeighbourProcessors[procPatchi],
913  cPp.name()
914  );
915  }
916  else
917  {
919  << "Sub patch ID set for non-cyclic patch type"
920  << exit(FatalError);
921  }
922 
923  curStart += size;
924 
925  nPatches++;
926  }
927  }
928 
929  // Add patches to the mesh
930  procMesh.addFvPatches(procPatches);
931 
932  // Create point zones
933  {
934  const meshPointZones& pz = completeMesh().pointZones();
935 
936  // Go through all the zoned points and find out if they
937  // belong to a zone. If so, add it to the zone as
938  // necessary
939  List<DynamicList<label>> zonePoints(pz.size());
940 
941  // Estimate size
942  forAll(zonePoints, zoneI)
943  {
944  zonePoints[zoneI].setCapacity(pz[zoneI].size()/nProcs());
945  }
946 
947  // Use the pointToZone map to find out the single zone (if any),
948  // use slow search only for shared points.
949  forAll(curPointLabels, pointi)
950  {
951  label curPoint = curPointLabels[pointi];
952 
953  label zoneI = pointToZone[curPoint];
954 
955  if (zoneI >= 0)
956  {
957  // Single zone.
958  zonePoints[zoneI].append(pointi);
959  }
960  else if (zoneI == -2)
961  {
962  // Multiple zones. Lookup.
963  forAll(pz, zoneI)
964  {
965  label index = pz[zoneI].whichPoint(curPoint);
966 
967  if (index != -1)
968  {
969  zonePoints[zoneI].append(pointi);
970  }
971  }
972  }
973  }
974 
975  procMesh.pointZones().clearAddressing();
976  procMesh.pointZones().setSize(zonePoints.size());
977  forAll(zonePoints, zoneI)
978  {
979  procMesh.pointZones().set
980  (
981  zoneI,
982  pz[zoneI].clone
983  (
984  procMesh.pointZones(),
985  zoneI,
986  zonePoints[zoneI].shrink()
987  )
988  );
989  }
990 
991  if (pz.size())
992  {
993  // Force writing on all processors
995  }
996  }
997 
998  // Create face zones
999  {
1000  const meshFaceZones& fz = completeMesh().faceZones();
1001 
1002  // Go through all the zoned face and find out if they
1003  // belong to a zone. If so, add it to the zone as
1004  // necessary
1005  List<DynamicList<label>> zoneFaces(fz.size());
1006  List<DynamicList<bool>> zoneFaceFlips(fz.size());
1007 
1008  // Estimate size
1009  forAll(zoneFaces, zoneI)
1010  {
1011  label procSize = fz[zoneI].size()/nProcs();
1012 
1013  zoneFaces[zoneI].setCapacity(procSize);
1014  zoneFaceFlips[zoneI].setCapacity(procSize);
1015  }
1016 
1017  // Go through all the zoned faces and find out if they
1018  // belong to a zone. If so, add it to the zone as
1019  // necessary
1020  forAll(curFaceLabels, facei)
1021  {
1022  // Remember to decrement the index by one (turning index)
1023  //
1024  label curF = mag(curFaceLabels[facei]) - 1;
1025 
1026  label zoneI = faceToZone[curF];
1027 
1028  if (zoneI >= 0)
1029  {
1030  // Single zone. Add the face
1031  zoneFaces[zoneI].append(facei);
1032 
1033  label index = fz[zoneI].whichFace(curF);
1034 
1035  bool flip = fz[zoneI].flipMap()[index];
1036 
1037  if (curFaceLabels[facei] < 0)
1038  {
1039  flip = !flip;
1040  }
1041 
1042  zoneFaceFlips[zoneI].append(flip);
1043  }
1044  else if (zoneI == -2)
1045  {
1046  // Multiple zones. Lookup.
1047  forAll(fz, zoneI)
1048  {
1049  label index = fz[zoneI].whichFace(curF);
1050 
1051  if (index != -1)
1052  {
1053  zoneFaces[zoneI].append(facei);
1054 
1055  bool flip = fz[zoneI].flipMap()[index];
1056 
1057  if (curFaceLabels[facei] < 0)
1058  {
1059  flip = !flip;
1060  }
1061 
1062  zoneFaceFlips[zoneI].append(flip);
1063  }
1064  }
1065  }
1066  }
1067 
1068  procMesh.faceZones().clearAddressing();
1069  procMesh.faceZones().setSize(zoneFaces.size());
1070  forAll(zoneFaces, zoneI)
1071  {
1072  procMesh.faceZones().set
1073  (
1074  zoneI,
1075  fz[zoneI].clone
1076  (
1077  zoneFaces[zoneI].shrink(), // addressing
1078  zoneFaceFlips[zoneI].shrink(), // flipmap
1079  zoneI,
1080  procMesh.faceZones()
1081  )
1082  );
1083  }
1084 
1085  if (fz.size())
1086  {
1087  // Force writing on all processors
1088  procMesh.faceZones().writeOpt() = IOobject::AUTO_WRITE;
1089  }
1090  }
1091 
1092  // Create cell zones
1093  {
1094  const meshCellZones& cz = completeMesh().cellZones();
1095 
1096  // Go through all the zoned cells and find out if they
1097  // belong to a zone. If so, add it to the zone as
1098  // necessary
1099  List<DynamicList<label>> zoneCells(cz.size());
1100 
1101  // Estimate size
1102  forAll(zoneCells, zoneI)
1103  {
1104  zoneCells[zoneI].setCapacity(cz[zoneI].size()/nProcs());
1105  }
1106 
1107  forAll(curCellLabels, celli)
1108  {
1109  label curCelli = curCellLabels[celli];
1110 
1111  label zoneI = cellToZone[curCelli];
1112 
1113  if (zoneI >= 0)
1114  {
1115  // Single zone.
1116  zoneCells[zoneI].append(celli);
1117  }
1118  else if (zoneI == -2)
1119  {
1120  // Multiple zones. Lookup.
1121  forAll(cz, zoneI)
1122  {
1123  label index = cz[zoneI].whichCell(curCelli);
1124 
1125  if (index != -1)
1126  {
1127  zoneCells[zoneI].append(celli);
1128  }
1129  }
1130  }
1131  }
1132 
1133  procMesh.cellZones().clearAddressing();
1134  procMesh.cellZones().setSize(zoneCells.size());
1135  forAll(zoneCells, zoneI)
1136  {
1137  procMesh.cellZones().set
1138  (
1139  zoneI,
1140  cz[zoneI].clone
1141  (
1142  zoneCells[zoneI].shrink(),
1143  zoneI,
1144  procMesh.cellZones()
1145  )
1146  );
1147  }
1148 
1149  if (cz.size())
1150  {
1151  // Force writing on all processors
1152  procMesh.cellZones().writeOpt() = IOobject::AUTO_WRITE;
1153  }
1154  }
1155 
1156  // Report processor and update global statistics
1157  {
1158  Info<< endl
1159  << "Processor " << proci << nl
1160  << " Number of cells = " << procMesh.nCells()
1161  << endl;
1162 
1163  maxProcCells = max(maxProcCells, procMesh.nCells());
1164 
1165  label nBoundaryFaces = 0;
1166  label nProcPatches = 0;
1167  label nProcFaces = 0;
1168 
1169  forAll(procMesh.boundaryMesh(), patchi)
1170  {
1171  if (isA<processorPolyPatch>(procMesh.boundaryMesh()[patchi]))
1172  {
1173  const processorPolyPatch& ppp =
1174  refCast<const processorPolyPatch>
1175  (
1176  procMesh.boundaryMesh()[patchi]
1177  );
1178 
1179  Info<< " Number of faces shared with processor "
1180  << ppp.neighbProcNo() << " = " << ppp.size() << endl;
1181 
1182  nProcPatches++;
1183  nProcFaces += ppp.size();
1184  }
1185  else
1186  {
1187  nBoundaryFaces += procMesh.boundaryMesh()[patchi].size();
1188  }
1189  }
1190 
1191  Info<< " Number of processor patches = " << nProcPatches << nl
1192  << " Number of processor faces = " << nProcFaces << nl
1193  << " Number of boundary faces = " << nBoundaryFaces << endl;
1194 
1195  totProcFaces += nProcFaces;
1196  totProcPatches += nProcPatches;
1197  maxProcPatches = max(maxProcPatches, nProcPatches);
1198  maxProcFaces = max(maxProcFaces, nProcFaces);
1199  }
1200  }
1201 
1202  // Determine the average number of processor elements
1203  scalar avgProcCells = scalar(completeMesh().nCells())/nProcs();
1204  scalar avgProcPatches = scalar(totProcPatches)/nProcs();
1205  scalar avgProcFaces = scalar(totProcFaces)/nProcs();
1206 
1207  // Prevent division by zero in the case of all faces on one processor
1208  if (totProcPatches == 0)
1209  {
1210  avgProcPatches = 1;
1211  }
1212  if (totProcFaces == 0)
1213  {
1214  avgProcFaces = 1;
1215  }
1216 
1217  Info<< nl
1218  << "Number of processor faces = " << totProcFaces/2 << nl
1219  << "Max number of cells = " << maxProcCells
1220  << " (" << 100.0*(maxProcCells-avgProcCells)/avgProcCells
1221  << "% above average " << avgProcCells << ")" << nl
1222  << "Max number of processor patches = " << maxProcPatches
1223  << " (" << 100.0*(maxProcPatches-avgProcPatches)/avgProcPatches
1224  << "% above average " << avgProcPatches << ")" << nl
1225  << "Max number of faces between processors = " << maxProcFaces
1226  << " (" << 100.0*(maxProcFaces-avgProcFaces)/avgProcFaces
1227  << "% above average " << avgProcFaces << ")" << nl
1228  << endl;
1229 
1230  // Unconform any non-conformal parts of the processor meshes
1231  unconform();
1232 }
1233 
1234 
1235 // ************************************************************************* //
label nPatches
Definition: readKivaGrid.H:402
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:453
void append(T *)
Append an element at the end of the list.
Definition: PtrListI.H:39
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
const word & name() const
Return name.
void sortedOrder(const UList< T > &, labelList &order)
Generate the (stable) sort order for the list.
bool set(const label) const
Is element set.
Definition: PtrListI.H:65
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
labelList pointLabels(nPoints, -1)
const meshCellZones & cellZones() const
Return cell zones.
Definition: polyMesh.H:501
error FatalError
void addFvPatches(const List< polyPatch *> &, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: fvMesh.C:643
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1243
labelListList invertOneToMany(const label len, const labelUList &)
Invert one-to-many map. Unmapped elements will be size 0.
Definition: ListOps.C:67
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
label nCells() const
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
void setPointsInstance(const fileName &)
Set the instance for the points files.
Definition: polyMeshIO.C:57
const cellList & cells() const
const PtrList< Time > & procTimes() const
Access the processor run times.
const cyclicPolyPatch & nbrPatch() const
UList< label > labelUList
Definition: UList.H:65
const meshPointZones & pointZones() const
Return point zones.
Definition: polyMesh.H:489
Neighbour processor patch.
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1211
Non-conformal processor cyclic poly patch. As nonConformalCyclicPolyPatch, but the neighbouring patch...
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
virtual autoPtr< polyPatch > clone(const polyBoundaryMesh &bm) const
Construct and return a clone, resetting the boundary mesh.
Definition: polyPatch.H:225
void decompose()
Decompose the complete mesh to create the processor meshes and.
A List obtained as a section of another List.
Definition: SubList.H:53
T clone(const T &t)
Definition: List.H:54
const pointField & points
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
const labelUList & faceCells() const
Return face-cell addressing.
Definition: polyPatch.C:340
label nPoints
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Cyclic plane patch.
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1237
static autoPtr< decompositionMethod > NewDecomposer(const dictionary &decompositionDict)
Return a reference to the selected decomposition method.
word timeName
Definition: getTimeIndex.H:3
List< label > labelList
A List of labels.
Definition: labelList.H:56
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1224
int neighbProcNo() const
Return neighbour processor number.
Non-conformal cyclic poly patch. As nonConformalCoupledPolyPatch, but the neighbouring patch is local...
static dictionary decomposeParDict(const Time &time)
Read and return the decomposeParDict.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:131
Foam::polyBoundaryMesh.
const polyPatch & origPatch() const
Original patch.
static const char nl
Definition: Ostream.H:260
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
void setSize(const label)
Reset size of List.
Definition: List.C:281
writeOption writeOpt() const
Definition: IOobject.H:375
label patchi
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:57
const meshFaceZones & faceZones() const
Return face zones.
Definition: polyMesh.H:495
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:306
void clearStorage()
Clear the list and delete storage.
Definition: DynamicListI.H:243
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
List< Key > toc() const
Return the table of contents.
Definition: HashTable.C:202
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
const fvMesh & completeMesh() const
Access the global mesh.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
void clearAddressing()
Clear addressing.
Definition: MeshZones.C:387
bool found
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:342
const Time & completeTime() const
Access the complete run time.
label nProcs() const
Return the number of processors in the decomposition.