Go to the documentation of this file.
53 #ifdef FOAM_USE_ZOLTAN
84 zeroGradientFvPatchScalarField::typeName
91 fld[cellI] = elems[cellI];
105 label diff = neighbour[faceI] - owner[faceI];
119 const bool calculateIntersect,
125 scalar& sumSqrIntersect
133 label own = owner[faceI];
134 label nei = neighbour[faceI];
138 cellBandwidth[nei] =
max(cellBandwidth[nei],
diff);
141 bandwidth =
max(cellBandwidth);
145 forAll(cellBandwidth, cellI)
147 profile += 1.0*cellBandwidth[cellI];
150 sumSqrIntersect = 0.0;
151 if (calculateIntersect)
155 for (
label colI = cellI-cellBandwidth[cellI]; colI <= cellI; colI++)
157 nIntersect[colI] += 1.0;
182 forAll(cellOrder, newCellI)
184 label oldCellI = cellOrder[newCellI];
193 label faceI = cFaces[i];
199 if (nbrCellI == newCellI)
204 if (newCellI < nbrCellI)
227 label index = order[i];
228 if (nbr[index] != -1)
230 oldToNewFace[cFaces[index]] = newFaceI++;
238 oldToNewFace[faceI] = faceI;
243 forAll(oldToNewFace, faceI)
245 if (oldToNewFace[faceI] == -1)
248 <<
"Did not determine new position" <<
" for face " << faceI
272 label prevRegion = -1;
274 forAll(cellOrder, newCellI)
276 label oldCellI = cellOrder[newCellI];
278 if (cellToRegion[oldCellI] != prevRegion)
280 prevRegion = cellToRegion[oldCellI];
289 label faceI = cFaces[i];
295 if (nbrCellI == newCellI)
300 if (cellToRegion[oldCellI] != cellToRegion[cellOrder[nbrCellI]])
305 else if (newCellI < nbrCellI)
329 oldToNewFace[cFaces[nbr.
indices()[i]]] = newFaceI++;
335 label nRegions =
max(cellToRegion)+1;
345 if (ownRegion != neiRegion)
348 min(ownRegion, neiRegion)*nRegions
349 +
max(ownRegion, neiRegion);
358 label key = sortKey[i];
370 oldToNewFace[sortKey.
indices()[i]] = newFaceI++;
377 oldToNewFace[faceI] = faceI;
382 forAll(oldToNewFace, faceI)
384 if (oldToNewFace[faceI] == -1)
387 <<
"Did not determine new position"
388 <<
" for face " << faceI
430 forAll(newNeighbour, faceI)
432 label own = newOwner[faceI];
433 label nei = newNeighbour[faceI];
437 newFaces[faceI].flip();
438 Swap(newOwner[faceI], newNeighbour[faceI]);
439 flipFaceFlux.
insert(faceI);
451 patchSizes[patchI] =
patches[patchI].size();
452 patchStarts[patchI] =
patches[patchI].start();
453 oldPatchNMeshPoints[patchI] =
patches[patchI].nPoints();
480 label oldFaceI = fZone[i];
481 newAddressing[i] = reverseFaceOrder[oldFaceI];
482 if (flipFaceFlux.
found(newAddressing[i]))
484 newFlipMap[i] = !fZone.
flipMap()[i];
488 newFlipMap[i] = fZone.
flipMap()[i];
561 Info<<
"Determining cell order:" <<
endl;
565 label nRegions =
max(cellToRegion)+1;
571 forAll(regionToCells, regionI)
573 Info<<
" region " << regionI <<
" starts at " << cellI <<
endl;
598 cellOrder[cellI++] = cellMap[subCellOrder[i]];
609 int main(
int argc,
char *argv[])
613 "Renumber mesh to minimise bandwidth"
623 "calculate the rms of the frontwidth"
629 "write cellMap, faceMap, pointMap in polyMesh/"
634 runTime.functionObjects().off();
639 #ifdef FOAM_USE_ZOLTAN
640 Info<<
"renumberMesh built with zoltan support." <<
nl <<
endl;
641 (void)zoltanRenumber::typeName;
662 scalar sumSqrIntersect;
686 <<
"Before renumbering :" <<
nl
687 <<
" band : " << band <<
nl
688 <<
" profile : " << profile <<
nl;
692 Info<<
" rms frontwidth : " << rmsFrontwidth <<
nl;
697 bool sortCoupledFaceCells =
false;
699 bool orderPoints =
false;
701 bool renumberSets =
true;
721 "sortCoupledFaceCells",
724 if (sortCoupledFaceCells)
726 Info<<
"Sorting cells on coupled boundaries to be last." <<
nl
733 Info<<
"Ordering cells into regions of size " << blockSize
734 <<
" (using decomposition);"
735 <<
" ordering faces into region-internal and region-external."
738 if (blockSize < 0 || blockSize >=
mesh.
nCells())
741 <<
"Block size " << blockSize
742 <<
" should be positive integer"
743 <<
" and less than the number of cells in the mesh."
751 Info<<
"Ordering points into internal and boundary points." <<
nl
755 renumberDict.
lookup(
"writeMaps") >> writeMaps;
758 Info<<
"Writing renumber maps (new to old) to polyMesh." <<
nl
766 Info<<
"Using default renumberMethod." <<
nl <<
endl;
771 Info<<
"Selecting renumberMethod " << renumberPtr().type() <<
nl <<
endl;
780 "cellProcAddressing",
793 "faceProcAddressing",
805 "pointProcAddressing",
817 "boundaryProcAddressing",
921 Info<<
"Reading pointSets:" <<
endl;
945 Info<<
"nBlocks = " << nBlocks <<
endl;
948 dictionary decomposeDict(renumberDictPtr().subDict(
"blockCoeffs"));
949 decomposeDict.
set(
"numberOfSubdomains", nBlocks);
961 decomposePtr().decompose
979 Info<<
nl <<
"Written decomposition as volScalarField to "
980 <<
"cellDist for use in postprocessing."
997 cellOrder = renumberPtr().renumber
1003 if (sortCoupledFaceCells)
1009 label nBndCells = 0;
1012 if (pbm[patchI].coupled())
1014 nBndCells += pbm[patchI].
size();
1025 if (pbm[patchI].coupled())
1027 const labelUList& faceCells = pbm[patchI].faceCells();
1030 label cellI = faceCells[i];
1032 if (reverseCellOrder[cellI] != -1)
1034 bndCells[nBndCells] = cellI;
1035 bndCellMap[nBndCells++] = reverseCellOrder[cellI];
1036 reverseCellOrder[cellI] = -1;
1042 bndCellMap.
setSize(nBndCells);
1054 label origCellI = bndCells[order[i]];
1055 newReverseCellOrder[origCellI] = --sortedI;
1058 Info<<
"Ordered all " << nBndCells <<
" cells with a coupled face"
1059 <<
" to the end of the cell list, starting at " << sortedI
1064 forAll(cellOrder, newCellI)
1066 label origCellI = cellOrder[newCellI];
1067 if (newReverseCellOrder[origCellI] == -1)
1069 newReverseCellOrder[origCellI] = sortedI++;
1113 pointOrderMap().pointMap()
1118 pointOrderMap().reversePointMap(),
1119 const_cast<labelList&
>(map().reversePointMap())
1134 Info<<
"Renumbering processor cell decomposition map "
1135 << cellProcAddressing.
name() <<
endl;
1148 Info<<
"Renumbering processor face decomposition map "
1160 label faceI = iter.key();
1165 if (masterFaceI == 0)
1178 Info<<
"Renumbering processor point decomposition map "
1179 << pointProcAddressing.
name() <<
endl;
1189 if (map().hasMotionPoints())
1198 scalar sumSqrIntersect;
1220 Info<<
"After renumbering :" <<
nl
1221 <<
" band : " << band <<
nl
1222 <<
" profile : " << profile <<
nl;
1227 Info<<
" rms frontwidth : " << rmsFrontwidth <<
nl;
1272 <<
" total : " << nTotPoints <<
nl
1273 <<
" internal: " << nTotIntPoints <<
nl
1274 <<
" boundary: " << nTotPoints-nTotIntPoints <<
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
1303 cellProcAddressing.
write();
1311 Info<<
"Deleting inconsistent processor cell decomposition"
1312 <<
" map " << fName <<
endl;
1330 Info<<
"Deleting inconsistent processor face decomposition"
1331 <<
" map " << fName <<
endl;
1337 if (pointProcAddressing.
headerOk())
1342 pointProcAddressing.
write();
1349 Info<<
"Deleting inconsistent processor point decomposition"
1350 <<
" map " << fName <<
endl;
1356 if (boundaryProcAddressing.
headerOk())
1361 boundaryProcAddressing.
write();
1368 Info<<
"Deleting inconsistent processor patch decomposition"
1369 <<
" map " << fName <<
endl;
1392 Info<<
nl <<
"Written current cellID and origCellID as volScalarField"
1393 <<
" for use in postprocessing."
1458 pointSets[i].updateMesh(map());
1460 pointSets[i].write();
void setLargeCellSubset(const labelList ®ion, const label currentRegion, const label patchID=-1, const bool syncCouples=true)
Set the subset from all cells with region == currentRegion.
vectorField pointField
pointField is a vectorField.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
fileName filePath() const
Return complete path + object name if the file exists.
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
T & last()
Return reference to the last element of the list.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
A class for handling words, derived from string.
A class for handling file names.
List< label > labelList
A List of labels.
ListType reorder(const labelUList &oldToNew, const ListType &)
Reorder the elements (indices, not values) of a list.
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
void append(T *)
Append an element at the end of the list.
label nInternalEdges() const
Internal edges using 0,1 or 2 boundary points.
void sort()
(stable) sort the list (if changed after construction time)
static void addNote(const string &)
Add extra notes for the usage information.
#define forAll(list, i)
Loop across all elements in list.
A class for managing temporary objects.
Post-processing mesh subset tool. Given the original mesh and the list of selected cells,...
autoPtr< mapPolyMesh > reorderMesh(polyMesh &mesh, const labelList &cellOrder, const labelList &faceOrder)
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.
labelListList invertOneToMany(const label len, const labelUList &)
Invert one-to-many map. Unmapped elements will be size 0.
static autoPtr< decompositionMethod > New(const dictionary &decompositionDict)
Return a reference to the selected decomposition method.
bool rm(const fileName &)
Remove a file, returning true if successful otherwise false.
label nInternal0Edges() const
Internal edges (i.e. not on boundary face) using.
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.
static const label labelMax
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
static bool & parRun()
Is this a parallel run?
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
const cellList & cells() const
Direct mesh changes based on v1.3 polyTopoChange syntax.
label nTotalCells() const
Return total number of cells in decomposed mesh.
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
const fileName & facesInstance() const
Return the current instance directory for faces.
const cellZoneMesh & cellZones() const
Return cell zone mesh.
virtual bool write() const
Write using setting from DB.
virtual void resetAddressing(const labelUList &, const boolList &)
Reset addressing and flip map (clearing demand-driven data)
virtual tmp< scalarField > movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Ostream & endl(Ostream &os)
Add newline and flush stream.
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
labelList getFaceOrder(const primitiveMesh &mesh, const labelList &cellOrder)
const fileName & instance() const
IOList< label > labelIOList
Label container classes.
A simple container for copying or transferring objects of type <T>.
Mesh consisting of general polyhedral cells.
static const pointMesh & New(const polyMesh &mesh)
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
const faceZoneMesh & faceZones() const
Return face zone mesh.
ListType renumber(const labelUList &oldToNew, const ListType &)
Renumber the values (not the indices) of a list.
void inplaceRenumber(const labelUList &oldToNew, ListType &)
Inplace renumber the values of a list.
virtual bool write() const
Write mesh using IO settings from time.
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
bool headerOk()
Read and check header info.
const word dictName("particleTrackDict")
const fileName & pointsInstance() const
Return the current instance directory for points.
void sortedOrder(const UList< T > &, labelList &order)
Generate the (stable) sort order for the list.
const labelList & cellMap() const
Return cell map.
virtual labelList renumber(const pointField &) const
Return the order in which cells need to be visited, i.e.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Pre-declare SubField and related Field type.
A subset of mesh faces organised as a primitive patch.
virtual const labelList & faceOwner() const
Return face owner.
virtual void updateMesh(const mapPolyMesh &mpm)
Update mesh corresponding to the given map.
const word & name() const
Return name.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
A list of keyword definitions, which are a keyword followed by any number of values (e....
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.
label nInternalFaces() const
#define forAllReverse(list, i)
Reverse loop across all elements in list.
A list that is sorted upon construction or when explicitly requested with the sort() method.
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))
Generic dimensioned Type class.
bool found(const Key &) const
Return true if hashedEntry is found in table.
label getBand(const labelList &owner, const labelList &neighbour)
Mesh data needed to do the Finite Volume discretisation.
label nInternalPoints() const
Points not on boundary.
errorManip< error > abort(error &err)
label nInternal1Edges() const
Internal edges using 0 or 1 boundary point.
A collection of cell labels.
IOobjectList lookupClass(const word &className) const
Return the list for all IOobjects of a given class.
List of IOobjects with searching and retrieving facilities.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Abstract base class for renumbering.
void setInstance(const fileName &)
Set the instance for mesh files.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
void setSize(const label)
Reset size of List.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
List< labelList > labelListList
A List of labelList.
const vectorField & cellCentres() const
virtual const faceList & faces() const
Return raw faces.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Xfer< T > xferMove(T &)
Construct by transferring the contents of the arg.
Helper routine to read fields.
labelList regionRenumber(const renumberMethod &method, const fvMesh &mesh, const labelList &cellToRegion)
tmp< volScalarField > createScalarField(const fvMesh &mesh, const word &name, const labelList &elems)
int main(int argc, char *argv[])
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
dimensionedScalar sqrt(const dimensionedScalar &ds)
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
label size() const
Return the number of elements in the PtrList.
bool insert(const Key &key)
Insert a new entry.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
const fvMesh & subMesh() const
Return reference to subset mesh.
bool optionFound(const word &opt) const
Return true if the named option is found.
const labelList & indices() const
Return the list of sorted indices. Updated every sort.
PtrList< labelIOList > & faceProcAddressing
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
void reset(T *=0)
If object pointer already set, delete object and set to given.
const Time & time() const
Return the top-level database.
A List with indirect addressing.
labelList invert(const label len, const labelUList &)
Invert one-to-one map. Unmapped elements will be -1.
static autoPtr< renumberMethod > New(const dictionary &renumberDict)
Return a reference to the selected renumbering method.
labelList getRegionFaceOrder(const primitiveMesh &mesh, const labelList &cellOrder, const labelList &cellToRegion)
void size(const label)
Override size to be inconsistent with allocated storage.
void clearAddressing()
Clear addressing.
const globalMeshData & globalData() const
Return parallel info.
Generic GeometricField class.
Cuthill-McKee renumbering.
Foam::argList args(argc, argv)
const boolList & flipMap() const
Return face flip map.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
A cell is defined as a list of faces with extra functionality.
IOobject dictIO(dictName, runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED, IOobject::NO_WRITE)
word name(const complex &)
Return a string representation of a complex.
virtual const labelList & faceNeighbour() const
Return face neighbour.
void set(entry *)
Assign a new entry, overwrite any existing entry.
Cell-face mesh analysis engine.