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-2015 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
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();
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
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.
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
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.
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 
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 
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
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 
627  (
628  "writeMaps",
629  "write cellMap, faceMap, pointMap in polyMesh/"
630  );
631 
632  #include "setRootCase.H"
633  #include "createTime.H"
634  runTime.functionObjects().off();
635 
636 
637  // Force linker to include zoltan symbols. This section is only needed since
638  // Zoltan is a static library
639  #ifdef FOAM_USE_ZOLTAN
640  Info<< "renumberMesh built with zoltan support." << nl << endl;
641  (void)zoltanRenumber::typeName;
642  #endif
643 
644 
645  // Get times list
646  instantList Times = runTime.times();
647 
648  // Set startTime and endTime depending on -time and -latestTime options
649  #include "checkTimeOptions.H"
650 
651  runTime.setTime(Times[startTime], startTime);
652 
653  #include "createNamedMesh.H"
654  const word oldInstance = mesh.pointsInstance();
655 
656  const bool readDict = args.optionFound("dict");
657  const bool doFrontWidth = args.optionFound("frontWidth");
658  const bool overwrite = args.optionFound("overwrite");
659 
660  label band;
661  scalar profile;
662  scalar sumSqrIntersect;
663  getBand
664  (
665  doFrontWidth,
666  mesh.nCells(),
667  mesh.faceOwner(),
669  band,
670  profile,
671  sumSqrIntersect
672  );
673 
674  reduce(band, maxOp<label>());
675  reduce(profile, sumOp<scalar>());
676  scalar rmsFrontwidth = Foam::sqrt
677  (
679  (
680  sumSqrIntersect,
681  sumOp<scalar>()
683  );
684 
685  Info<< "Mesh size: " << mesh.globalData().nTotalCells() << nl
686  << "Before renumbering :" << nl
687  << " band : " << band << nl
688  << " profile : " << profile << nl;
689 
690  if (doFrontWidth)
691  {
692  Info<< " rms frontwidth : " << rmsFrontwidth << nl;
693  }
694 
695  Info<< endl;
696 
697  bool sortCoupledFaceCells = false;
698  bool writeMaps = args.optionFound("writeMaps");
699  bool orderPoints = false;
700  label blockSize = 0;
701  bool renumberSets = true;
702 
703  // Construct renumberMethod
704  autoPtr<IOdictionary> renumberDictPtr;
705  autoPtr<renumberMethod> renumberPtr;
706 
707  if (readDict)
708  {
709  const word dictName("renumberMeshDict");
710  #include "setSystemMeshDictionaryIO.H"
711 
712  Info<< "Renumber according to " << dictName << nl << endl;
713 
714  renumberDictPtr.reset(new IOdictionary(dictIO));
715  const IOdictionary& renumberDict = renumberDictPtr();
716 
717  renumberPtr = renumberMethod::New(renumberDict);
718 
719  sortCoupledFaceCells = renumberDict.lookupOrDefault
720  (
721  "sortCoupledFaceCells",
722  false
723  );
724  if (sortCoupledFaceCells)
725  {
726  Info<< "Sorting cells on coupled boundaries to be last." << nl
727  << endl;
728  }
729 
730  blockSize = renumberDict.lookupOrDefault("blockSize", 0);
731  if (blockSize > 0)
732  {
733  Info<< "Ordering cells into regions of size " << blockSize
734  << " (using decomposition);"
735  << " ordering faces into region-internal and region-external."
736  << nl << endl;
737 
738  if (blockSize < 0 || blockSize >= mesh.nCells())
739  {
741  << "Block size " << blockSize
742  << " should be positive integer"
743  << " and less than the number of cells in the mesh."
744  << exit(FatalError);
745  }
746  }
747 
748  orderPoints = renumberDict.lookupOrDefault("orderPoints", false);
749  if (orderPoints)
750  {
751  Info<< "Ordering points into internal and boundary points." << nl
752  << endl;
753  }
754 
755  renumberDict.lookup("writeMaps") >> writeMaps;
756  if (writeMaps)
757  {
758  Info<< "Writing renumber maps (new to old) to polyMesh." << nl
759  << endl;
760  }
761 
762  renumberSets = renumberDict.lookupOrDefault("renumberSets", true);
763  }
764  else
765  {
766  Info<< "Using default renumberMethod." << nl << endl;
767  dictionary renumberDict;
768  renumberPtr.reset(new CuthillMcKeeRenumber(renumberDict));
769  }
770 
771  Info<< "Selecting renumberMethod " << renumberPtr().type() << nl << endl;
772 
773 
774 
775  // Read parallel reconstruct maps
776  labelIOList cellProcAddressing
777  (
778  IOobject
779  (
780  "cellProcAddressing",
783  mesh,
785  ),
786  labelList(0)
787  );
788 
790  (
791  IOobject
792  (
793  "faceProcAddressing",
796  mesh,
798  ),
799  labelList(0)
800  );
801  labelIOList pointProcAddressing
802  (
803  IOobject
804  (
805  "pointProcAddressing",
808  mesh,
810  ),
811  labelList(0)
812  );
813  labelIOList boundaryProcAddressing
814  (
815  IOobject
816  (
817  "boundaryProcAddressing",
820  mesh,
822  ),
823  labelList(0)
824  );
825 
826 
827  // Read objects in time directory
828  IOobjectList objects(mesh, runTime.timeName());
829 
830 
831  // Read vol fields.
832 
834  ReadFields(mesh, objects, vsFlds);
835 
837  ReadFields(mesh, objects, vvFlds);
838 
840  ReadFields(mesh, objects, vstFlds);
841 
842  PtrList<volSymmTensorField> vsymtFlds;
843  ReadFields(mesh, objects, vsymtFlds);
844 
846  ReadFields(mesh, objects, vtFlds);
847 
848 
849  // Read surface fields.
850 
852  ReadFields(mesh, objects, ssFlds);
853 
855  ReadFields(mesh, objects, svFlds);
856 
858  ReadFields(mesh, objects, sstFlds);
859 
861  ReadFields(mesh, objects, ssymtFlds);
862 
864  ReadFields(mesh, objects, stFlds);
865 
866 
867  // Read point fields.
868 
870  ReadFields(pointMesh::New(mesh), objects, psFlds);
871 
873  ReadFields(pointMesh::New(mesh), objects, pvFlds);
874 
876  ReadFields(pointMesh::New(mesh), objects, pstFlds);
877 
879  ReadFields(pointMesh::New(mesh), objects, psymtFlds);
880 
882  ReadFields(pointMesh::New(mesh), objects, ptFlds);
883 
884 
885  // Read sets
888  PtrList<pointSet> pointSets;
889  if (renumberSets)
890  {
891  // Read sets
892  IOobjectList objects(mesh, mesh.facesInstance(), "polyMesh/sets");
893  {
894  IOobjectList cSets(objects.lookupClass(cellSet::typeName));
895  if (cSets.size())
896  {
897  Info<< "Reading cellSets:" << endl;
898  forAllConstIter(IOobjectList, cSets, iter)
899  {
900  cellSets.append(new cellSet(*iter()));
901  Info<< " " << cellSets.last().name() << endl;
902  }
903  }
904  }
905  {
906  IOobjectList fSets(objects.lookupClass(faceSet::typeName));
907  if (fSets.size())
908  {
909  Info<< "Reading faceSets:" << endl;
910  forAllConstIter(IOobjectList, fSets, iter)
911  {
912  faceSets.append(new faceSet(*iter()));
913  Info<< " " << faceSets.last().name() << endl;
914  }
915  }
916  }
917  {
918  IOobjectList pSets(objects.lookupClass(pointSet::typeName));
919  if (pSets.size())
920  {
921  Info<< "Reading pointSets:" << endl;
922  forAllConstIter(IOobjectList, pSets, iter)
923  {
924  pointSets.append(new pointSet(*iter()));
925  Info<< " " << pointSets.last().name() << endl;
926  }
927  }
928  }
929  }
930 
931 
932  Info<< endl;
933 
934  // From renumbering:
935  // - from new cell/face back to original cell/face
936  labelList cellOrder;
937  labelList faceOrder;
938 
939  if (blockSize > 0)
940  {
941  // Renumbering in two phases. Should be done in one so mapping of
942  // fields is done correctly!
943 
944  label nBlocks = mesh.nCells()/blockSize;
945  Info<< "nBlocks = " << nBlocks << endl;
946 
947  // Read decompositionMethod dictionary
948  dictionary decomposeDict(renumberDictPtr().subDict("blockCoeffs"));
949  decomposeDict.set("numberOfSubdomains", nBlocks);
950 
951  bool oldParRun = UPstream::parRun();
952  UPstream::parRun() = false;
953 
955  (
956  decomposeDict
957  );
958 
959  labelList cellToRegion
960  (
961  decomposePtr().decompose
962  (
963  mesh,
964  mesh.cellCentres()
965  )
966  );
967 
968  // Restore state
969  UPstream::parRun() = oldParRun;
970 
971  // For debugging: write out region
973  (
974  mesh,
975  "cellDist",
976  cellToRegion
977  )().write();
978 
979  Info<< nl << "Written decomposition as volScalarField to "
980  << "cellDist for use in postprocessing."
981  << nl << endl;
982 
983 
984  cellOrder = regionRenumber(renumberPtr(), mesh, cellToRegion);
985 
986  // Determine new to old face order with new cell numbering
987  faceOrder = getRegionFaceOrder
988  (
989  mesh,
990  cellOrder,
991  cellToRegion
992  );
993  }
994  else
995  {
996  // Detemines sorted back to original cell ordering
997  cellOrder = renumberPtr().renumber
998  (
999  mesh,
1000  mesh.cellCentres()
1001  );
1002 
1003  if (sortCoupledFaceCells)
1004  {
1005  // Change order so all coupled patch faceCells are at the end.
1006  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
1007 
1008  // Collect all boundary cells on coupled patches
1009  label nBndCells = 0;
1010  forAll(pbm, patchI)
1011  {
1012  if (pbm[patchI].coupled())
1013  {
1014  nBndCells += pbm[patchI].size();
1015  }
1016  }
1017 
1018  labelList reverseCellOrder = invert(mesh.nCells(), cellOrder);
1019 
1020  labelList bndCells(nBndCells);
1021  labelList bndCellMap(nBndCells);
1022  nBndCells = 0;
1023  forAll(pbm, patchI)
1024  {
1025  if (pbm[patchI].coupled())
1026  {
1027  const labelUList& faceCells = pbm[patchI].faceCells();
1028  forAll(faceCells, i)
1029  {
1030  label cellI = faceCells[i];
1031 
1032  if (reverseCellOrder[cellI] != -1)
1033  {
1034  bndCells[nBndCells] = cellI;
1035  bndCellMap[nBndCells++] = reverseCellOrder[cellI];
1036  reverseCellOrder[cellI] = -1;
1037  }
1038  }
1039  }
1040  }
1041  bndCells.setSize(nBndCells);
1042  bndCellMap.setSize(nBndCells);
1043 
1044  // Sort
1045  labelList order;
1046  sortedOrder(bndCellMap, order);
1047 
1048  // Redo newReverseCellOrder
1049  labelList newReverseCellOrder(mesh.nCells(), -1);
1050 
1051  label sortedI = mesh.nCells();
1052  forAllReverse(order, i)
1053  {
1054  label origCellI = bndCells[order[i]];
1055  newReverseCellOrder[origCellI] = --sortedI;
1056  }
1057 
1058  Info<< "Ordered all " << nBndCells << " cells with a coupled face"
1059  << " to the end of the cell list, starting at " << sortedI
1060  << endl;
1061 
1062  // Compact
1063  sortedI = 0;
1064  forAll(cellOrder, newCellI)
1065  {
1066  label origCellI = cellOrder[newCellI];
1067  if (newReverseCellOrder[origCellI] == -1)
1068  {
1069  newReverseCellOrder[origCellI] = sortedI++;
1070  }
1071  }
1072 
1073  // Update sorted back to original (unsorted) map
1074  cellOrder = invert(mesh.nCells(), newReverseCellOrder);
1075  }
1076 
1077 
1078  // Determine new to old face order with new cell numbering
1079  faceOrder = getFaceOrder
1080  (
1081  mesh,
1082  cellOrder // New to old cell
1083  );
1084  }
1085 
1086 
1087  if (!overwrite)
1088  {
1089  runTime++;
1090  }
1091 
1092 
1093  // Change the mesh.
1094  autoPtr<mapPolyMesh> map = reorderMesh(mesh, cellOrder, faceOrder);
1095 
1096 
1097  if (orderPoints)
1098  {
1099  polyTopoChange meshMod(mesh);
1100  autoPtr<mapPolyMesh> pointOrderMap = meshMod.changeMesh
1101  (
1102  mesh,
1103  false, // inflate
1104  true, // syncParallel
1105  false, // orderCells
1106  orderPoints // orderPoints
1107  );
1108 
1109  // Combine point reordering into map.
1110  const_cast<labelList&>(map().pointMap()) = UIndirectList<label>
1111  (
1112  map().pointMap(),
1113  pointOrderMap().pointMap()
1114  )();
1115 
1117  (
1118  pointOrderMap().reversePointMap(),
1119  const_cast<labelList&>(map().reversePointMap())
1120  );
1121  }
1122 
1123 
1124  // Update fields
1125  mesh.updateMesh(map);
1126 
1127  // Update proc maps
1128  if
1129  (
1130  cellProcAddressing.headerOk()
1131  && cellProcAddressing.size() == mesh.nCells()
1132  )
1133  {
1134  Info<< "Renumbering processor cell decomposition map "
1135  << cellProcAddressing.name() << endl;
1136 
1137  cellProcAddressing = labelList
1138  (
1139  UIndirectList<label>(cellProcAddressing, map().cellMap())
1140  );
1141  }
1142  if
1143  (
1144  faceProcAddressing.headerOk()
1145  && faceProcAddressing.size() == mesh.nFaces()
1146  )
1147  {
1148  Info<< "Renumbering processor face decomposition map "
1149  << faceProcAddressing.name() << endl;
1150 
1152  (
1154  );
1155 
1156  // Detect any flips.
1157  const labelHashSet& fff = map().flipFaceFlux();
1158  forAllConstIter(labelHashSet, fff, iter)
1159  {
1160  label faceI = iter.key();
1161  label masterFaceI = faceProcAddressing[faceI];
1162 
1163  faceProcAddressing[faceI] = -masterFaceI;
1164 
1165  if (masterFaceI == 0)
1166  {
1168  << " masterFaceI:" << masterFaceI << exit(FatalError);
1169  }
1170  }
1171  }
1172  if
1173  (
1174  pointProcAddressing.headerOk()
1175  && pointProcAddressing.size() == mesh.nPoints()
1176  )
1177  {
1178  Info<< "Renumbering processor point decomposition map "
1179  << pointProcAddressing.name() << endl;
1180 
1181  pointProcAddressing = labelList
1182  (
1183  UIndirectList<label>(pointProcAddressing, map().pointMap())
1184  );
1185  }
1186 
1187 
1188  // Move mesh (since morphing might not do this)
1189  if (map().hasMotionPoints())
1190  {
1191  mesh.movePoints(map().preMotionPoints());
1192  }
1193 
1194 
1195  {
1196  label band;
1197  scalar profile;
1198  scalar sumSqrIntersect;
1199  getBand
1200  (
1201  doFrontWidth,
1202  mesh.nCells(),
1203  mesh.faceOwner(),
1204  mesh.faceNeighbour(),
1205  band,
1206  profile,
1207  sumSqrIntersect
1208  );
1209  reduce(band, maxOp<label>());
1210  reduce(profile, sumOp<scalar>());
1211  scalar rmsFrontwidth = Foam::sqrt
1212  (
1213  returnReduce
1214  (
1215  sumSqrIntersect,
1216  sumOp<scalar>()
1218  );
1219 
1220  Info<< "After renumbering :" << nl
1221  << " band : " << band << nl
1222  << " profile : " << profile << nl;
1223 
1224  if (doFrontWidth)
1225  {
1226 
1227  Info<< " rms frontwidth : " << rmsFrontwidth << nl;
1228  }
1229 
1230  Info<< endl;
1231  }
1232 
1233  if (orderPoints)
1234  {
1235  // Force edge calculation (since only reason that points would need to
1236  // be sorted)
1237  (void)mesh.edges();
1238 
1239  label nTotPoints = returnReduce
1240  (
1241  mesh.nPoints(),
1242  sumOp<label>()
1243  );
1244  label nTotIntPoints = returnReduce
1245  (
1247  sumOp<label>()
1248  );
1249 
1250  label nTotEdges = returnReduce
1251  (
1252  mesh.nEdges(),
1253  sumOp<label>()
1254  );
1255  label nTotIntEdges = returnReduce
1256  (
1257  mesh.nInternalEdges(),
1258  sumOp<label>()
1259  );
1260  label nTotInt0Edges = returnReduce
1261  (
1263  sumOp<label>()
1264  );
1265  label nTotInt1Edges = returnReduce
1266  (
1268  sumOp<label>()
1269  );
1270 
1271  Info<< "Points:" << nl
1272  << " total : " << nTotPoints << nl
1273  << " internal: " << nTotIntPoints << nl
1274  << " boundary: " << nTotPoints-nTotIntPoints << nl
1275  << "Edges:" << nl
1276  << " total : " << nTotEdges << nl
1277  << " internal: " << nTotIntEdges << nl
1278  << " internal using 0 boundary points: "
1279  << nTotInt0Edges << nl
1280  << " internal using 1 boundary points: "
1281  << nTotInt1Edges-nTotInt0Edges << nl
1282  << " internal using 2 boundary points: "
1283  << nTotIntEdges-nTotInt1Edges << nl
1284  << " boundary: " << nTotEdges-nTotIntEdges << nl
1285  << endl;
1286  }
1287 
1288 
1289  if (overwrite)
1290  {
1291  mesh.setInstance(oldInstance);
1292  }
1293 
1294  Info<< "Writing mesh to " << mesh.facesInstance() << endl;
1295 
1296  mesh.write();
1297 
1298  if (cellProcAddressing.headerOk())
1299  {
1300  cellProcAddressing.instance() = mesh.facesInstance();
1301  if (cellProcAddressing.size() == mesh.nCells())
1302  {
1303  cellProcAddressing.write();
1304  }
1305  else
1306  {
1307  // procAddressing file no longer valid. Delete it.
1308  const fileName fName(cellProcAddressing.filePath());
1309  if (fName.size())
1310  {
1311  Info<< "Deleting inconsistent processor cell decomposition"
1312  << " map " << fName << endl;
1313  rm(fName);
1314  }
1315  }
1316  }
1317 
1318  if (faceProcAddressing.headerOk())
1319  {
1320  faceProcAddressing.instance() = mesh.facesInstance();
1321  if (faceProcAddressing.size() == mesh.nFaces())
1322  {
1323  faceProcAddressing.write();
1324  }
1325  else
1326  {
1327  const fileName fName(faceProcAddressing.filePath());
1328  if (fName.size())
1329  {
1330  Info<< "Deleting inconsistent processor face decomposition"
1331  << " map " << fName << endl;
1332  rm(fName);
1333  }
1334  }
1335  }
1336 
1337  if (pointProcAddressing.headerOk())
1338  {
1339  pointProcAddressing.instance() = mesh.facesInstance();
1340  if (pointProcAddressing.size() == mesh.nPoints())
1341  {
1342  pointProcAddressing.write();
1343  }
1344  else
1345  {
1346  const fileName fName(pointProcAddressing.filePath());
1347  if (fName.size())
1348  {
1349  Info<< "Deleting inconsistent processor point decomposition"
1350  << " map " << fName << endl;
1351  rm(fName);
1352  }
1353  }
1354  }
1355 
1356  if (boundaryProcAddressing.headerOk())
1357  {
1358  boundaryProcAddressing.instance() = mesh.facesInstance();
1359  if (boundaryProcAddressing.size() == mesh.boundaryMesh().size())
1360  {
1361  boundaryProcAddressing.write();
1362  }
1363  else
1364  {
1365  const fileName fName(boundaryProcAddressing.filePath());
1366  if (fName.size())
1367  {
1368  Info<< "Deleting inconsistent processor patch decomposition"
1369  << " map " << fName << endl;
1370  rm(fName);
1371  }
1372  }
1373  }
1374 
1375  if (writeMaps)
1376  {
1377  // For debugging: write out region
1379  (
1380  mesh,
1381  "origCellID",
1382  map().cellMap()
1383  )().write();
1384 
1386  (
1387  mesh,
1388  "cellID",
1389  identity(mesh.nCells())
1390  )().write();
1391 
1392  Info<< nl << "Written current cellID and origCellID as volScalarField"
1393  << " for use in postprocessing."
1394  << nl << endl;
1395 
1396  labelIOList
1397  (
1398  IOobject
1399  (
1400  "cellMap",
1401  mesh.facesInstance(),
1403  mesh,
1406  false
1407  ),
1408  map().cellMap()
1409  ).write();
1410 
1411  labelIOList
1412  (
1413  IOobject
1414  (
1415  "faceMap",
1416  mesh.facesInstance(),
1418  mesh,
1421  false
1422  ),
1423  map().faceMap()
1424  ).write();
1425 
1426  labelIOList
1427  (
1428  IOobject
1429  (
1430  "pointMap",
1431  mesh.facesInstance(),
1433  mesh,
1436  false
1437  ),
1438  map().pointMap()
1439  ).write();
1440  }
1441 
1442  if (renumberSets)
1443  {
1444  forAll(cellSets, i)
1445  {
1446  cellSets[i].updateMesh(map());
1447  cellSets[i].instance() = mesh.facesInstance();
1448  cellSets[i].write();
1449  }
1450  forAll(faceSets, i)
1451  {
1452  faceSets[i].updateMesh(map());
1453  faceSets[i].instance() = mesh.facesInstance();
1454  faceSets[i].write();
1455  }
1456  forAll(pointSets, i)
1457  {
1458  pointSets[i].updateMesh(map());
1459  pointSets[i].instance() = mesh.facesInstance();
1460  pointSets[i].write();
1461  }
1462  }
1463 
1464  Info<< "\nEnd\n" << endl;
1465 
1466  return 0;
1467 }
1468 
1469 
1470 // ************************************************************************* //
Foam::fvMeshSubset::setLargeCellSubset
void setLargeCellSubset(const labelList &region, const label currentRegion, const label patchID=-1, const bool syncCouples=true)
Set the subset from all cells with region == currentRegion.
Definition: fvMeshSubset.C:795
volFields.H
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
Foam::maxOp
Definition: ops.H:172
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::IOobject::filePath
fileName filePath() const
Return complete path + object name if the file exists.
Definition: IOobject.C:336
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
Foam::PtrList::last
T & last()
Return reference to the last element of the list.
Definition: PtrListI.H:60
Foam::faceMap
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Definition: blockMeshMergeFast.C:90
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::labelList
List< label > labelList
A List of labels.
Definition: labelList.H:56
Foam::reorder
ListType reorder(const labelUList &oldToNew, const ListType &)
Reorder the elements (indices, not values) of a list.
Definition: ListOpsTemplates.C:74
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:86
Foam::PtrList::append
void append(T *)
Append an element at the end of the list.
Foam::primitiveMesh::nInternalEdges
label nInternalEdges() const
Internal edges using 0,1 or 2 boundary points.
Definition: primitiveMeshI.H:63
Foam::SortableList::sort
void sort()
(stable) sort the list (if changed after construction time)
Definition: SortableList.C:95
Foam::polyBoundaryMesh
Foam::polyBoundaryMesh.
Definition: polyBoundaryMesh.H:60
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:117
Foam::argList::addNote
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:139
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::fvMeshSubset
Post-processing mesh subset tool. Given the original mesh and the list of selected cells,...
Definition: fvMeshSubset.H:73
reorderMesh
autoPtr< mapPolyMesh > reorderMesh(polyMesh &mesh, const labelList &cellOrder, const labelList &faceOrder)
Definition: renumberMesh.C:398
Foam::polyTopoChange::changeMesh
autoPtr< mapPolyMesh > changeMesh(polyMesh &mesh, const bool inflate, const bool syncParallel=true, const bool orderCells=false, const bool orderPoints=false)
Inplace changes mesh without change of patches.
Definition: polyTopoChange.C:3039
Foam::invertOneToMany
labelListList invertOneToMany(const label len, const labelUList &)
Invert one-to-many map. Unmapped elements will be size 0.
Definition: ListOps.C:67
Foam::decompositionMethod::New
static autoPtr< decompositionMethod > New(const dictionary &decompositionDict)
Return a reference to the selected decomposition method.
Definition: decompositionMethod.C:179
Foam::rm
bool rm(const fileName &)
Remove a file, returning true if successful otherwise false.
Definition: POSIX.C:954
fvMeshSubset.H
Foam::primitiveMesh::nInternal0Edges
label nInternal0Edges() const
Internal edges (i.e. not on boundary face) using.
Definition: primitiveMeshI.H:47
addOverwriteOption.H
Foam::ReadFields
wordList ReadFields(const typename GeoMesh::Mesh &mesh, const IOobjectList &objects, PtrList< GeometricField< Type, PatchField, GeoMesh > > &fields, const bool syncPar=true, const bool readOldTime=false)
Helper routine to read Geometric fields.
Foam::labelMax
static const label labelMax
Definition: label.H:62
Foam::polyMesh::meshSubDir
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:309
polyTopoChange.H
Foam::UPstream::parRun
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:377
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
Foam::primitiveMesh::cells
const cellList & cells() const
Definition: primitiveMeshCells.C:136
Foam::polyTopoChange
Direct mesh changes based on v1.3 polyTopoChange syntax.
Definition: polyTopoChange.H:97
Foam::globalMeshData::nTotalCells
label nTotalCells() const
Return total number of cells in decomposed mesh.
Definition: globalMeshData.H:394
Foam::argList::addBoolOption
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:98
Foam::dictionary::lookupOrDefault
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Definition: dictionaryTemplates.C:33
Foam::polyMesh::facesInstance
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:778
Foam::polyMesh::cellZones
const cellZoneMesh & cellZones() const
Return cell zone mesh.
Definition: polyMesh.H:469
Foam::regIOobject::write
virtual bool write() const
Write using setting from DB.
Definition: regIOobjectWrite.C:126
Foam::faceSet
A list of face labels.
Definition: faceSet.H:48
Foam::faceZone::resetAddressing
virtual void resetAddressing(const labelUList &, const boolList &)
Reset addressing and flip map (clearing demand-driven data)
Definition: faceZone.C:408
IOobjectList.H
Foam::fvMesh::movePoints
virtual tmp< scalarField > movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: fvMesh.C:726
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
surfaceFields.H
Foam::surfaceFields.
decompositionMethod.H
Foam::primitiveMesh::edges
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
Definition: primitiveMeshEdges.C:353
getFaceOrder
labelList getFaceOrder(const primitiveMesh &mesh, const labelList &cellOrder)
Definition: renumberMesh.C:165
Foam::IOobject::instance
const fileName & instance() const
Definition: IOobject.H:350
Foam::HashSet< label, Hash< label > >
Foam::labelIOList
IOList< label > labelIOList
Label container classes.
Definition: labelIOList.H:42
Foam::primitiveMesh::nEdges
label nEdges() const
Definition: primitiveMeshI.H:41
addTimeOptions.H
setSystemMeshDictionaryIO.H
Foam::Xfer
A simple container for copying or transferring objects of type <T>.
Definition: Xfer.H:85
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::MeshObject< polyMesh, UpdateableMeshObject, pointMesh >::New
static const pointMesh & New(const polyMesh &mesh)
Definition: MeshObject.C:44
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Foam::polyMesh::faceZones
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:463
addDictOption.H
Foam::faceSets
Definition: faceSets.H:44
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
Foam::renumber
ListType renumber(const labelUList &oldToNew, const ListType &)
Renumber the values (not the indices) of a list.
Definition: ListOpsTemplates.C:32
Foam::inplaceRenumber
void inplaceRenumber(const labelUList &oldToNew, ListType &)
Inplace renumber the values of a list.
Definition: ListOpsTemplates.C:57
Foam::primitiveMesh::nCells
label nCells() const
Definition: primitiveMeshI.H:64
Foam::fvMesh::write
virtual bool write() const
Write mesh using IO settings from time.
Definition: fvMesh.C:873
Foam::diff
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition: triad.C:407
Foam::IOobject::headerOk
bool headerOk()
Read and check header info.
Definition: IOobject.C:439
dictName
const word dictName("particleTrackDict")
Foam::polyMesh::pointsInstance
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:772
SortableList.H
Foam::sortedOrder
void sortedOrder(const UList< T > &, labelList &order)
Generate the (stable) sort order for the list.
Definition: ListOpsTemplates.C:179
Foam::fvMeshSubset::cellMap
const labelList & cellMap() const
Return cell map.
Definition: fvMeshSubset.C:1542
Foam::renumberMethod::renumber
virtual labelList renumber(const pointField &) const
Return the order in which cells need to be visited, i.e.
Definition: renumberMethod.H:114
Foam::primitiveMesh::nPoints
label nPoints() const
Definition: primitiveMeshI.H:35
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:43
renumberMethod.H
Foam::label
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
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::faceZone
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:64
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
argList.H
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1017
faceSet.H
addRegionOption.H
Foam::fvMesh::updateMesh
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
Definition: fvMesh.C:802
Foam::ZoneMesh< faceZone, polyMesh >
Foam::IOobject::name
const word & name() const
Return name.
Definition: IOobject.H:273
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
Foam::cellSets
Definition: cellSets.H:44
Foam::identity
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::polyMesh::resetPrimitives
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
Foam::primitiveMesh::nInternalFaces
label nInternalFaces() const
Definition: primitiveMeshI.H:52
createNamedMesh.H
forAllReverse
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: UList.H:418
Foam::SortableList
A list that is sorted upon construction or when explicitly requested with the sort() method.
Definition: List.H:65
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
fld
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::HashTable::found
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:109
getBand
label getBand(const labelList &owner, const labelList &neighbour)
Definition: renumberMesh.C:96
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::primitiveMesh::nInternalPoints
label nInternalPoints() const
Points not on boundary.
Definition: primitiveMeshI.H:35
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::primitiveMesh::nInternal1Edges
label nInternal1Edges() const
Internal edges using 0 or 1 boundary point.
Definition: primitiveMeshI.H:55
Foam::cellSet
A collection of cell labels.
Definition: cellSet.H:48
Foam::IOobjectList::lookupClass
IOobjectList lookupClass(const word &className) const
Return the list for all IOobjects of a given class.
Definition: IOobjectList.C:199
checkTimeOptions.H
Foam::IOobjectList
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::renumberMethod
Abstract base class for renumbering.
Definition: renumberMethod.H:48
Foam::polyMesh::setInstance
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:32
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
setRootCase.H
Foam::primitiveMesh::cellCentres
const vectorField & cellCentres() const
Definition: primitiveMeshCellCentresAndVols.C:211
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1004
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::xferMove
Xfer< T > xferMove(T &)
Construct by transferring the contents of the arg.
Foam::primitiveMesh::nFaces
label nFaces() const
Definition: primitiveMeshI.H:58
ReadFields.H
Helper routine to read fields.
Foam::sumOp
Definition: ops.H:162
Foam::pointSet
A set of point labels.
Definition: pointSet.H:48
regionRenumber
labelList regionRenumber(const renumberMethod &method, const fvMesh &mesh, const labelList &cellToRegion)
Definition: renumberMesh.C:552
createScalarField
tmp< volScalarField > createScalarField(const fvMesh &mesh, const word &name, const labelList &elems)
Definition: renumberMesh.C:60
main
int main(int argc, char *argv[])
Definition: renumberMesh.C:606
Foam::Time::timeName
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
Foam::primitiveMesh::isInternalFace
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
Definition: primitiveMeshI.H:70
Foam::UList
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
Foam::IOList< label >
Foam::PtrList::size
label size() const
Return the number of elements in the PtrList.
Definition: PtrListI.H:32
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:116
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:333
Foam::fvMeshSubset::subMesh
const fvMesh & subMesh() const
Return reference to subset mesh.
Definition: fvMeshSubset.C:1472
createTime.H
Foam::argList::optionFound
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
Foam::SortableList::indices
const labelList & indices() const
Return the list of sorted indices. Updated every sort.
Definition: SortableList.H:90
faceProcAddressing
PtrList< labelIOList > & faceProcAddressing
Definition: checkFaceAddressingComp.H:9
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
Foam::autoPtr::reset
void reset(T *=0)
If object pointer already set, delete object and set to given.
Definition: autoPtrI.H:114
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
zoltanRenumber.H
Foam::UIndirectList
A List with indirect addressing.
Definition: fvMatrix.H:106
Foam::invert
labelList invert(const label len, const labelUList &)
Invert one-to-one map. Unmapped elements will be -1.
Definition: ListOps.C:37
patches
patches[0]
Definition: createSingleCellMesh.H:36
Foam::renumberMethod::New
static autoPtr< renumberMethod > New(const dictionary &renumberDict)
Return a reference to the selected renumbering method.
Definition: renumberMethod.C:49
getRegionFaceOrder
labelList getRegionFaceOrder(const primitiveMesh &mesh, const labelList &cellOrder, const labelList &cellToRegion)
Definition: renumberMesh.C:257
startTime
Foam::label startTime
Definition: checkTimeOptions.H:5
Foam::sort
void sort(UList< T > &)
Definition: UList.C:107
cellSet.H
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::ZoneMesh::clearAddressing
void clearAddressing()
Clear addressing.
Definition: ZoneMesh.C:394
Foam::polyMesh::globalData
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1143
Foam::Swap
void Swap(T &a, T &b)
Definition: Swap.H:43
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::CuthillMcKeeRenumber
Cuthill-McKee renumbering.
Definition: CuthillMcKeeRenumber.H:47
args
Foam::argList args(argc, argv)
Foam::IOobject::READ_IF_PRESENT
@ READ_IF_PRESENT
Definition: IOobject.H:110
Foam::faceZone::flipMap
const boolList & flipMap() const
Return face flip map.
Definition: faceZone.H:249
zeroGradientFvPatchFields.H
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
write
Tcoeff write()
dictIO
IOobject dictIO(dictName, runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED, IOobject::NO_WRITE)
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::polyMesh::faceNeighbour
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1023
CuthillMcKeeRenumber.H
Foam::dictionary::set
void set(entry *)
Assign a new entry, overwrite any existing entry.
Definition: dictionary.C:856
pointSet.H
Foam::primitiveMesh
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:79