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