Go to the documentation of this file.
44 Info <<
"Renumbering the mesh" <<
endl;
58 label cellInOrder = 0;
69 nextCell.
append(currentCell);
78 while( nextCell.
size() > 0 )
82 if( !visited[currentCell] )
84 visited[currentCell] =
true;
87 newOrder[cellInOrder] = currentCell;
93 const label nei = cellCells(currentCell, nI);
118 newCells[cellI].
transfer(oldCells[newOrder[cellI]]);
120 reverseOrder[newOrder[cellI]] = cellI;
130 const cell&
c = newCells[cellI];
134 --reverseFaceOrder[
c[faceI]];
139 label nMarkedFaces = 0;
149 const cell&
c = newCells[cellI];
154 label nNeighbours(0);
158 if( reverseFaceOrder[
c[faceI]] == -2 )
161 if( cellI == reverseOrder[oldOwner[
c[faceI]]] )
163 neiCells[faceI] = reverseOrder[oldNeighbour[
c[faceI]]];
165 else if( cellI == reverseOrder[oldNeighbour[
c[faceI]]] )
167 neiCells[faceI] = reverseOrder[oldOwner[
c[faceI]]];
179 for(
label i=0;i<nNeighbours;++i)
187 if( (neiCells[ncI] > -1) && (neiCells[ncI] < minNei) )
190 minNei = neiCells[ncI];
197 reverseFaceOrder[
c[nextNei]] = nMarkedFaces;
200 neiCells[nextNei] = -1;
208 "void polyMeshGenModifier::renumberedMesh() const"
209 ) <<
"Error in internal face insertion"
216 forAll(reverseFaceOrder, faceI)
218 if( reverseFaceOrder[faceI] < 0 )
220 reverseFaceOrder[faceI] = nMarkedFaces;
231 faceOrder[reverseFaceOrder[faceI]] = faceI;
237 cell&
c = newCells[cellI];
241 c[fI] = reverseFaceOrder[
c[fI]];
250 newFaces[faceI].
transfer(oldFaces[faceOrder[faceI]]);
255 forAll(oldNeighbour, faceI)
257 const label oldFaceI = faceOrder[faceI];
258 if( oldNeighbour[oldFaceI] < 0 )
263 reverseOrder[oldNeighbour[oldFaceI]]
264 < reverseOrder[oldOwner[oldFaceI]]
267 newFaces[faceI] = newFaces[faceI].reverseFace();
273 oldCells[cellI].
transfer(newCells[cellI]);
275 oldFaces[faceI].
transfer(newFaces[faceI]);
282 Info <<
"Finished renumbering the mesh" <<
endl;
const labelList & neighbour() const
const labelList & owner() const
owner and neighbour cells for faces
const polyMeshGenAddressing & addressingData() const
addressing which may be needed
void append(const T &e)
Append an element at the end of the list.
#define forAll(list, i)
Loop across all elements in list.
Template functions to aid in the implementation of demand driven data.
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Ostream & endl(Ostream &os)
Add newline and flush stream.
void clearOut() const
clear all pointer data
label size() const
Size of the active part of the list.
polyMeshGen & mesh_
reference to the mesh
void transfer(cellList &)
void updateFaceSubsets(const ListType &)
cellListPMG & cellsAccess()
access to cells
void updateCellSubsets(const ListType &)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
void clearOut()
clear out unnecessary data (pointFacesPtr_);
const cellListPMG & cells() const
access to cells
label size() const
Returns the number of rows.
#define forAllRow(graph, rowI, index)
errorManip< error > abort(error &err)
void renumberMesh()
reorder the cells and faces to reduce the matrix bandwidth
label size() const
return the number of used elements
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
const dimensionedScalar c
Speed of light in a vacuum.
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
faceListPMG & facesAccess()
access to mesh faces
void size(const label)
Override size to be inconsistent with allocated storage.
A cell is defined as a list of faces with extra functionality.
void transfer(faceList &)
const VRWGraph & cellCells() const
label size() const
return the number of used elements