A deep-copy description of an OpenFOAM volume mesh in data structures suitable for VTK UnstructuredGrid, including the possibility of decomposing polyhedral cells into primitive cell types. More...
#include <foamVtuCells.H>


Public Member Functions | |
| vtuCells (const vtuCells &)=delete | |
| No copy construct. | |
| void | operator= (const vtuCells &)=delete |
| No copy assignment. | |
| vtuCells (const enum contentType output=contentType::XML, const bool decompose=false) | |
| Default construct (XML format, no polyhedral decomposition). | |
| vtuCells (const polyMesh &mesh, const enum contentType output=contentType::XML, const bool decompose=false) | |
| Construct from components, create output information immediately. | |
| vtuCells (const vtk::outputOptions opts, const bool decompose=false) | |
| Construct from components. | |
| vtuCells (const polyMesh &mesh, const vtk::outputOptions opts, const bool decompose=false) | |
| Construct from components, create output information immediately. | |
| enum contentType | content () const noexcept |
| The output content type. | |
| bool | decomposeRequested () const noexcept |
| Query the polyhedral decompose requested flag. | |
| bool | empty () const noexcept |
| True if no cellTypes are populated. | |
| label | size () const noexcept |
| The size of populated cellTypes (the number of cells). | |
| bool | useCellMap () const noexcept |
| The cellMap is non-identity for a SUBSET_MESH or when there are decomposed cells. | |
| bool | usePointMap () const noexcept |
| The pointMap is available and non-identity [FUTURE USE]. | |
| bool | is_hdf () const noexcept |
| Test for hdf content type. | |
| void | clear () |
| Reset all sizes to zero. | |
| void | reset (const polyMesh &mesh) |
| Create the geometry using the previously requested output and decomposition types. | |
| void | reset (const polyMesh &mesh, const labelUList &subsetCellsIds) |
| Create the geometry for a mesh subset, using previously requested output and decomposition types. | |
| void | reset (const polyMesh &mesh, const enum contentType output, const bool decompose) |
| Respecify requested output and decomposition type prior to creating the geometry. | |
| void | resetShapes (const UList< cellShape > &shapes) |
| Reset sizing using primitive shapes only (ADVANCED USAGE). | |
| void | renumberCells (const labelUList &mapping) |
| Renumber cell ids to account for subset meshes. | |
| void | renumberPoints (const labelUList &mapping) |
| Renumber point ids to account for subset meshes. | |
| void | addPointCellLabels (const labelUList &cellIds) |
| Define which additional cell-centres are to be used (ADVANCED!). | |
| const List< uint8_t > & | cellTypes () const noexcept |
The cell "types" (legacy: "CELL_TYPES"). | |
| const labelList & | vertLabels () const noexcept |
The cell "connectivity" (legacy: "CELLS"). | |
| const labelList & | vertOffsets () const noexcept |
The "offsets" into the connectivity list vertLabels(). | |
| const labelList & | faceLabels () const noexcept |
The polyhedral "faces" or "FaceConnectivity". | |
| const labelList & | faceOffsets () const noexcept |
The "facesoffset" or the "FaceOffsets". | |
| const labelList & | polyFaceIds () const |
| Face ids per polyhedral cell (HDF only). | |
| const labelList & | polyFaceOffsets () const noexcept |
| Offsets into the polyhedral face ids (HDF only). | |
| const labelList & | addPointCellLabels () const noexcept |
| Additional point addressing (from added point to original cell). | |
| const labelList & | cellMap () const noexcept |
| The mesh cell ids for all cells (regular and decomposed). | |
| const labelUList & | pointMap () const noexcept |
| The mesh point ids [FUTURE USE]. | |
| refPtr< labelList > | vertLabels (label pointOffset) const |
| The "connectivity" (legacy: "CELLS") adjusted by the specified (global) point offset. | |
| refPtr< labelList > | vertOffsets (label beginOffset, bool syncPar=false) const |
| Connectivity (vertices) offsets for each cell adjusted by the specified (global) cell offset. | |
| refPtr< labelList > | faceLabels (label pointOffset) const |
| The "faces" stream/connectivity (non-legacy only) adjusted by the specified (global) point offset. | |
| refPtr< labelList > | faceOffsets (label beginOffset, bool syncPar=false) const |
| The offsets into "faces" (non-legacy only) adjusted by the specified (global) face offset. | |
| refPtr< labelList > | polyFaceIds (label beginOffset) const |
| The "PolyhedronToFaces" ids (HDF only) adjusted by the specified (global) face offset. | |
| refPtr< labelList > | polyFaceOffsets (label beginOffset, bool syncPar=false) const |
| The "PolyhedronOffsets" (HDF only), which are the offsets into the polyhedral face ids, adjusted by the specified (global) cell offset. | |
| template<class PointType = Foam::point> | |
| tmp< Field< PointType > > | points (const polyMesh &mesh) const |
| Return the mesh points, possibly with cell centres for addPointCellLabels(). | |
| Public Member Functions inherited from vtuSizing | |
| vtuSizing () noexcept | |
| Default construct. | |
| vtuSizing (const polyMesh &mesh) | |
| Construct sizing by analyzing the mesh. | |
| vtuSizing (const polyMesh &mesh, const bool decompose) | |
| Construct sizing by analyzing the mesh. | |
| void | reset (const polyMesh &mesh, const bool decompose=false) |
| Reset sizing by analyzing the mesh. | |
| void | reset (const polyMesh &mesh, const labelUList &subsetCellsIds, const bool decompose=false) |
| Reset sizing by analyzing a subset of the mesh. | |
| void | resetShapes (const UList< cellShape > &shapes) |
| Reset sizing using primitive shapes only (ADVANCED USAGE). | |
| void | clear () noexcept |
| Reset all sizes to zero. | |
| bool | decompose () const noexcept |
| Query the decompose flag (normally off). | |
| bool | manifold () const noexcept |
| Manifold mesh cells detected? Globally consistent quantity. | |
| selectionModeType | selectionMode () const noexcept |
| Query how the mesh cells have been selected or defined. | |
| bool | hasPolyCells () const noexcept |
| Has polyhedral cells? | |
| label | nCells () const noexcept |
| Number of cells for the mesh. | |
| label | nPoints () const noexcept |
| Number of points for the mesh. | |
| label | nVertLabels () const noexcept |
| Number of vertex labels for the mesh. | |
| label | nFaceLabels () const noexcept |
| Number of polyhedral face labels for the mesh, without any prefixing information. | |
| label | nFaceLabels (contentType output) const |
| Number of polyhedral face labels for the mesh with format-dependent prefixing information. | |
| label | nCellsPoly () const noexcept |
| Number of output polyhedral cells for the mesh. | |
| label | nFacesPoly () const noexcept |
| Number of (non-unique) faces used by polyhedral cells. | |
| label | nVertPoly () const noexcept |
| Number of vertex labels for polyhedral cells of the mesh. | |
| label | nAddCells () const noexcept |
| Number of additional (decomposed) cells for the mesh. | |
| label | nAddPoints () const noexcept |
| Number of additional (decomposed) points for the mesh. | |
| label | nAddVerts () const noexcept |
| Number of additional (decomposed) vertices for the mesh. | |
| label | nFieldCells () const noexcept |
| Number of field cells = nCells + nAddCells. | |
| label | nFieldPoints () const noexcept |
| Number of field points = nPoints + nAddPoints. | |
| void | setNumPoints (label n) noexcept |
| Alter number of mesh points (ADVANCED USAGE). | |
| void | setNumAddPoints (label n) noexcept |
| Alter number of additional (cell-centre) points (ADVANCED USAGE). | |
| label | sizeOf (const enum contentType output, const enum slotType slot) const noexcept |
| Return the required size for the storage slot. | |
| template<vtuSizing::slotType slot> | |
| label | sizeOf (contentType output) const noexcept |
| Return the required size for the storage slot. | |
| void | populateLegacy (const polyMesh &mesh, UList< uint8_t > &cellTypes, labelUList &connectivity, foamVtkMeshMaps &maps) const |
| Populate lists for Legacy output. | |
| void | populateXml (const polyMesh &mesh, UList< uint8_t > &cellTypes, labelUList &connectivity, labelUList &offsets, labelUList &faces, labelUList &facesOffsets, foamVtkMeshMaps &maps) const |
| Populate lists for XML output. | |
| void | populateHdf (const polyMesh &mesh, UList< uint8_t > &cellTypes, labelUList &connectivity, labelUList &offsets, labelUList &faces, labelUList &facesOffsets, labelUList &polyFaceIds, labelUList &polyFaceOffsets, foamVtkMeshMaps &maps) const |
| Populate lists for HDF output. | |
| void | populateInternal (const polyMesh &mesh, UList< uint8_t > &cellTypes, UList< int > &connectivity, UList< int > &offsets, UList< int > &faces, UList< int > &facesOffsets, foamVtkMeshMaps &maps, const enum contentType output) const |
| Populate lists for Internal VTK format. | |
| void | populateInternal (const polyMesh &mesh, UList< uint8_t > &cellTypes, UList< int > &connectivity, UList< int > &offsets, UList< int > &faces, UList< int > &facesOffsets, labelUList &cellMap, labelUList &addPointsIds, const enum contentType output) const |
| Populate lists for Internal VTK format. | |
| void | populateInternal (const polyMesh &mesh, UList< uint8_t > &cellTypes, UList< long > &connectivity, UList< long > &offsets, UList< long > &faces, UList< long > &facesOffsets, foamVtkMeshMaps &maps, const enum contentType output) const |
| Populate lists for Internal VTK format. | |
| void | populateInternal (const polyMesh &mesh, UList< uint8_t > &cellTypes, UList< long > &connectivity, UList< long > &offsets, UList< long > &faces, UList< long > &facesOffsets, labelUList &cellMap, labelUList &addPointsIds, const enum contentType output) const |
| Populate lists for Internal VTK format. | |
| void | populateInternal (const polyMesh &mesh, UList< uint8_t > &cellTypes, UList< long long > &connectivity, UList< long long > &offsets, UList< long long > &faces, UList< long long > &facesOffsets, foamVtkMeshMaps &maps, const enum contentType output) const |
| Populate lists for Internal VTK format. | |
| void | populateInternal (const polyMesh &mesh, UList< uint8_t > &cellTypes, UList< long long > &connectivity, UList< long long > &offsets, UList< long long > &faces, UList< long long > &facesOffsets, labelUList &cellMap, labelUList &addPointsIds, const enum contentType output) const |
| Populate lists for Internal VTK format. | |
| void | info (Ostream &os) const |
| Report some information. | |
| bool | operator== (const vtuSizing &rhs) const |
| Test equality. | |
| bool | operator!= (const vtuSizing &rhs) const |
| Test inequality. | |
| label | sizeLegacy () const |
| The calculated size for legacy storage. | |
| label | sizeLegacy (slotType slot) const |
| The calculated size for legacy storage of the specified slot. | |
| label | sizeXml (slotType slot) const |
| The calculated size for xml storage of the specified slot. | |
| label | sizeInternal1 (slotType slot) const |
| The calculated size for vtk-internal storage of the specified slot. | |
| label | sizeInternal2 (slotType slot) const |
| The calculated size for vtk-internal storage of the specified slot. | |
Additional Inherited Members | |
| Public Types inherited from vtuSizing | |
| enum class | contentType : unsigned char { LEGACY , XML , INTERNAL1 , INTERNAL2 , HDF } |
| Types of content that the storage may represent. More... | |
| enum class | slotType : unsigned char { CELLS , CELLS_OFFSETS , FACES , FACES_OFFSETS , POLY_FACEIDS , POLY_FACEIDS_OFFSETS } |
| The possible storage 'slots' that can be used. More... | |
| enum class | selectionModeType : unsigned char { FULL_MESH , SUBSET_MESH , SHAPE_MESH } |
| How the mesh cells have been selected or defined. More... | |
| Static Public Member Functions inherited from vtuSizing | |
| static void | renumberVertLabels (labelUList &connectivity, const label pointOffset, const contentType output) |
| Renumber vertex labels by (global) point offset. | |
| static void | renumberFaceLabels (labelUList &faceLabels, const label pointOffset, const contentType output) |
| Renumber faces stream labels by (global) point offset. | |
| static void | renumberFaceOffsets (labelUList &faceOffsets, const label beginOffset, const contentType output) |
| Renumber face offsets with a specified (global) begin offset. | |
| static labelList | copyVertLabels (const labelUList &connectivity, const label pointOffset, const contentType output) |
| Copy vertex labels with a (global) point offset. | |
| static labelList | copyFaceLabels (const labelUList &faceLabels, const label pointOffset, const contentType output) |
| Copy faces stream labels with a (global) point offset. | |
| static labelList | copyFaceOffsets (const labelUList &faceOffsets, const label beginOffset, const contentType output) |
| Copy face offsets with a specified (global) begin offset. | |
| static void | renumberVertLabelsLegacy (labelUList &connectivity, const label pointOffset) |
| Renumber vertex labels by (global) point offset - legacy format. | |
| static void | renumberVertLabelsXml (labelUList &connectivity, const label pointOffset) |
| Renumber vertex labels by (global) point offset - XML format. | |
| static void | renumberFaceLabelsXml (labelUList &faceLabels, const label pointOffset) |
| Renumber faces stream labels by (global) point offset - XML format. | |
| static void | renumberFaceOffsetsXml (labelUList &faceOffsets, const label beginOffset) |
| Renumber face offsets with a specified begin offset - XML format. | |
| static labelList | copyVertLabelsLegacy (const labelUList &connectivity, const label pointOffset) |
| Copy vertex labels with a (global) point offset - legacy format. | |
| static labelList | copyVertLabelsXml (const labelUList &connectivity, const label pointOffset) |
| Copy vertex labels with a (global) point offset - XML format. | |
| static labelList | copyFaceLabelsXml (const labelUList &faceLabels, const label pointOffset) |
| Copy faces stream labels with a (global) point offset - XML format. | |
| static labelList | copyFaceOffsetsXml (const labelUList &faceOffsets, const label beginOffset) |
| Copy face offsets with an offset from previous - XML format. | |
| Protected Member Functions inherited from vtuSizing | |
| void | populateShapesLegacy (const UList< cellShape > &shapes, UList< uint8_t > &cellTypes, labelUList &connectivity, foamVtkMeshMaps &maps) const |
| Reset list for primitive shapes only (ADVANCED USAGE). | |
| void | populateShapesXml (const UList< cellShape > &shapes, UList< uint8_t > &cellTypes, labelUList &connectivity, labelUList &offsets, labelUList &faces, labelUList &facesOffsets, foamVtkMeshMaps &maps) const |
| Reset list for primitive shapes only (ADVANCED USAGE). | |
| Static Protected Member Functions inherited from vtuSizing | |
| static labelList | dummyFaceOffsets (const label numCells, const enum contentType output, label beginOffset=0) |
| Return a list of dummy faceOffsets. | |
A deep-copy description of an OpenFOAM volume mesh in data structures suitable for VTK UnstructuredGrid, including the possibility of decomposing polyhedral cells into primitive cell types.
Knowledge of the vtkUnstructuredGrid and the corresponding .vtu xml file-format aids in understanding this class. The class can be used for the VTK xml format, legacy format, as well as a VTK internal representation. The internal representation is somewhat related to the xml format, but not entirely.
Definition at line 73 of file foamVtuCells.H.
|
delete |
No copy construct.
References vtuCells().
Referenced by operator=(), and vtuCells().


|
explicit |
Default construct (XML format, no polyhedral decomposition).
References vtuSizing::decompose(), and vtuSizing::XML.

|
explicit |
Construct from components, create output information immediately.
References vtuSizing::decompose(), mesh, and vtuSizing::XML.

|
explicit |
Construct from components.
Optionally with polyhedral decomposition.
References vtuSizing::decompose().

| vtuCells | ( | const polyMesh & | mesh, |
| const vtk::outputOptions | opts, | ||
| const bool | decompose = false ) |
Construct from components, create output information immediately.
References vtuSizing::decompose(), decomposeRequested(), empty(), is_hdf(), mesh, Foam::noexcept, size(), useCellMap(), and usePointMap().

|
delete |
|
inlinenoexcept |
|
inlinenoexcept |
Query the polyhedral decompose requested flag.
Definition at line 38 of file foamVtuCellsI.H.
References Foam::noexcept.
Referenced by vtuCells().

|
inlinenoexcept |
True if no cellTypes are populated.
Definition at line 44 of file foamVtuCellsI.H.
References Foam::noexcept.
Referenced by vtuCells().

|
inlinenoexcept |
The size of populated cellTypes (the number of cells).
Definition at line 50 of file foamVtuCellsI.H.
References Foam::noexcept.
Referenced by vtuCells().

|
inlinenoexcept |
The cellMap is non-identity for a SUBSET_MESH or when there are decomposed cells.
Definition at line 105 of file foamVtuCellsI.H.
References vtuSizing::nAddCells(), Foam::noexcept, vtuSizing::selectionMode(), and vtuSizing::SUBSET_MESH.
Referenced by vtuCells().


|
inlinenoexcept |
The pointMap is available and non-identity [FUTURE USE].
Definition at line 122 of file foamVtuCellsI.H.
References Foam::noexcept.
Referenced by vtuCells().

|
inlinenoexcept |
Test for hdf content type.
Definition at line 32 of file foamVtuCellsI.H.
References vtuSizing::HDF, and Foam::noexcept.
Referenced by vtuCells().

| void clear | ( | ) |
| void reset | ( | const polyMesh & | mesh | ) |
| void reset | ( | const polyMesh & | mesh, |
| const labelUList & | subsetCellsIds ) |
| void reset | ( | const polyMesh & | mesh, |
| const enum contentType | output, | ||
| const bool | decompose ) |
Respecify requested output and decomposition type prior to creating the geometry.
References vtuSizing::decompose(), mesh, and reset().

Reset sizing using primitive shapes only (ADVANCED USAGE).
Effectively removes any polyhedrals!
References resetShapes().
Referenced by isoSurfaceTopo::isoSurfaceTopo(), and resetShapes().


| void renumberCells | ( | const labelUList & | mapping | ) |
Renumber cell ids to account for subset meshes.
References renumberCells().
Referenced by renumberCells().


| void renumberPoints | ( | const labelUList & | mapping | ) |
Renumber point ids to account for subset meshes.
References renumberPoints().
Referenced by renumberPoints().


| void addPointCellLabels | ( | const labelUList & | cellIds | ) |
Define which additional cell-centres are to be used (ADVANCED!).
References addPointCellLabels(), cellTypes(), faceLabels(), faceOffsets(), Foam::noexcept, vertLabels(), and vertOffsets().
Referenced by addPointCellLabels(), isoSurfaceTopo::isoSurfaceTopo(), and polyFaceIds().


|
inlinenoexcept |
The cell "types" (legacy: "CELL_TYPES").
Definition at line 57 of file foamVtuCellsI.H.
References Foam::noexcept.
Referenced by addPointCellLabels().

|
inlinenoexcept |
The cell "connectivity" (legacy: "CELLS").
Definition at line 64 of file foamVtuCellsI.H.
References Foam::noexcept.
Referenced by addPointCellLabels(), and vertLabels().

|
inlinenoexcept |
The "offsets" into the connectivity list vertLabels().
The content depends on the format:
Definition at line 71 of file foamVtuCellsI.H.
References Foam::noexcept.
Referenced by addPointCellLabels(), and vertOffsets().

|
inlinenoexcept |
The polyhedral "faces" or "FaceConnectivity".
The content depends on the format:
Definition at line 78 of file foamVtuCellsI.H.
References Foam::noexcept.
Referenced by addPointCellLabels(), and faceLabels().

|
inlinenoexcept |
The "facesoffset" or the "FaceOffsets".
The content depends on the format:
end offset into faceLabels(),begin offset into faceLabels(), with -1 for primitive typesDefinition at line 85 of file foamVtuCellsI.H.
References Foam::noexcept.
Referenced by addPointCellLabels(), and faceOffsets().

| const labelList & polyFaceIds | ( | ) | const |
Face ids per polyhedral cell (HDF only).
Possibly demand-driven content
References addPointCellLabels(), cellMap(), Foam::noexcept, pointMap(), polyFaceIds(), and polyFaceOffsets().
Referenced by polyFaceIds(), and polyFaceIds().


|
inlinenoexcept |
Offsets into the polyhedral face ids (HDF only).
Definition at line 92 of file foamVtuCellsI.H.
References Foam::noexcept.
Referenced by polyFaceIds(), and polyFaceOffsets().

|
inlinenoexcept |
Additional point addressing (from added point to original cell).
Definition at line 99 of file foamVtuCellsI.H.
References Foam::noexcept.
|
inlinenoexcept |
The mesh cell ids for all cells (regular and decomposed).
For a volume field these would be the cell ids to extract from the field.
Definition at line 116 of file foamVtuCellsI.H.
References Foam::noexcept.
Referenced by polyFaceIds().

|
inlinenoexcept |
The mesh point ids [FUTURE USE].
Definition at line 129 of file foamVtuCellsI.H.
References Foam::noexcept.
Referenced by polyFaceIds().

The "connectivity" (legacy: "CELLS") adjusted by the specified (global) point offset.
References vertLabels().

Connectivity (vertices) offsets for each cell adjusted by the specified (global) cell offset.
| beginOffset | For formats with begin/end, can specify -1 for an automatic calculation (in parallel) |
| syncPar | For begin/end offsets for parallel concatenation |
References vertOffsets().

The "faces" stream/connectivity (non-legacy only) adjusted by the specified (global) point offset.
References faceLabels().

The offsets into "faces" (non-legacy only) adjusted by the specified (global) face offset.
Returns a dummy set of offsets if there are no "faces"
| beginOffset | For formats with begin/end, can specify -1 for an automatic calculation (in parallel) |
| syncPar | For begin/end offsets for parallel concatenation |
References faceOffsets().

The "PolyhedronToFaces" ids (HDF only) adjusted by the specified (global) face offset.
References polyFaceIds().

The "PolyhedronOffsets" (HDF only), which are the offsets into the polyhedral face ids, adjusted by the specified (global) cell offset.
| beginOffset | For formats with begin/end, can specify -1 for an automatic calculation (in parallel) |
| syncPar | For begin/end offsets for parallel concatenation |
References polyFaceOffsets().

| tmp< Field< PointType > > points | ( | const polyMesh & | mesh | ) | const |
Return the mesh points, possibly with cell centres for addPointCellLabels().
References mesh, and points().
Referenced by points().

