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