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 
626 
627  #include "setRootCase.H"
628  #include "createTime.H"
629  runTime.functionObjects().off();
630 
631 
632  // Force linker to include zoltan symbols. This section is only needed since
633  // Zoltan is a static library
634  #ifdef FOAM_USE_ZOLTAN
635  Info<< "renumberMesh built with zoltan support." << nl << endl;
636  (void)zoltanRenumber::typeName;
637  #endif
638 
639 
640  // Get times list
641  instantList Times = runTime.times();
642 
643  // Set startTime and endTime depending on -time and -latestTime options
644  #include "checkTimeOptions.H"
645 
646  runTime.setTime(Times[startTime], startTime);
647 
648  #include "createNamedMesh.H"
649  const word oldInstance = mesh.pointsInstance();
650 
651  const bool readDict = args.optionFound("dict");
652  const bool doFrontWidth = args.optionFound("frontWidth");
653  const bool overwrite = args.optionFound("overwrite");
654 
655  label band;
656  scalar profile;
657  scalar sumSqrIntersect;
658  getBand
659  (
660  doFrontWidth,
661  mesh.nCells(),
662  mesh.faceOwner(),
664  band,
665  profile,
666  sumSqrIntersect
667  );
668 
669  reduce(band, maxOp<label>());
670  reduce(profile, sumOp<scalar>());
671  scalar rmsFrontwidth = Foam::sqrt
672  (
674  (
675  sumSqrIntersect,
676  sumOp<scalar>()
678  );
679 
680  Info<< "Mesh size: " << mesh.globalData().nTotalCells() << nl
681  << "Before renumbering :" << nl
682  << " band : " << band << nl
683  << " profile : " << profile << nl;
684 
685  if (doFrontWidth)
686  {
687  Info<< " rms frontwidth : " << rmsFrontwidth << nl;
688  }
689 
690  Info<< endl;
691 
692  bool sortCoupledFaceCells = false;
693  bool writeMaps = false;
694  bool orderPoints = false;
695  label blockSize = 0;
696  bool renumberSets = true;
697 
698  // Construct renumberMethod
699  autoPtr<IOdictionary> renumberDictPtr;
700  autoPtr<renumberMethod> renumberPtr;
701 
702  if (readDict)
703  {
704  const word dictName("renumberMeshDict");
705  #include "setSystemMeshDictionaryIO.H"
706 
707  Info<< "Renumber according to " << dictName << nl << endl;
708 
709  renumberDictPtr.reset(new IOdictionary(dictIO));
710  const IOdictionary& renumberDict = renumberDictPtr();
711 
712  renumberPtr = renumberMethod::New(renumberDict);
713 
714  sortCoupledFaceCells = renumberDict.lookupOrDefault
715  (
716  "sortCoupledFaceCells",
717  false
718  );
719  if (sortCoupledFaceCells)
720  {
721  Info<< "Sorting cells on coupled boundaries to be last." << nl
722  << endl;
723  }
724 
725  blockSize = renumberDict.lookupOrDefault("blockSize", 0);
726  if (blockSize > 0)
727  {
728  Info<< "Ordering cells into regions of size " << blockSize
729  << " (using decomposition);"
730  << " ordering faces into region-internal and region-external."
731  << nl << endl;
732 
733  if (blockSize < 0 || blockSize >= mesh.nCells())
734  {
736  << "Block size " << blockSize
737  << " should be positive integer"
738  << " and less than the number of cells in the mesh."
739  << exit(FatalError);
740  }
741  }
742 
743  orderPoints = renumberDict.lookupOrDefault("orderPoints", false);
744  if (orderPoints)
745  {
746  Info<< "Ordering points into internal and boundary points." << nl
747  << endl;
748  }
749 
750  renumberDict.lookup("writeMaps") >> writeMaps;
751  if (writeMaps)
752  {
753  Info<< "Writing renumber maps (new to old) to polyMesh." << nl
754  << endl;
755  }
756 
757  renumberSets = renumberDict.lookupOrDefault("renumberSets", true);
758  }
759  else
760  {
761  Info<< "Using default renumberMethod." << nl << endl;
762  dictionary renumberDict;
763  renumberPtr.reset(new CuthillMcKeeRenumber(renumberDict));
764  }
765 
766  Info<< "Selecting renumberMethod " << renumberPtr().type() << nl << endl;
767 
768 
769 
770  // Read parallel reconstruct maps
771  labelIOList cellProcAddressing
772  (
773  IOobject
774  (
775  "cellProcAddressing",
778  mesh,
780  ),
781  labelList(0)
782  );
783 
785  (
786  IOobject
787  (
788  "faceProcAddressing",
791  mesh,
793  ),
794  labelList(0)
795  );
796  labelIOList pointProcAddressing
797  (
798  IOobject
799  (
800  "pointProcAddressing",
803  mesh,
805  ),
806  labelList(0)
807  );
808  labelIOList boundaryProcAddressing
809  (
810  IOobject
811  (
812  "boundaryProcAddressing",
815  mesh,
817  ),
818  labelList(0)
819  );
820 
821 
822  // Read objects in time directory
823  IOobjectList objects(mesh, runTime.timeName());
824 
825 
826  // Read vol fields.
827 
829  ReadFields(mesh, objects, vsFlds);
830 
832  ReadFields(mesh, objects, vvFlds);
833 
835  ReadFields(mesh, objects, vstFlds);
836 
837  PtrList<volSymmTensorField> vsymtFlds;
838  ReadFields(mesh, objects, vsymtFlds);
839 
841  ReadFields(mesh, objects, vtFlds);
842 
843 
844  // Read surface fields.
845 
847  ReadFields(mesh, objects, ssFlds);
848 
850  ReadFields(mesh, objects, svFlds);
851 
853  ReadFields(mesh, objects, sstFlds);
854 
856  ReadFields(mesh, objects, ssymtFlds);
857 
859  ReadFields(mesh, objects, stFlds);
860 
861 
862  // Read point fields.
863 
865  ReadFields(pointMesh::New(mesh), objects, psFlds);
866 
868  ReadFields(pointMesh::New(mesh), objects, pvFlds);
869 
871  ReadFields(pointMesh::New(mesh), objects, pstFlds);
872 
874  ReadFields(pointMesh::New(mesh), objects, psymtFlds);
875 
877  ReadFields(pointMesh::New(mesh), objects, ptFlds);
878 
879 
880  // Read sets
883  PtrList<pointSet> pointSets;
884  if (renumberSets)
885  {
886  // Read sets
887  IOobjectList objects(mesh, mesh.facesInstance(), "polyMesh/sets");
888  {
889  IOobjectList cSets(objects.lookupClass(cellSet::typeName));
890  if (cSets.size())
891  {
892  Info<< "Reading cellSets:" << endl;
893  forAllConstIter(IOobjectList, cSets, iter)
894  {
895  cellSets.append(new cellSet(*iter()));
896  Info<< " " << cellSets.last().name() << endl;
897  }
898  }
899  }
900  {
901  IOobjectList fSets(objects.lookupClass(faceSet::typeName));
902  if (fSets.size())
903  {
904  Info<< "Reading faceSets:" << endl;
905  forAllConstIter(IOobjectList, fSets, iter)
906  {
907  faceSets.append(new faceSet(*iter()));
908  Info<< " " << faceSets.last().name() << endl;
909  }
910  }
911  }
912  {
913  IOobjectList pSets(objects.lookupClass(pointSet::typeName));
914  if (pSets.size())
915  {
916  Info<< "Reading pointSets:" << endl;
917  forAllConstIter(IOobjectList, pSets, iter)
918  {
919  pointSets.append(new pointSet(*iter()));
920  Info<< " " << pointSets.last().name() << endl;
921  }
922  }
923  }
924  }
925 
926 
927  Info<< endl;
928 
929  // From renumbering:
930  // - from new cell/face back to original cell/face
931  labelList cellOrder;
932  labelList faceOrder;
933 
934  if (blockSize > 0)
935  {
936  // Renumbering in two phases. Should be done in one so mapping of
937  // fields is done correctly!
938 
939  label nBlocks = mesh.nCells()/blockSize;
940  Info<< "nBlocks = " << nBlocks << endl;
941 
942  // Read decompositionMethod dictionary
943  dictionary decomposeDict(renumberDictPtr().subDict("blockCoeffs"));
944  decomposeDict.set("numberOfSubdomains", nBlocks);
945 
946  bool oldParRun = UPstream::parRun();
947  UPstream::parRun() = false;
948 
950  (
951  decomposeDict
952  );
953 
954  labelList cellToRegion
955  (
956  decomposePtr().decompose
957  (
958  mesh,
959  mesh.cellCentres()
960  )
961  );
962 
963  // Restore state
964  UPstream::parRun() = oldParRun;
965 
966  // For debugging: write out region
968  (
969  mesh,
970  "cellDist",
971  cellToRegion
972  )().write();
973 
974  Info<< nl << "Written decomposition as volScalarField to "
975  << "cellDist for use in postprocessing."
976  << nl << endl;
977 
978 
979  cellOrder = regionRenumber(renumberPtr(), mesh, cellToRegion);
980 
981  // Determine new to old face order with new cell numbering
982  faceOrder = getRegionFaceOrder
983  (
984  mesh,
985  cellOrder,
986  cellToRegion
987  );
988  }
989  else
990  {
991  // Detemines sorted back to original cell ordering
992  cellOrder = renumberPtr().renumber
993  (
994  mesh,
995  mesh.cellCentres()
996  );
997 
998  if (sortCoupledFaceCells)
999  {
1000  // Change order so all coupled patch faceCells are at the end.
1001  const polyBoundaryMesh& pbm = mesh.boundaryMesh();
1002 
1003  // Collect all boundary cells on coupled patches
1004  label nBndCells = 0;
1005  forAll(pbm, patchI)
1006  {
1007  if (pbm[patchI].coupled())
1008  {
1009  nBndCells += pbm[patchI].size();
1010  }
1011  }
1012 
1013  labelList reverseCellOrder = invert(mesh.nCells(), cellOrder);
1014 
1015  labelList bndCells(nBndCells);
1016  labelList bndCellMap(nBndCells);
1017  nBndCells = 0;
1018  forAll(pbm, patchI)
1019  {
1020  if (pbm[patchI].coupled())
1021  {
1022  const labelUList& faceCells = pbm[patchI].faceCells();
1023  forAll(faceCells, i)
1024  {
1025  label cellI = faceCells[i];
1026 
1027  if (reverseCellOrder[cellI] != -1)
1028  {
1029  bndCells[nBndCells] = cellI;
1030  bndCellMap[nBndCells++] = reverseCellOrder[cellI];
1031  reverseCellOrder[cellI] = -1;
1032  }
1033  }
1034  }
1035  }
1036  bndCells.setSize(nBndCells);
1037  bndCellMap.setSize(nBndCells);
1038 
1039  // Sort
1040  labelList order;
1041  sortedOrder(bndCellMap, order);
1042 
1043  // Redo newReverseCellOrder
1044  labelList newReverseCellOrder(mesh.nCells(), -1);
1045 
1046  label sortedI = mesh.nCells();
1047  forAllReverse(order, i)
1048  {
1049  label origCellI = bndCells[order[i]];
1050  newReverseCellOrder[origCellI] = --sortedI;
1051  }
1052 
1053  Info<< "Ordered all " << nBndCells << " cells with a coupled face"
1054  << " to the end of the cell list, starting at " << sortedI
1055  << endl;
1056 
1057  // Compact
1058  sortedI = 0;
1059  forAll(cellOrder, newCellI)
1060  {
1061  label origCellI = cellOrder[newCellI];
1062  if (newReverseCellOrder[origCellI] == -1)
1063  {
1064  newReverseCellOrder[origCellI] = sortedI++;
1065  }
1066  }
1067 
1068  // Update sorted back to original (unsorted) map
1069  cellOrder = invert(mesh.nCells(), newReverseCellOrder);
1070  }
1071 
1072 
1073  // Determine new to old face order with new cell numbering
1074  faceOrder = getFaceOrder
1075  (
1076  mesh,
1077  cellOrder // New to old cell
1078  );
1079  }
1080 
1081 
1082  if (!overwrite)
1083  {
1084  runTime++;
1085  }
1086 
1087 
1088  // Change the mesh.
1089  autoPtr<mapPolyMesh> map = reorderMesh(mesh, cellOrder, faceOrder);
1090 
1091 
1092  if (orderPoints)
1093  {
1094  polyTopoChange meshMod(mesh);
1095  autoPtr<mapPolyMesh> pointOrderMap = meshMod.changeMesh
1096  (
1097  mesh,
1098  false, // inflate
1099  true, // syncParallel
1100  false, // orderCells
1101  orderPoints // orderPoints
1102  );
1103 
1104  // Combine point reordering into map.
1105  const_cast<labelList&>(map().pointMap()) = UIndirectList<label>
1106  (
1107  map().pointMap(),
1108  pointOrderMap().pointMap()
1109  )();
1110 
1112  (
1113  pointOrderMap().reversePointMap(),
1114  const_cast<labelList&>(map().reversePointMap())
1115  );
1116  }
1117 
1118 
1119  // Update fields
1120  mesh.updateMesh(map);
1121 
1122  // Update proc maps
1123  if
1124  (
1125  cellProcAddressing.headerOk()
1126  && cellProcAddressing.size() == mesh.nCells()
1127  )
1128  {
1129  Info<< "Renumbering processor cell decomposition map "
1130  << cellProcAddressing.name() << endl;
1131 
1132  cellProcAddressing = labelList
1133  (
1134  UIndirectList<label>(cellProcAddressing, map().cellMap())
1135  );
1136  }
1137  if
1138  (
1139  faceProcAddressing.headerOk()
1140  && faceProcAddressing.size() == mesh.nFaces()
1141  )
1142  {
1143  Info<< "Renumbering processor face decomposition map "
1144  << faceProcAddressing.name() << endl;
1145 
1147  (
1149  );
1150 
1151  // Detect any flips.
1152  const labelHashSet& fff = map().flipFaceFlux();
1153  forAllConstIter(labelHashSet, fff, iter)
1154  {
1155  label faceI = iter.key();
1156  label masterFaceI = faceProcAddressing[faceI];
1157 
1158  faceProcAddressing[faceI] = -masterFaceI;
1159 
1160  if (masterFaceI == 0)
1161  {
1163  << " masterFaceI:" << masterFaceI << exit(FatalError);
1164  }
1165  }
1166  }
1167  if
1168  (
1169  pointProcAddressing.headerOk()
1170  && pointProcAddressing.size() == mesh.nPoints()
1171  )
1172  {
1173  Info<< "Renumbering processor point decomposition map "
1174  << pointProcAddressing.name() << endl;
1175 
1176  pointProcAddressing = labelList
1177  (
1178  UIndirectList<label>(pointProcAddressing, map().pointMap())
1179  );
1180  }
1181 
1182 
1183  // Move mesh (since morphing might not do this)
1184  if (map().hasMotionPoints())
1185  {
1186  mesh.movePoints(map().preMotionPoints());
1187  }
1188 
1189 
1190  {
1191  label band;
1192  scalar profile;
1193  scalar sumSqrIntersect;
1194  getBand
1195  (
1196  doFrontWidth,
1197  mesh.nCells(),
1198  mesh.faceOwner(),
1199  mesh.faceNeighbour(),
1200  band,
1201  profile,
1202  sumSqrIntersect
1203  );
1204  reduce(band, maxOp<label>());
1205  reduce(profile, sumOp<scalar>());
1206  scalar rmsFrontwidth = Foam::sqrt
1207  (
1208  returnReduce
1209  (
1210  sumSqrIntersect,
1211  sumOp<scalar>()
1213  );
1214 
1215  Info<< "After renumbering :" << nl
1216  << " band : " << band << nl
1217  << " profile : " << profile << nl;
1218 
1219  if (doFrontWidth)
1220  {
1221 
1222  Info<< " rms frontwidth : " << rmsFrontwidth << nl;
1223  }
1224 
1225  Info<< endl;
1226  }
1227 
1228  if (orderPoints)
1229  {
1230  // Force edge calculation (since only reason that points would need to
1231  // be sorted)
1232  (void)mesh.edges();
1233 
1234  label nTotPoints = returnReduce
1235  (
1236  mesh.nPoints(),
1237  sumOp<label>()
1238  );
1239  label nTotIntPoints = returnReduce
1240  (
1242  sumOp<label>()
1243  );
1244 
1245  label nTotEdges = returnReduce
1246  (
1247  mesh.nEdges(),
1248  sumOp<label>()
1249  );
1250  label nTotIntEdges = returnReduce
1251  (
1252  mesh.nInternalEdges(),
1253  sumOp<label>()
1254  );
1255  label nTotInt0Edges = returnReduce
1256  (
1258  sumOp<label>()
1259  );
1260  label nTotInt1Edges = returnReduce
1261  (
1263  sumOp<label>()
1264  );
1265 
1266  Info<< "Points:" << nl
1267  << " total : " << nTotPoints << nl
1268  << " internal: " << nTotIntPoints << nl
1269  << " boundary: " << nTotPoints-nTotIntPoints << nl
1270  << "Edges:" << nl
1271  << " total : " << nTotEdges << nl
1272  << " internal: " << nTotIntEdges << nl
1273  << " internal using 0 boundary points: "
1274  << nTotInt0Edges << nl
1275  << " internal using 1 boundary points: "
1276  << nTotInt1Edges-nTotInt0Edges << nl
1277  << " internal using 2 boundary points: "
1278  << nTotIntEdges-nTotInt1Edges << nl
1279  << " boundary: " << nTotEdges-nTotIntEdges << nl
1280  << endl;
1281  }
1282 
1283 
1284  if (overwrite)
1285  {
1286  mesh.setInstance(oldInstance);
1287  }
1288 
1289  Info<< "Writing mesh to " << mesh.facesInstance() << endl;
1290 
1291  mesh.write();
1292 
1293  if (cellProcAddressing.headerOk())
1294  {
1295  cellProcAddressing.instance() = mesh.facesInstance();
1296  if (cellProcAddressing.size() == mesh.nCells())
1297  {
1298  cellProcAddressing.write();
1299  }
1300  else
1301  {
1302  // procAddressing file no longer valid. Delete it.
1303  const fileName fName(cellProcAddressing.filePath());
1304  if (fName.size())
1305  {
1306  Info<< "Deleting inconsistent processor cell decomposition"
1307  << " map " << fName << endl;
1308  rm(fName);
1309  }
1310  }
1311  }
1312 
1313  if (faceProcAddressing.headerOk())
1314  {
1315  faceProcAddressing.instance() = mesh.facesInstance();
1316  if (faceProcAddressing.size() == mesh.nFaces())
1317  {
1318  faceProcAddressing.write();
1319  }
1320  else
1321  {
1322  const fileName fName(faceProcAddressing.filePath());
1323  if (fName.size())
1324  {
1325  Info<< "Deleting inconsistent processor face decomposition"
1326  << " map " << fName << endl;
1327  rm(fName);
1328  }
1329  }
1330  }
1331 
1332  if (pointProcAddressing.headerOk())
1333  {
1334  pointProcAddressing.instance() = mesh.facesInstance();
1335  if (pointProcAddressing.size() == mesh.nPoints())
1336  {
1337  pointProcAddressing.write();
1338  }
1339  else
1340  {
1341  const fileName fName(pointProcAddressing.filePath());
1342  if (fName.size())
1343  {
1344  Info<< "Deleting inconsistent processor point decomposition"
1345  << " map " << fName << endl;
1346  rm(fName);
1347  }
1348  }
1349  }
1350 
1351  if (boundaryProcAddressing.headerOk())
1352  {
1353  boundaryProcAddressing.instance() = mesh.facesInstance();
1354  if (boundaryProcAddressing.size() == mesh.boundaryMesh().size())
1355  {
1356  boundaryProcAddressing.write();
1357  }
1358  else
1359  {
1360  const fileName fName(boundaryProcAddressing.filePath());
1361  if (fName.size())
1362  {
1363  Info<< "Deleting inconsistent processor patch decomposition"
1364  << " map " << fName << endl;
1365  rm(fName);
1366  }
1367  }
1368  }
1369 
1370  if (writeMaps)
1371  {
1372  // For debugging: write out region
1374  (
1375  mesh,
1376  "origCellID",
1377  map().cellMap()
1378  )().write();
1379 
1381  (
1382  mesh,
1383  "cellID",
1384  identity(mesh.nCells())
1385  )().write();
1386 
1387  Info<< nl << "Written current cellID and origCellID as volScalarField"
1388  << " for use in postprocessing."
1389  << nl << endl;
1390 
1391  labelIOList
1392  (
1393  IOobject
1394  (
1395  "cellMap",
1396  mesh.facesInstance(),
1398  mesh,
1401  false
1402  ),
1403  map().cellMap()
1404  ).write();
1405 
1406  labelIOList
1407  (
1408  IOobject
1409  (
1410  "faceMap",
1411  mesh.facesInstance(),
1413  mesh,
1416  false
1417  ),
1418  map().faceMap()
1419  ).write();
1420 
1421  labelIOList
1422  (
1423  IOobject
1424  (
1425  "pointMap",
1426  mesh.facesInstance(),
1428  mesh,
1431  false
1432  ),
1433  map().pointMap()
1434  ).write();
1435  }
1436 
1437  if (renumberSets)
1438  {
1439  forAll(cellSets, i)
1440  {
1441  cellSets[i].updateMesh(map());
1442  cellSets[i].instance() = mesh.facesInstance();
1443  cellSets[i].write();
1444  }
1445  forAll(faceSets, i)
1446  {
1447  faceSets[i].updateMesh(map());
1448  faceSets[i].instance() = mesh.facesInstance();
1449  faceSets[i].write();
1450  }
1451  forAll(pointSets, i)
1452  {
1453  pointSets[i].updateMesh(map());
1454  pointSets[i].instance() = mesh.facesInstance();
1455  pointSets[i].write();
1456  }
1457  }
1458 
1459  Info<< "\nEnd\n" << endl;
1460 
1461  return 0;
1462 }
1463 
1464 
1465 // ************************************************************************* //
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::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::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::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::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
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::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::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
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
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::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:333
createTime.H
Foam::argList::optionFound
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
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
pointSet.H
Foam::primitiveMesh
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:79