domainDecomposition.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2019-2025 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 "domainDecomposition.H"
30 #include "dictionary.H"
31 #include "labelIOList.H"
32 #include "processorPolyPatch.H"
34 #include "fvMesh.H"
35 #include "OSspecific.H"
36 #include "Map.H"
37 #include "DynamicList.H"
38 #include "fvFieldDecomposer.H"
39 #include "IOobjectList.H"
40 #include "PtrDynList.H"
41 #include "cellSet.H"
42 #include "faceSet.H"
43 #include "pointSet.H"
44 #include "decompositionModel.H"
45 #include "hexRef8Data.H"
46 
47 // For handling pointMeshes with additional patches
48 #include "pointMesh.H"
49 #include "meshPointPatch.H"
50 #include "processorPointPatch.H"
51 #include "DynamicField.H"
52 
53 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
54 
55 void Foam::domainDecomposition::mark
56 (
57  const labelList& zoneElems,
58  const label zoneI,
59  labelList& elementToZone
60 )
61 {
62  for (const label pointi : zoneElems)
63  {
64  if (elementToZone[pointi] == -1)
65  {
66  // First occurrence
67  elementToZone[pointi] = zoneI;
68  }
69  else if (elementToZone[pointi] >= 0)
70  {
71  // Multiple zones
72  elementToZone[pointi] = -2;
73  }
74  }
75 }
76 
77 
78 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
79 
81 (
82  const IOobject& io,
83  const fileName& decompDictFile
84 )
85 :
86  fvMesh(io),
87  facesInstancePointsPtr_
88  (
89  pointsInstance() != facesInstance()
90  ? new pointIOField
91  (
92  IOobject
93  (
94  "points",
95  facesInstance(),
96  polyMesh::meshSubDir,
97  *this,
98  IOobject::MUST_READ,
99  IOobject::NO_WRITE,
100  IOobject::NO_REGISTER
101  )
102  )
103  : nullptr
104  ),
105  decompDictFile_(decompDictFile),
106  nProcs_
107  (
108  decompositionMethod::nDomains
109  (
110  decompositionModel::New
111  (
112  *this,
113  decompDictFile
114  )
115  )
116  ),
117  distributed_(false),
118  cellToProc_(nCells()),
119  procPointAddressing_(nProcs_),
120  procFaceAddressing_(nProcs_),
121  procCellAddressing_(nProcs_),
122  procPatchSize_(nProcs_),
123  procPatchStartIndex_(nProcs_),
124  procNeighbourProcessors_(nProcs_),
125  procProcessorPatchSize_(nProcs_),
126  procProcessorPatchStartIndex_(nProcs_),
127  procProcessorPatchSubPatchIDs_(nProcs_),
128  procProcessorPatchSubPatchStarts_(nProcs_)
129 {
130  updateParameters(this->model());
131 }
132 
133 
134 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
135 
137 {
138  return decompositionModel::New(*this, decompDictFile_);
139 }
140 
141 
143 (
144  const dictionary& params
145 )
146 {
147  params.readIfPresent("distributed", distributed_);
148 }
149 
150 
151 bool Foam::domainDecomposition::writeDecomposition(const bool decomposeSets)
152 {
153  Info<< "\nConstructing processor meshes" << endl;
154 
155  // Mark point/faces/cells that are in zones.
156  // -1 : not in zone
157  // -2 : in multiple zones
158  // >= 0 : in single given zone
159  // This will give direct lookup of elements that are in a single zone
160  // and we'll only have to revert back to searching through all zones
161  // for the duplicate elements
162 
163  // Point zones
164  labelList pointToZone(points().size(), -1);
165 
166  forAll(pointZones(), zonei)
167  {
168  mark(pointZones()[zonei], zonei, pointToZone);
169  }
170 
171  // Face zones
172  labelList faceToZone(faces().size(), -1);
173 
174  forAll(faceZones(), zonei)
175  {
176  mark(faceZones()[zonei], zonei, faceToZone);
177  }
178 
179  // Cell zones
180  labelList cellToZone(nCells(), -1);
181 
182  forAll(cellZones(), zonei)
183  {
184  mark(cellZones()[zonei], zonei, cellToZone);
185  }
186 
187 
188  PtrDynList<const cellSet> cellSets;
189  PtrDynList<const faceSet> faceSets;
190  PtrDynList<const pointSet> pointSets;
191  if (decomposeSets)
192  {
193  // Read sets
194  IOobjectList objects(*this, facesInstance(), "polyMesh/sets");
195  for (const IOobject& io : objects.csorted<cellSet>())
196  {
197  cellSets.emplace_back(io);
198  }
199  for (const IOobject& io : objects.csorted<faceSet>())
200  {
201  faceSets.emplace_back(io);
202  }
203  for (const IOobject& io : objects.csorted<pointSet>())
204  {
205  pointSets.emplace_back(io);
206  }
207  }
208 
209 
210  // Load refinement data (if any)
211  hexRef8Data baseMeshData
212  (
213  IOobject
214  (
215  "dummy",
216  facesInstance(),
218  *this,
222  )
223  );
224 
225 
226  label maxProcCells = 0;
227  label maxProcFaces = 0;
228  label totProcFaces = 0;
229  label maxProcPatches = 0;
230  label totProcPatches = 0;
231 
232  // Write out the meshes
233  for (label proci = 0; proci < nProcs_; proci++)
234  {
235  // Create processor points
236  const labelList& curPointLabels = procPointAddressing_[proci];
237 
238  const pointField& meshPoints = points();
239 
240  labelList pointLookup(nPoints(), -1);
241 
242  pointField procPoints(curPointLabels.size());
243 
244  forAll(curPointLabels, pointi)
245  {
246  procPoints[pointi] = meshPoints[curPointLabels[pointi]];
247 
248  pointLookup[curPointLabels[pointi]] = pointi;
249  }
250 
251  // Create processor faces
252  const labelList& curFaceLabels = procFaceAddressing_[proci];
253 
254  const faceList& meshFaces = faces();
255 
256  labelList faceLookup(nFaces(), -1);
257 
258  faceList procFaces(curFaceLabels.size());
259 
260  forAll(curFaceLabels, facei)
261  {
262  // Mark the original face as used
263  // Remember to decrement the index by one (turning index)
264  label curF = mag(curFaceLabels[facei]) - 1;
265 
266  faceLookup[curF] = facei;
267 
268  // get the original face
269  labelList origFaceLabels;
270 
271  if (curFaceLabels[facei] >= 0)
272  {
273  // face not turned
274  origFaceLabels = meshFaces[curF];
275  }
276  else
277  {
278  origFaceLabels = meshFaces[curF].reverseFace();
279  }
280 
281  // translate face labels into local point list
282  face& procFaceLabels = procFaces[facei];
283 
284  procFaceLabels.setSize(origFaceLabels.size());
285 
286  forAll(origFaceLabels, pointi)
287  {
288  procFaceLabels[pointi] = pointLookup[origFaceLabels[pointi]];
289  }
290  }
291 
292  // Create processor cells
293  const labelList& curCellLabels = procCellAddressing_[proci];
294 
295  const cellList& meshCells = cells();
296 
297  cellList procCells(curCellLabels.size());
298 
299  forAll(curCellLabels, celli)
300  {
301  const labelList& origCellLabels = meshCells[curCellLabels[celli]];
302 
303  cell& curCell = procCells[celli];
304 
305  curCell.setSize(origCellLabels.size());
306 
307  forAll(origCellLabels, cellFacei)
308  {
309  curCell[cellFacei] = faceLookup[origCellLabels[cellFacei]];
310  }
311  }
312 
313  // Create processor mesh without a boundary
314 
315  // create a database
316  Time processorDb
317  (
319  time().rootPath(),
320  time().caseName()/("processor" + Foam::name(proci)),
321  false, // No function objects
322  false // No extra controlDict libs
323  );
324  processorDb.setTime(time());
325 
326  // create the mesh. Two situations:
327  // - points and faces come from the same time ('instance'). The mesh
328  // will get constructed in the same instance.
329  // - points come from a different time (moving mesh cases).
330  // It will read the points belonging to the faces instance and
331  // construct the procMesh with it which then gets handled as above.
332  // (so with 'old' geometry).
333  // Only at writing time will it additionally write the current
334  // points.
335 
336  autoPtr<polyMesh> procMeshPtr;
337 
338  if (facesInstancePointsPtr_)
339  {
340  // Construct mesh from facesInstance.
341  pointField facesInstancePoints
342  (
343  facesInstancePointsPtr_(),
344  curPointLabels
345  );
346 
347  procMeshPtr = autoPtr<polyMesh>::New
348  (
349  IOobject
350  (
351  this->polyMesh::name(), // region of undecomposed mesh
352  facesInstance(),
353  processorDb,
356  ),
357  std::move(facesInstancePoints),
358  std::move(procFaces),
359  std::move(procCells)
360  );
361  }
362  else
363  {
364  procMeshPtr = autoPtr<polyMesh>::New
365  (
366  IOobject
367  (
368  this->polyMesh::name(), // region of undecomposed mesh
369  facesInstance(),
370  processorDb,
373  ),
374  std::move(procPoints),
375  std::move(procFaces),
376  std::move(procCells)
377  );
378  }
379  polyMesh& procMesh = procMeshPtr();
380 
381 
382  // Create processor boundary patches
383  const labelList& curPatchSizes = procPatchSize_[proci];
384 
385  const labelList& curPatchStarts = procPatchStartIndex_[proci];
386 
387  const labelList& curNeighbourProcessors =
388  procNeighbourProcessors_[proci];
389 
390  const labelList& curProcessorPatchSizes =
391  procProcessorPatchSize_[proci];
392 
393  const labelList& curProcessorPatchStarts =
394  procProcessorPatchStartIndex_[proci];
395 
396  const labelListList& curSubPatchIDs =
397  procProcessorPatchSubPatchIDs_[proci];
398 
399  const labelListList& curSubStarts =
400  procProcessorPatchSubPatchStarts_[proci];
401 
402  const polyPatchList& meshPatches = boundaryMesh();
403 
404 
405  // Count the number of inter-proc patches
406  label nInterProcPatches = 0;
407  forAll(curSubPatchIDs, procPatchi)
408  {
409  nInterProcPatches += curSubPatchIDs[procPatchi].size();
410  }
411 
412  polyPatchList procPatches
413  (
414  curPatchSizes.size() + nInterProcPatches
415  );
416 
417  label nPatches = 0;
418 
419  forAll(curPatchSizes, patchi)
420  {
421  // Get the face labels consistent with the field mapping
422  // (reuse the patch field mappers)
423  const polyPatch& meshPatch = meshPatches[patchi];
424 
425  fvFieldDecomposer::patchFieldDecomposer patchMapper
426  (
427  SubList<label>
428  (
429  curFaceLabels,
430  curPatchSizes[patchi],
431  curPatchStarts[patchi]
432  ),
433  meshPatch.start()
434  );
435 
436  // Map existing patches
437  procPatches.set
438  (
439  nPatches,
440  meshPatch.clone
441  (
442  procMesh.boundaryMesh(),
443  nPatches,
444  patchMapper.directAddressing(),
445  curPatchStarts[patchi]
446  )
447  );
448 
449  nPatches++;
450  }
451 
452  forAll(curProcessorPatchSizes, procPatchi)
453  {
454  const labelList& subPatchID = curSubPatchIDs[procPatchi];
455  const labelList& subStarts = curSubStarts[procPatchi];
456 
457  label curStart = curProcessorPatchStarts[procPatchi];
458 
459  forAll(subPatchID, i)
460  {
461  label size =
462  (
463  i < subPatchID.size()-1
464  ? subStarts[i+1] - subStarts[i]
465  : curProcessorPatchSizes[procPatchi] - subStarts[i]
466  );
467 
468  if (subPatchID[i] == -1)
469  {
470  // From internal faces
471  procPatches.set
472  (
473  nPatches,
474  new processorPolyPatch
475  (
476  size,
477  curStart,
478  nPatches,
479  procMesh.boundaryMesh(),
480  proci,
481  curNeighbourProcessors[procPatchi]
482  )
483  );
484  }
485  else
486  {
487  const coupledPolyPatch& pcPatch
488  = refCast<const coupledPolyPatch>
489  (
490  boundaryMesh()[subPatchID[i]]
491  );
492 
493  procPatches.set
494  (
495  nPatches,
496  new processorCyclicPolyPatch
497  (
498  size,
499  curStart,
500  nPatches,
501  procMesh.boundaryMesh(),
502  proci,
503  curNeighbourProcessors[procPatchi],
504  pcPatch.name(),
505  pcPatch.transform()
506  )
507  );
508  }
509 
510  curStart += size;
511  ++nPatches;
512  }
513  }
514 
515  // Add boundary patches
516  procMesh.addPatches(procPatches);
517 
518  // Create and add zones
519 
520  // Point zones
521  {
522  const pointZoneMesh& pz = pointZones();
523 
524  // Go through all the zoned points and find out if they
525  // belong to a zone. If so, add it to the zone as
526  // necessary
527  List<DynamicList<label>> zonePoints(pz.size());
528 
529  // Estimate size
530  forAll(zonePoints, zonei)
531  {
532  zonePoints[zonei].setCapacity(pz[zonei].size()/nProcs_);
533  }
534 
535  // Use the pointToZone map to find out the single zone (if any),
536  // use slow search only for shared points.
537  forAll(curPointLabels, pointi)
538  {
539  label curPoint = curPointLabels[pointi];
540 
541  label zonei = pointToZone[curPoint];
542 
543  if (zonei >= 0)
544  {
545  // Single zone.
546  zonePoints[zonei].append(pointi);
547  }
548  else if (zonei == -2)
549  {
550  // Multiple zones. Lookup.
551  forAll(pz, zonei)
552  {
553  label index = pz[zonei].whichPoint(curPoint);
554 
555  if (index != -1)
556  {
557  zonePoints[zonei].append(pointi);
558  }
559  }
560  }
561  }
562 
563  procMesh.pointZones().clearAddressing();
564  procMesh.pointZones().setSize(zonePoints.size());
565  forAll(zonePoints, zonei)
566  {
567  procMesh.pointZones().set
568  (
569  zonei,
570  pz[zonei].clone
571  (
572  procMesh.pointZones(),
573  zonei,
574  zonePoints[zonei].shrink()
575  )
576  );
577  }
578 
579  if (pz.size())
580  {
581  // Force writing on all processors
582  procMesh.pointZones().writeOpt(IOobject::AUTO_WRITE);
583  }
584  }
585 
586  // Face zones
587  {
588  const faceZoneMesh& fz = faceZones();
589 
590  // Go through all the zoned face and find out if they
591  // belong to a zone. If so, add it to the zone as
592  // necessary
593  List<DynamicList<label>> zoneFaces(fz.size());
594  List<DynamicList<bool>> zoneFaceFlips(fz.size());
595 
596  // Estimate size
597  forAll(zoneFaces, zonei)
598  {
599  label procSize = fz[zonei].size() / nProcs_;
600 
601  zoneFaces[zonei].setCapacity(procSize);
602  zoneFaceFlips[zonei].setCapacity(procSize);
603  }
604 
605  // Go through all the zoned faces and find out if they
606  // belong to a zone. If so, add it to the zone as
607  // necessary
608  forAll(curFaceLabels, facei)
609  {
610  // Remember to decrement the index by one (turning index)
611  //
612  label curF = mag(curFaceLabels[facei]) - 1;
613 
614  label zonei = faceToZone[curF];
615 
616  if (zonei >= 0)
617  {
618  // Single zone. Add the face
619  zoneFaces[zonei].append(facei);
620 
621  label index = fz[zonei].whichFace(curF);
622 
623  bool flip = fz[zonei].flipMap()[index];
624 
625  if (curFaceLabels[facei] < 0)
626  {
627  flip = !flip;
628  }
629 
630  zoneFaceFlips[zonei].append(flip);
631  }
632  else if (zonei == -2)
633  {
634  // Multiple zones. Lookup.
635  forAll(fz, zonei)
636  {
637  label index = fz[zonei].whichFace(curF);
638 
639  if (index != -1)
640  {
641  zoneFaces[zonei].append(facei);
642 
643  bool flip = fz[zonei].flipMap()[index];
644 
645  if (curFaceLabels[facei] < 0)
646  {
647  flip = !flip;
648  }
649 
650  zoneFaceFlips[zonei].append(flip);
651  }
652  }
653  }
654  }
655 
656  procMesh.faceZones().clearAddressing();
657  procMesh.faceZones().setSize(zoneFaces.size());
658  forAll(zoneFaces, zonei)
659  {
660  procMesh.faceZones().set
661  (
662  zonei,
663  fz[zonei].clone
664  (
665  zoneFaces[zonei].shrink(), // addressing
666  zoneFaceFlips[zonei].shrink(), // flipmap
667  zonei,
668  procMesh.faceZones()
669  )
670  );
671  }
672 
673  if (fz.size())
674  {
675  // Force writing on all processors
676  procMesh.faceZones().writeOpt(IOobject::AUTO_WRITE);
677  }
678  }
679 
680  // Cell zones
681  {
682  const cellZoneMesh& cz = cellZones();
683 
684  // Go through all the zoned cells and find out if they
685  // belong to a zone. If so, add it to the zone as
686  // necessary
687  List<DynamicList<label>> zoneCells(cz.size());
688 
689  // Estimate size
690  forAll(zoneCells, zonei)
691  {
692  zoneCells[zonei].setCapacity(cz[zonei].size()/nProcs_);
693  }
694 
695  forAll(curCellLabels, celli)
696  {
697  label curCelli = curCellLabels[celli];
698 
699  label zonei = cellToZone[curCelli];
700 
701  if (zonei >= 0)
702  {
703  // Single zone.
704  zoneCells[zonei].append(celli);
705  }
706  else if (zonei == -2)
707  {
708  // Multiple zones. Lookup.
709  forAll(cz, zonei)
710  {
711  label index = cz[zonei].whichCell(curCelli);
712 
713  if (index != -1)
714  {
715  zoneCells[zonei].append(celli);
716  }
717  }
718  }
719  }
720 
721  procMesh.cellZones().clearAddressing();
722  procMesh.cellZones().setSize(zoneCells.size());
723  forAll(zoneCells, zonei)
724  {
725  procMesh.cellZones().set
726  (
727  zonei,
728  cz[zonei].clone
729  (
730  zoneCells[zonei].shrink(),
731  zonei,
732  procMesh.cellZones()
733  )
734  );
735  }
736 
737  if (cz.size())
738  {
739  // Force writing on all processors
740  procMesh.cellZones().writeOpt(IOobject::AUTO_WRITE);
741  }
742  }
743 
744  // More precision (for points data)
746 
747  procMesh.write();
748 
749  // Add pointMesh if it was available
750  const auto* pMeshPtr =
751  thisDb().cfindObject<pointMesh>(pointMesh::typeName);
752  if (pMeshPtr)
753  {
754  const auto& pMesh = *pMeshPtr;
755  const auto& pMeshBoundary = pMesh.boundary();
756 
757 
758  // 1. Generate pointBoundaryMesh from polyBoundaryMesh (so ignoring
759  // any additional patches
760  const auto& procPointMesh = pointMesh::New(procMesh);
761 
762  pointBoundaryMesh& procBoundary =
763  const_cast<pointBoundaryMesh&>(procPointMesh.boundary());
764 
765  // Keep track if it differs from the polyBoundaryMesh since then
766  // we need to write the boundary file.
767  bool differsFromPoly = false;
768 
769  // 2. Explicitly add subsetted meshPointPatches
770  forAll(pMeshBoundary, patchi)
771  {
772  const auto* mppPtr = isA<meshPointPatch>(pMeshBoundary[patchi]);
773  if (mppPtr && (procBoundary.findPatchID(mppPtr->name()) == -1))
774  {
775  const auto& mpp = *mppPtr;
776 
777  DynamicList<label> procMeshPoints(mpp.size());
778  DynamicField<vector> procNormals(mpp.size());
779  forAll(mpp.meshPoints(), i)
780  {
781  const label pointi = mpp.meshPoints()[i];
782  const label procPointi = pointLookup[pointi];
783  if (procPointi != -1)
784  {
785  procMeshPoints.append(procPointi);
786  procNormals.append(mpp.pointNormals()[i]);
787  }
788  }
789 
790  procBoundary.push_back
791  (
792  new meshPointPatch
793  (
794  mpp.name(),
795  procMeshPoints,
796  procNormals,
797  procBoundary.size(),
798  procBoundary,
799  meshPointPatch::typeName
800  )
801  );
802  differsFromPoly = true;
803  }
804  }
805 
806  // 3. Shuffle new patches before any processor patches
807  labelList oldToNew(procBoundary.size());
808  label newPatchi = 0;
809  forAll(procBoundary, patchi)
810  {
811  if (!isA<processorPointPatch>(procBoundary[patchi]))
812  {
813  oldToNew[patchi] = newPatchi;
814 
815  if (newPatchi != patchi)
816  {
817  differsFromPoly = true;
818  }
819 
820  newPatchi++;
821  }
822  }
823 
824  // decomposed-to-undecomposed patch numbering
825  labelList boundaryProcAddressing(identity(newPatchi));
826  boundaryProcAddressing.setSize(procBoundary.size(), -1);
827 
828  forAll(procBoundary, patchi)
829  {
830  if (isA<processorPointPatch>(procBoundary[patchi]))
831  {
832  oldToNew[patchi] = newPatchi++;
833  }
834  }
835  procBoundary.reorder(oldToNew, true);
836 
837 
838  if (differsFromPoly)
839  {
840  // Write pointMesh/boundary
841  procBoundary.write();
842 
843  // Write pointMesh/boundaryProcAddressing
844  IOobject ioAddr
845  (
846  "boundaryProcAddressing",
847  procMesh.facesInstance(),
849  procPointMesh.thisDb()
850  );
851  IOList<label>::writeContents(ioAddr, boundaryProcAddressing);
852  }
853  }
854 
855  // Write points if pointsInstance differing from facesInstance
856  if (facesInstancePointsPtr_)
857  {
858  pointIOField pointsInstancePoints
859  (
860  IOobject
861  (
862  "points",
863  pointsInstance(),
865  procMesh,
869  ),
870  std::move(procPoints)
871  );
872  pointsInstancePoints.write();
873  }
874 
875 
876  // Decompose any sets
877  if (decomposeSets)
878  {
879  forAll(cellSets, i)
880  {
881  const cellSet& cs = cellSets[i];
882  cellSet set(procMesh, cs.name(), cs.size()/nProcs_);
883  forAll(curCellLabels, i)
884  {
885  if (cs.found(curCellLabels[i]))
886  {
887  set.insert(i);
888  }
889  }
890  set.write();
891  }
892  forAll(faceSets, i)
893  {
894  const faceSet& cs = faceSets[i];
895  faceSet set(procMesh, cs.name(), cs.size()/nProcs_);
896  forAll(curFaceLabels, i)
897  {
898  if (cs.found(mag(curFaceLabels[i])-1))
899  {
900  set.insert(i);
901  }
902  }
903  set.write();
904  }
905  forAll(pointSets, i)
906  {
907  const pointSet& cs = pointSets[i];
908  pointSet set(procMesh, cs.name(), cs.size()/nProcs_);
909  forAll(curPointLabels, i)
910  {
911  if (cs.found(curPointLabels[i]))
912  {
913  set.insert(i);
914  }
915  }
916  set.write();
917  }
918  }
919 
920 
921  // Optional hexRef8 data
922  hexRef8Data
923  (
924  IOobject
925  (
926  "dummy",
927  facesInstance(),
929  procMesh,
933  ),
934  baseMeshData,
935  procCellAddressing_[proci],
936  procPointAddressing_[proci]
937  ).write();
938 
939 
940  // Statistics
941  Info<< nl << "Processor " << proci;
942 
943  if (procMesh.nCells())
944  {
945  Info<< nl << " ";
946  }
947  else
948  {
949  Info<< ": ";
950  }
951 
952  Info<< "Number of cells = " << procMesh.nCells() << nl;
953 
954  if (procMesh.nCells())
955  {
956  Info<< " Number of points = " << procMesh.nPoints() << nl;
957  }
958 
959  maxProcCells = max(maxProcCells, procMesh.nCells());
960 
961  label nBoundaryFaces = 0;
962  label nProcPatches = 0;
963  label nProcFaces = 0;
964 
965  for (const polyPatch& pp : procMesh.boundaryMesh())
966  {
967  const auto* cpp = isA<processorPolyPatch>(pp);
968 
969  if (cpp)
970  {
971  const auto& procPatch = *cpp;
972 
973  Info<< " Number of faces shared with processor "
974  << procPatch.neighbProcNo() << " = "
975  << procPatch.size() << nl;
976 
977  nProcFaces += procPatch.size();
978  ++nProcPatches;
979  }
980  else
981  {
982  nBoundaryFaces += pp.size();
983  }
984  }
985 
986  if (procMesh.nCells() && (nBoundaryFaces || nProcFaces))
987  {
988  Info<< " Number of processor patches = " << nProcPatches << nl
989  << " Number of processor faces = " << nProcFaces << nl
990  << " Number of boundary faces = " << nBoundaryFaces << nl;
991  }
992 
993  totProcFaces += nProcFaces;
994  totProcPatches += nProcPatches;
995  maxProcFaces = max(maxProcFaces, nProcFaces);
996  maxProcPatches = max(maxProcPatches, nProcPatches);
997 
998  // Write the addressing information
999 
1000  IOobject ioAddr
1001  (
1002  "procAddressing",
1003  procMesh.facesInstance(),
1005  procMesh.thisDb(),
1009  );
1010 
1011  // pointProcAddressing
1012  ioAddr.rename("pointProcAddressing");
1013  IOList<label>::writeContents(ioAddr, procPointAddressing_[proci]);
1014 
1015  // faceProcAddressing
1016  ioAddr.rename("faceProcAddressing");
1017  IOList<label>::writeContents(ioAddr, procFaceAddressing_[proci]);
1018 
1019  // cellProcAddressing
1020  ioAddr.rename("cellProcAddressing");
1021  IOList<label>::writeContents(ioAddr, procCellAddressing_[proci]);
1022 
1023  // Write patch map for backwards compatibility.
1024  // (= identity map for original patches, -1 for processor patches)
1025  label nMeshPatches = curPatchSizes.size();
1026  labelList procBoundaryAddr(identity(nMeshPatches));
1027  procBoundaryAddr.resize(nMeshPatches+nProcPatches, -1);
1028 
1029  // boundaryProcAddressing
1030  ioAddr.rename("boundaryProcAddressing");
1031  IOList<label>::writeContents(ioAddr, procBoundaryAddr);
1032  }
1033 
1034 
1035  // Summary stats
1036  Info<< nl
1037  << "Number of processor faces = " << (totProcFaces/2) << nl
1038  << "Max number of cells = " << maxProcCells;
1039 
1040  if (maxProcCells != nCells())
1041  {
1042  scalar avgValue = scalar(nCells())/nProcs_;
1043 
1044  Info<< " (" << 100.0*(maxProcCells-avgValue)/avgValue
1045  << "% above average " << avgValue << ')';
1046  }
1047  Info<< nl;
1048 
1049  Info<< "Max number of processor patches = " << maxProcPatches;
1050  if (totProcPatches)
1051  {
1052  scalar avgValue = scalar(totProcPatches)/nProcs_;
1053 
1054  Info<< " (" << 100.0*(maxProcPatches-avgValue)/avgValue
1055  << "% above average " << avgValue << ')';
1056  }
1057  Info<< nl;
1058 
1059  Info<< "Max number of faces between processors = " << maxProcFaces;
1060  if (totProcFaces)
1061  {
1062  scalar avgValue = scalar(totProcFaces)/nProcs_;
1063 
1064  Info<< " (" << 100.0*(maxProcFaces-avgValue)/avgValue
1065  << "% above average " << avgValue << ')';
1066  }
1067  Info<< nl << endl;
1068 
1069  return true;
1070 }
1071 
1072 
1073 // ************************************************************************* //
label nPatches
Definition: readKivaGrid.H:394
bool writeDecomposition(const bool decomposeSets)
Write decomposition.
ZoneMesh< faceZone, polyMesh > faceZoneMesh
A ZoneMesh with the type faceZone.
List< cell > cellList
List of cell.
Definition: cellListFwd.H:39
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
vectorIOField pointIOField
pointIOField is a vectorIOField.
Definition: pointIOField.H:38
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
static FOAM_NO_DANGLING_REFERENCE const pointMesh & New(const polyMesh &mesh, Args &&... args)
Get existing or create MeshObject registered with typeName.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:411
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:529
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Ignore writing from objectRegistry::writeObject()
static const decompositionModel & New(const polyMesh &mesh, const fileName &decompDictFile="", const dictionary *fallback=nullptr)
Read and register on mesh, optionally with alternative decomposeParDict path/name or with fallback co...
static void writeContents(const IOobject &io, const UList< T > &content)
Write contents. The IOobject is never registered.
Definition: IOList.C:141
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:286
List< face > faceList
List of faces.
Definition: faceListFwd.H:39
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
PtrList< polyPatch > polyPatchList
Store lists of polyPatch as a PtrList.
Definition: polyPatch.H:56
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
void setSize(const label n)
Alias for resize()
Definition: List.H:325
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const cellShapeList & cells
const pointField & points
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
Definition: labelLists.C:44
label nPoints
Reading is optional [identical to LAZY_READ].
ZoneMesh< pointZone, polyMesh > pointZoneMesh
A ZoneMesh with the type pointZone.
MeshObject wrapper of decompositionMethod.
static unsigned int minPrecision(unsigned int prec) noexcept
Set the minimum default precision.
Definition: IOstream.H:440
static word controlDictName
The default control dictionary name (normally "controlDict")
Definition: Time.H:268
void updateParameters(const dictionary &params)
Update flags based on the decomposition model settings.
const label nProcPatches
static word meshSubDir
Return the mesh sub-directory name (usually "pointMesh")
Definition: pointMesh.H:107
Nothing to be read.
Automatically write from objectRegistry::writeObject()
messageStream Info
Information stream (stdout output on master, null elsewhere)
List< label > labelList
A List of labels.
Definition: List.H:61
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
domainDecomposition(const IOobject &io, const fileName &decompDictFile="")
Construct from components.
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition: autoPtr.H:178
const decompositionModel & model() const
Return decomposition model used.
Do not request registration (bool: false)
ZoneMesh< cellZone, polyMesh > cellZoneMesh
A ZoneMesh with the type cellZone.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())