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-2024 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 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48 
49 void Foam::domainDecomposition::mark
50 (
51  const labelList& zoneElems,
52  const label zoneI,
53  labelList& elementToZone
54 )
55 {
56  for (const label pointi : zoneElems)
57  {
58  if (elementToZone[pointi] == -1)
59  {
60  // First occurrence
61  elementToZone[pointi] = zoneI;
62  }
63  else if (elementToZone[pointi] >= 0)
64  {
65  // Multiple zones
66  elementToZone[pointi] = -2;
67  }
68  }
69 }
70 
71 
72 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
73 
75 (
76  const IOobject& io,
77  const fileName& decompDictFile
78 )
79 :
80  fvMesh(io),
81  facesInstancePointsPtr_
82  (
83  pointsInstance() != facesInstance()
84  ? new pointIOField
85  (
86  IOobject
87  (
88  "points",
89  facesInstance(),
90  polyMesh::meshSubDir,
91  *this,
92  IOobject::MUST_READ,
93  IOobject::NO_WRITE,
94  IOobject::NO_REGISTER
95  )
96  )
97  : nullptr
98  ),
99  decompDictFile_(decompDictFile),
100  nProcs_
101  (
102  decompositionMethod::nDomains
103  (
104  decompositionModel::New
105  (
106  *this,
107  decompDictFile
108  )
109  )
110  ),
111  distributed_(false),
112  cellToProc_(nCells()),
113  procPointAddressing_(nProcs_),
114  procFaceAddressing_(nProcs_),
115  procCellAddressing_(nProcs_),
116  procPatchSize_(nProcs_),
117  procPatchStartIndex_(nProcs_),
118  procNeighbourProcessors_(nProcs_),
119  procProcessorPatchSize_(nProcs_),
120  procProcessorPatchStartIndex_(nProcs_),
121  procProcessorPatchSubPatchIDs_(nProcs_),
122  procProcessorPatchSubPatchStarts_(nProcs_)
123 {
124  updateParameters(this->model());
125 }
126 
127 
128 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
129 
131 {
132  return decompositionModel::New(*this, decompDictFile_);
133 }
134 
135 
137 (
138  const dictionary& params
139 )
140 {
141  params.readIfPresent("distributed", distributed_);
142 }
143 
144 
145 bool Foam::domainDecomposition::writeDecomposition(const bool decomposeSets)
146 {
147  Info<< "\nConstructing processor meshes" << endl;
148 
149  // Mark point/faces/cells that are in zones.
150  // -1 : not in zone
151  // -2 : in multiple zones
152  // >= 0 : in single given zone
153  // This will give direct lookup of elements that are in a single zone
154  // and we'll only have to revert back to searching through all zones
155  // for the duplicate elements
156 
157  // Point zones
158  labelList pointToZone(points().size(), -1);
159 
160  forAll(pointZones(), zonei)
161  {
162  mark(pointZones()[zonei], zonei, pointToZone);
163  }
164 
165  // Face zones
166  labelList faceToZone(faces().size(), -1);
167 
168  forAll(faceZones(), zonei)
169  {
170  mark(faceZones()[zonei], zonei, faceToZone);
171  }
172 
173  // Cell zones
174  labelList cellToZone(nCells(), -1);
175 
176  forAll(cellZones(), zonei)
177  {
178  mark(cellZones()[zonei], zonei, cellToZone);
179  }
180 
181 
182  PtrDynList<const cellSet> cellSets;
183  PtrDynList<const faceSet> faceSets;
184  PtrDynList<const pointSet> pointSets;
185  if (decomposeSets)
186  {
187  // Read sets
188  IOobjectList objects(*this, facesInstance(), "polyMesh/sets");
189  for (const IOobject& io : objects.csorted<cellSet>())
190  {
191  cellSets.emplace_back(io);
192  }
193  for (const IOobject& io : objects.csorted<faceSet>())
194  {
195  faceSets.emplace_back(io);
196  }
197  for (const IOobject& io : objects.csorted<pointSet>())
198  {
199  pointSets.emplace_back(io);
200  }
201  }
202 
203 
204  // Load refinement data (if any)
205  hexRef8Data baseMeshData
206  (
207  IOobject
208  (
209  "dummy",
210  facesInstance(),
212  *this,
216  )
217  );
218 
219 
220  label maxProcCells = 0;
221  label maxProcFaces = 0;
222  label totProcFaces = 0;
223  label maxProcPatches = 0;
224  label totProcPatches = 0;
225 
226  // Write out the meshes
227  for (label proci = 0; proci < nProcs_; proci++)
228  {
229  // Create processor points
230  const labelList& curPointLabels = procPointAddressing_[proci];
231 
232  const pointField& meshPoints = points();
233 
234  labelList pointLookup(nPoints(), -1);
235 
236  pointField procPoints(curPointLabels.size());
237 
238  forAll(curPointLabels, pointi)
239  {
240  procPoints[pointi] = meshPoints[curPointLabels[pointi]];
241 
242  pointLookup[curPointLabels[pointi]] = pointi;
243  }
244 
245  // Create processor faces
246  const labelList& curFaceLabels = procFaceAddressing_[proci];
247 
248  const faceList& meshFaces = faces();
249 
250  labelList faceLookup(nFaces(), -1);
251 
252  faceList procFaces(curFaceLabels.size());
253 
254  forAll(curFaceLabels, facei)
255  {
256  // Mark the original face as used
257  // Remember to decrement the index by one (turning index)
258  label curF = mag(curFaceLabels[facei]) - 1;
259 
260  faceLookup[curF] = facei;
261 
262  // get the original face
263  labelList origFaceLabels;
264 
265  if (curFaceLabels[facei] >= 0)
266  {
267  // face not turned
268  origFaceLabels = meshFaces[curF];
269  }
270  else
271  {
272  origFaceLabels = meshFaces[curF].reverseFace();
273  }
274 
275  // translate face labels into local point list
276  face& procFaceLabels = procFaces[facei];
277 
278  procFaceLabels.setSize(origFaceLabels.size());
279 
280  forAll(origFaceLabels, pointi)
281  {
282  procFaceLabels[pointi] = pointLookup[origFaceLabels[pointi]];
283  }
284  }
285 
286  // Create processor cells
287  const labelList& curCellLabels = procCellAddressing_[proci];
288 
289  const cellList& meshCells = cells();
290 
291  cellList procCells(curCellLabels.size());
292 
293  forAll(curCellLabels, celli)
294  {
295  const labelList& origCellLabels = meshCells[curCellLabels[celli]];
296 
297  cell& curCell = procCells[celli];
298 
299  curCell.setSize(origCellLabels.size());
300 
301  forAll(origCellLabels, cellFacei)
302  {
303  curCell[cellFacei] = faceLookup[origCellLabels[cellFacei]];
304  }
305  }
306 
307  // Create processor mesh without a boundary
308 
309  // create a database
310  Time processorDb
311  (
313  time().rootPath(),
314  time().caseName()/("processor" + Foam::name(proci)),
315  false, // No function objects
316  false // No extra controlDict libs
317  );
318  processorDb.setTime(time());
319 
320  // create the mesh. Two situations:
321  // - points and faces come from the same time ('instance'). The mesh
322  // will get constructed in the same instance.
323  // - points come from a different time (moving mesh cases).
324  // It will read the points belonging to the faces instance and
325  // construct the procMesh with it which then gets handled as above.
326  // (so with 'old' geometry).
327  // Only at writing time will it additionally write the current
328  // points.
329 
330  autoPtr<polyMesh> procMeshPtr;
331 
332  if (facesInstancePointsPtr_)
333  {
334  // Construct mesh from facesInstance.
335  pointField facesInstancePoints
336  (
337  facesInstancePointsPtr_(),
338  curPointLabels
339  );
340 
341  procMeshPtr = autoPtr<polyMesh>::New
342  (
343  IOobject
344  (
345  this->polyMesh::name(), // region of undecomposed mesh
346  facesInstance(),
347  processorDb,
350  ),
351  std::move(facesInstancePoints),
352  std::move(procFaces),
353  std::move(procCells)
354  );
355  }
356  else
357  {
358  procMeshPtr = autoPtr<polyMesh>::New
359  (
360  IOobject
361  (
362  this->polyMesh::name(), // region of undecomposed mesh
363  facesInstance(),
364  processorDb,
367  ),
368  std::move(procPoints),
369  std::move(procFaces),
370  std::move(procCells)
371  );
372  }
373  polyMesh& procMesh = procMeshPtr();
374 
375 
376  // Create processor boundary patches
377  const labelList& curPatchSizes = procPatchSize_[proci];
378 
379  const labelList& curPatchStarts = procPatchStartIndex_[proci];
380 
381  const labelList& curNeighbourProcessors =
382  procNeighbourProcessors_[proci];
383 
384  const labelList& curProcessorPatchSizes =
385  procProcessorPatchSize_[proci];
386 
387  const labelList& curProcessorPatchStarts =
388  procProcessorPatchStartIndex_[proci];
389 
390  const labelListList& curSubPatchIDs =
391  procProcessorPatchSubPatchIDs_[proci];
392 
393  const labelListList& curSubStarts =
394  procProcessorPatchSubPatchStarts_[proci];
395 
396  const polyPatchList& meshPatches = boundaryMesh();
397 
398 
399  // Count the number of inter-proc patches
400  label nInterProcPatches = 0;
401  forAll(curSubPatchIDs, procPatchi)
402  {
403  nInterProcPatches += curSubPatchIDs[procPatchi].size();
404  }
405 
406  polyPatchList procPatches
407  (
408  curPatchSizes.size() + nInterProcPatches
409  );
410 
411  label nPatches = 0;
412 
413  forAll(curPatchSizes, patchi)
414  {
415  // Get the face labels consistent with the field mapping
416  // (reuse the patch field mappers)
417  const polyPatch& meshPatch = meshPatches[patchi];
418 
419  fvFieldDecomposer::patchFieldDecomposer patchMapper
420  (
421  SubList<label>
422  (
423  curFaceLabels,
424  curPatchSizes[patchi],
425  curPatchStarts[patchi]
426  ),
427  meshPatch.start()
428  );
429 
430  // Map existing patches
431  procPatches.set
432  (
433  nPatches,
434  meshPatch.clone
435  (
436  procMesh.boundaryMesh(),
437  nPatches,
438  patchMapper.directAddressing(),
439  curPatchStarts[patchi]
440  )
441  );
442 
443  nPatches++;
444  }
445 
446  forAll(curProcessorPatchSizes, procPatchi)
447  {
448  const labelList& subPatchID = curSubPatchIDs[procPatchi];
449  const labelList& subStarts = curSubStarts[procPatchi];
450 
451  label curStart = curProcessorPatchStarts[procPatchi];
452 
453  forAll(subPatchID, i)
454  {
455  label size =
456  (
457  i < subPatchID.size()-1
458  ? subStarts[i+1] - subStarts[i]
459  : curProcessorPatchSizes[procPatchi] - subStarts[i]
460  );
461 
462  if (subPatchID[i] == -1)
463  {
464  // From internal faces
465  procPatches.set
466  (
467  nPatches,
468  new processorPolyPatch
469  (
470  size,
471  curStart,
472  nPatches,
473  procMesh.boundaryMesh(),
474  proci,
475  curNeighbourProcessors[procPatchi]
476  )
477  );
478  }
479  else
480  {
481  const coupledPolyPatch& pcPatch
482  = refCast<const coupledPolyPatch>
483  (
484  boundaryMesh()[subPatchID[i]]
485  );
486 
487  procPatches.set
488  (
489  nPatches,
490  new processorCyclicPolyPatch
491  (
492  size,
493  curStart,
494  nPatches,
495  procMesh.boundaryMesh(),
496  proci,
497  curNeighbourProcessors[procPatchi],
498  pcPatch.name(),
499  pcPatch.transform()
500  )
501  );
502  }
503 
504  curStart += size;
505  ++nPatches;
506  }
507  }
508 
509  // Add boundary patches
510  procMesh.addPatches(procPatches);
511 
512  // Create and add zones
513 
514  // Point zones
515  {
516  const pointZoneMesh& pz = pointZones();
517 
518  // Go through all the zoned points and find out if they
519  // belong to a zone. If so, add it to the zone as
520  // necessary
521  List<DynamicList<label>> zonePoints(pz.size());
522 
523  // Estimate size
524  forAll(zonePoints, zonei)
525  {
526  zonePoints[zonei].setCapacity(pz[zonei].size()/nProcs_);
527  }
528 
529  // Use the pointToZone map to find out the single zone (if any),
530  // use slow search only for shared points.
531  forAll(curPointLabels, pointi)
532  {
533  label curPoint = curPointLabels[pointi];
534 
535  label zonei = pointToZone[curPoint];
536 
537  if (zonei >= 0)
538  {
539  // Single zone.
540  zonePoints[zonei].append(pointi);
541  }
542  else if (zonei == -2)
543  {
544  // Multiple zones. Lookup.
545  forAll(pz, zonei)
546  {
547  label index = pz[zonei].whichPoint(curPoint);
548 
549  if (index != -1)
550  {
551  zonePoints[zonei].append(pointi);
552  }
553  }
554  }
555  }
556 
557  procMesh.pointZones().clearAddressing();
558  procMesh.pointZones().setSize(zonePoints.size());
559  forAll(zonePoints, zonei)
560  {
561  procMesh.pointZones().set
562  (
563  zonei,
564  pz[zonei].clone
565  (
566  procMesh.pointZones(),
567  zonei,
568  zonePoints[zonei].shrink()
569  )
570  );
571  }
572 
573  if (pz.size())
574  {
575  // Force writing on all processors
576  procMesh.pointZones().writeOpt(IOobject::AUTO_WRITE);
577  }
578  }
579 
580  // Face zones
581  {
582  const faceZoneMesh& fz = faceZones();
583 
584  // Go through all the zoned face and find out if they
585  // belong to a zone. If so, add it to the zone as
586  // necessary
587  List<DynamicList<label>> zoneFaces(fz.size());
588  List<DynamicList<bool>> zoneFaceFlips(fz.size());
589 
590  // Estimate size
591  forAll(zoneFaces, zonei)
592  {
593  label procSize = fz[zonei].size() / nProcs_;
594 
595  zoneFaces[zonei].setCapacity(procSize);
596  zoneFaceFlips[zonei].setCapacity(procSize);
597  }
598 
599  // Go through all the zoned faces and find out if they
600  // belong to a zone. If so, add it to the zone as
601  // necessary
602  forAll(curFaceLabels, facei)
603  {
604  // Remember to decrement the index by one (turning index)
605  //
606  label curF = mag(curFaceLabels[facei]) - 1;
607 
608  label zonei = faceToZone[curF];
609 
610  if (zonei >= 0)
611  {
612  // Single zone. Add the face
613  zoneFaces[zonei].append(facei);
614 
615  label index = fz[zonei].whichFace(curF);
616 
617  bool flip = fz[zonei].flipMap()[index];
618 
619  if (curFaceLabels[facei] < 0)
620  {
621  flip = !flip;
622  }
623 
624  zoneFaceFlips[zonei].append(flip);
625  }
626  else if (zonei == -2)
627  {
628  // Multiple zones. Lookup.
629  forAll(fz, zonei)
630  {
631  label index = fz[zonei].whichFace(curF);
632 
633  if (index != -1)
634  {
635  zoneFaces[zonei].append(facei);
636 
637  bool flip = fz[zonei].flipMap()[index];
638 
639  if (curFaceLabels[facei] < 0)
640  {
641  flip = !flip;
642  }
643 
644  zoneFaceFlips[zonei].append(flip);
645  }
646  }
647  }
648  }
649 
650  procMesh.faceZones().clearAddressing();
651  procMesh.faceZones().setSize(zoneFaces.size());
652  forAll(zoneFaces, zonei)
653  {
654  procMesh.faceZones().set
655  (
656  zonei,
657  fz[zonei].clone
658  (
659  zoneFaces[zonei].shrink(), // addressing
660  zoneFaceFlips[zonei].shrink(), // flipmap
661  zonei,
662  procMesh.faceZones()
663  )
664  );
665  }
666 
667  if (fz.size())
668  {
669  // Force writing on all processors
670  procMesh.faceZones().writeOpt(IOobject::AUTO_WRITE);
671  }
672  }
673 
674  // Cell zones
675  {
676  const cellZoneMesh& cz = cellZones();
677 
678  // Go through all the zoned cells and find out if they
679  // belong to a zone. If so, add it to the zone as
680  // necessary
681  List<DynamicList<label>> zoneCells(cz.size());
682 
683  // Estimate size
684  forAll(zoneCells, zonei)
685  {
686  zoneCells[zonei].setCapacity(cz[zonei].size()/nProcs_);
687  }
688 
689  forAll(curCellLabels, celli)
690  {
691  label curCelli = curCellLabels[celli];
692 
693  label zonei = cellToZone[curCelli];
694 
695  if (zonei >= 0)
696  {
697  // Single zone.
698  zoneCells[zonei].append(celli);
699  }
700  else if (zonei == -2)
701  {
702  // Multiple zones. Lookup.
703  forAll(cz, zonei)
704  {
705  label index = cz[zonei].whichCell(curCelli);
706 
707  if (index != -1)
708  {
709  zoneCells[zonei].append(celli);
710  }
711  }
712  }
713  }
714 
715  procMesh.cellZones().clearAddressing();
716  procMesh.cellZones().setSize(zoneCells.size());
717  forAll(zoneCells, zonei)
718  {
719  procMesh.cellZones().set
720  (
721  zonei,
722  cz[zonei].clone
723  (
724  zoneCells[zonei].shrink(),
725  zonei,
726  procMesh.cellZones()
727  )
728  );
729  }
730 
731  if (cz.size())
732  {
733  // Force writing on all processors
734  procMesh.cellZones().writeOpt(IOobject::AUTO_WRITE);
735  }
736  }
737 
738  // More precision (for points data)
740 
741  procMesh.write();
742 
743  // Write points if pointsInstance differing from facesInstance
744  if (facesInstancePointsPtr_)
745  {
746  pointIOField pointsInstancePoints
747  (
748  IOobject
749  (
750  "points",
751  pointsInstance(),
753  procMesh,
757  ),
758  std::move(procPoints)
759  );
760  pointsInstancePoints.write();
761  }
762 
763 
764  // Decompose any sets
765  if (decomposeSets)
766  {
767  forAll(cellSets, i)
768  {
769  const cellSet& cs = cellSets[i];
770  cellSet set(procMesh, cs.name(), cs.size()/nProcs_);
771  forAll(curCellLabels, i)
772  {
773  if (cs.found(curCellLabels[i]))
774  {
775  set.insert(i);
776  }
777  }
778  set.write();
779  }
780  forAll(faceSets, i)
781  {
782  const faceSet& cs = faceSets[i];
783  faceSet set(procMesh, cs.name(), cs.size()/nProcs_);
784  forAll(curFaceLabels, i)
785  {
786  if (cs.found(mag(curFaceLabels[i])-1))
787  {
788  set.insert(i);
789  }
790  }
791  set.write();
792  }
793  forAll(pointSets, i)
794  {
795  const pointSet& cs = pointSets[i];
796  pointSet set(procMesh, cs.name(), cs.size()/nProcs_);
797  forAll(curPointLabels, i)
798  {
799  if (cs.found(curPointLabels[i]))
800  {
801  set.insert(i);
802  }
803  }
804  set.write();
805  }
806  }
807 
808 
809  // Optional hexRef8 data
810  hexRef8Data
811  (
812  IOobject
813  (
814  "dummy",
815  facesInstance(),
817  procMesh,
821  ),
822  baseMeshData,
823  procCellAddressing_[proci],
824  procPointAddressing_[proci]
825  ).write();
826 
827 
828  // Statistics
829  Info<< nl << "Processor " << proci;
830 
831  if (procMesh.nCells())
832  {
833  Info<< nl << " ";
834  }
835  else
836  {
837  Info<< ": ";
838  }
839 
840  Info<< "Number of cells = " << procMesh.nCells() << nl;
841 
842  if (procMesh.nCells())
843  {
844  Info<< " Number of points = " << procMesh.nPoints() << nl;
845  }
846 
847  maxProcCells = max(maxProcCells, procMesh.nCells());
848 
849  label nBoundaryFaces = 0;
850  label nProcPatches = 0;
851  label nProcFaces = 0;
852 
853  for (const polyPatch& pp : procMesh.boundaryMesh())
854  {
855  const auto* cpp = isA<processorPolyPatch>(pp);
856 
857  if (cpp)
858  {
859  const auto& procPatch = *cpp;
860 
861  Info<< " Number of faces shared with processor "
862  << procPatch.neighbProcNo() << " = "
863  << procPatch.size() << nl;
864 
865  nProcFaces += procPatch.size();
866  ++nProcPatches;
867  }
868  else
869  {
870  nBoundaryFaces += pp.size();
871  }
872  }
873 
874  if (procMesh.nCells() && (nBoundaryFaces || nProcFaces))
875  {
876  Info<< " Number of processor patches = " << nProcPatches << nl
877  << " Number of processor faces = " << nProcFaces << nl
878  << " Number of boundary faces = " << nBoundaryFaces << nl;
879  }
880 
881  totProcFaces += nProcFaces;
882  totProcPatches += nProcPatches;
883  maxProcFaces = max(maxProcFaces, nProcFaces);
884  maxProcPatches = max(maxProcPatches, nProcPatches);
885 
886  // Write the addressing information
887 
888  IOobject ioAddr
889  (
890  "procAddressing",
891  procMesh.facesInstance(),
893  procMesh.thisDb(),
897  );
898 
899  // pointProcAddressing
900  ioAddr.rename("pointProcAddressing");
901  IOListRef<label>(ioAddr, procPointAddressing_[proci]).write();
902 
903  // faceProcAddressing
904  ioAddr.rename("faceProcAddressing");
905  IOListRef<label>(ioAddr, procFaceAddressing_[proci]).write();
906 
907  // cellProcAddressing
908  ioAddr.rename("cellProcAddressing");
909  IOListRef<label>(ioAddr, procCellAddressing_[proci]).write();
910 
911  // Write patch map for backwards compatibility.
912  // (= identity map for original patches, -1 for processor patches)
913  label nMeshPatches = curPatchSizes.size();
914  labelList procBoundaryAddr(identity(nMeshPatches));
915  procBoundaryAddr.resize(nMeshPatches+nProcPatches, -1);
916 
917  // boundaryProcAddressing
918  ioAddr.rename("boundaryProcAddressing");
919  IOListRef<label>(ioAddr, procBoundaryAddr).write();
920  }
921 
922 
923  // Summary stats
924  Info<< nl
925  << "Number of processor faces = " << (totProcFaces/2) << nl
926  << "Max number of cells = " << maxProcCells;
927 
928  if (maxProcCells != nCells())
929  {
930  scalar avgValue = scalar(nCells())/nProcs_;
931 
932  Info<< " (" << 100.0*(maxProcCells-avgValue)/avgValue
933  << "% above average " << avgValue << ')';
934  }
935  Info<< nl;
936 
937  Info<< "Max number of processor patches = " << maxProcPatches;
938  if (totProcPatches)
939  {
940  scalar avgValue = scalar(totProcPatches)/nProcs_;
941 
942  Info<< " (" << 100.0*(maxProcPatches-avgValue)/avgValue
943  << "% above average " << avgValue << ')';
944  }
945  Info<< nl;
946 
947  Info<< "Max number of faces between processors = " << maxProcFaces;
948  if (totProcFaces)
949  {
950  scalar avgValue = scalar(totProcFaces)/nProcs_;
951 
952  Info<< " (" << 100.0*(maxProcFaces-avgValue)/avgValue
953  << "% above average " << avgValue << ')';
954  }
955  Info<< nl << endl;
956 
957  return true;
958 }
959 
960 
961 // ************************************************************************* //
label nPatches
Definition: readKivaGrid.H:396
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
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:531
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...
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
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
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
void setSize(const label n)
Alias for resize()
Definition: List.H:320
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
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:62
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())