VTK  9.3.20240319
Public Member Functions | Static Public Member Functions | Protected Member Functions | Static Protected Member Functions | Protected Attributes | List of all members
vtkUnstructuredGrid Class Reference

dataset represents arbitrary combinations of all possible cell types More...

#include <vtkUnstructuredGrid.h>

Inheritance diagram for vtkUnstructuredGrid:
[legend]
Collaboration diagram for vtkUnstructuredGrid:
[legend]

Public Member Functions

int GetDataObjectType () override
 Standard vtkDataSet API methods. More...
 
bool AllocateEstimate (vtkIdType numCells, vtkIdType maxCellSize)
 Pre-allocate memory in internal data structures. More...
 
bool AllocateExact (vtkIdType numCells, vtkIdType connectivitySize)
 Pre-allocate memory in internal data structures. More...
 
void Allocate (vtkIdType numCells=1000, int vtkNotUsed(extSize)=1000) override
 Method allocates initial storage for the cell connectivity. More...
 
int GetCellType (vtkIdType cellId) override
 Get the type of the cell with the given cellId. More...
 
vtkIdType GetCellSize (vtkIdType cellId) override
 Get the size of the cell with given cellId. More...
 
void GetCellTypes (vtkCellTypes *types) override
 Get a list of types of cells in a dataset. More...
 
vtkUnsignedCharArrayGetDistinctCellTypesArray ()
 Get a list of types of cells in a dataset. More...
 
void GetCellPoints (vtkIdType cellId, vtkIdType &npts, vtkIdType const *&pts)
 A higher-performing variant of the virtual vtkDataSet::GetCellPoints() for unstructured grids. More...
 
void GetCellPoints (vtkIdType cellId, vtkIdType &npts, vtkIdType const *&pts, vtkIdList *ptIds) override
 A higher-performing variant of the virtual vtkDataSet::GetCellPoints() for unstructured grids. More...
 
vtkUnsignedCharArrayGetCellTypesArray ()
 Get the array of all cell types in the grid. More...
 
void Squeeze () override
 Squeeze all arrays in the grid to conserve memory. More...
 
void Initialize () override
 Reset the grid to an empty state and free any memory. More...
 
int GetMaxCellSize () override
 Get the size, in number of points, of the largest cell. More...
 
int GetMaxSpatialDimension () override
 Get the maximum spatial dimensionality of the data which is the maximum dimension of all cells. More...
 
void BuildLinks ()
 Build topological links from points to lists of cells that use each point. More...
 
vtkAbstractCellLinksGetCellLinks ()
 Get the cell links. More...
 
void GetFaceStream (vtkIdType cellId, vtkIdList *ptIds)
 Get the face stream of a polyhedron cell in the following format: (numCellFaces, numFace0Pts, id1, id2, id3, numFace1Pts,id1, id2, id3, ...). More...
 
void GetFaceStream (vtkIdType cellId, vtkIdType &nfaces, vtkIdType const *&ptIds)
 Get the number of faces and the face stream of a polyhedral cell. More...
 
vtkCellArrayGetCells ()
 Return the unstructured grid connectivity array. More...
 
virtual int GetGhostLevel ()
 Get the ghost level. More...
 
unsigned long GetActualMemorySize () override
 Return the actual size of the data in kibibytes (1024 bytes). More...
 
void GetIdsOfCellsOfType (int type, vtkIdTypeArray *array) override
 Fill vtkIdTypeArray container with list of cell Ids. More...
 
int IsHomogeneous () override
 Returns whether cells are all of the same type. More...
 
void RemoveGhostCells ()
 This method will remove any cell that is marked as ghost (has the vtkDataSetAttributes::DUPLICATECELL or the vtkDataSetAttributes::HIDDENCELL bit set). More...
 
vtkIdTypeGetFaces (vtkIdType cellId)
 Special support for polyhedron. More...
 
void GetPolyhedronFaces (vtkIdType cellId, vtkCellArray *faces)
 Special support for polyhedron. More...
 
int InitializeFacesRepresentation (vtkIdType numPrevCells)
 Special function used by vtkUnstructuredGridReader. More...
 
virtual vtkMTimeType GetMeshMTime ()
 Return the mesh (geometry/topology) modification time. More...
 
vtkIdTypeArrayGetCellLocationsArray ()
 Get the array of all the starting indices of cell definitions in the cell array. More...
 
void Reset ()
 Standard vtkDataSet methods; see vtkDataSet.h for documentation. More...
 
void CopyStructure (vtkDataSet *ds) override
 Standard vtkDataSet methods; see vtkDataSet.h for documentation. More...
 
vtkIdType GetNumberOfCells () override
 Standard vtkDataSet methods; see vtkDataSet.h for documentation. More...
 
vtkCellGetCell (vtkIdType cellId) override
 Standard vtkDataSet methods; see vtkDataSet.h for documentation. More...
 
void GetCell (vtkIdType cellId, vtkGenericCell *cell) override
 Standard vtkDataSet methods; see vtkDataSet.h for documentation. More...
 
void GetCellBounds (vtkIdType cellId, double bounds[6]) override
 Standard vtkDataSet methods; see vtkDataSet.h for documentation. More...
 
void GetCellPoints (vtkIdType cellId, vtkIdList *ptIds) override
 Standard vtkDataSet methods; see vtkDataSet.h for documentation. More...
 
void GetPointCells (vtkIdType ptId, vtkIdList *cellIds) override
 Standard vtkDataSet methods; see vtkDataSet.h for documentation. More...
 
vtkCellIteratorNewCellIterator () override
 Standard vtkDataSet methods; see vtkDataSet.h for documentation. More...
 
virtual vtkCellGetCell (vtkIdType cellId)=0
 Get cell with cellId such that: 0 <= cellId < NumberOfCells. More...
 
virtual vtkCellGetCell (int vtkNotUsed(i), int vtkNotUsed(j), int vtkNotUsed(k))
 Standard vtkDataSet methods; see vtkDataSet.h for documentation. More...
 
virtual void GetCell (vtkIdType cellId, vtkGenericCell *cell)=0
 Get cell with cellId such that: 0 <= cellId < NumberOfCells. More...
 
void GetPointCells (vtkIdType ptId, vtkIdType &ncells, vtkIdType *&cells)
 Special (efficient) operation to return the list of cells using the specified point ptId. More...
 
 vtkSetSmartPointerMacro (Links, vtkAbstractCellLinks)
 Set/Get the links that you created possibly without using BuildLinks. More...
 
 vtkGetSmartPointerMacro (Links, vtkAbstractCellLinks)
 Set/Get the links that you created possibly without using BuildLinks. More...
 
void SetCells (int type, vtkCellArray *cells)
 Provide cell information to define the dataset. More...
 
void SetCells (int *types, vtkCellArray *cells)
 Provide cell information to define the dataset. More...
 
void SetCells (vtkUnsignedCharArray *cellTypes, vtkCellArray *cells)
 Provide cell information to define the dataset. More...
 
void SetPolyhedralCells (vtkUnsignedCharArray *cellTypes, vtkCellArray *cells, vtkCellArray *faceLocations, vtkCellArray *faces)
 Provide cell information to define the dataset. More...
 
void SetCells (vtkUnsignedCharArray *cellTypes, vtkCellArray *cells, vtkIdTypeArray *faceLocations, vtkIdTypeArray *faces)
 Provide cell information to define the dataset. More...
 
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 specified (e.g., cellId). More...
 
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 specified (e.g., cellId). More...
 
bool IsCellBoundary (vtkIdType cellId, vtkIdType npts, const vtkIdType *ptIds, vtkIdType &neighborCellId)
 A topological inquiry to determine whether a topological entity (e.g., point, edge, or face) defined by the point ids (ptIds of length npts) is a boundary entity of a specified cell (indicated by cellId). More...
 
bool IsCellBoundary (vtkIdType cellId, vtkIdType npts, const vtkIdType *ptIds)
 A topological inquiry to determine whether a topological entity (e.g., point, edge, or face) defined by the point ids (ptIds of length npts) is a boundary entity of a specified cell (indicated by cellId). More...
 
bool IsCellBoundary (vtkIdType cellId, vtkIdType npts, const vtkIdType *ptIds, vtkIdList *vtkNotUsed(cellIds))
 A topological inquiry to determine whether a topological entity (e.g., point, edge, or face) defined by the point ids (ptIds of length npts) is a boundary entity of a specified cell (indicated by cellId). More...
 
vtkIdType InsertNextLinkedCell (int type, int npts, const vtkIdType pts[])
 Use these methods only if the dataset has been specified as Editable. More...
 
void RemoveReferenceToCell (vtkIdType ptId, vtkIdType cellId)
 Use these methods only if the dataset has been specified as Editable. More...
 
void AddReferenceToCell (vtkIdType ptId, vtkIdType cellId)
 Use these methods only if the dataset has been specified as Editable. More...
 
void ResizeCellList (vtkIdType ptId, int size)
 Use these methods only if the dataset has been specified as Editable. More...
 
virtual int GetPiece ()
 Set / Get the piece and the number of pieces. More...
 
virtual int GetNumberOfPieces ()
 Set / Get the piece and the number of pieces. More...
 
void ShallowCopy (vtkDataObject *src) override
 Shallow and Deep copy. More...
 
void DeepCopy (vtkDataObject *src) override
 Shallow and Deep copy. More...
 
vtkIdTypeArrayGetFaces ()
 Get pointer to faces and facelocations. More...
 
vtkIdTypeArrayGetFaceLocations ()
 Get pointer to faces and facelocations. More...
 
vtkCellArrayGetPolyhedronFaces ()
 Get pointer to faces and facelocations for polyhedron cells. More...
 
vtkCellArrayGetPolyhedronFaceLocations ()
 Get pointer to faces and facelocations. More...
 
void SetCells (vtkUnsignedCharArray *cellTypes, vtkIdTypeArray *cellLocations, vtkCellArray *cells)
 Special methods specific to vtkUnstructuredGrid for defining the cells composing the dataset. More...
 
void SetCells (vtkUnsignedCharArray *cellTypes, vtkIdTypeArray *cellLocations, vtkCellArray *cells, vtkIdTypeArray *faceLocations, vtkIdTypeArray *faces)
 Special methods specific to vtkUnstructuredGrid for defining the cells composing the dataset. More...
 
- Public Member Functions inherited from vtkUnstructuredGridBase
 vtkAbstractTypeMacro (vtkUnstructuredGridBase, vtkPointSet)
 
void PrintSelf (ostream &os, vtkIndent indent) override
 Standard methods for type information and printing. More...
 
int GetDataObjectType () override
 Standard vtkDataSet API methods. More...
 
virtual void Allocate (vtkIdType numCells=1000, int extSize=1000)=0
 Allocate memory for the number of cells indicated. More...
 
void DeepCopy (vtkDataObject *src) override
 Shallow and Deep copy. More...
 
vtkIdType InsertNextCell (int type, vtkIdType npts, const vtkIdType ptIds[])
 Insert/create cell in object by type and list of point ids defining cell topology. More...
 
vtkIdType InsertNextCell (int type, vtkIdList *ptIds)
 Insert/create cell in object by a list of point ids defining cell topology. More...
 
vtkIdType InsertNextCell (int type, vtkIdType npts, const vtkIdType ptIds[], vtkIdType nfaces, const vtkIdType faces[])
 
vtkIdType InsertNextCell (int type, vtkIdType npts, const vtkIdType ptIds[], vtkCellArray *faces)
 
void ReplaceCell (vtkIdType cellId, int npts, const vtkIdType pts[])
 Replace the points defining cell "cellId" with a new set of points. More...
 
- Public Member Functions inherited from vtkPointSet
double * GetPoint (vtkIdType ptId) override
 See vtkDataSet for additional information. More...
 
void BuildCellLocator ()
 Build the cell locator. More...
 
vtkMTimeType GetMTime () override
 Get MTime which also considers its vtkPoints MTime. More...
 
void ComputeBounds () override
 Compute the (X, Y, Z) bounds of the data. More...
 
virtual void SetEditable (bool)
 Specify whether this dataset is editable after creation. More...
 
virtual bool GetEditable ()
 Specify whether this dataset is editable after creation. More...
 
virtual void EditableOn ()
 Specify whether this dataset is editable after creation. More...
 
virtual void EditableOff ()
 Specify whether this dataset is editable after creation. More...
 
vtkIdType GetNumberOfPoints () override
 See vtkDataSet for additional information. More...
 
void GetPoint (vtkIdType ptId, double x[3]) override
 See vtkDataSet for additional information. More...
 
vtkIdType FindPoint (double x[3]) override
 See vtkDataSet for additional information. More...
 
vtkIdType FindCell (double x[3], vtkCell *cell, vtkIdType cellId, double tol2, int &subId, double pcoords[3], double *weights) override
 See vtkDataSet for additional information. More...
 
vtkIdType FindCell (double x[3], vtkCell *cell, vtkGenericCell *gencell, vtkIdType cellId, double tol2, int &subId, double pcoords[3], double *weights) override
 See vtkDataSet for additional information. More...
 
virtual void GetCellPoints (vtkIdType cellId, vtkIdList *ptIds)=0
 This method resets parameter idList, as there is no cell in a vtkPointSet. More...
 
virtual void GetCellPoints (vtkIdType cellId, vtkIdType &npts, vtkIdType const *&pts, vtkIdList *ptIds)
 This method resets parameter idList, as there is no cell in a vtkPointSet. More...
 
void BuildPointLocator ()
 Build the internal point locator . More...
 
void BuildLocator ()
 Build the internal point locator . More...
 
virtual void SetPointLocator (vtkAbstractPointLocator *)
 Set / get an instance of vtkAbstractPointLocator which is used to support the FindPoint() and FindCell() methods. More...
 
virtual vtkAbstractPointLocatorGetPointLocator ()
 Set / get an instance of vtkAbstractPointLocator which is used to support the FindPoint() and FindCell() methods. More...
 
virtual void SetCellLocator (vtkAbstractCellLocator *)
 Set / get an instance of vtkAbstractCellLocator which may be used when a vtkCellLocatorStrategy is used during a FindCell() operation. More...
 
virtual vtkAbstractCellLocatorGetCellLocator ()
 Set / get an instance of vtkAbstractCellLocator which may be used when a vtkCellLocatorStrategy is used during a FindCell() operation. More...
 
virtual void SetPoints (vtkPoints *)
 Specify point array to define point coordinates. More...
 
vtkPointsGetPoints () override
 Specify point array to define point coordinates. More...
 
bool UsesGarbageCollector () const override
 Overwritten to handle the data/locator loop. More...
 
vtkPointSetNewInstance () const
 Standard methods for type information and printing. More...
 
- Public Member Functions inherited from vtkDataSet
vtkDataSetNewInstance () const
 
virtual void CopyAttributes (vtkDataSet *ds)
 Copy the attributes associated with the specified dataset to this instance of vtkDataSet. More...
 
virtual vtkCellGetCell (int vtkNotUsed(i), int vtkNotUsed(j), int vtkNotUsed(k))
 
void SetCellOrderAndRationalWeights (vtkIdType cellId, vtkGenericCell *cell)
 
int GetCellNumberOfFaces (vtkIdType cellId, unsigned char &cellType, vtkGenericCell *cell)
 Get the number of faces of a cell. More...
 
virtual vtkCellFindAndGetCell (double x[3], vtkCell *cell, vtkIdType cellId, double tol2, int &subId, double pcoords[3], double *weights)
 Locate the cell that contains a point and return the cell. More...
 
vtkCellDataGetCellData ()
 Return a pointer to this dataset's cell data. More...
 
vtkPointDataGetPointData ()
 Return a pointer to this dataset's point data. More...
 
double * GetBounds ()
 Return a pointer to the geometry bounding box in the form (xmin,xmax, ymin,ymax, zmin,zmax). More...
 
void GetBounds (double bounds[6])
 Return a pointer to the geometry bounding box in the form (xmin,xmax, ymin,ymax, zmin,zmax). More...
 
double * GetCenter ()
 Get the center of the bounding box. More...
 
void GetCenter (double center[3])
 Get the center of the bounding box. More...
 
double GetLength ()
 Return the length of the diagonal of the bounding box. More...
 
double GetLength2 ()
 Return the squared length of the diagonal of the bounding box. More...
 
virtual void GetScalarRange (double range[2])
 Convenience method to get the range of the first component (and only the first component) of any scalars in the data set. More...
 
double * GetScalarRange ()
 Convenience method to get the range of the first component (and only the first component) of any scalars in the data set. More...
 
int CheckAttributes ()
 This method checks to see if the cell and point attributes match the geometry. More...
 
vtkFieldDataGetAttributesAsFieldData (int type) override
 Returns the attributes of the data object as a vtkFieldData. More...
 
vtkIdType GetNumberOfElements (int type) override
 Get the number of elements for a specific attribute type (POINT, CELL, etc.). More...
 
bool HasAnyGhostCells ()
 Returns 1 if there are any ghost cells 0 otherwise. More...
 
bool HasAnyGhostPoints ()
 Returns 1 if there are any ghost points 0 otherwise. More...
 
virtual bool HasAnyBlankCells ()
 Returns 1 if there are any blanking cells 0 otherwise. More...
 
virtual bool HasAnyBlankPoints ()
 Returns 1 if there are any blanking points 0 otherwise. More...
 
vtkUnsignedCharArrayGetPointGhostArray ()
 Gets the array that defines the ghost type of each point. More...
 
void UpdatePointGhostArrayCache ()
 Updates the pointer to the point ghost array. More...
 
vtkUnsignedCharArrayAllocatePointGhostArray ()
 Allocate ghost array for points. More...
 
vtkUnsignedCharArrayGetCellGhostArray ()
 Get the array that defines the ghost type of each cell. More...
 
void UpdateCellGhostArrayCache ()
 Updates the pointer to the cell ghost array. More...
 
vtkUnsignedCharArrayAllocateCellGhostArray ()
 Allocate ghost array for cells. More...
 
vtkUnsignedCharArrayGetGhostArray (int type) override
 Returns the ghost array for the given type (point or cell). More...
 
bool SupportsGhostArray (int type) override
 Returns true for POINT or CELL, false otherwise. More...
 
vtkIdType FindPoint (double x, double y, double z)
 Locate the closest point to the global coordinate x. More...
 
virtual void GenerateGhostArray (int zeroExt[6])
 Normally called by pipeline executives or algorithms only. More...
 
virtual void GenerateGhostArray (int zeroExt[6], bool cellOnly)
 Normally called by pipeline executives or algorithms only. More...
 
- Public Member Functions inherited from vtkDataObject
vtkDataObjectNewInstance () const
 
void ReleaseData ()
 Release data back to system to conserve memory resource. More...
 
vtkMTimeType GetUpdateTime ()
 Used by Threaded ports to determine if they should initiate an asynchronous update (still in development). More...
 
virtual void CopyInformationFromPipeline (vtkInformation *vtkNotUsed(info))
 Copy from the pipeline information to the data object's own information. More...
 
virtual void CopyInformationToPipeline (vtkInformation *vtkNotUsed(info))
 Copy information from this data object to the pipeline information. More...
 
void DataHasBeenGenerated ()
 This method is called by the source when it executes to generate data. More...
 
virtual void PrepareForNewData ()
 make the output data ready for new data to be inserted. More...
 
virtual int GetExtentType ()
 The ExtentType will be left as VTK_PIECES_EXTENT for data objects such as vtkPolyData and vtkUnstructuredGrid. More...
 
virtual void Crop (const int *updateExtent)
 This method crops the data object (if necessary) so that the extent matches the update extent. More...
 
virtual vtkDataSetAttributesGetAttributes (int type)
 Returns the attributes of the data object of the specified attribute type. More...
 
virtual int GetAttributeTypeForArray (vtkAbstractArray *arr)
 Retrieves the attribute type that an array came from. More...
 
virtual vtkInformationGetInformation ()
 Set/Get the information object associated with this data object. More...
 
virtual void SetInformation (vtkInformation *)
 Set/Get the information object associated with this data object. More...
 
virtual vtkTypeBool GetDataReleased ()
 Get the flag indicating the data has been released. More...
 
virtual void SetFieldData (vtkFieldData *)
 Assign or retrieve a general field data to this data object. More...
 
virtual vtkFieldDataGetFieldData ()
 Assign or retrieve a general field data to this data object. More...
 
void GlobalReleaseDataFlagOn ()
 Turn on/off flag to control whether every object releases its data after being used by a filter. More...
 
void GlobalReleaseDataFlagOff ()
 Turn on/off flag to control whether every object releases its data after being used by a filter. More...
 
- Public Member Functions inherited from vtkObject
 vtkBaseTypeMacro (vtkObject, vtkObjectBase)
 
virtual void DebugOn ()
 Turn debugging output on. More...
 
virtual void DebugOff ()
 Turn debugging output off. More...
 
bool GetDebug ()
 Get the value of the debug flag. More...
 
void SetDebug (bool debugFlag)
 Set the value of the debug flag. More...
 
virtual void Modified ()
 Update the modification time for this object. More...
 
void RemoveObserver (unsigned long tag)
 
void RemoveObservers (unsigned long event)
 
void RemoveObservers (const char *event)
 
void RemoveAllObservers ()
 
vtkTypeBool HasObserver (unsigned long event)
 
vtkTypeBool HasObserver (const char *event)
 
vtkTypeBool InvokeEvent (unsigned long event)
 
vtkTypeBool InvokeEvent (const char *event)
 
std::string GetObjectDescription () const override
 The object description printed in messages and PrintSelf output. More...
 
unsigned long AddObserver (unsigned long event, vtkCommand *, float priority=0.0f)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
unsigned long AddObserver (const char *event, vtkCommand *, float priority=0.0f)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
vtkCommandGetCommand (unsigned long tag)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
void RemoveObserver (vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
void RemoveObservers (unsigned long event, vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
void RemoveObservers (const char *event, vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
vtkTypeBool HasObserver (unsigned long event, vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
vtkTypeBool HasObserver (const char *event, vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
template<class U , class T >
unsigned long AddObserver (unsigned long event, U observer, void(T::*callback)(), float priority=0.0f)
 Overloads to AddObserver that allow developers to add class member functions as callbacks for events. More...
 
template<class U , class T >
unsigned long AddObserver (unsigned long event, U observer, void(T::*callback)(vtkObject *, unsigned long, void *), float priority=0.0f)
 Overloads to AddObserver that allow developers to add class member functions as callbacks for events. More...
 
template<class U , class T >
unsigned long AddObserver (unsigned long event, U observer, bool(T::*callback)(vtkObject *, unsigned long, void *), float priority=0.0f)
 Allow user to set the AbortFlagOn() with the return value of the callback method. More...
 
vtkTypeBool InvokeEvent (unsigned long event, void *callData)
 This method invokes an event and return whether the event was aborted or not. More...
 
vtkTypeBool InvokeEvent (const char *event, void *callData)
 This method invokes an event and return whether the event was aborted or not. More...
 
virtual void SetObjectName (const std::string &objectName)
 Set/get the name of this object for reporting purposes. More...
 
virtual std::string GetObjectName () const
 Set/get the name of this object for reporting purposes. More...
 
- Public Member Functions inherited from vtkObjectBase
const char * GetClassName () const
 Return the class name as a string. More...
 
virtual vtkIdType GetNumberOfGenerationsFromBase (const char *name)
 Given the name of a base class of this class type, return the distance of inheritance between this class type and the named class (how many generations of inheritance are there between this class and the named class). More...
 
virtual void Delete ()
 Delete a VTK object. More...
 
virtual void FastDelete ()
 Delete a reference to this object. More...
 
void InitializeObjectBase ()
 
void Print (ostream &os)
 Print an object to an ostream. More...
 
void Register (vtkObjectBase *o)
 Increase the reference count (mark as used by another object). More...
 
virtual void UnRegister (vtkObjectBase *o)
 Decrease the reference count (release by another object). More...
 
int GetReferenceCount ()
 Return the current reference count of this object. More...
 
void SetReferenceCount (int)
 Sets the reference count. More...
 
bool GetIsInMemkind () const
 A local state flag that remembers whether this object lives in the normal or extended memory space. More...
 
virtual void PrintHeader (ostream &os, vtkIndent indent)
 Methods invoked by print to print information about the object including superclasses. More...
 
virtual void PrintTrailer (ostream &os, vtkIndent indent)
 Methods invoked by print to print information about the object including superclasses. More...
 

Static Public Member Functions

static vtkUnstructuredGridNew ()
 Standard instantiation method. More...
 
static vtkUnstructuredGridExtendedNew ()
 
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, i, j, k, nFace1Pts, i, j, k, ...] into three components: (1) an integer indicating the number of faces (2) a standard vtkCellArray storing point ids [nCell0Pts, i, j, k] and (3) an vtkIdTypeArray storing face connectivity in format [nFace0Pts, i, j, k, nFace1Pts, i, j, k, ...] Note: input is assumed to contain only one polyhedron cell. More...
 
static void DecomposeAPolyhedronCell (const vtkIdType *polyhedronCellStream, vtkIdType &nCellpts, vtkIdType &nCellfaces, vtkCellArray *cellArray, vtkIdTypeArray *faces)
 
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, i, j, k, nFace1Pts, i, j, k, ...] into three components: (1) an integer indicating the number of faces (2) a standard vtkCellArray storing point ids [nCell0Pts, i, j, k] and (3) an vtkIdTypeArray storing face connectivity in format [nFace0Pts, i, j, k, nFace1Pts, i, j, k, ...] Note: input is assumed to contain only one polyhedron cell. More...
 
static void ConvertFaceStreamPointIds (vtkIdList *faceStream, vtkIdType *idMap)
 Convert pid in a face stream into idMap[pid]. More...
 
static void ConvertFaceStreamPointIds (vtkIdType nfaces, vtkIdType *faceStream, vtkIdType *idMap)
 Convert pid in a face stream into idMap[pid]. More...
 
static void ConvertFaceStreamPointIds (vtkCellArray *faces, vtkIdType *idMap)
 Convert pid in a face stream into idMap[pid]. More...
 
static vtkUnstructuredGridGetData (vtkInformation *info)
 Retrieve an instance of this class from an information object. More...
 
static vtkUnstructuredGridGetData (vtkInformationVector *v, int i=0)
 Retrieve an instance of this class from an information object. More...
 
- Static Public Member Functions inherited from vtkUnstructuredGridBase
static vtkUnstructuredGridBaseGetData (vtkInformation *info)
 Retrieve an instance of this class from an information object. More...
 
static vtkUnstructuredGridBaseGetData (vtkInformationVector *v, int i=0)
 Retrieve an instance of this class from an information object. More...
 
- Static Public Member Functions inherited from vtkPointSet
static vtkPointSetNew ()
 Standard instantiation method. More...
 
static vtkPointSetExtendedNew ()
 
static vtkPointSetGetData (vtkInformation *info)
 Retrieve an instance of this class from an information object. More...
 
static vtkPointSetGetData (vtkInformationVector *v, int i=0)
 Retrieve an instance of this class from an information object. More...
 
static vtkTypeBool IsTypeOf (const char *type)
 Standard methods for type information and printing. More...
 
static vtkPointSetSafeDownCast (vtkObjectBase *o)
 Standard methods for type information and printing. More...
 
- Static Public Member Functions inherited from vtkDataSet
static vtkTypeBool IsTypeOf (const char *type)
 
static vtkDataSetSafeDownCast (vtkObjectBase *o)
 
static vtkDataSetGetData (vtkInformation *info)
 Retrieve an instance of this class from an information object. More...
 
static vtkDataSetGetData (vtkInformationVector *v, int i=0)
 Retrieve an instance of this class from an information object. More...
 
- Static Public Member Functions inherited from vtkDataObject
static vtkDataObjectNew ()
 
static vtkTypeBool IsTypeOf (const char *type)
 
static vtkDataObjectSafeDownCast (vtkObjectBase *o)
 
static vtkInformationGetActiveFieldInformation (vtkInformation *info, int fieldAssociation, int attributeType)
 Return the information object within the input information object's field data corresponding to the specified association (FIELD_ASSOCIATION_POINTS or FIELD_ASSOCIATION_CELLS) and attribute (SCALARS, VECTORS, NORMALS, TCOORDS, or TENSORS) More...
 
static vtkInformationGetNamedFieldInformation (vtkInformation *info, int fieldAssociation, const char *name)
 Return the information object within the input information object's field data corresponding to the specified association (FIELD_ASSOCIATION_POINTS or FIELD_ASSOCIATION_CELLS) and name. More...
 
static void RemoveNamedFieldInformation (vtkInformation *info, int fieldAssociation, const char *name)
 Remove the info associated with an array. More...
 
static vtkInformationSetActiveAttribute (vtkInformation *info, int fieldAssociation, const char *attributeName, int attributeType)
 Set the named array to be the active field for the specified type (SCALARS, VECTORS, NORMALS, TCOORDS, or TENSORS) and association (FIELD_ASSOCIATION_POINTS or FIELD_ASSOCIATION_CELLS). More...
 
static void SetActiveAttributeInfo (vtkInformation *info, int fieldAssociation, int attributeType, const char *name, int arrayType, int numComponents, int numTuples)
 Set the name, array type, number of components, and number of tuples within the passed information object for the active attribute of type attributeType (in specified association, FIELD_ASSOCIATION_POINTS or FIELD_ASSOCIATION_CELLS). More...
 
static void SetPointDataActiveScalarInfo (vtkInformation *info, int arrayType, int numComponents)
 Convenience version of previous method for use (primarily) by the Imaging filters. More...
 
static const char * GetAssociationTypeAsString (int associationType)
 Given an integer association type, this static method returns a string type for the attribute (i.e. More...
 
static int GetAssociationTypeFromString (const char *associationName)
 Given a string association name, this static method returns an integer association type for the attribute (i.e. More...
 
static vtkInformationStringKeyDATA_TYPE_NAME ()
 
static vtkInformationDataObjectKeyDATA_OBJECT ()
 
static vtkInformationIntegerKeyDATA_EXTENT_TYPE ()
 
static vtkInformationIntegerPointerKeyDATA_EXTENT ()
 
static vtkInformationIntegerVectorKeyALL_PIECES_EXTENT ()
 
static vtkInformationIntegerKeyDATA_PIECE_NUMBER ()
 
static vtkInformationIntegerKeyDATA_NUMBER_OF_PIECES ()
 
static vtkInformationIntegerKeyDATA_NUMBER_OF_GHOST_LEVELS ()
 
static vtkInformationDoubleKeyDATA_TIME_STEP ()
 
static vtkInformationInformationVectorKeyPOINT_DATA_VECTOR ()
 
static vtkInformationInformationVectorKeyCELL_DATA_VECTOR ()
 
static vtkInformationInformationVectorKeyVERTEX_DATA_VECTOR ()
 
static vtkInformationInformationVectorKeyEDGE_DATA_VECTOR ()
 
static vtkInformationIntegerKeyFIELD_ARRAY_TYPE ()
 
static vtkInformationIntegerKeyFIELD_ASSOCIATION ()
 
static vtkInformationIntegerKeyFIELD_ATTRIBUTE_TYPE ()
 
static vtkInformationIntegerKeyFIELD_ACTIVE_ATTRIBUTE ()
 
static vtkInformationIntegerKeyFIELD_NUMBER_OF_COMPONENTS ()
 
static vtkInformationIntegerKeyFIELD_NUMBER_OF_TUPLES ()
 
static vtkInformationIntegerKeyFIELD_OPERATION ()
 
static vtkInformationDoubleVectorKeyFIELD_RANGE ()
 
static vtkInformationIntegerVectorKeyPIECE_EXTENT ()
 
static vtkInformationStringKeyFIELD_NAME ()
 
static vtkInformationDoubleVectorKeyORIGIN ()
 
static vtkInformationDoubleVectorKeySPACING ()
 
static vtkInformationDoubleVectorKeyDIRECTION ()
 
static vtkInformationDoubleVectorKeyBOUNDING_BOX ()
 
static vtkInformationDataObjectKeySIL ()
 
static vtkDataObjectGetData (vtkInformation *info)
 Retrieve an instance of this class from an information object. More...
 
static vtkDataObjectGetData (vtkInformationVector *v, int i=0)
 Retrieve an instance of this class from an information object. More...
 
static void SetGlobalReleaseDataFlag (vtkTypeBool val)
 Turn on/off flag to control whether every object releases its data after being used by a filter. More...
 
static vtkTypeBool GetGlobalReleaseDataFlag ()
 Turn on/off flag to control whether every object releases its data after being used by a filter. More...
 
- Static Public Member Functions inherited from vtkObject
static vtkObjectNew ()
 Create an object with Debug turned off, modified time initialized to zero, and reference counting on. More...
 
static void BreakOnError ()
 This method is called when vtkErrorMacro executes. More...
 
static void SetGlobalWarningDisplay (vtkTypeBool val)
 This is a global flag that controls whether any debug, warning or error messages are displayed. More...
 
static void GlobalWarningDisplayOn ()
 This is a global flag that controls whether any debug, warning or error messages are displayed. More...
 
static void GlobalWarningDisplayOff ()
 This is a global flag that controls whether any debug, warning or error messages are displayed. More...
 
static vtkTypeBool GetGlobalWarningDisplay ()
 This is a global flag that controls whether any debug, warning or error messages are displayed. More...
 
- Static Public Member Functions inherited from vtkObjectBase
static vtkTypeBool IsTypeOf (const char *name)
 Return 1 if this class type is the same type of (or a subclass of) the named class. More...
 
static vtkIdType GetNumberOfGenerationsFromBaseType (const char *name)
 Given a the name of a base class of this class type, return the distance of inheritance between this class type and the named class (how many generations of inheritance are there between this class and the named class). More...
 
static vtkObjectBaseNew ()
 Create an object with Debug turned off, modified time initialized to zero, and reference counting on. More...
 
static void SetMemkindDirectory (const char *directoryname)
 The name of a directory, ideally mounted -o dax, to memory map an extended memory space within. More...
 
static bool GetUsingMemkind ()
 A global state flag that controls whether vtkObjects are constructed in the usual way (the default) or within the extended memory space. More...
 

Protected Member Functions

 vtkUnstructuredGrid ()
 
 ~vtkUnstructuredGrid () override
 
void ReportReferences (vtkGarbageCollector *) override
 
vtkIdType InternalInsertNextCell (int type, vtkIdType npts, const vtkIdType ptIds[]) override
 
vtkIdType InternalInsertNextCell (int type, vtkIdList *ptIds) override
 
vtkIdType InternalInsertNextCell (int type, vtkIdType npts, const vtkIdType ptIds[], vtkIdType nfaces, const vtkIdType faces[])
 
vtkIdType InternalInsertNextCell (int type, vtkIdType npts, const vtkIdType pts[], vtkCellArray *faces) override
 
void InternalReplaceCell (vtkIdType cellId, int npts, const vtkIdType pts[]) override
 
- Protected Member Functions inherited from vtkUnstructuredGridBase
 vtkUnstructuredGridBase ()
 
 ~vtkUnstructuredGridBase () override
 
- Protected Member Functions inherited from vtkPointSet
 vtkPointSet ()
 
 ~vtkPointSet () override
 
- Protected Member Functions inherited from vtkDataSet
 vtkDataSet ()
 
 ~vtkDataSet () override
 
vtkMTimeType GetGhostCellsTime ()
 Return the MTime of the ghost cells array. More...
 
virtual void ComputeScalarRange ()
 Compute the range of the scalars and cache it into ScalarRange only if the cache became invalid (ScalarRangeComputeTime). More...
 
- Protected Member Functions inherited from vtkDataObject
 vtkDataObject ()
 
 ~vtkDataObject () override
 
- Protected Member Functions inherited from vtkObject
 vtkObject ()
 
 ~vtkObject () override
 
void RegisterInternal (vtkObjectBase *, vtkTypeBool check) override
 
void UnRegisterInternal (vtkObjectBase *, vtkTypeBool check) override
 
void InternalGrabFocus (vtkCommand *mouseEvents, vtkCommand *keypressEvents=nullptr)
 These methods allow a command to exclusively grab all events. More...
 
void InternalReleaseFocus ()
 These methods allow a command to exclusively grab all events. More...
 
- Protected Member Functions inherited from vtkObjectBase
 vtkObjectBase ()
 
virtual ~vtkObjectBase ()
 
 vtkObjectBase (const vtkObjectBase &)
 
void operator= (const vtkObjectBase &)
 

Static Protected Member Functions

static void DecomposeAPolyhedronCell (const vtkIdType *cellStream, vtkIdType &numCellPts, vtkIdType &nCellFaces, vtkCellArray *cellArray, vtkCellArray *faces)
 A static method for converting an input polyhedron cell stream of format [nFace0Pts, i, j, k, nFace1Pts, i, j, k, ...] into three components: (1) an integer indicating the number of faces (2) a standard vtkCellArray storing point ids [nCell0Pts, i, j, k] and (3) an vtkIdTypeArray storing face connectivity in format [nFace0Pts, i, j, k, nFace1Pts, i, j, k, ...] Note: input is assumed to contain only one polyhedron cell. More...
 
static void DecomposeAPolyhedronCell (vtkIdType nCellFaces, const vtkIdType *cellStream, vtkIdType &numCellPts, vtkCellArray *cellArray, vtkCellArray *facesArray)
 
static void DecomposeAPolyhedronCell (vtkCellArray *polyhedronCell, vtkIdType &numCellPts, vtkIdType &nCellfaces, vtkCellArray *cellArray, vtkCellArray *faces)
 
static void DecomposeAPolyhedronCell (const vtkIdType *cellStream, vtkIdType &numCellPts, vtkIdType &nCellFaces, vtkCellArray *cellArray, vtkCellArray *faces, vtkCellArray *faceLocations)
 
static void DecomposeAPolyhedronCell (vtkIdType nCellFaces, const vtkIdType *cellStream, vtkIdType &numCellPts, vtkCellArray *cellArray, vtkCellArray *faces, vtkCellArray *faceLocations)
 
static int CopyPolyhedronToFaceStream (vtkCellArray *faceArray, vtkCellArray *faceLocationArray, vtkIdTypeArray *faceStream, vtkIdTypeArray *faceLocation)
 Backward compatibility function to convert new polyhedron storage to legacy. More...
 
- Static Protected Member Functions inherited from vtkObjectBase
static vtkMallocingFunction GetCurrentMallocFunction ()
 
static vtkReallocingFunction GetCurrentReallocFunction ()
 
static vtkFreeingFunction GetCurrentFreeFunction ()
 
static vtkFreeingFunction GetAlternateFreeFunction ()
 

Protected Attributes

vtkSmartPointer< vtkCellArrayConnectivity
 
vtkSmartPointer< vtkAbstractCellLinksLinks
 
vtkSmartPointer< vtkUnsignedCharArrayTypes
 
vtkSmartPointer< vtkCellTypesDistinctCellTypes
 
vtkMTimeType DistinctCellTypesUpdateMTime
 
vtkSmartPointer< vtkCellArrayFaces
 Special support for polyhedra/cells with explicit face representations. More...
 
vtkSmartPointer< vtkCellArrayFaceLocations
 
vtkSmartPointer< vtkIdTypeArrayCellLocations
 
vtkSmartPointer< vtkIdTypeArrayLegacyFaces
 Legacy support – stores the old-style Faces && FaceLocations Special support for polyhedra/cells with explicit face representations. More...
 
vtkSmartPointer< vtkIdTypeArrayLegacyFaceLocations
 
vtkSmartPointer< vtkIdListLegacyPointIdsBuffer
 Legacy backward compatibility for GetFaceStream This member should be removed simultaneously with the deprecated GetFaceStream. More...
 
- Protected Attributes inherited from vtkPointSet
bool Editable
 
vtkPointsPoints
 
vtkAbstractPointLocatorPointLocator
 
vtkAbstractCellLocatorCellLocator
 
- Protected Attributes inherited from vtkDataSet
vtkNew< vtkGenericCellGenericCell
 
vtkCellDataCellData
 
vtkPointDataPointData
 
vtkCallbackCommandDataObserver
 
vtkTimeStamp ComputeTime
 
double Bounds [6]
 
double Center [3]
 
double ScalarRange [2]
 
vtkTimeStamp ScalarRangeComputeTime
 
vtkUnsignedCharArrayPointGhostArray
 These arrays pointers are caches used to avoid a string comparison (when getting ghost arrays using GetArray(name)) More...
 
vtkUnsignedCharArrayCellGhostArray
 These arrays pointers are caches used to avoid a string comparison (when getting ghost arrays using GetArray(name)) More...
 
bool PointGhostArrayCached
 These arrays pointers are caches used to avoid a string comparison (when getting ghost arrays using GetArray(name)) More...
 
bool CellGhostArrayCached
 These arrays pointers are caches used to avoid a string comparison (when getting ghost arrays using GetArray(name)) More...
 
- Protected Attributes inherited from vtkDataObject
vtkFieldDataFieldData
 
vtkTypeBool DataReleased
 
vtkTimeStamp UpdateTime
 
vtkInformationInformation
 
- Protected Attributes inherited from vtkObject
bool Debug
 
vtkTimeStamp MTime
 
vtkSubjectHelper * SubjectHelper
 
std::string ObjectName
 
- Protected Attributes inherited from vtkObjectBase
std::atomic< int32_t > ReferenceCount
 
vtkWeakPointerBase ** WeakPointers
 
typedef vtkUnstructuredGridBase Superclass
 Standard methods for type information and printing. More...
 
static vtkTypeBool IsTypeOf (const char *type)
 Standard methods for type information and printing. More...
 
static vtkUnstructuredGridSafeDownCast (vtkObjectBase *o)
 Standard methods for type information and printing. More...
 
virtual vtkTypeBool IsA (const char *type)
 Standard methods for type information and printing. More...
 
vtkUnstructuredGridNewInstance () const
 Standard methods for type information and printing. More...
 
void PrintSelf (ostream &os, vtkIndent indent) override
 Standard methods for type information and printing. More...
 
virtual vtkObjectBaseNewInstanceInternal () const
 Standard methods for type information and printing. More...
 

Additional Inherited Members

- Public Types inherited from vtkPointSet
typedef vtkDataSet Superclass
 Standard methods for type information and printing. More...
 
- Public Types inherited from vtkDataSet
enum  FieldDataType { DATA_OBJECT_FIELD = 0 , POINT_DATA_FIELD = 1 , CELL_DATA_FIELD = 2 }
 
typedef vtkDataObject Superclass
 
- Public Types inherited from vtkDataObject
enum  FieldAssociations {
  FIELD_ASSOCIATION_POINTS , FIELD_ASSOCIATION_CELLS , FIELD_ASSOCIATION_NONE , FIELD_ASSOCIATION_POINTS_THEN_CELLS ,
  FIELD_ASSOCIATION_VERTICES , FIELD_ASSOCIATION_EDGES , FIELD_ASSOCIATION_ROWS , NUMBER_OF_ASSOCIATIONS
}
 Possible values for the FIELD_ASSOCIATION information entry. More...
 
enum  AttributeTypes {
  POINT , CELL , FIELD , POINT_THEN_CELL ,
  VERTEX , EDGE , ROW , NUMBER_OF_ATTRIBUTE_TYPES
}
 Possible attribute types. More...
 
enum  FieldOperations { FIELD_OPERATION_PRESERVED , FIELD_OPERATION_REINTERPOLATED , FIELD_OPERATION_MODIFIED , FIELD_OPERATION_REMOVED }
 Possible values for the FIELD_OPERATION information entry. More...
 
typedef vtkObject Superclass
 

Detailed Description

dataset represents arbitrary combinations of all possible cell types

vtkUnstructuredGrid is a data object that is a concrete implementation of vtkDataSet. vtkUnstructuredGrid represents any combinations of any cell types. This includes 0D (e.g., points), 1D (e.g., lines, polylines), 2D (e.g., triangles, polygons), and 3D (e.g., hexahedron, tetrahedron, polyhedron, etc.). vtkUnstructuredGrid provides random access to cells, as well as topological information (such as lists of cells using each point).

Examples:
vtkUnstructuredGrid (Examples)
Online Examples:

Tests:
vtkUnstructuredGrid (Tests)

Definition at line 149 of file vtkUnstructuredGrid.h.

Member Typedef Documentation

◆ Superclass

Standard methods for type information and printing.

Definition at line 162 of file vtkUnstructuredGrid.h.

Constructor & Destructor Documentation

◆ vtkUnstructuredGrid()

vtkUnstructuredGrid::vtkUnstructuredGrid ( )
protected

◆ ~vtkUnstructuredGrid()

vtkUnstructuredGrid::~vtkUnstructuredGrid ( )
overrideprotected

Member Function Documentation

◆ New()

static vtkUnstructuredGrid* vtkUnstructuredGrid::New ( )
static

Standard instantiation method.

◆ ExtendedNew()

static vtkUnstructuredGrid* vtkUnstructuredGrid::ExtendedNew ( )
static

◆ IsTypeOf()

static vtkTypeBool vtkUnstructuredGrid::IsTypeOf ( const char *  type)
static

Standard methods for type information and printing.

◆ IsA()

virtual vtkTypeBool vtkUnstructuredGrid::IsA ( const char *  type)
virtual

Standard methods for type information and printing.

Reimplemented from vtkPointSet.

◆ SafeDownCast()

static vtkUnstructuredGrid* vtkUnstructuredGrid::SafeDownCast ( vtkObjectBase o)
static

Standard methods for type information and printing.

◆ NewInstanceInternal()

virtual vtkObjectBase* vtkUnstructuredGrid::NewInstanceInternal ( ) const
protectedvirtual

Standard methods for type information and printing.

Reimplemented from vtkPointSet.

◆ NewInstance()

vtkUnstructuredGrid* vtkUnstructuredGrid::NewInstance ( ) const

Standard methods for type information and printing.

◆ PrintSelf()

void vtkUnstructuredGrid::PrintSelf ( ostream &  os,
vtkIndent  indent 
)
overridevirtual

Standard methods for type information and printing.

Reimplemented from vtkPointSet.

◆ GetDataObjectType()

int vtkUnstructuredGrid::GetDataObjectType ( )
inlineoverridevirtual

Standard vtkDataSet API methods.

See vtkDataSet for more information.

Reimplemented from vtkPointSet.

Definition at line 169 of file vtkUnstructuredGrid.h.

◆ AllocateEstimate()

bool vtkUnstructuredGrid::AllocateEstimate ( vtkIdType  numCells,
vtkIdType  maxCellSize 
)
inline

Pre-allocate memory in internal data structures.

Does not change the number of cells, only the array capacities. Existing data is NOT preserved.

Parameters
numCellsThe number of expected cells in the dataset.
maxCellSizeThe number of points per cell to allocate memory for.
Returns
True if allocation succeeds.
See also
Squeeze();

Definition at line 180 of file vtkUnstructuredGrid.h.

◆ AllocateExact()

bool vtkUnstructuredGrid::AllocateExact ( vtkIdType  numCells,
vtkIdType  connectivitySize 
)

Pre-allocate memory in internal data structures.

Does not change the number of cells, only the array capacities. Existing data is NOT preserved.

Parameters
numCellsThe number of expected cells in the dataset.
connectivitySizeThe total number of pointIds stored for all cells.
Returns
True if allocation succeeds.
See also
Squeeze();

◆ Allocate()

void vtkUnstructuredGrid::Allocate ( vtkIdType  numCells = 1000,
int   vtkNotUsedextSize = 1000 
)
inlineoverride

Method allocates initial storage for the cell connectivity.

Use this method before the method InsertNextCell(). The array capacity is doubled when the inserting a cell exceeds the current capacity. extSize is no longer used.

Note
Prefer AllocateExact or AllocateEstimate, which give more control over how allocations are distributed.

Definition at line 205 of file vtkUnstructuredGrid.h.

◆ Reset()

void vtkUnstructuredGrid::Reset ( )

Standard vtkDataSet methods; see vtkDataSet.h for documentation.

◆ CopyStructure()

void vtkUnstructuredGrid::CopyStructure ( vtkDataSet ds)
overridevirtual

Standard vtkDataSet methods; see vtkDataSet.h for documentation.

Reimplemented from vtkPointSet.

◆ GetNumberOfCells()

vtkIdType vtkUnstructuredGrid::GetNumberOfCells ( )
overridevirtual

Standard vtkDataSet methods; see vtkDataSet.h for documentation.

Reimplemented from vtkPointSet.

◆ GetCell() [1/5]

vtkCell* vtkUnstructuredGrid::GetCell ( vtkIdType  cellId)
overridevirtual

Standard vtkDataSet methods; see vtkDataSet.h for documentation.

Reimplemented from vtkPointSet.

◆ GetCell() [2/5]

void vtkUnstructuredGrid::GetCell ( vtkIdType  cellId,
vtkGenericCell cell 
)
overridevirtual

Standard vtkDataSet methods; see vtkDataSet.h for documentation.

Reimplemented from vtkPointSet.

◆ GetCellBounds()

void vtkUnstructuredGrid::GetCellBounds ( vtkIdType  cellId,
double  bounds[6] 
)
overridevirtual

Standard vtkDataSet methods; see vtkDataSet.h for documentation.

Reimplemented from vtkDataSet.

◆ GetCellPoints() [1/3]

void vtkUnstructuredGrid::GetCellPoints ( vtkIdType  cellId,
vtkIdList ptIds 
)
overridevirtual

Standard vtkDataSet methods; see vtkDataSet.h for documentation.

Reimplemented from vtkPointSet.

◆ GetPointCells() [1/2]

void vtkUnstructuredGrid::GetPointCells ( vtkIdType  ptId,
vtkIdList cellIds 
)
overridevirtual

Standard vtkDataSet methods; see vtkDataSet.h for documentation.

Reimplemented from vtkPointSet.

◆ NewCellIterator()

vtkCellIterator* vtkUnstructuredGrid::NewCellIterator ( )
overridevirtual

Standard vtkDataSet methods; see vtkDataSet.h for documentation.

Reimplemented from vtkPointSet.

◆ GetCellType()

int vtkUnstructuredGrid::GetCellType ( vtkIdType  cellId)
overridevirtual

Get the type of the cell with the given cellId.

Reimplemented from vtkPointSet.

◆ GetCellSize()

vtkIdType vtkUnstructuredGrid::GetCellSize ( vtkIdType  cellId)
overridevirtual

Get the size of the cell with given cellId.

Reimplemented from vtkPointSet.

◆ GetCellTypes()

void vtkUnstructuredGrid::GetCellTypes ( vtkCellTypes types)
overridevirtual

Get a list of types of cells in a dataset.

The list consists of an array of types (not necessarily in any order), with a single entry per type. For example a dataset with 5 triangles, 3 lines, and 100 hexahedra would result in a list of three entries, corresponding to the types VTK_TRIANGLE, VTK_LINE, and VTK_HEXAHEDRON. This override implements an optimization that recomputes cell types only when the types of cells may have changed.

THIS METHOD IS THREAD SAFE IF FIRST CALLED FROM A SINGLE THREAD AND THE DATASET IS NOT MODIFIED

Reimplemented from vtkDataSet.

◆ GetDistinctCellTypesArray()

vtkUnsignedCharArray* vtkUnstructuredGrid::GetDistinctCellTypesArray ( )

Get a list of types of cells in a dataset.

The list consists of an array of types (not necessarily in any order), with a single entry per type. For example a dataset with 5 triangles, 3 lines, and 100 hexahedra would result in a list of three entries, corresponding to the types VTK_TRIANGLE, VTK_LINE, and VTK_HEXAHEDRON. This override implements an optimization that recomputes cell types only when the types of cells may have changed. This method never returns nullptr.

THIS METHOD IS THREAD SAFE IF FIRST CALLED FROM A SINGLE THREAD AND THE DATASET IS NOT MODIFIED

◆ GetCellPoints() [2/3]

void vtkUnstructuredGrid::GetCellPoints ( vtkIdType  cellId,
vtkIdType npts,
vtkIdType const *&  pts 
)
inline

A higher-performing variant of the virtual vtkDataSet::GetCellPoints() for unstructured grids.

Given a cellId, return the number of defining points and the list of points defining the cell.

Warning
Subsequent calls to this method may invalidate previous call results.

The pts pointer must not be modified.

Note: This method MAY NOT be thread-safe. (See GetCellAtId at vtkCellArray)

Definition at line 275 of file vtkUnstructuredGrid.h.

◆ GetCellPoints() [3/3]

void vtkUnstructuredGrid::GetCellPoints ( vtkIdType  cellId,
vtkIdType npts,
vtkIdType const *&  pts,
vtkIdList ptIds 
)
inlineoverridevirtual

A higher-performing variant of the virtual vtkDataSet::GetCellPoints() for unstructured grids.

Given a cellId, return the number of defining points and the list of points defining the cell.

This function MAY use ptIds, which is an object that is created by each thread, to guarantee thread safety.

Warning
Subsequent calls to this method may invalidate previous call results.

The pts pointer must not be modified.

Note: This method is thread-safe.

Reimplemented from vtkDataSet.

Definition at line 295 of file vtkUnstructuredGrid.h.

◆ GetPointCells() [2/2]

void vtkUnstructuredGrid::GetPointCells ( vtkIdType  ptId,
vtkIdType ncells,
vtkIdType *&  cells 
)

Special (efficient) operation to return the list of cells using the specified point ptId.

Use carefully (i.e., make sure that BuildLinks() has been called).

◆ GetCellTypesArray()

vtkUnsignedCharArray* vtkUnstructuredGrid::GetCellTypesArray ( )

Get the array of all cell types in the grid.

Each single-component tuple in the array at an index that corresponds to the type of the cell with the same index. To get an array of only the distinct cell types in the dataset, use GetCellTypes().

◆ Squeeze()

void vtkUnstructuredGrid::Squeeze ( )
overridevirtual

Squeeze all arrays in the grid to conserve memory.

Reimplemented from vtkPointSet.

◆ Initialize()

void vtkUnstructuredGrid::Initialize ( )
overridevirtual

Reset the grid to an empty state and free any memory.

Reimplemented from vtkPointSet.

◆ GetMaxCellSize()

int vtkUnstructuredGrid::GetMaxCellSize ( )
overridevirtual

Get the size, in number of points, of the largest cell.

Reimplemented from vtkPointSet.

◆ GetMaxSpatialDimension()

int vtkUnstructuredGrid::GetMaxSpatialDimension ( )
overridevirtual

Get the maximum spatial dimensionality of the data which is the maximum dimension of all cells.

Reimplemented from vtkDataSet.

◆ BuildLinks()

void vtkUnstructuredGrid::BuildLinks ( )

Build topological links from points to lists of cells that use each point.

See vtkAbstractCellLinks for more information.

◆ vtkSetSmartPointerMacro()

vtkUnstructuredGrid::vtkSetSmartPointerMacro ( Links  ,
vtkAbstractCellLinks   
)

Set/Get the links that you created possibly without using BuildLinks.

◆ vtkGetSmartPointerMacro()

vtkUnstructuredGrid::vtkGetSmartPointerMacro ( Links  ,
vtkAbstractCellLinks   
)

Set/Get the links that you created possibly without using BuildLinks.

◆ GetCellLinks()

vtkAbstractCellLinks* vtkUnstructuredGrid::GetCellLinks ( )

Get the cell links.

The cell links will be one of nullptr=0; vtkCellLinks=1; vtkStaticCellLinksTemplate<VTK_UNSIGNED_SHORT>=2; vtkStaticCellLinksTemplate<VTK_UNSIGNED_INT>=3; vtkStaticCellLinksTemplate<VTK_ID_TYPE>=4. (See enum types defined in vtkAbstractCellLinks.)

◆ GetFaceStream() [1/2]

void vtkUnstructuredGrid::GetFaceStream ( vtkIdType  cellId,
vtkIdList ptIds 
)

Get the face stream of a polyhedron cell in the following format: (numCellFaces, numFace0Pts, id1, id2, id3, numFace1Pts,id1, id2, id3, ...).

If the requested cell is not a polyhedron, then the standard GetCellPoints is called to return a list of unique point ids (id1, id2, id3, ...). This function is threadsafe.

◆ GetFaceStream() [2/2]

void vtkUnstructuredGrid::GetFaceStream ( vtkIdType  cellId,
vtkIdType nfaces,
vtkIdType const *&  ptIds 
)

Get the number of faces and the face stream of a polyhedral cell.

The output ptIds has the following format: (numFace0Pts, id1, id2, id3, numFace1Pts,id1, id2, id3, ...). If the requested cell is not a polyhedron, then the standard GetCellPoints is called to return the number of points and a list of unique point ids (id1, id2, id3, ...). This function is NOT THREADSAFE.

◆ SetCells() [1/6]

void vtkUnstructuredGrid::SetCells ( int  type,
vtkCellArray cells 
)

Provide cell information to define the dataset.

Cells like vtkPolyhedron require points plus a list of faces. To handle vtkPolyhedron, use SetPolyhedralCells() SetPolyhedralCells requires a faces vtkCellArray that will describe the faces use by the polyhedral cells. SetPolyhedralCells also requires a faceLocations vtkCellArray to fully describe a polyhedron cell The faceLocations is a collection of face ids pointing to the faces vtkCellArray.

◆ SetCells() [2/6]

void vtkUnstructuredGrid::SetCells ( int *  types,
vtkCellArray cells 
)

Provide cell information to define the dataset.

Cells like vtkPolyhedron require points plus a list of faces. To handle vtkPolyhedron, use SetPolyhedralCells() SetPolyhedralCells requires a faces vtkCellArray that will describe the faces use by the polyhedral cells. SetPolyhedralCells also requires a faceLocations vtkCellArray to fully describe a polyhedron cell The faceLocations is a collection of face ids pointing to the faces vtkCellArray.

◆ SetCells() [3/6]

void vtkUnstructuredGrid::SetCells ( vtkUnsignedCharArray cellTypes,
vtkCellArray cells 
)

Provide cell information to define the dataset.

Cells like vtkPolyhedron require points plus a list of faces. To handle vtkPolyhedron, use SetPolyhedralCells() SetPolyhedralCells requires a faces vtkCellArray that will describe the faces use by the polyhedral cells. SetPolyhedralCells also requires a faceLocations vtkCellArray to fully describe a polyhedron cell The faceLocations is a collection of face ids pointing to the faces vtkCellArray.

◆ SetPolyhedralCells()

void vtkUnstructuredGrid::SetPolyhedralCells ( vtkUnsignedCharArray cellTypes,
vtkCellArray cells,
vtkCellArray faceLocations,
vtkCellArray faces 
)

Provide cell information to define the dataset.

Cells like vtkPolyhedron require points plus a list of faces. To handle vtkPolyhedron, use SetPolyhedralCells() SetPolyhedralCells requires a faces vtkCellArray that will describe the faces use by the polyhedral cells. SetPolyhedralCells also requires a faceLocations vtkCellArray to fully describe a polyhedron cell The faceLocations is a collection of face ids pointing to the faces vtkCellArray.

◆ SetCells() [4/6]

void vtkUnstructuredGrid::SetCells ( vtkUnsignedCharArray cellTypes,
vtkCellArray cells,
vtkIdTypeArray faceLocations,
vtkIdTypeArray faces 
)

Provide cell information to define the dataset.

Cells like vtkPolyhedron require points plus a list of faces. To handle vtkPolyhedron, use SetPolyhedralCells() SetPolyhedralCells requires a faces vtkCellArray that will describe the faces use by the polyhedral cells. SetPolyhedralCells also requires a faceLocations vtkCellArray to fully describe a polyhedron cell The faceLocations is a collection of face ids pointing to the faces vtkCellArray.

◆ GetCells()

vtkCellArray* vtkUnstructuredGrid::GetCells ( )
inline

Return the unstructured grid connectivity array.

Definition at line 409 of file vtkUnstructuredGrid.h.

◆ GetCellNeighbors() [1/2]

void vtkUnstructuredGrid::GetCellNeighbors ( vtkIdType  cellId,
vtkIdList ptIds,
vtkIdList cellIds 
)
inlineoverridevirtual

A topological inquiry to retrieve all of the cells using list of points exclusive of the current cell specified (e.g., cellId).

THIS METHOD IS THREAD SAFE IF FIRST CALLED FROM A SINGLE THREAD AND THE DATASET IS NOT MODIFIED.

Reimplemented from vtkDataSet.

Definition at line 418 of file vtkUnstructuredGrid.h.

◆ GetCellNeighbors() [2/2]

void vtkUnstructuredGrid::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 specified (e.g., cellId).

THIS METHOD IS THREAD SAFE IF FIRST CALLED FROM A SINGLE THREAD AND THE DATASET IS NOT MODIFIED.

◆ IsCellBoundary() [1/3]

bool vtkUnstructuredGrid::IsCellBoundary ( vtkIdType  cellId,
vtkIdType  npts,
const vtkIdType ptIds,
vtkIdType neighborCellId 
)

A topological inquiry to determine whether a topological entity (e.g., point, edge, or face) defined by the point ids (ptIds of length npts) is a boundary entity of a specified cell (indicated by cellId).

A boundary entity is a topological feature used by exactly one cell. This method is related to GetCellNeighbors() except that it simply indicates whether a topological feature is boundary - hence the method is faster. THIS METHOD IS THREAD SAFE IF FIRST CALLED FROM A SINGLE THREAD AND THE DATASET IS NOT MODIFIED.

◆ IsCellBoundary() [2/3]

bool vtkUnstructuredGrid::IsCellBoundary ( vtkIdType  cellId,
vtkIdType  npts,
const vtkIdType ptIds 
)
inline

A topological inquiry to determine whether a topological entity (e.g., point, edge, or face) defined by the point ids (ptIds of length npts) is a boundary entity of a specified cell (indicated by cellId).

A boundary entity is a topological feature used by exactly one cell. This method is related to GetCellNeighbors() except that it simply indicates whether a topological feature is boundary - hence the method is faster. THIS METHOD IS THREAD SAFE IF FIRST CALLED FROM A SINGLE THREAD AND THE DATASET IS NOT MODIFIED.

Definition at line 439 of file vtkUnstructuredGrid.h.

◆ IsCellBoundary() [3/3]

bool vtkUnstructuredGrid::IsCellBoundary ( vtkIdType  cellId,
vtkIdType  npts,
const vtkIdType ptIds,
vtkIdList vtkNotUsedcellIds 
)
inline

A topological inquiry to determine whether a topological entity (e.g., point, edge, or face) defined by the point ids (ptIds of length npts) is a boundary entity of a specified cell (indicated by cellId).

A boundary entity is a topological feature used by exactly one cell. This method is related to GetCellNeighbors() except that it simply indicates whether a topological feature is boundary - hence the method is faster. THIS METHOD IS THREAD SAFE IF FIRST CALLED FROM A SINGLE THREAD AND THE DATASET IS NOT MODIFIED.

Definition at line 445 of file vtkUnstructuredGrid.h.

◆ InsertNextLinkedCell()

vtkIdType vtkUnstructuredGrid::InsertNextLinkedCell ( int  type,
int  npts,
const vtkIdType  pts[] 
)

Use these methods only if the dataset has been specified as Editable.

See vtkPointSet for more information.

◆ RemoveReferenceToCell()

void vtkUnstructuredGrid::RemoveReferenceToCell ( vtkIdType  ptId,
vtkIdType  cellId 
)

Use these methods only if the dataset has been specified as Editable.

See vtkPointSet for more information.

◆ AddReferenceToCell()

void vtkUnstructuredGrid::AddReferenceToCell ( vtkIdType  ptId,
vtkIdType  cellId 
)

Use these methods only if the dataset has been specified as Editable.

See vtkPointSet for more information.

◆ ResizeCellList()

void vtkUnstructuredGrid::ResizeCellList ( vtkIdType  ptId,
int  size 
)

Use these methods only if the dataset has been specified as Editable.

See vtkPointSet for more information.

◆ GetPiece()

virtual int vtkUnstructuredGrid::GetPiece ( )
virtual

Set / Get the piece and the number of pieces.

Similar to extent in 3D.

◆ GetNumberOfPieces()

virtual int vtkUnstructuredGrid::GetNumberOfPieces ( )
virtual

Set / Get the piece and the number of pieces.

Similar to extent in 3D.

◆ GetGhostLevel()

virtual int vtkUnstructuredGrid::GetGhostLevel ( )
virtual

Get the ghost level.

◆ GetActualMemorySize()

unsigned long vtkUnstructuredGrid::GetActualMemorySize ( )
overridevirtual

Return the actual size of the data in kibibytes (1024 bytes).

This number is valid only after the pipeline has updated. The memory size returned is guaranteed to be greater than or equal to the memory required to represent the data (e.g., extra space in arrays, etc. are not included in the return value). THIS METHOD IS THREAD SAFE.

Reimplemented from vtkPointSet.

◆ ShallowCopy()

void vtkUnstructuredGrid::ShallowCopy ( vtkDataObject src)
overridevirtual

Shallow and Deep copy.

Reimplemented from vtkPointSet.

◆ DeepCopy()

void vtkUnstructuredGrid::DeepCopy ( vtkDataObject src)
overridevirtual

Shallow and Deep copy.

Reimplemented from vtkPointSet.

◆ GetIdsOfCellsOfType()

void vtkUnstructuredGrid::GetIdsOfCellsOfType ( int  type,
vtkIdTypeArray array 
)
overridevirtual

Fill vtkIdTypeArray container with list of cell Ids.

This method traverses all cells and, for a particular cell type, inserts the cell Id into the container.

Implements vtkUnstructuredGridBase.

◆ IsHomogeneous()

int vtkUnstructuredGrid::IsHomogeneous ( )
overridevirtual

Returns whether cells are all of the same type.

Implements vtkUnstructuredGridBase.

◆ RemoveGhostCells()

void vtkUnstructuredGrid::RemoveGhostCells ( )

This method will remove any cell that is marked as ghost (has the vtkDataSetAttributes::DUPLICATECELL or the vtkDataSetAttributes::HIDDENCELL bit set).

◆ GetData() [1/2]

static vtkUnstructuredGrid* vtkUnstructuredGrid::GetData ( vtkInformation info)
static

Retrieve an instance of this class from an information object.

◆ GetData() [2/2]

static vtkUnstructuredGrid* vtkUnstructuredGrid::GetData ( vtkInformationVector v,
int  i = 0 
)
static

Retrieve an instance of this class from an information object.

◆ GetFaces() [1/2]

vtkIdType* vtkUnstructuredGrid::GetFaces ( vtkIdType  cellId)

Special support for polyhedron.

Return nullptr for all other cell types.

◆ GetPolyhedronFaces() [1/2]

void vtkUnstructuredGrid::GetPolyhedronFaces ( vtkIdType  cellId,
vtkCellArray faces 
)

Special support for polyhedron.

Do not handle all other cell types.

◆ GetFaces() [2/2]

vtkIdTypeArray* vtkUnstructuredGrid::GetFaces ( )

Get pointer to faces and facelocations.

Support for polyhedron cells. Use an internal cache to handle legacy layout

◆ GetFaceLocations()

vtkIdTypeArray* vtkUnstructuredGrid::GetFaceLocations ( )

Get pointer to faces and facelocations.

Support for polyhedron cells. Use an internal cache to handle legacy layout

◆ GetPolyhedronFaces() [2/2]

vtkCellArray* vtkUnstructuredGrid::GetPolyhedronFaces ( )

Get pointer to faces and facelocations for polyhedron cells.

This is a direct access to internal vtkCellArray structures without any copy.

◆ GetPolyhedronFaceLocations()

vtkCellArray* vtkUnstructuredGrid::GetPolyhedronFaceLocations ( )

Get pointer to faces and facelocations.

Support for polyhedron cells. Use an internal cache to handle legacy layout

◆ InitializeFacesRepresentation()

int vtkUnstructuredGrid::InitializeFacesRepresentation ( vtkIdType  numPrevCells)

Special function used by vtkUnstructuredGridReader.

By default vtkUnstructuredGrid does not contain face information, which is only used by polyhedron cells. If so far no polyhedron cells have been added, Faces and FaceLocations pointers will be nullptr. In this case, need to initialize the arrays and assign values to the previous non-polyhedron cells.

◆ GetMeshMTime()

virtual vtkMTimeType vtkUnstructuredGrid::GetMeshMTime ( )
virtual

Return the mesh (geometry/topology) modification time.

This time is different from the usual MTime which also takes into account the modification of data arrays. This function can be used to track the changes on the mesh separately from the data arrays (eg. static mesh over time with transient data).

◆ DecomposeAPolyhedronCell() [1/8]

static void vtkUnstructuredGrid::DecomposeAPolyhedronCell ( vtkCellArray polyhedronCellArray,
vtkIdType nCellpts,
vtkIdType nCellfaces,
vtkCellArray cellArray,
vtkIdTypeArray faces 
)
static

A static method for converting a polyhedron vtkCellArray of format [nCellFaces, nFace0Pts, i, j, k, nFace1Pts, i, j, k, ...] into three components: (1) an integer indicating the number of faces (2) a standard vtkCellArray storing point ids [nCell0Pts, i, j, k] and (3) an vtkIdTypeArray storing face connectivity in format [nFace0Pts, i, j, k, nFace1Pts, i, j, k, ...] Note: input is assumed to contain only one polyhedron cell.

Outputs (2) and (3) will be stacked at the end of the input cellArray and faces. The original data in the input will not be touched.

◆ DecomposeAPolyhedronCell() [2/8]

static void vtkUnstructuredGrid::DecomposeAPolyhedronCell ( const vtkIdType polyhedronCellStream,
vtkIdType nCellpts,
vtkIdType nCellfaces,
vtkCellArray cellArray,
vtkIdTypeArray faces 
)
static

◆ DecomposeAPolyhedronCell() [3/8]

static void vtkUnstructuredGrid::DecomposeAPolyhedronCell ( vtkIdType  nCellFaces,
const vtkIdType inFaceStream,
vtkIdType nCellpts,
vtkCellArray cellArray,
vtkIdTypeArray faces 
)
static

A static method for converting an input polyhedron cell stream of format [nFace0Pts, i, j, k, nFace1Pts, i, j, k, ...] into three components: (1) an integer indicating the number of faces (2) a standard vtkCellArray storing point ids [nCell0Pts, i, j, k] and (3) an vtkIdTypeArray storing face connectivity in format [nFace0Pts, i, j, k, nFace1Pts, i, j, k, ...] Note: input is assumed to contain only one polyhedron cell.

Outputs (2) and (3) will be stacked at the end of the input cellArray and faces. The original data in the input will not be touched.

◆ ConvertFaceStreamPointIds() [1/3]

static void vtkUnstructuredGrid::ConvertFaceStreamPointIds ( vtkIdList faceStream,
vtkIdType idMap 
)
static

Convert pid in a face stream into idMap[pid].

The face stream is of format [nCellFaces, nFace0Pts, i, j, k, nFace1Pts, i, j, k, ...]. The user is responsible to make sure all the Ids in faceStream do not exceed the range of idMap.

◆ ConvertFaceStreamPointIds() [2/3]

static void vtkUnstructuredGrid::ConvertFaceStreamPointIds ( vtkIdType  nfaces,
vtkIdType faceStream,
vtkIdType idMap 
)
static

Convert pid in a face stream into idMap[pid].

The face stream is of format [nFace0Pts, i, j, k, nFace1Pts, i, j, k, ...]. The user is responsible to make sure all the Ids in faceStream do not exceed the range of idMap.

◆ ConvertFaceStreamPointIds() [3/3]

static void vtkUnstructuredGrid::ConvertFaceStreamPointIds ( vtkCellArray faces,
vtkIdType idMap 
)
static

Convert pid in a face stream into idMap[pid].

The face stream is of format vtkCellArray. The user is responsible to make sure all the Ids in faceStream do not exceed the range of idMap.

◆ GetCellLocationsArray()

vtkIdTypeArray* vtkUnstructuredGrid::GetCellLocationsArray ( )

Get the array of all the starting indices of cell definitions in the cell array.

Warning
vtkCellArray supports random access now. This array is no longer used.

◆ SetCells() [5/6]

void vtkUnstructuredGrid::SetCells ( vtkUnsignedCharArray cellTypes,
vtkIdTypeArray cellLocations,
vtkCellArray cells 
)

Special methods specific to vtkUnstructuredGrid for defining the cells composing the dataset.

Most cells require just arrays of cellTypes, cellLocations and cellConnectivities which implicitly define the set of points in each cell and their ordering. In those cases the cellConnectivities are of the format (numFace0Pts, id1, id2, id3, numFace1Pts, id1, id2, id3...). However, some cells like vtkPolyhedron require points plus a list of faces. To handle vtkPolyhedron, SetCells() support a special input cellConnectivities format (numCellFaces, numFace0Pts, id1, id2, id3, numFace1Pts,id1, id2, id3, ...) The functions use vtkPolyhedron::DecomposeAPolyhedronCell() to convert polyhedron cells into standard format.

Warning
The cellLocations array is no longer used; this information is stored in vtkCellArray. Use the other SetCells overloads.

◆ SetCells() [6/6]

void vtkUnstructuredGrid::SetCells ( vtkUnsignedCharArray cellTypes,
vtkIdTypeArray cellLocations,
vtkCellArray cells,
vtkIdTypeArray faceLocations,
vtkIdTypeArray faces 
)

Special methods specific to vtkUnstructuredGrid for defining the cells composing the dataset.

Most cells require just arrays of cellTypes, cellLocations and cellConnectivities which implicitly define the set of points in each cell and their ordering. In those cases the cellConnectivities are of the format (numFace0Pts, id1, id2, id3, numFace1Pts, id1, id2, id3...). However, some cells like vtkPolyhedron require points plus a list of faces. To handle vtkPolyhedron, SetCells() support a special input cellConnectivities format (numCellFaces, numFace0Pts, id1, id2, id3, numFace1Pts,id1, id2, id3, ...) The functions use vtkPolyhedron::DecomposeAPolyhedronCell() to convert polyhedron cells into standard format.

Warning
The cellLocations array is no longer used; this information is stored in vtkCellArray. Use the other SetCells overloads.

◆ ReportReferences()

void vtkUnstructuredGrid::ReportReferences ( vtkGarbageCollector )
overrideprotectedvirtual

Reimplemented from vtkPointSet.

◆ InternalInsertNextCell() [1/4]

vtkIdType vtkUnstructuredGrid::InternalInsertNextCell ( int  type,
vtkIdType  npts,
const vtkIdType  ptIds[] 
)
overrideprotectedvirtual

◆ InternalInsertNextCell() [2/4]

vtkIdType vtkUnstructuredGrid::InternalInsertNextCell ( int  type,
vtkIdList ptIds 
)
overrideprotectedvirtual

◆ InternalInsertNextCell() [3/4]

vtkIdType vtkUnstructuredGrid::InternalInsertNextCell ( int  type,
vtkIdType  npts,
const vtkIdType  ptIds[],
vtkIdType  nfaces,
const vtkIdType  faces[] 
)
protected

◆ InternalInsertNextCell() [4/4]

vtkIdType vtkUnstructuredGrid::InternalInsertNextCell ( int  type,
vtkIdType  npts,
const vtkIdType  pts[],
vtkCellArray faces 
)
overrideprotectedvirtual

◆ InternalReplaceCell()

void vtkUnstructuredGrid::InternalReplaceCell ( vtkIdType  cellId,
int  npts,
const vtkIdType  pts[] 
)
overrideprotectedvirtual

◆ DecomposeAPolyhedronCell() [4/8]

static void vtkUnstructuredGrid::DecomposeAPolyhedronCell ( const vtkIdType cellStream,
vtkIdType numCellPts,
vtkIdType nCellFaces,
vtkCellArray cellArray,
vtkCellArray faces 
)
staticprotected

A static method for converting an input polyhedron cell stream of format [nFace0Pts, i, j, k, nFace1Pts, i, j, k, ...] into three components: (1) an integer indicating the number of faces (2) a standard vtkCellArray storing point ids [nCell0Pts, i, j, k] and (3) an vtkIdTypeArray storing face connectivity in format [nFace0Pts, i, j, k, nFace1Pts, i, j, k, ...] Note: input is assumed to contain only one polyhedron cell.

Outputs (2) and (3) will be stacked at the end of the input cellArray and faces. The original data in the input will not be touched.

◆ DecomposeAPolyhedronCell() [5/8]

static void vtkUnstructuredGrid::DecomposeAPolyhedronCell ( vtkIdType  nCellFaces,
const vtkIdType cellStream,
vtkIdType numCellPts,
vtkCellArray cellArray,
vtkCellArray facesArray 
)
staticprotected

◆ DecomposeAPolyhedronCell() [6/8]

static void vtkUnstructuredGrid::DecomposeAPolyhedronCell ( vtkCellArray polyhedronCell,
vtkIdType numCellPts,
vtkIdType nCellfaces,
vtkCellArray cellArray,
vtkCellArray faces 
)
staticprotected

◆ DecomposeAPolyhedronCell() [7/8]

static void vtkUnstructuredGrid::DecomposeAPolyhedronCell ( const vtkIdType cellStream,
vtkIdType numCellPts,
vtkIdType nCellFaces,
vtkCellArray cellArray,
vtkCellArray faces,
vtkCellArray faceLocations 
)
staticprotected

◆ DecomposeAPolyhedronCell() [8/8]

static void vtkUnstructuredGrid::DecomposeAPolyhedronCell ( vtkIdType  nCellFaces,
const vtkIdType cellStream,
vtkIdType numCellPts,
vtkCellArray cellArray,
vtkCellArray faces,
vtkCellArray faceLocations 
)
staticprotected

◆ CopyPolyhedronToFaceStream()

static int vtkUnstructuredGrid::CopyPolyhedronToFaceStream ( vtkCellArray faceArray,
vtkCellArray faceLocationArray,
vtkIdTypeArray faceStream,
vtkIdTypeArray faceLocation 
)
staticprotected

Backward compatibility function to convert new polyhedron storage to legacy.

◆ GetCell() [3/5]

virtual vtkCell* vtkDataSet::GetCell

Get cell with cellId such that: 0 <= cellId < NumberOfCells.

The returned vtkCell is an object owned by this instance, hence the return value must not be deleted by the caller.

Warning
Repeat calls to this function for different face ids will change the data stored in the internal member object whose pointer is returned by this function.
THIS METHOD IS NOT THREAD SAFE. For a thread-safe version, please use void GetCell(vtkIdType cellId, vtkGenericCell* cell).

◆ GetCell() [4/5]

virtual vtkCell* vtkDataSet::GetCell
inline

Standard vtkDataSet methods; see vtkDataSet.h for documentation.

Definition at line 238 of file vtkDataSet.h.

◆ GetCell() [5/5]

virtual void vtkDataSet::GetCell

Get cell with cellId such that: 0 <= cellId < NumberOfCells.

This is a thread-safe alternative to the previous GetCell() method. THIS METHOD IS THREAD SAFE IF FIRST CALLED FROM A SINGLE THREAD AND THE DATASET IS NOT MODIFIED

Member Data Documentation

◆ Connectivity

vtkSmartPointer<vtkCellArray> vtkUnstructuredGrid::Connectivity
protected

Definition at line 674 of file vtkUnstructuredGrid.h.

◆ Links

vtkSmartPointer<vtkAbstractCellLinks> vtkUnstructuredGrid::Links
protected

Definition at line 675 of file vtkUnstructuredGrid.h.

◆ Types

vtkSmartPointer<vtkUnsignedCharArray> vtkUnstructuredGrid::Types
protected

Definition at line 676 of file vtkUnstructuredGrid.h.

◆ DistinctCellTypes

vtkSmartPointer<vtkCellTypes> vtkUnstructuredGrid::DistinctCellTypes
protected

Definition at line 679 of file vtkUnstructuredGrid.h.

◆ DistinctCellTypesUpdateMTime

vtkMTimeType vtkUnstructuredGrid::DistinctCellTypesUpdateMTime
protected

Definition at line 683 of file vtkUnstructuredGrid.h.

◆ Faces

vtkSmartPointer<vtkCellArray> vtkUnstructuredGrid::Faces
protected

Special support for polyhedra/cells with explicit face representations.

The Faces class represents polygonal faces using a vtkCellArray structure. The FaceLocations store a polyhedron as a list of faces defined in Faces using a vtkCellArray structure.

Definition at line 691 of file vtkUnstructuredGrid.h.

◆ FaceLocations

vtkSmartPointer<vtkCellArray> vtkUnstructuredGrid::FaceLocations
protected

Definition at line 692 of file vtkUnstructuredGrid.h.

◆ CellLocations

vtkSmartPointer<vtkIdTypeArray> vtkUnstructuredGrid::CellLocations
protected

Definition at line 695 of file vtkUnstructuredGrid.h.

◆ LegacyFaces

vtkSmartPointer<vtkIdTypeArray> vtkUnstructuredGrid::LegacyFaces
protected

Legacy support – stores the old-style Faces && FaceLocations Special support for polyhedra/cells with explicit face representations.

The Faces class represents polygonal faces using a modified vtkCellArray structure. Each cell face list begins with the total number of faces in the cell, followed by a vtkCellArray data organization (n,i,j,k,n,i,j,k,...).

Warning
The Faces and FaceLocations arrays are no longer used; this information is stored in vtkCellArrays ElementFaces and ElementFaceLocations. Use SetPolyhedralCells.

Definition at line 717 of file vtkUnstructuredGrid.h.

◆ LegacyFaceLocations

vtkSmartPointer<vtkIdTypeArray> vtkUnstructuredGrid::LegacyFaceLocations
protected

Definition at line 718 of file vtkUnstructuredGrid.h.

◆ LegacyPointIdsBuffer

vtkSmartPointer<vtkIdList> vtkUnstructuredGrid::LegacyPointIdsBuffer
protected

Legacy backward compatibility for GetFaceStream This member should be removed simultaneously with the deprecated GetFaceStream.

Definition at line 725 of file vtkUnstructuredGrid.h.


The documentation for this class was generated from the following file: