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-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "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  {
190  IOobjectList sets(objects.lookupClass<cellSet>());
191  forAllConstIters(sets, iter)
192  {
193  cellSets.append(new cellSet(*(iter.val())));
194  }
195  }
196  {
197  IOobjectList sets(objects.lookupClass<faceSet>());
198  forAllConstIters(sets, iter)
199  {
200  faceSets.append(new faceSet(*(iter.val())));
201  }
202  }
203  {
204  IOobjectList sets(objects.lookupClass<pointSet>());
205  forAllConstIters(sets, iter)
206  {
207  pointSets.append(new pointSet(*(iter.val())));
208  }
209  }
210  }
211 
212 
213  // Load refinement data (if any)
214  hexRef8Data baseMeshData
215  (
216  IOobject
217  (
218  "dummy",
219  facesInstance(),
221  *this,
225  )
226  );
227 
228 
229  label maxProcCells = 0;
230  label maxProcFaces = 0;
231  label totProcFaces = 0;
232  label maxProcPatches = 0;
233  label totProcPatches = 0;
234 
235  // Write out the meshes
236  for (label proci = 0; proci < nProcs_; proci++)
237  {
238  // Create processor points
239  const labelList& curPointLabels = procPointAddressing_[proci];
240 
241  const pointField& meshPoints = points();
242 
243  labelList pointLookup(nPoints(), -1);
244 
245  pointField procPoints(curPointLabels.size());
246 
247  forAll(curPointLabels, pointi)
248  {
249  procPoints[pointi] = meshPoints[curPointLabels[pointi]];
250 
251  pointLookup[curPointLabels[pointi]] = pointi;
252  }
253 
254  // Create processor faces
255  const labelList& curFaceLabels = procFaceAddressing_[proci];
256 
257  const faceList& meshFaces = faces();
258 
259  labelList faceLookup(nFaces(), -1);
260 
261  faceList procFaces(curFaceLabels.size());
262 
263  forAll(curFaceLabels, facei)
264  {
265  // Mark the original face as used
266  // Remember to decrement the index by one (turning index)
267  label curF = mag(curFaceLabels[facei]) - 1;
268 
269  faceLookup[curF] = facei;
270 
271  // get the original face
272  labelList origFaceLabels;
273 
274  if (curFaceLabels[facei] >= 0)
275  {
276  // face not turned
277  origFaceLabels = meshFaces[curF];
278  }
279  else
280  {
281  origFaceLabels = meshFaces[curF].reverseFace();
282  }
283 
284  // translate face labels into local point list
285  face& procFaceLabels = procFaces[facei];
286 
287  procFaceLabels.setSize(origFaceLabels.size());
288 
289  forAll(origFaceLabels, pointi)
290  {
291  procFaceLabels[pointi] = pointLookup[origFaceLabels[pointi]];
292  }
293  }
294 
295  // Create processor cells
296  const labelList& curCellLabels = procCellAddressing_[proci];
297 
298  const cellList& meshCells = cells();
299 
300  cellList procCells(curCellLabels.size());
301 
302  forAll(curCellLabels, celli)
303  {
304  const labelList& origCellLabels = meshCells[curCellLabels[celli]];
305 
306  cell& curCell = procCells[celli];
307 
308  curCell.setSize(origCellLabels.size());
309 
310  forAll(origCellLabels, cellFacei)
311  {
312  curCell[cellFacei] = faceLookup[origCellLabels[cellFacei]];
313  }
314  }
315 
316  // Create processor mesh without a boundary
317 
318  // create a database
319  Time processorDb
320  (
322  time().rootPath(),
323  time().caseName()/("processor" + Foam::name(proci)),
324  false, // No function objects
325  false // No extra controlDict libs
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_)
342  {
343  // Construct mesh from facesInstance.
344  pointField facesInstancePoints
345  (
346  facesInstancePointsPtr_(),
347  curPointLabels
348  );
349 
350  procMeshPtr = autoPtr<polyMesh>::New
351  (
352  IOobject
353  (
354  this->polyMesh::name(), // region of undecomposed mesh
355  facesInstance(),
356  processorDb,
359  ),
360  std::move(facesInstancePoints),
361  std::move(procFaces),
362  std::move(procCells)
363  );
364  }
365  else
366  {
367  procMeshPtr = autoPtr<polyMesh>::New
368  (
369  IOobject
370  (
371  this->polyMesh::name(), // region of undecomposed mesh
372  facesInstance(),
373  processorDb,
376  ),
377  std::move(procPoints),
378  std::move(procFaces),
379  std::move(procCells)
380  );
381  }
382  polyMesh& procMesh = procMeshPtr();
383 
384 
385  // Create processor boundary patches
386  const labelList& curPatchSizes = procPatchSize_[proci];
387 
388  const labelList& curPatchStarts = procPatchStartIndex_[proci];
389 
390  const labelList& curNeighbourProcessors =
391  procNeighbourProcessors_[proci];
392 
393  const labelList& curProcessorPatchSizes =
394  procProcessorPatchSize_[proci];
395 
396  const labelList& curProcessorPatchStarts =
397  procProcessorPatchStartIndex_[proci];
398 
399  const labelListList& curSubPatchIDs =
400  procProcessorPatchSubPatchIDs_[proci];
401 
402  const labelListList& curSubStarts =
403  procProcessorPatchSubPatchStarts_[proci];
404 
405  const polyPatchList& meshPatches = boundaryMesh();
406 
407 
408  // Count the number of inter-proc patches
409  label nInterProcPatches = 0;
410  forAll(curSubPatchIDs, procPatchi)
411  {
412  nInterProcPatches += curSubPatchIDs[procPatchi].size();
413  }
414 
415  polyPatchList procPatches
416  (
417  curPatchSizes.size() + nInterProcPatches
418  );
419 
420  label nPatches = 0;
421 
422  forAll(curPatchSizes, patchi)
423  {
424  // Get the face labels consistent with the field mapping
425  // (reuse the patch field mappers)
426  const polyPatch& meshPatch = meshPatches[patchi];
427 
428  fvFieldDecomposer::patchFieldDecomposer patchMapper
429  (
430  SubList<label>
431  (
432  curFaceLabels,
433  curPatchSizes[patchi],
434  curPatchStarts[patchi]
435  ),
436  meshPatch.start()
437  );
438 
439  // Map existing patches
440  procPatches.set
441  (
442  nPatches,
443  meshPatch.clone
444  (
445  procMesh.boundaryMesh(),
446  nPatches,
447  patchMapper.directAddressing(),
448  curPatchStarts[patchi]
449  )
450  );
451 
452  nPatches++;
453  }
454 
455  forAll(curProcessorPatchSizes, procPatchi)
456  {
457  const labelList& subPatchID = curSubPatchIDs[procPatchi];
458  const labelList& subStarts = curSubStarts[procPatchi];
459 
460  label curStart = curProcessorPatchStarts[procPatchi];
461 
462  forAll(subPatchID, i)
463  {
464  label size =
465  (
466  i < subPatchID.size()-1
467  ? subStarts[i+1] - subStarts[i]
468  : curProcessorPatchSizes[procPatchi] - subStarts[i]
469  );
470 
471  if (subPatchID[i] == -1)
472  {
473  // From internal faces
474  procPatches.set
475  (
476  nPatches,
477  new processorPolyPatch
478  (
479  size,
480  curStart,
481  nPatches,
482  procMesh.boundaryMesh(),
483  proci,
484  curNeighbourProcessors[procPatchi]
485  )
486  );
487  }
488  else
489  {
490  const coupledPolyPatch& pcPatch
491  = refCast<const coupledPolyPatch>
492  (
493  boundaryMesh()[subPatchID[i]]
494  );
495 
496  procPatches.set
497  (
498  nPatches,
499  new processorCyclicPolyPatch
500  (
501  size,
502  curStart,
503  nPatches,
504  procMesh.boundaryMesh(),
505  proci,
506  curNeighbourProcessors[procPatchi],
507  pcPatch.name(),
508  pcPatch.transform()
509  )
510  );
511  }
512 
513  curStart += size;
514  ++nPatches;
515  }
516  }
517 
518  // Add boundary patches
519  procMesh.addPatches(procPatches);
520 
521  // Create and add zones
522 
523  // Point zones
524  {
525  const pointZoneMesh& pz = pointZones();
526 
527  // Go through all the zoned points and find out if they
528  // belong to a zone. If so, add it to the zone as
529  // necessary
530  List<DynamicList<label>> zonePoints(pz.size());
531 
532  // Estimate size
533  forAll(zonePoints, zonei)
534  {
535  zonePoints[zonei].setCapacity(pz[zonei].size()/nProcs_);
536  }
537 
538  // Use the pointToZone map to find out the single zone (if any),
539  // use slow search only for shared points.
540  forAll(curPointLabels, pointi)
541  {
542  label curPoint = curPointLabels[pointi];
543 
544  label zonei = pointToZone[curPoint];
545 
546  if (zonei >= 0)
547  {
548  // Single zone.
549  zonePoints[zonei].append(pointi);
550  }
551  else if (zonei == -2)
552  {
553  // Multiple zones. Lookup.
554  forAll(pz, zonei)
555  {
556  label index = pz[zonei].whichPoint(curPoint);
557 
558  if (index != -1)
559  {
560  zonePoints[zonei].append(pointi);
561  }
562  }
563  }
564  }
565 
566  procMesh.pointZones().clearAddressing();
567  procMesh.pointZones().setSize(zonePoints.size());
568  forAll(zonePoints, zonei)
569  {
570  procMesh.pointZones().set
571  (
572  zonei,
573  pz[zonei].clone
574  (
575  procMesh.pointZones(),
576  zonei,
577  zonePoints[zonei].shrink()
578  )
579  );
580  }
581 
582  if (pz.size())
583  {
584  // Force writing on all processors
585  procMesh.pointZones().writeOpt(IOobject::AUTO_WRITE);
586  }
587  }
588 
589  // Face zones
590  {
591  const faceZoneMesh& fz = faceZones();
592 
593  // Go through all the zoned face and find out if they
594  // belong to a zone. If so, add it to the zone as
595  // necessary
596  List<DynamicList<label>> zoneFaces(fz.size());
597  List<DynamicList<bool>> zoneFaceFlips(fz.size());
598 
599  // Estimate size
600  forAll(zoneFaces, zonei)
601  {
602  label procSize = fz[zonei].size() / nProcs_;
603 
604  zoneFaces[zonei].setCapacity(procSize);
605  zoneFaceFlips[zonei].setCapacity(procSize);
606  }
607 
608  // Go through all the zoned faces and find out if they
609  // belong to a zone. If so, add it to the zone as
610  // necessary
611  forAll(curFaceLabels, facei)
612  {
613  // Remember to decrement the index by one (turning index)
614  //
615  label curF = mag(curFaceLabels[facei]) - 1;
616 
617  label zonei = faceToZone[curF];
618 
619  if (zonei >= 0)
620  {
621  // Single zone. Add the face
622  zoneFaces[zonei].append(facei);
623 
624  label index = fz[zonei].whichFace(curF);
625 
626  bool flip = fz[zonei].flipMap()[index];
627 
628  if (curFaceLabels[facei] < 0)
629  {
630  flip = !flip;
631  }
632 
633  zoneFaceFlips[zonei].append(flip);
634  }
635  else if (zonei == -2)
636  {
637  // Multiple zones. Lookup.
638  forAll(fz, zonei)
639  {
640  label index = fz[zonei].whichFace(curF);
641 
642  if (index != -1)
643  {
644  zoneFaces[zonei].append(facei);
645 
646  bool flip = fz[zonei].flipMap()[index];
647 
648  if (curFaceLabels[facei] < 0)
649  {
650  flip = !flip;
651  }
652 
653  zoneFaceFlips[zonei].append(flip);
654  }
655  }
656  }
657  }
658 
659  procMesh.faceZones().clearAddressing();
660  procMesh.faceZones().setSize(zoneFaces.size());
661  forAll(zoneFaces, zonei)
662  {
663  procMesh.faceZones().set
664  (
665  zonei,
666  fz[zonei].clone
667  (
668  zoneFaces[zonei].shrink(), // addressing
669  zoneFaceFlips[zonei].shrink(), // flipmap
670  zonei,
671  procMesh.faceZones()
672  )
673  );
674  }
675 
676  if (fz.size())
677  {
678  // Force writing on all processors
679  procMesh.faceZones().writeOpt(IOobject::AUTO_WRITE);
680  }
681  }
682 
683  // Cell zones
684  {
685  const cellZoneMesh& cz = cellZones();
686 
687  // Go through all the zoned cells and find out if they
688  // belong to a zone. If so, add it to the zone as
689  // necessary
690  List<DynamicList<label>> zoneCells(cz.size());
691 
692  // Estimate size
693  forAll(zoneCells, zonei)
694  {
695  zoneCells[zonei].setCapacity(cz[zonei].size()/nProcs_);
696  }
697 
698  forAll(curCellLabels, celli)
699  {
700  label curCelli = curCellLabels[celli];
701 
702  label zonei = cellToZone[curCelli];
703 
704  if (zonei >= 0)
705  {
706  // Single zone.
707  zoneCells[zonei].append(celli);
708  }
709  else if (zonei == -2)
710  {
711  // Multiple zones. Lookup.
712  forAll(cz, zonei)
713  {
714  label index = cz[zonei].whichCell(curCelli);
715 
716  if (index != -1)
717  {
718  zoneCells[zonei].append(celli);
719  }
720  }
721  }
722  }
723 
724  procMesh.cellZones().clearAddressing();
725  procMesh.cellZones().setSize(zoneCells.size());
726  forAll(zoneCells, zonei)
727  {
728  procMesh.cellZones().set
729  (
730  zonei,
731  cz[zonei].clone
732  (
733  zoneCells[zonei].shrink(),
734  zonei,
735  procMesh.cellZones()
736  )
737  );
738  }
739 
740  if (cz.size())
741  {
742  // Force writing on all processors
743  procMesh.cellZones().writeOpt(IOobject::AUTO_WRITE);
744  }
745  }
746 
747  // Set the precision of the points data to be min 10
749 
750  procMesh.write();
751 
752  // Write points if pointsInstance differing from facesInstance
753  if (facesInstancePointsPtr_)
754  {
755  pointIOField pointsInstancePoints
756  (
757  IOobject
758  (
759  "points",
760  pointsInstance(),
762  procMesh,
766  ),
767  std::move(procPoints)
768  );
769  pointsInstancePoints.write();
770  }
771 
772 
773  // Decompose any sets
774  if (decomposeSets)
775  {
776  forAll(cellSets, i)
777  {
778  const cellSet& cs = cellSets[i];
779  cellSet set(procMesh, cs.name(), cs.size()/nProcs_);
780  forAll(curCellLabels, i)
781  {
782  if (cs.found(curCellLabels[i]))
783  {
784  set.insert(i);
785  }
786  }
787  set.write();
788  }
789  forAll(faceSets, i)
790  {
791  const faceSet& cs = faceSets[i];
792  faceSet set(procMesh, cs.name(), cs.size()/nProcs_);
793  forAll(curFaceLabels, i)
794  {
795  if (cs.found(mag(curFaceLabels[i])-1))
796  {
797  set.insert(i);
798  }
799  }
800  set.write();
801  }
802  forAll(pointSets, i)
803  {
804  const pointSet& cs = pointSets[i];
805  pointSet set(procMesh, cs.name(), cs.size()/nProcs_);
806  forAll(curPointLabels, i)
807  {
808  if (cs.found(curPointLabels[i]))
809  {
810  set.insert(i);
811  }
812  }
813  set.write();
814  }
815  }
816 
817 
818  // Optional hexRef8 data
819  hexRef8Data
820  (
821  IOobject
822  (
823  "dummy",
824  facesInstance(),
826  procMesh,
830  ),
831  baseMeshData,
832  procCellAddressing_[proci],
833  procPointAddressing_[proci]
834  ).write();
835 
836 
837  // Statistics
838  Info<< nl << "Processor " << proci;
839 
840  if (procMesh.nCells())
841  {
842  Info<< nl << " ";
843  }
844  else
845  {
846  Info<< ": ";
847  }
848 
849  Info<< "Number of cells = " << procMesh.nCells() << nl;
850 
851  if (procMesh.nCells())
852  {
853  Info<< " Number of points = " << procMesh.nPoints() << nl;
854  }
855 
856  maxProcCells = max(maxProcCells, procMesh.nCells());
857 
858  label nBoundaryFaces = 0;
859  label nProcPatches = 0;
860  label nProcFaces = 0;
861 
862  for (const polyPatch& pp : procMesh.boundaryMesh())
863  {
864  const auto* cpp = isA<processorPolyPatch>(pp);
865 
866  if (cpp)
867  {
868  const auto& procPatch = *cpp;
869 
870  Info<< " Number of faces shared with processor "
871  << procPatch.neighbProcNo() << " = "
872  << procPatch.size() << nl;
873 
874  nProcFaces += procPatch.size();
875  ++nProcPatches;
876  }
877  else
878  {
879  nBoundaryFaces += pp.size();
880  }
881  }
882 
883  if (procMesh.nCells() && (nBoundaryFaces || nProcFaces))
884  {
885  Info<< " Number of processor patches = " << nProcPatches << nl
886  << " Number of processor faces = " << nProcFaces << nl
887  << " Number of boundary faces = " << nBoundaryFaces << nl;
888  }
889 
890  totProcFaces += nProcFaces;
891  totProcPatches += nProcPatches;
892  maxProcFaces = max(maxProcFaces, nProcFaces);
893  maxProcPatches = max(maxProcPatches, nProcPatches);
894 
895  // Write the addressing information
896 
897  IOobject ioAddr
898  (
899  "procAddressing",
900  procMesh.facesInstance(),
902  procMesh.thisDb(),
906  );
907 
908  // pointProcAddressing
909  ioAddr.rename("pointProcAddressing");
910  IOListRef<label>(ioAddr, procPointAddressing_[proci]).write();
911 
912  // faceProcAddressing
913  ioAddr.rename("faceProcAddressing");
914  IOListRef<label>(ioAddr, procFaceAddressing_[proci]).write();
915 
916  // cellProcAddressing
917  ioAddr.rename("cellProcAddressing");
918  IOListRef<label>(ioAddr, procCellAddressing_[proci]).write();
919 
920  // Write patch map for backwards compatibility.
921  // (= identity map for original patches, -1 for processor patches)
922  label nMeshPatches = curPatchSizes.size();
923  labelList procBoundaryAddr(identity(nMeshPatches));
924  procBoundaryAddr.resize(nMeshPatches+nProcPatches, -1);
925 
926  // boundaryProcAddressing
927  ioAddr.rename("boundaryProcAddressing");
928  IOListRef<label>(ioAddr, procBoundaryAddr).write();
929  }
930 
931 
932  // Summary stats
933  Info<< nl
934  << "Number of processor faces = " << (totProcFaces/2) << nl
935  << "Max number of cells = " << maxProcCells;
936 
937  if (maxProcCells != nCells())
938  {
939  scalar avgValue = scalar(nCells())/nProcs_;
940 
941  Info<< " (" << 100.0*(maxProcCells-avgValue)/avgValue
942  << "% above average " << avgValue << ')';
943  }
944  Info<< nl;
945 
946  Info<< "Max number of processor patches = " << maxProcPatches;
947  if (totProcPatches)
948  {
949  scalar avgValue = scalar(totProcPatches)/nProcs_;
950 
951  Info<< " (" << 100.0*(maxProcPatches-avgValue)/avgValue
952  << "% above average " << avgValue << ')';
953  }
954  Info<< nl;
955 
956  Info<< "Max number of faces between processors = " << maxProcFaces;
957  if (totProcFaces)
958  {
959  scalar avgValue = scalar(totProcFaces)/nProcs_;
960 
961  Info<< " (" << 100.0*(maxProcFaces-avgValue)/avgValue
962  << "% above average " << avgValue << ')';
963  }
964  Info<< nl << endl;
965 
966  return true;
967 }
968 
969 
970 // ************************************************************************* //
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:162
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:49
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:402
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
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.
static unsigned int defaultPrecision() noexcept
Return the default precision.
Definition: IOstream.H:423
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:414
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:289
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
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)
Definition: labelList.C:31
label nPoints
Reading is optional [identical to LAZY_READ].
ZoneMesh< pointZone, polyMesh > pointZoneMesh
A ZoneMesh with the type pointZone.
MeshObject wrapper of decompositionMethod.
static word controlDictName
The default control dictionary name (normally "controlDict")
Definition: Time.H:280
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())
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28