renumberMesh.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-2016 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 Application
25  renumberMesh
26 
27 Description
28  Renumbers the cell list in order to reduce the bandwidth, reading and
29  renumbering all fields from all the time directories.
30 
31  By default uses bandCompression (CuthillMcKee) but will
32  read system/renumberMeshDict if -dict option is present
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #include "argList.H"
37 #include "IOobjectList.H"
38 #include "fvMesh.H"
39 #include "polyTopoChange.H"
40 #include "ReadFields.H"
41 #include "volFields.H"
42 #include "surfaceFields.H"
43 #include "SortableList.H"
44 #include "decompositionMethod.H"
45 #include "renumberMethod.H"
47 #include "CuthillMcKeeRenumber.H"
48 #include "fvMeshSubset.H"
49 #include "cellSet.H"
50 #include "faceSet.H"
51 #include "pointSet.H"
52 
53 #ifdef FOAM_USE_ZOLTAN
54  #include "zoltanRenumber.H"
55 #endif
56 
57 
58 using namespace Foam;
59 
60 
61 // Create named field from labelList for postprocessing
62 tmp<volScalarField> createScalarField
63 (
64  const fvMesh& mesh,
65  const word& name,
66  const labelList& elems
67 )
68 {
70  (
71  new volScalarField
72  (
73  IOobject
74  (
75  name,
76  mesh.time().timeName(),
77  mesh,
80  false
81  ),
82  mesh,
83  dimensionedScalar("zero", dimless, 0),
84  zeroGradientFvPatchScalarField::typeName
85  )
86  );
87  volScalarField& fld = tfld.ref();
88 
89  forAll(fld, celli)
90  {
91  fld[celli] = elems[celli];
92  }
93 
94  return tfld;
95 }
96 
97 
98 // Calculate band of matrix
99 label getBand(const labelList& owner, const labelList& neighbour)
100 {
101  label band = 0;
102 
103  forAll(neighbour, facei)
104  {
105  label diff = neighbour[facei] - owner[facei];
106 
107  if (diff > band)
108  {
109  band = diff;
110  }
111  }
112  return band;
113 }
114 
115 
116 // Calculate band of matrix
117 void getBand
118 (
119  const bool calculateIntersect,
120  const label nCells,
121  const labelList& owner,
122  const labelList& neighbour,
123  label& bandwidth,
124  scalar& profile, // scalar to avoid overflow
125  scalar& sumSqrIntersect // scalar to avoid overflow
126 )
127 {
128  labelList cellBandwidth(nCells, 0);
129  scalarField nIntersect(nCells, 0.0);
130 
131  forAll(neighbour, facei)
132  {
133  label own = owner[facei];
134  label nei = neighbour[facei];
135 
136  // Note: mag not necessary for correct (upper-triangular) ordering.
137  label diff = nei-own;
138  cellBandwidth[nei] = max(cellBandwidth[nei], diff);
139  }
140 
141  bandwidth = max(cellBandwidth);
142 
143  // Do not use field algebra because of conversion label to scalar
144  profile = 0.0;
145  forAll(cellBandwidth, celli)
146  {
147  profile += 1.0*cellBandwidth[celli];
148  }
149 
150  sumSqrIntersect = 0.0;
151  if (calculateIntersect)
152  {
153  forAll(nIntersect, celli)
154  {
155  for (label colI = celli-cellBandwidth[celli]; colI <= celli; colI++)
156  {
157  nIntersect[colI] += 1.0;
158  }
159  }
160 
161  sumSqrIntersect = sum(Foam::sqr(nIntersect));
162  }
163 }
164 
165 
166 // Determine upper-triangular face order
167 labelList getFaceOrder
168 (
169  const primitiveMesh& mesh,
170  const labelList& cellOrder // New to old cell
171 )
172 {
173  labelList reverseCellOrder(invert(cellOrder.size(), cellOrder));
174 
175  labelList oldToNewFace(mesh.nFaces(), -1);
176 
177  label newFacei = 0;
178 
179  labelList nbr;
180  labelList order;
181 
182  forAll(cellOrder, newCelli)
183  {
184  label oldCelli = cellOrder[newCelli];
185 
186  const cell& cFaces = mesh.cells()[oldCelli];
187 
188  // Neighbouring cells
189  nbr.setSize(cFaces.size());
190 
191  forAll(cFaces, i)
192  {
193  label facei = cFaces[i];
194 
195  if (mesh.isInternalFace(facei))
196  {
197  // Internal face. Get cell on other side.
198  label nbrCelli = reverseCellOrder[mesh.faceNeighbour()[facei]];
199  if (nbrCelli == newCelli)
200  {
201  nbrCelli = reverseCellOrder[mesh.faceOwner()[facei]];
202  }
203 
204  if (newCelli < nbrCelli)
205  {
206  // Celli is master
207  nbr[i] = nbrCelli;
208  }
209  else
210  {
211  // nbrCell is master. Let it handle this face.
212  nbr[i] = -1;
213  }
214  }
215  else
216  {
217  // External face. Do later.
218  nbr[i] = -1;
219  }
220  }
221 
222  order.setSize(nbr.size());
223  sortedOrder(nbr, order);
224 
225  forAll(order, i)
226  {
227  label index = order[i];
228  if (nbr[index] != -1)
229  {
230  oldToNewFace[cFaces[index]] = newFacei++;
231  }
232  }
233  }
234 
235  // Leave patch faces intact.
236  for (label facei = newFacei; facei < mesh.nFaces(); facei++)
237  {
238  oldToNewFace[facei] = facei;
239  }
240 
241 
242  // Check done all faces.
243  forAll(oldToNewFace, facei)
244  {
245  if (oldToNewFace[facei] == -1)
246  {
248  << "Did not determine new position" << " for face " << facei
249  << abort(FatalError);
250  }
251  }
252 
253  return invert(mesh.nFaces(), oldToNewFace);
254 }
255 
256 
257 // Determine face order such that inside region faces are sorted
258 // upper-triangular but inbetween region faces are handled like boundary faces.
259 labelList getRegionFaceOrder
260 (
261  const primitiveMesh& mesh,
262  const labelList& cellOrder, // New to old cell
263  const labelList& cellToRegion // Old cell to region
264 )
265 {
266  labelList reverseCellOrder(invert(cellOrder.size(), cellOrder));
267 
268  labelList oldToNewFace(mesh.nFaces(), -1);
269 
270  label newFacei = 0;
271 
272  label prevRegion = -1;
273 
274  forAll(cellOrder, newCelli)
275  {
276  label oldCelli = cellOrder[newCelli];
277 
278  if (cellToRegion[oldCelli] != prevRegion)
279  {
280  prevRegion = cellToRegion[oldCelli];
281  }
282 
283  const cell& cFaces = mesh.cells()[oldCelli];
284 
285  SortableList<label> nbr(cFaces.size());
286 
287  forAll(cFaces, i)
288  {
289  label facei = cFaces[i];
290 
291  if (mesh.isInternalFace(facei))
292  {
293  // Internal face. Get cell on other side.
294  label nbrCelli = reverseCellOrder[mesh.faceNeighbour()[facei]];
295  if (nbrCelli == newCelli)
296  {
297  nbrCelli = reverseCellOrder[mesh.faceOwner()[facei]];
298  }
299 
300  if (cellToRegion[oldCelli] != cellToRegion[cellOrder[nbrCelli]])
301  {
302  // Treat like external face. Do later.
303  nbr[i] = -1;
304  }
305  else if (newCelli < nbrCelli)
306  {
307  // Celli is master
308  nbr[i] = nbrCelli;
309  }
310  else
311  {
312  // nbrCell is master. Let it handle this face.
313  nbr[i] = -1;
314  }
315  }
316  else
317  {
318  // External face. Do later.
319  nbr[i] = -1;
320  }
321  }
322 
323  nbr.sort();
324 
325  forAll(nbr, i)
326  {
327  if (nbr[i] != -1)
328  {
329  oldToNewFace[cFaces[nbr.indices()[i]]] = newFacei++;
330  }
331  }
332  }
333 
334  // Do region interfaces
335  label nRegions = max(cellToRegion)+1;
336  {
337  // Sort in increasing region
338  SortableList<label> sortKey(mesh.nFaces(), labelMax);
339 
340  for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
341  {
342  label ownRegion = cellToRegion[mesh.faceOwner()[facei]];
343  label neiRegion = cellToRegion[mesh.faceNeighbour()[facei]];
344 
345  if (ownRegion != neiRegion)
346  {
347  sortKey[facei] =
348  min(ownRegion, neiRegion)*nRegions
349  +max(ownRegion, neiRegion);
350  }
351  }
352  sortKey.sort();
353 
354  // Extract.
355  label prevKey = -1;
356  forAll(sortKey, i)
357  {
358  label key = sortKey[i];
359 
360  if (key == labelMax)
361  {
362  break;
363  }
364 
365  if (prevKey != key)
366  {
367  prevKey = key;
368  }
369 
370  oldToNewFace[sortKey.indices()[i]] = newFacei++;
371  }
372  }
373 
374  // Leave patch faces intact.
375  for (label facei = newFacei; facei < mesh.nFaces(); facei++)
376  {
377  oldToNewFace[facei] = facei;
378  }
379 
380 
381  // Check done all faces.
382  forAll(oldToNewFace, facei)
383  {
384  if (oldToNewFace[facei] == -1)
385  {
387  << "Did not determine new position"
388  << " for face " << facei
389  << abort(FatalError);
390  }
391  }
392 
393  return invert(mesh.nFaces(), oldToNewFace);
394 }
395 
396 
397 // cellOrder: old cell for every new cell
398 // faceOrder: old face for every new face. Ordering of boundary faces not
399 // changed.
400 autoPtr<mapPolyMesh> reorderMesh
401 (
402  polyMesh& mesh,
403  const labelList& cellOrder,
404  const labelList& faceOrder
405 )
406 {
407  labelList reverseCellOrder(invert(cellOrder.size(), cellOrder));
408  labelList reverseFaceOrder(invert(faceOrder.size(), faceOrder));
409 
410  faceList newFaces(reorder(reverseFaceOrder, mesh.faces()));
411  labelList newOwner
412  (
413  renumber
414  (
415  reverseCellOrder,
416  reorder(reverseFaceOrder, mesh.faceOwner())
417  )
418  );
419  labelList newNeighbour
420  (
421  renumber
422  (
423  reverseCellOrder,
424  reorder(reverseFaceOrder, mesh.faceNeighbour())
425  )
426  );
427 
428  // Check if any faces need swapping.
429  labelHashSet flipFaceFlux(newOwner.size());
430  forAll(newNeighbour, facei)
431  {
432  label own = newOwner[facei];
433  label nei = newNeighbour[facei];
434 
435  if (nei < own)
436  {
437  newFaces[facei].flip();
438  Swap(newOwner[facei], newNeighbour[facei]);
439  flipFaceFlux.insert(facei);
440  }
441  }
442 
443  const polyBoundaryMesh& patches = mesh.boundaryMesh();
444  labelList patchSizes(patches.size());
445  labelList patchStarts(patches.size());
446  labelList oldPatchNMeshPoints(patches.size());
447  labelListList patchPointMap(patches.size());
448 
449  forAll(patches, patchi)
450  {
451  patchSizes[patchi] = patches[patchi].size();
452  patchStarts[patchi] = patches[patchi].start();
453  oldPatchNMeshPoints[patchi] = patches[patchi].nPoints();
454  patchPointMap[patchi] = identity(patches[patchi].nPoints());
455  }
456 
457  mesh.resetPrimitives
458  (
460  xferMove(newFaces),
461  xferMove(newOwner),
462  xferMove(newNeighbour),
463  patchSizes,
464  patchStarts,
465  true
466  );
467 
468 
469  // Re-do the faceZones
470  {
471  faceZoneMesh& faceZones = mesh.faceZones();
472  faceZones.clearAddressing();
473  forAll(faceZones, zoneI)
474  {
475  faceZone& fZone = faceZones[zoneI];
476  labelList newAddressing(fZone.size());
477  boolList newFlipMap(fZone.size());
478  forAll(fZone, i)
479  {
480  label oldFacei = fZone[i];
481  newAddressing[i] = reverseFaceOrder[oldFacei];
482  if (flipFaceFlux.found(newAddressing[i]))
483  {
484  newFlipMap[i] = !fZone.flipMap()[i];
485  }
486  else
487  {
488  newFlipMap[i] = fZone.flipMap()[i];
489  }
490  }
491  labelList newToOld;
492  sortedOrder(newAddressing, newToOld);
493  fZone.resetAddressing
494  (
495  UIndirectList<label>(newAddressing, newToOld)(),
496  UIndirectList<bool>(newFlipMap, newToOld)()
497  );
498  }
499  }
500  // Re-do the cellZones
501  {
502  cellZoneMesh& cellZones = mesh.cellZones();
503  cellZones.clearAddressing();
504  forAll(cellZones, zoneI)
505  {
506  cellZones[zoneI] = UIndirectList<label>
507  (
508  reverseCellOrder,
509  cellZones[zoneI]
510  )();
511  Foam::sort(cellZones[zoneI]);
512  }
513  }
514 
515 
516  return autoPtr<mapPolyMesh>
517  (
518  new mapPolyMesh
519  (
520  mesh, // const polyMesh& mesh,
521  mesh.nPoints(), // nOldPoints,
522  mesh.nFaces(), // nOldFaces,
523  mesh.nCells(), // nOldCells,
524  identity(mesh.nPoints()), // pointMap,
525  List<objectMap>(0), // pointsFromPoints,
526  faceOrder, // faceMap,
527  List<objectMap>(0), // facesFromPoints,
528  List<objectMap>(0), // facesFromEdges,
529  List<objectMap>(0), // facesFromFaces,
530  cellOrder, // cellMap,
531  List<objectMap>(0), // cellsFromPoints,
532  List<objectMap>(0), // cellsFromEdges,
533  List<objectMap>(0), // cellsFromFaces,
534  List<objectMap>(0), // cellsFromCells,
535  identity(mesh.nPoints()), // reversePointMap,
536  reverseFaceOrder, // reverseFaceMap,
537  reverseCellOrder, // reverseCellMap,
538  flipFaceFlux, // flipFaceFlux,
539  patchPointMap, // patchPointMap,
540  labelListList(0), // pointZoneMap,
541  labelListList(0), // faceZonePointMap,
542  labelListList(0), // faceZoneFaceMap,
543  labelListList(0), // cellZoneMap,
544  pointField(0), // preMotionPoints,
545  patchStarts, // oldPatchStarts,
546  oldPatchNMeshPoints, // oldPatchNMeshPoints
547  autoPtr<scalarField>() // oldCellVolumes
548  )
549  );
550 }
551 
552 
553 // Return new to old cell numbering
554 labelList regionRenumber
555 (
556  const renumberMethod& method,
557  const fvMesh& mesh,
558  const labelList& cellToRegion
559 )
560 {
561  Info<< "Determining cell order:" << endl;
562 
563  labelList cellOrder(cellToRegion.size());
564 
565  label nRegions = max(cellToRegion)+1;
566 
567  labelListList regionToCells(invertOneToMany(nRegions, cellToRegion));
568 
569  label celli = 0;
570 
571  forAll(regionToCells, regionI)
572  {
573  Info<< " region " << regionI << " starts at " << celli << endl;
574 
575  // Make sure no parallel comms
576  bool oldParRun = UPstream::parRun();
577  UPstream::parRun() = false;
578 
579  // Per region do a reordering.
580  fvMeshSubset subsetter(mesh);
581  subsetter.setLargeCellSubset(cellToRegion, regionI);
582 
583  const fvMesh& subMesh = subsetter.subMesh();
584 
585  labelList subCellOrder = method.renumber
586  (
587  subMesh,
588  subMesh.cellCentres()
589  );
590 
591  // Restore state
592  UPstream::parRun() = oldParRun;
593 
594  const labelList& cellMap = subsetter.cellMap();
595 
596  forAll(subCellOrder, i)
597  {
598  cellOrder[celli++] = cellMap[subCellOrder[i]];
599  }
600  }
601  Info<< endl;
602 
603  return cellOrder;
604 }
605 
606 
607 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
608 
609 int main(int argc, char *argv[])
610 {
612  (
613  "Renumber mesh to minimise bandwidth"
614  );
615 
616  #include "addRegionOption.H"
617  #include "addOverwriteOption.H"
618  #include "addTimeOptions.H"
619  #include "addDictOption.H"
621  (
622  "frontWidth",
623  "calculate the rms of the frontwidth"
624  );
625 
626 
627  #include "setRootCase.H"
628  #include "createTime.H"
629  runTime.functionObjects().off();
630 
631 
632  // Force linker to include zoltan symbols. This section is only needed since
633  // Zoltan is a static library
634  #ifdef FOAM_USE_ZOLTAN
635  Info<< "renumberMesh built with zoltan support." << nl << endl;
636  (void)zoltanRenumber::typeName;
637  #endif
638 
639 
640  // Get times list
641  instantList Times = runTime.times();
642 
643  // Set startTime and endTime depending on -time and -latestTime options
644  #include "checkTimeOptions.H"
645 
646  runTime.setTime(Times[startTime], startTime);
647 
648  #include "createNamedMesh.H"
649  const word oldInstance = mesh.pointsInstance();
650 
651  const bool readDict = args.optionFound("dict");
652  const bool doFrontWidth = args.optionFound("frontWidth");
653  const bool overwrite = args.optionFound("overwrite");
654 
655  label band;
656  scalar profile;
657  scalar sumSqrIntersect;
658  getBand
659  (
660  doFrontWidth,
661  mesh.nCells(),
662  mesh.faceOwner(),
663  mesh.faceNeighbour(),
664  band,
665  profile,
666  sumSqrIntersect
667  );
668 
669  reduce(band, maxOp<label>());
670  reduce(profile, sumOp<scalar>());
671  scalar rmsFrontwidth = Foam::sqrt
672  (
674  (
675  sumSqrIntersect,
676  sumOp<scalar>()
677  )/mesh.globalData().nTotalCells()
678  );
679 
680  Info<< "Mesh size: " << mesh.globalData().nTotalCells() << nl
681  << "Before renumbering :" << nl
682  << " band : " << band << nl
683  << " profile : " << profile << nl;
684 
685  if (doFrontWidth)
686  {
687  Info<< " rms frontwidth : " << rmsFrontwidth << nl;
688  }
689 
690  Info<< endl;
691 
692  bool sortCoupledFaceCells = false;
693  bool writeMaps = false;
694  bool orderPoints = false;
695  label blockSize = 0;
696  bool renumberSets = true;
697 
698  // Construct renumberMethod
699  autoPtr<IOdictionary> renumberDictPtr;
700  autoPtr<renumberMethod> renumberPtr;
701 
702  if (readDict)
703  {
704  const word dictName("renumberMeshDict");
705  #include "setSystemMeshDictionaryIO.H"
706 
707  Info<< "Renumber according to " << dictName << nl << endl;
708 
709  renumberDictPtr.reset(new IOdictionary(dictIO));
710  const IOdictionary& renumberDict = renumberDictPtr();
711 
712  renumberPtr = renumberMethod::New(renumberDict);
713 
714  sortCoupledFaceCells = renumberDict.lookupOrDefault
715  (
716  "sortCoupledFaceCells",
717  false
718  );
719  if (sortCoupledFaceCells)
720  {
721  Info<< "Sorting cells on coupled boundaries to be last." << nl
722  << endl;
723  }
724 
725  blockSize = renumberDict.lookupOrDefault("blockSize", 0);
726  if (blockSize > 0)
727  {
728  Info<< "Ordering cells into regions of size " << blockSize
729  << " (using decomposition);"
730  << " ordering faces into region-internal and region-external."
731  << nl << endl;
732 
733  if (blockSize < 0 || blockSize >= mesh.nCells())
734  {
736  << "Block size " << blockSize
737  << " should be positive integer"
738  << " and less than the number of cells in the mesh."
739  << exit(FatalError);
740  }
741  }
742 
743  orderPoints = renumberDict.lookupOrDefault("orderPoints", false);
744  if (orderPoints)
745  {
746  Info<< "Ordering points into internal and boundary points." << nl
747  << endl;
748  }
749 
750  renumberDict.lookup("writeMaps") >> writeMaps;
751  if (writeMaps)
752  {
753  Info<< "Writing renumber maps (new to old) to polyMesh." << nl
754  << endl;
755  }
756 
757  renumberSets = renumberDict.lookupOrDefault("renumberSets", true);
758  }
759  else
760  {
761  Info<< "Using default renumberMethod." << nl << endl;
762  dictionary renumberDict;
763  renumberPtr.reset(new CuthillMcKeeRenumber(renumberDict));
764  }
765 
766  Info<< "Selecting renumberMethod " << renumberPtr().type() << nl << endl;
767 
768 
769 
770  // Read parallel reconstruct maps
771  labelIOList cellProcAddressing
772  (
773  IOobject
774  (
775  "cellProcAddressing",
776  mesh.facesInstance(),
778  mesh,
780  ),
781  labelList(0)
782  );
783 
785  (
786  IOobject
787  (
788  "faceProcAddressing",
789  mesh.facesInstance(),
791  mesh,
793  ),
794  labelList(0)
795  );
796  labelIOList pointProcAddressing
797  (
798  IOobject
799  (
800  "pointProcAddressing",
801  mesh.pointsInstance(),
803  mesh,
805  ),
806  labelList(0)
807  );
808  labelIOList boundaryProcAddressing
809  (
810  IOobject
811  (
812  "boundaryProcAddressing",
813  mesh.pointsInstance(),
815  mesh,
817  ),
818  labelList(0)
819  );
820 
821 
822  // Read objects in time directory
823  IOobjectList objects(mesh, runTime.timeName());
824 
825 
826  // Read vol fields.
827 
829  ReadFields(mesh, objects, vsFlds);
830 
832  ReadFields(mesh, objects, vvFlds);
833 
835  ReadFields(mesh, objects, vstFlds);
836 
837  PtrList<volSymmTensorField> vsymtFlds;
838  ReadFields(mesh, objects, vsymtFlds);
839 
841  ReadFields(mesh, objects, vtFlds);
842 
843 
844  // Read surface fields.
845 
847  ReadFields(mesh, objects, ssFlds);
848 
850  ReadFields(mesh, objects, svFlds);
851 
853  ReadFields(mesh, objects, sstFlds);
854 
856  ReadFields(mesh, objects, ssymtFlds);
857 
859  ReadFields(mesh, objects, stFlds);
860 
861 
862  // Read point fields.
863 
865  ReadFields(pointMesh::New(mesh), objects, psFlds);
866 
868  ReadFields(pointMesh::New(mesh), objects, pvFlds);
869 
871  ReadFields(pointMesh::New(mesh), objects, pstFlds);
872 
874  ReadFields(pointMesh::New(mesh), objects, psymtFlds);
875 
877  ReadFields(pointMesh::New(mesh), objects, ptFlds);
878 
879 
880  Info<< endl;
881 
882  // From renumbering:
883  // - from new cell/face back to original cell/face
884  labelList cellOrder;
885  labelList faceOrder;
886 
887  if (blockSize > 0)
888  {
889  // Renumbering in two phases. Should be done in one so mapping of
890  // fields is done correctly!
891 
892  label nBlocks = mesh.nCells()/blockSize;
893  Info<< "nBlocks = " << nBlocks << endl;
894 
895  // Read decompositionMethod dictionary
896  dictionary decomposeDict(renumberDictPtr().subDict("blockCoeffs"));
897  decomposeDict.set("numberOfSubdomains", nBlocks);
898 
899  bool oldParRun = UPstream::parRun();
900  UPstream::parRun() = false;
901 
903  (
904  decomposeDict
905  );
906 
907  labelList cellToRegion
908  (
909  decomposePtr().decompose
910  (
911  mesh,
912  mesh.cellCentres()
913  )
914  );
915 
916  // Restore state
917  UPstream::parRun() = oldParRun;
918 
919  // For debugging: write out region
920  createScalarField
921  (
922  mesh,
923  "cellDist",
924  cellToRegion
925  )().write();
926 
927  Info<< nl << "Written decomposition as volScalarField to "
928  << "cellDist for use in postprocessing."
929  << nl << endl;
930 
931 
932  cellOrder = regionRenumber(renumberPtr(), mesh, cellToRegion);
933 
934  // Determine new to old face order with new cell numbering
935  faceOrder = getRegionFaceOrder
936  (
937  mesh,
938  cellOrder,
939  cellToRegion
940  );
941  }
942  else
943  {
944  // Detemines sorted back to original cell ordering
945  cellOrder = renumberPtr().renumber
946  (
947  mesh,
948  mesh.cellCentres()
949  );
950 
951  if (sortCoupledFaceCells)
952  {
953  // Change order so all coupled patch faceCells are at the end.
954  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
955 
956  // Collect all boundary cells on coupled patches
957  label nBndCells = 0;
958  forAll(pbm, patchi)
959  {
960  if (pbm[patchi].coupled())
961  {
962  nBndCells += pbm[patchi].size();
963  }
964  }
965 
966  labelList reverseCellOrder = invert(mesh.nCells(), cellOrder);
967 
968  labelList bndCells(nBndCells);
969  labelList bndCellMap(nBndCells);
970  nBndCells = 0;
971  forAll(pbm, patchi)
972  {
973  if (pbm[patchi].coupled())
974  {
975  const labelUList& faceCells = pbm[patchi].faceCells();
976  forAll(faceCells, i)
977  {
978  label celli = faceCells[i];
979 
980  if (reverseCellOrder[celli] != -1)
981  {
982  bndCells[nBndCells] = celli;
983  bndCellMap[nBndCells++] = reverseCellOrder[celli];
984  reverseCellOrder[celli] = -1;
985  }
986  }
987  }
988  }
989  bndCells.setSize(nBndCells);
990  bndCellMap.setSize(nBndCells);
991 
992  // Sort
993  labelList order;
994  sortedOrder(bndCellMap, order);
995 
996  // Redo newReverseCellOrder
997  labelList newReverseCellOrder(mesh.nCells(), -1);
998 
999  label sortedI = mesh.nCells();
1000  forAllReverse(order, i)
1001  {
1002  label origCelli = bndCells[order[i]];
1003  newReverseCellOrder[origCelli] = --sortedI;
1004  }
1005 
1006  Info<< "Ordered all " << nBndCells << " cells with a coupled face"
1007  << " to the end of the cell list, starting at " << sortedI
1008  << endl;
1009 
1010  // Compact
1011  sortedI = 0;
1012  forAll(cellOrder, newCelli)
1013  {
1014  label origCelli = cellOrder[newCelli];
1015  if (newReverseCellOrder[origCelli] == -1)
1016  {
1017  newReverseCellOrder[origCelli] = sortedI++;
1018  }
1019  }
1020 
1021  // Update sorted back to original (unsorted) map
1022  cellOrder = invert(mesh.nCells(), newReverseCellOrder);
1023  }
1024 
1025 
1026  // Determine new to old face order with new cell numbering
1027  faceOrder = getFaceOrder
1028  (
1029  mesh,
1030  cellOrder // New to old cell
1031  );
1032  }
1033 
1034 
1035  if (!overwrite)
1036  {
1037  runTime++;
1038  }
1039 
1040 
1041  // Change the mesh.
1042  autoPtr<mapPolyMesh> map = reorderMesh(mesh, cellOrder, faceOrder);
1043 
1044 
1045  if (orderPoints)
1046  {
1047  polyTopoChange meshMod(mesh);
1048  autoPtr<mapPolyMesh> pointOrderMap = meshMod.changeMesh
1049  (
1050  mesh,
1051  false, // inflate
1052  true, // syncParallel
1053  false, // orderCells
1054  orderPoints // orderPoints
1055  );
1056 
1057  // Combine point reordering into map.
1058  const_cast<labelList&>(map().pointMap()) = UIndirectList<label>
1059  (
1060  map().pointMap(),
1061  pointOrderMap().pointMap()
1062  )();
1063 
1065  (
1066  pointOrderMap().reversePointMap(),
1067  const_cast<labelList&>(map().reversePointMap())
1068  );
1069  }
1070 
1071 
1072  // Update fields
1073  mesh.updateMesh(map);
1074 
1075  // Update proc maps
1076  if
1077  (
1078  cellProcAddressing.headerOk()
1079  && cellProcAddressing.size() == mesh.nCells()
1080  )
1081  {
1082  Info<< "Renumbering processor cell decomposition map "
1083  << cellProcAddressing.name() << endl;
1084 
1085  cellProcAddressing = labelList
1086  (
1087  UIndirectList<label>(cellProcAddressing, map().cellMap())
1088  );
1089  }
1090  if
1091  (
1092  faceProcAddressing.headerOk()
1093  && faceProcAddressing.size() == mesh.nFaces()
1094  )
1095  {
1096  Info<< "Renumbering processor face decomposition map "
1097  << faceProcAddressing.name() << endl;
1098 
1100  (
1102  );
1103 
1104  // Detect any flips.
1105  const labelHashSet& fff = map().flipFaceFlux();
1106  forAllConstIter(labelHashSet, fff, iter)
1107  {
1108  label facei = iter.key();
1109  label masterFacei = faceProcAddressing[facei];
1110 
1111  faceProcAddressing[facei] = -masterFacei;
1112 
1113  if (masterFacei == 0)
1114  {
1116  << " masterFacei:" << masterFacei << exit(FatalError);
1117  }
1118  }
1119  }
1120  if
1121  (
1122  pointProcAddressing.headerOk()
1123  && pointProcAddressing.size() == mesh.nPoints()
1124  )
1125  {
1126  Info<< "Renumbering processor point decomposition map "
1127  << pointProcAddressing.name() << endl;
1128 
1129  pointProcAddressing = labelList
1130  (
1131  UIndirectList<label>(pointProcAddressing, map().pointMap())
1132  );
1133  }
1134 
1135 
1136  // Move mesh (since morphing might not do this)
1137  if (map().hasMotionPoints())
1138  {
1139  mesh.movePoints(map().preMotionPoints());
1140  }
1141 
1142 
1143  {
1144  label band;
1145  scalar profile;
1146  scalar sumSqrIntersect;
1147  getBand
1148  (
1149  doFrontWidth,
1150  mesh.nCells(),
1151  mesh.faceOwner(),
1152  mesh.faceNeighbour(),
1153  band,
1154  profile,
1155  sumSqrIntersect
1156  );
1157  reduce(band, maxOp<label>());
1158  reduce(profile, sumOp<scalar>());
1159  scalar rmsFrontwidth = Foam::sqrt
1160  (
1161  returnReduce
1162  (
1163  sumSqrIntersect,
1164  sumOp<scalar>()
1165  )/mesh.globalData().nTotalCells()
1166  );
1167 
1168  Info<< "After renumbering :" << nl
1169  << " band : " << band << nl
1170  << " profile : " << profile << nl;
1171 
1172  if (doFrontWidth)
1173  {
1174 
1175  Info<< " rms frontwidth : " << rmsFrontwidth << nl;
1176  }
1177 
1178  Info<< endl;
1179  }
1180 
1181  if (orderPoints)
1182  {
1183  // Force edge calculation (since only reason that points would need to
1184  // be sorted)
1185  (void)mesh.edges();
1186 
1187  label nTotPoints = returnReduce
1188  (
1189  mesh.nPoints(),
1190  sumOp<label>()
1191  );
1192  label nTotIntPoints = returnReduce
1193  (
1194  mesh.nInternalPoints(),
1195  sumOp<label>()
1196  );
1197 
1198  label nTotEdges = returnReduce
1199  (
1200  mesh.nEdges(),
1201  sumOp<label>()
1202  );
1203  label nTotIntEdges = returnReduce
1204  (
1205  mesh.nInternalEdges(),
1206  sumOp<label>()
1207  );
1208  label nTotInt0Edges = returnReduce
1209  (
1210  mesh.nInternal0Edges(),
1211  sumOp<label>()
1212  );
1213  label nTotInt1Edges = returnReduce
1214  (
1215  mesh.nInternal1Edges(),
1216  sumOp<label>()
1217  );
1218 
1219  Info<< "Points:" << nl
1220  << " total : " << nTotPoints << nl
1221  << " internal: " << nTotIntPoints << nl
1222  << " boundary: " << nTotPoints-nTotIntPoints << nl
1223  << "Edges:" << nl
1224  << " total : " << nTotEdges << nl
1225  << " internal: " << nTotIntEdges << nl
1226  << " internal using 0 boundary points: "
1227  << nTotInt0Edges << nl
1228  << " internal using 1 boundary points: "
1229  << nTotInt1Edges-nTotInt0Edges << nl
1230  << " internal using 2 boundary points: "
1231  << nTotIntEdges-nTotInt1Edges << nl
1232  << " boundary: " << nTotEdges-nTotIntEdges << nl
1233  << endl;
1234  }
1235 
1236 
1237  if (overwrite)
1238  {
1239  mesh.setInstance(oldInstance);
1240  }
1241 
1242  Info<< "Writing mesh to " << mesh.facesInstance() << endl;
1243 
1244  mesh.write();
1245 
1246  if (cellProcAddressing.headerOk())
1247  {
1248  cellProcAddressing.instance() = mesh.facesInstance();
1249  if (cellProcAddressing.size() == mesh.nCells())
1250  {
1251  cellProcAddressing.write();
1252  }
1253  else
1254  {
1255  // procAddressing file no longer valid. Delete it.
1256  const fileName fName(cellProcAddressing.filePath());
1257  if (fName.size())
1258  {
1259  Info<< "Deleting inconsistent processor cell decomposition"
1260  << " map " << fName << endl;
1261  rm(fName);
1262  }
1263  }
1264  }
1265 
1266  if (faceProcAddressing.headerOk())
1267  {
1268  faceProcAddressing.instance() = mesh.facesInstance();
1269  if (faceProcAddressing.size() == mesh.nFaces())
1270  {
1271  faceProcAddressing.write();
1272  }
1273  else
1274  {
1275  const fileName fName(faceProcAddressing.filePath());
1276  if (fName.size())
1277  {
1278  Info<< "Deleting inconsistent processor face decomposition"
1279  << " map " << fName << endl;
1280  rm(fName);
1281  }
1282  }
1283  }
1284 
1285  if (pointProcAddressing.headerOk())
1286  {
1287  pointProcAddressing.instance() = mesh.facesInstance();
1288  if (pointProcAddressing.size() == mesh.nPoints())
1289  {
1290  pointProcAddressing.write();
1291  }
1292  else
1293  {
1294  const fileName fName(pointProcAddressing.filePath());
1295  if (fName.size())
1296  {
1297  Info<< "Deleting inconsistent processor point decomposition"
1298  << " map " << fName << endl;
1299  rm(fName);
1300  }
1301  }
1302  }
1303 
1304  if (boundaryProcAddressing.headerOk())
1305  {
1306  boundaryProcAddressing.instance() = mesh.facesInstance();
1307  if (boundaryProcAddressing.size() == mesh.boundaryMesh().size())
1308  {
1309  boundaryProcAddressing.write();
1310  }
1311  else
1312  {
1313  const fileName fName(boundaryProcAddressing.filePath());
1314  if (fName.size())
1315  {
1316  Info<< "Deleting inconsistent processor patch decomposition"
1317  << " map " << fName << endl;
1318  rm(fName);
1319  }
1320  }
1321  }
1322 
1323  if (writeMaps)
1324  {
1325  // For debugging: write out region
1326  createScalarField
1327  (
1328  mesh,
1329  "origCellID",
1330  map().cellMap()
1331  )().write();
1332 
1333  createScalarField
1334  (
1335  mesh,
1336  "cellID",
1337  identity(mesh.nCells())
1338  )().write();
1339 
1340  Info<< nl << "Written current cellID and origCellID as volScalarField"
1341  << " for use in postprocessing."
1342  << nl << endl;
1343 
1344  labelIOList
1345  (
1346  IOobject
1347  (
1348  "cellMap",
1349  mesh.facesInstance(),
1351  mesh,
1354  false
1355  ),
1356  map().cellMap()
1357  ).write();
1358 
1359  labelIOList
1360  (
1361  IOobject
1362  (
1363  "faceMap",
1364  mesh.facesInstance(),
1366  mesh,
1369  false
1370  ),
1371  map().faceMap()
1372  ).write();
1373 
1374  labelIOList
1375  (
1376  IOobject
1377  (
1378  "pointMap",
1379  mesh.facesInstance(),
1381  mesh,
1384  false
1385  ),
1386  map().pointMap()
1387  ).write();
1388  }
1389 
1390 
1391  // Renumber sets if required
1392  if (renumberSets)
1393  {
1394  Info<< endl;
1395 
1396  // Read sets
1397  IOobjectList objects(mesh, mesh.facesInstance(), "polyMesh/sets");
1398 
1399  {
1400  IOobjectList cSets(objects.lookupClass(cellSet::typeName));
1401  if (cSets.size())
1402  {
1403  Info<< "Renumbering cellSets:" << endl;
1404  forAllConstIter(IOobjectList, cSets, iter)
1405  {
1406  cellSet cs(*iter());
1407  Info<< " " << cs.name() << endl;
1408  cs.updateMesh(map());
1409  cs.instance() = mesh.facesInstance();
1410  cs.write();
1411  }
1412  }
1413  }
1414 
1415  {
1416  IOobjectList fSets(objects.lookupClass(faceSet::typeName));
1417  if (fSets.size())
1418  {
1419  Info<< "Renumbering faceSets:" << endl;
1420  forAllConstIter(IOobjectList, fSets, iter)
1421  {
1422  faceSet fs(*iter());
1423  Info<< " " << fs.name() << endl;
1424  fs.updateMesh(map());
1425  fs.instance() = mesh.facesInstance();
1426  fs.write();
1427  }
1428  }
1429  }
1430 
1431  {
1432  IOobjectList pSets(objects.lookupClass(pointSet::typeName));
1433  if (pSets.size())
1434  {
1435  Info<< "Renumbering pointSets:" << endl;
1436  forAllConstIter(IOobjectList, pSets, iter)
1437  {
1438  pointSet ps(*iter());
1439  Info<< " " << ps.name() << endl;
1440  ps.updateMesh(map());
1441  ps.instance() = mesh.facesInstance();
1442  ps.write();
1443  }
1444  }
1445  }
1446  }
1447 
1448  Info<< "End\n" << endl;
1449 
1450  return 0;
1451 }
1452 
1453 
1454 // ************************************************************************* //
Foam::surfaceFields.
ListType renumber(const labelUList &oldToNew, const ListType &)
Renumber the values (not the indices) of a list.
A simple container for copying or transferring objects of type <T>.
Definition: Xfer.H:85
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition: triad.C:394
wordList ReadFields(const Mesh &mesh, const IOobjectList &objects, PtrList< GeoField > &fields, const bool syncPar=true)
Read all fields of the specified type.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
void clearAddressing()
Clear addressing.
Definition: ZoneMesh.C:387
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
virtual void resetAddressing(const labelUList &, const boolList &)
Reset addressing and flip map (clearing demand-driven data)
Definition: faceZone.C:377
void sortedOrder(const UList< T > &, labelList &order)
Generate the (stable) sort order for the list.
A class for handling file names.
Definition: fileName.H:69
A list of face labels.
Definition: faceSet.H:48
void inplaceRenumber(const labelUList &oldToNew, ListType &)
Inplace renumber the values of a list.
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const boolList & flipMap() const
Return face flip map.
Definition: faceZone.H:253
A set of point labels.
Definition: pointSet.H:48
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:778
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
A list that is sorted upon construction or when explicitly requested with the sort() method...
Definition: List.H:68
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:74
PtrList< labelIOList > & faceProcAddressing
dimensionedSymmTensor sqr(const dimensionedVector &dv)
labelListList invertOneToMany(const label len, const labelUList &)
Invert one-to-many map. Unmapped elements will be size 0.
Definition: ListOps.C:67
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:309
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
static autoPtr< renumberMethod > New(const dictionary &renumberDict)
Return a reference to the selected renumbering method.
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: UList.H:440
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
Field reading functions for post-processing utilities.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:715
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
const cellList & cells() const
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
label nCells() const
static const pointMesh & New(const polyMesh &mesh)
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
virtual const labelList & faceNeighbour() const =0
Face face-neighbour addressing.
Abstract base class for renumbering.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
labelList invert(const label len, const labelUList &)
Invert one-to-one map. Unmapped elements will be -1.
Definition: ListOps.C:37
void reset(T *=0)
If object pointer already set, delete object and set to given.
Definition: autoPtrI.H:114
dynamicFvMesh & mesh
void Swap(T &a, T &b)
Definition: Swap.H:43
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
void resetPrimitives(const Xfer< pointField > &points, const Xfer< faceList > &faces, const Xfer< labelList > &owner, const Xfer< labelList > &neighbour, const labelList &patchSizes, const labelList &patchStarts, const bool validBoundary=true)
Reset mesh primitive data. Assumes all patch info correct.
Definition: polyMesh.C:642
Xfer< T > xferMove(T &)
Construct by transferring the contents of the arg.
A class for handling words, derived from string.
Definition: word.H:59
label nPoints
Cuthill-McKee renumbering.
void sort(UList< T > &)
Definition: UList.C:115
const vectorField & cellCentres() const
static const label labelMax
Definition: label.H:62
List< label > labelList
A List of labels.
Definition: labelList.H:56
virtual labelList renumber(const pointField &) const
Return the order in which cells need to be visited, i.e.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
errorManip< error > abort(error &err)
Definition: errorManip.H:131
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
static autoPtr< decompositionMethod > New(const dictionary &decompositionDict)
Return a reference to the selected decomposition method.
Foam::polyBoundaryMesh.
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1140
ListType reorder(const labelUList &oldToNew, const ListType &)
Reorder the elements (indices, not values) of a list.
static const char nl
Definition: Ostream.H:262
Post-processing mesh subset tool. Given the original mesh and the list of selected cells...
Definition: fvMeshSubset.H:73
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
IOobject dictIO(dictName, runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED, IOobject::NO_WRITE)
Internal & ref()
Return a reference to the dimensioned internal field.
label nFaces() const
void setSize(const label)
Reset size of List.
Definition: List.C:295
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:393
const cellZoneMesh & cellZones() const
Return cell zone mesh.
Definition: polyMesh.H:469
label patchi
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1023
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
word dictName("noiseDict")
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:62
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:772
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
A collection of cell labels.
Definition: cellSet.H:48
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
A List with indirect addressing.
Definition: fvMatrix.H:106
Direct mesh changes based on v1.3 polyTopoChange syntax.
label nPoints() const
messageStream Info
virtual const labelList & faceOwner() const =0
Face face-owner addresing.
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:463
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:53
virtual Ostream & write(const token &)=0
Write next token to stream.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:83
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:64
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:124
A class for managing temporary objects.
Definition: PtrList.H:54
bool rm(const fileName &)
Remove a file, returning true if successful otherwise false.
Definition: POSIX.C:819
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1017
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1004
Foam::argList args(argc, argv)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
label nTotalCells() const
Return total number of cells in decomposed mesh.
label nInternalFaces() const
runTime write()
IOList< label > labelIOList
Label container classes.
Definition: labelIOList.H:42
Namespace for OpenFOAM.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:451