VTK  9.5.20250731
vtkUnstructuredGrid.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-License-Identifier: BSD-3-Clause
131#ifndef vtkUnstructuredGrid_h
132#define vtkUnstructuredGrid_h
133
134#include "vtkAbstractCellLinks.h" // For vtkAbstractCellLinks
135#include "vtkCellArray.h" // inline GetCellPoints()
136#include "vtkCommonDataModelModule.h" // For export macro
137#include "vtkDeprecation.h" // VTK_DEPRECATED_IN_9_6_0()
138#include "vtkIdTypeArray.h" // inline GetCellPoints()
139#include "vtkSmartPointer.h" // for smart pointer
141#include "vtkWrappingHints.h" // For VTK_MARSHALMANUAL
142
143VTK_ABI_NAMESPACE_BEGIN
144class vtkCellArray;
145class vtkIdList;
146class vtkIdTypeArray;
148class vtkIdTypeArray;
149
150class VTKCOMMONDATAMODEL_EXPORT VTK_MARSHALMANUAL vtkUnstructuredGrid
152{
153public:
159
161
165 void PrintSelf(ostream& os, vtkIndent indent) override;
167
171 int GetDataObjectType() VTK_FUTURE_CONST override { return VTK_UNSTRUCTURED_GRID; }
172
182 bool AllocateEstimate(vtkIdType numCells, vtkIdType maxCellSize)
183 {
184 return this->AllocateExact(numCells, numCells * maxCellSize);
185 }
186
196 bool AllocateExact(vtkIdType numCells, vtkIdType connectivitySize);
197
207 void Allocate(vtkIdType numCells = 1000, int vtkNotUsed(extSize) = 1000) override
208 {
209 this->AllocateExact(numCells, numCells);
210 }
211
213
216 void Reset();
217 void CopyStructure(vtkDataSet* ds) override;
220 vtkCell* GetCell(vtkIdType cellId) override;
221 void GetCell(vtkIdType cellId, vtkGenericCell* cell) override;
222 void GetCellBounds(vtkIdType cellId, double bounds[6]) override;
223 void GetCellPoints(vtkIdType cellId, vtkIdList* ptIds) override;
224 void GetPointCells(vtkIdType ptId, vtkIdList* cellIds) override;
227
231 int GetCellType(vtkIdType cellId) override;
232
237
249 void GetCellTypes(vtkCellTypes* types) override;
250
264
277 void GetCellPoints(vtkIdType cellId, vtkIdType& npts, vtkIdType const*& pts)
278 {
279 this->Connectivity->GetCellAtId(cellId, npts, pts);
280 }
281
298 vtkIdType cellId, vtkIdType& npts, vtkIdType const*& pts, vtkIdList* ptIds) override
299 {
300 this->Connectivity->GetCellAtId(cellId, npts, pts, ptIds);
301 }
302
304
309 void GetPointCells(vtkIdType ptId, vtkIdType& ncells, vtkIdType*& cells)
310 VTK_SIZEHINT(cells, ncells);
312
320
324 void Squeeze() override;
325
329 void Initialize() override;
330
334 int GetMaxCellSize() override;
335
337
344
350
352
355 vtkSetSmartPointerMacro(Links, vtkAbstractCellLinks);
356 vtkGetSmartPointerMacro(Links, vtkAbstractCellLinks);
358
366 void GetFaceStream(vtkIdType cellId, vtkIdList* ptIds);
367
369
379 void SetCells(int type, vtkCellArray* cells);
380 void SetCells(int* types, vtkCellArray* cells);
383 vtkCellArray* faceLocations, vtkCellArray* faces);
384
385 // This was incorrectly deprecated in v9.4, so it should only been fully removed when removing
386 // VTK_DEPRECATED_IN_9_6_0
387 VTK_DEPRECATED_IN_9_5_0("This function is deprecated, use SetPolyhedralCells")
388 void SetCells(vtkUnsignedCharArray* cellTypes, vtkCellArray* cells, vtkIdTypeArray* faceLocations,
389 vtkIdTypeArray* faces);
391
395 vtkCellArray* GetCells() { return this->Connectivity; }
396
398
404 void GetCellNeighbors(vtkIdType cellId, vtkIdList* ptIds, vtkIdList* cellIds) override
405 {
406 this->GetCellNeighbors(cellId, ptIds->GetNumberOfIds(), ptIds->GetPointer(0), cellIds);
407 }
409 vtkIdType cellId, vtkIdType npts, const vtkIdType* ptIds, vtkIdList* cellIds);
411
421 vtkIdType cellId, unsigned char& cellType, vtkGenericCell* cell) override;
422
424
435 vtkIdType cellId, vtkIdType npts, const vtkIdType* ptIds, vtkIdType& neighborCellId);
436 bool IsCellBoundary(vtkIdType cellId, vtkIdType npts, const vtkIdType* ptIds)
437 {
438 vtkIdType neighborCellId;
439 return this->IsCellBoundary(cellId, npts, ptIds, neighborCellId);
440 }
442
444
448 vtkIdType InsertNextLinkedCell(int type, int npts, const vtkIdType pts[]) VTK_SIZEHINT(pts, npts);
451 void ResizeCellList(vtkIdType ptId, int size);
453
455
458 virtual int GetPiece();
459 virtual int GetNumberOfPieces();
461
465 virtual int GetGhostLevel();
466
475 unsigned long GetActualMemorySize() override;
476
478
481 void ShallowCopy(vtkDataObject* src) override;
482 void DeepCopy(vtkDataObject* src) override;
484
490 void GetIdsOfCellsOfType(int type, vtkIdTypeArray* array) override;
491
495 int IsHomogeneous() override;
496
503
505
511
516
518
525
534
543
556 static void DecomposeAPolyhedronCell(vtkCellArray* polyhedronCellArray, vtkIdType& nCellpts,
557 vtkIdType& nCellfaces, vtkCellArray* cellArray, vtkIdTypeArray* faces);
558
559 static void DecomposeAPolyhedronCell(const vtkIdType* polyhedronCellStream, vtkIdType& nCellpts,
560 vtkIdType& nCellfaces, vtkCellArray* cellArray, vtkIdTypeArray* faces);
561
574 static void DecomposeAPolyhedronCell(vtkIdType nCellFaces, const vtkIdType* inFaceStream,
575 vtkIdType& nCellpts, vtkCellArray* cellArray, vtkIdTypeArray* faces);
576
583 static void ConvertFaceStreamPointIds(vtkIdList* faceStream, vtkIdType* idMap);
584
590 static void ConvertFaceStreamPointIds(vtkIdType nfaces, vtkIdType* faceStream, vtkIdType* idMap);
591
598
599 //====================== Begin Legacy Methods ================================
600
608 VTK_DEPRECATED_IN_9_6_0("CellLocations is not longer used")
609 vtkIdTypeArray* GetCellLocationsArray();
610
612
628 VTK_DEPRECATED_IN_9_6_0("CellLocations is not longer used, use other SetCells methods")
629 void SetCells(
630 vtkUnsignedCharArray* cellTypes, vtkIdTypeArray* cellLocations, vtkCellArray* cells);
631 VTK_DEPRECATED_IN_9_6_0("This function is deprecated, use SetPolyhedralCells")
632 void SetCells(vtkUnsignedCharArray* cellTypes, vtkIdTypeArray* cellLocations, vtkCellArray* cells,
633 vtkIdTypeArray* faceLocations, vtkIdTypeArray* faces);
635
636 //====================== End Legacy Methods ==================================
637
638protected:
641
642 void ReportReferences(vtkGarbageCollector*) override;
643
644 // Points derived from vtkPointSet.
645 // Attribute data (i.e., point and cell data (i.e., scalars, vectors, normals, tcoords)
646 // derived from vtkDataSet.
647
648 // The heart of the data representation. The points are managed by the
649 // superclass vtkPointSet. A cell is defined by its connectivity (i.e., the
650 // point ids that define the cell) and the cell type, represented by the
651 // Connectivity and Types arrays.
652 // Finally, when certain topological information is needed (e.g.,
653 // all the cells that use a point), the cell links array is built.
657
658 // Set of all cell types present in the grid. All entries are unique.
659 vtkSmartPointer<vtkCellTypes> DistinctCellTypes;
660
661 // The DistinctCellTypes is cached, so we keep track of the last time it was
662 // updated so we can compare it to the modified time of the Types array.
663 vtkMTimeType DistinctCellTypesUpdateMTime;
664
673
674 // VTK_DEPRECATED_IN_9_6_0()
675 // Legacy support -- stores the old-style cell array locations.
677
678 vtkIdType InternalInsertNextCell(int type, vtkIdType npts, const vtkIdType ptIds[]) override;
679 vtkIdType InternalInsertNextCell(int type, vtkIdList* ptIds) override;
680 vtkIdType InternalInsertNextCell(
681 int type, vtkIdType npts, const vtkIdType ptIds[], vtkIdType nfaces, const vtkIdType faces[]);
682 vtkIdType InternalInsertNextCell(
683 int type, vtkIdType npts, const vtkIdType pts[], vtkCellArray* faces) override;
684 void InternalReplaceCell(vtkIdType cellId, int npts, const vtkIdType pts[]) override;
685
698 static void DecomposeAPolyhedronCell(const vtkIdType* cellStream, vtkIdType& numCellPts,
699 vtkIdType& nCellFaces, vtkCellArray* cellArray, vtkCellArray* faces);
700
701 static void DecomposeAPolyhedronCell(vtkIdType nCellFaces, const vtkIdType* cellStream,
702 vtkIdType& numCellPts, vtkCellArray* cellArray, vtkCellArray* facesArray);
703
704 static void DecomposeAPolyhedronCell(vtkCellArray* polyhedronCell, vtkIdType& numCellPts,
705 vtkIdType& nCellfaces, vtkCellArray* cellArray, vtkCellArray* faces);
706
707 static void DecomposeAPolyhedronCell(const vtkIdType* cellStream, vtkIdType& numCellPts,
708 vtkIdType& nCellFaces, vtkCellArray* cellArray, vtkCellArray* faces,
709 vtkCellArray* faceLocations);
710
711 static void DecomposeAPolyhedronCell(vtkIdType nCellFaces, const vtkIdType* cellStream,
712 vtkIdType& numCellPts, vtkCellArray* cellArray, vtkCellArray* faces,
713 vtkCellArray* faceLocations);
714
715private:
716 // Hide these from the user and the compiler.
718 void operator=(const vtkUnstructuredGrid&) = delete;
719
720 void Cleanup();
721};
722
723VTK_ABI_NAMESPACE_END
724#endif
object to represent cell connectivity
Efficient cell iterator for vtkDataSet topologies.
object provides direct access to cells in vtkCellArray and type information
abstract class to specify cell behavior
Definition vtkCell.h:129
general representation of visualization data
abstract class to specify dataset behavior
Definition vtkDataSet.h:165
virtual vtkCell * GetCell(vtkIdType cellId)=0
Get cell with cellId such that: 0 <= cellId < NumberOfCells.
virtual void GetCellNeighbors(vtkIdType cellId, vtkIdList *ptIds, vtkIdList *cellIds)
Topological inquiry to get all cells using list of points exclusive of cell specified (e....
Detect and break reference loops.
provides thread-safe access to cells
list of point or cell ids
Definition vtkIdList.h:133
vtkIdType GetNumberOfIds() const noexcept
Return the number of id's in the list.
Definition vtkIdList.h:159
vtkIdType * GetPointer(vtkIdType i)
Get a pointer to a particular data index.
Definition vtkIdList.h:225
dynamic, self-adjusting array of vtkIdType
a simple class to control print indentation
Definition vtkIndent.h:108
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
Hold a reference to a vtkObjectBase instance.
dynamic, self-adjusting array of unsigned char
dataset represents arbitrary combinations of all possible cell types.
dataset represents arbitrary combinations of all possible cell types
void RemoveReferenceToCell(vtkIdType ptId, vtkIdType cellId)
Use these methods only if the dataset has been specified as Editable.
void GetCellTypes(vtkCellTypes *types) override
Get a list of types of cells in a dataset.
void GetCell(vtkIdType cellId, vtkGenericCell *cell) override
Standard vtkDataSet methods; see vtkDataSet.h for documentation.
int GetDataObjectType() VTK_FUTURE_CONST override
Standard vtkDataSet API methods.
void Allocate(vtkIdType numCells=1000, int extSize=1000) override
Method allocates initial storage for the cell connectivity.
void GetCellPoints(vtkIdType cellId, vtkIdType &npts, vtkIdType const *&pts)
A higher-performing variant of the virtual vtkDataSet::GetCellPoints() for unstructured grids.
void SetCells(int *types, vtkCellArray *cells)
Provide cell information to define the dataset.
vtkCell * GetCell(vtkIdType cellId) override
Standard vtkDataSet methods; see vtkDataSet.h for documentation.
void GetCellPoints(vtkIdType cellId, vtkIdList *ptIds) override
Standard vtkDataSet methods; see vtkDataSet.h for documentation.
vtkCellArray * GetPolyhedronFaces()
Get pointer to faces and facelocations for polyhedron cells.
int GetMaxSpatialDimension() override
Get the maximum/minimum spatial dimensionality of the data which is the maximum/minimum dimension of ...
void Squeeze() override
Squeeze all arrays in the grid to conserve memory.
vtkIdType GetCellSize(vtkIdType cellId) override
Get the size of the cell with given cellId.
static void DecomposeAPolyhedronCell(vtkIdType nCellFaces, const vtkIdType *inFaceStream, vtkIdType &nCellpts, vtkCellArray *cellArray, vtkIdTypeArray *faces)
A static method for converting an input polyhedron cell stream of format [nFace0Pts,...
int IsHomogeneous() override
Returns whether cells are all of the same type.
void SetPolyhedralCells(vtkUnsignedCharArray *cellTypes, vtkCellArray *cells, vtkCellArray *faceLocations, vtkCellArray *faces)
Provide cell information to define the dataset.
static void ConvertFaceStreamPointIds(vtkIdList *faceStream, vtkIdType *idMap)
Convert pid in a face stream into idMap[pid].
static vtkUnstructuredGrid * GetData(vtkInformation *info)
Retrieve an instance of this class from an information object.
int GetCellNumberOfFaces(vtkIdType cellId, unsigned char &cellType, vtkGenericCell *cell) override
Get the number of faces of a cell.
void GetCellNeighbors(vtkIdType cellId, vtkIdType npts, const vtkIdType *ptIds, vtkIdList *cellIds)
A topological inquiry to retrieve all of the cells using list of points exclusive of the current cell...
void GetPolyhedronFaces(vtkIdType cellId, vtkCellArray *faces)
Special support for polyhedron.
void ResizeCellList(vtkIdType ptId, int size)
Use these methods only if the dataset has been specified as Editable.
bool AllocateEstimate(vtkIdType numCells, vtkIdType maxCellSize)
Pre-allocate memory in internal data structures.
int InitializeFacesRepresentation(vtkIdType numPrevCells)
Special function used by vtkUnstructuredGridReader.
bool IsCellBoundary(vtkIdType cellId, vtkIdType npts, const vtkIdType *ptIds, vtkIdType &neighborCellId)
A topological inquiry to determine whether a topological entity (e.g., point, edge,...
void GetCellBounds(vtkIdType cellId, double bounds[6]) override
Standard vtkDataSet methods; see vtkDataSet.h for documentation.
void SetCells(int type, vtkCellArray *cells)
Provide cell information to define the dataset.
vtkCellArray * GetPolyhedronFaceLocations()
Get pointer to faces and facelocations for polyhedron cells.
vtkUnsignedCharArray * GetDistinctCellTypesArray()
Get a list of types of cells in a dataset.
static void DecomposeAPolyhedronCell(const vtkIdType *polyhedronCellStream, vtkIdType &nCellpts, vtkIdType &nCellfaces, vtkCellArray *cellArray, vtkIdTypeArray *faces)
bool IsCellBoundary(vtkIdType cellId, vtkIdType npts, const vtkIdType *ptIds)
A topological inquiry to determine whether a topological entity (e.g., point, edge,...
vtkMTimeType GetMeshMTime() override
Return the mesh (geometry/topology) modification time.
vtkIdType InsertNextLinkedCell(int type, int npts, const vtkIdType pts[])
Use these methods only if the dataset has been specified as Editable.
void GetFaceStream(vtkIdType cellId, vtkIdList *ptIds)
Get the face stream of a polyhedron cell in the following format: (numCellFaces, numFace0Pts,...
vtkIdType GetNumberOfCells() override
Standard vtkDataSet methods; see vtkDataSet.h for documentation.
bool AllocateExact(vtkIdType numCells, vtkIdType connectivitySize)
Pre-allocate memory in internal data structures.
void Reset()
Standard vtkDataSet methods; see vtkDataSet.h for documentation.
void GetPointCells(vtkIdType ptId, vtkIdList *cellIds) override
Standard vtkDataSet methods; see vtkDataSet.h for documentation.
vtkCellIterator * NewCellIterator() override
Standard vtkDataSet methods; see vtkDataSet.h for documentation.
void GetPointCells(vtkIdType ptId, vtkIdType &ncells, vtkIdType *&cells)
Special (efficient) operation to return the list of cells using the specified point ptId.
void SetCells(vtkUnsignedCharArray *cellTypes, vtkCellArray *cells)
Provide cell information to define the dataset.
void RemoveGhostCells()
This method will remove any cell that is marked as ghost (has the vtkDataSetAttributes::DUPLICATECELL...
void GetCellNeighbors(vtkIdType cellId, vtkIdList *ptIds, vtkIdList *cellIds) override
A topological inquiry to retrieve all of the cells using list of points exclusive of the current cell...
static vtkUnstructuredGrid * GetData(vtkInformationVector *v, int i=0)
Retrieve an instance of this class from an information object.
void CopyStructure(vtkDataSet *ds) override
Standard vtkDataSet methods; see vtkDataSet.h for documentation.
void Initialize() override
Reset the grid to an empty state and free any memory.
int GetMaxCellSize() override
Get the size, in number of points, of the largest cell.
static void DecomposeAPolyhedronCell(vtkCellArray *polyhedronCellArray, vtkIdType &nCellpts, vtkIdType &nCellfaces, vtkCellArray *cellArray, vtkIdTypeArray *faces)
A static method for converting a polyhedron vtkCellArray of format [nCellFaces, nFace0Pts,...
static void ConvertFaceStreamPointIds(vtkCellArray *faces, vtkIdType *idMap)
Convert pid in a face stream into idMap[pid].
int GetMinSpatialDimension() override
Get the maximum/minimum spatial dimensionality of the data which is the maximum/minimum dimension of ...
virtual int GetGhostLevel()
Get the ghost level.
int GetCellType(vtkIdType cellId) override
Get the type of the cell with the given cellId.
void PrintSelf(ostream &os, vtkIndent indent) override
Standard methods for type information and printing.
static void ConvertFaceStreamPointIds(vtkIdType nfaces, vtkIdType *faceStream, vtkIdType *idMap)
Convert pid in a face stream into idMap[pid].
void BuildLinks()
Build topological links from points to lists of cells that use each point.
void GetIdsOfCellsOfType(int type, vtkIdTypeArray *array) override
Fill vtkIdTypeArray container with list of cell Ids.
virtual int GetPiece()
Set / Get the piece and the number of pieces.
unsigned long GetActualMemorySize() override
Return the actual size of the data in kibibytes (1024 bytes).
virtual int GetNumberOfPieces()
Set / Get the piece and the number of pieces.
void AddReferenceToCell(vtkIdType ptId, vtkIdType cellId)
Use these methods only if the dataset has been specified as Editable.
static vtkUnstructuredGrid * ExtendedNew()
static vtkUnstructuredGrid * New()
Standard instantiation method.
vtkUnsignedCharArray * GetCellTypesArray()
Get the array of all cell types in the grid.
void GetCellPoints(vtkIdType cellId, vtkIdType &npts, vtkIdType const *&pts, vtkIdList *ptIds) override
A higher-performing variant of the virtual vtkDataSet::GetCellPoints() for unstructured grids.
void ShallowCopy(vtkDataObject *src) override
Shallow and Deep copy.
void DeepCopy(vtkDataObject *src) override
Shallow and Deep copy.
#define VTK_DEPRECATED_IN_9_6_0(reason)
#define VTK_DEPRECATED_IN_9_5_0(reason)
int vtkIdType
Definition vtkType.h:332
vtkTypeUInt32 vtkMTimeType
Definition vtkType.h:287
@ VTK_UNSTRUCTURED_GRID
Definition vtkType.h:80
#define VTK_SIZEHINT(...)
#define VTK_MARSHALMANUAL