VTK  9.1.0
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
vtkCellIterator Class Referenceabstract

Efficient cell iterator for vtkDataSet topologies. More...

#include <vtkCellIterator.h>

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

Public Member Functions

void PrintSelf (ostream &os, vtkIndent indent) override
 Methods invoked by print to print information about the object including superclasses.
 
 vtkAbstractTypeMacro (vtkCellIterator, vtkObject)
 
void InitTraversal ()
 Reset to the first cell.
 
void GoToNextCell ()
 Increment to next cell.
 
virtual bool IsDoneWithTraversal ()=0
 Returns false while the iterator is valid.
 
int GetCellType ()
 Get the current cell type (e.g.
 
int GetCellDimension ()
 Get the current cell dimension (0, 1, 2, or 3).
 
virtual vtkIdType GetCellId ()=0
 Get the id of the current cell.
 
vtkIdListGetPointIds ()
 Get the ids of the points in the current cell.
 
vtkPointsGetPoints ()
 Get the points in the current cell.
 
vtkIdListGetFaces ()
 Get the faces for a polyhedral cell.
 
void GetCell (vtkGenericCell *cell)
 Write the current full cell information into the argument.
 
vtkIdType GetNumberOfPoints ()
 Return the number of points in the current cell.
 
vtkIdType GetNumberOfFaces ()
 Return the number of faces in the current cell.
 
- Public Member Functions inherited from vtkObject
 vtkBaseTypeMacro (vtkObject, vtkObjectBase)
 
virtual void DebugOn ()
 Turn debugging output on.
 
virtual void DebugOff ()
 Turn debugging output off.
 
bool GetDebug ()
 Get the value of the debug flag.
 
void SetDebug (bool debugFlag)
 Set the value of the debug flag.
 
virtual void Modified ()
 Update the modification time for this object.
 
virtual vtkMTimeType GetMTime ()
 Return this object's modified time.
 
void PrintSelf (ostream &os, vtkIndent indent) override
 Methods invoked by print to print information about the object including superclasses.
 
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)
 
int InvokeEvent (unsigned long event)
 
int InvokeEvent (const char *event)
 
unsigned long AddObserver (unsigned long event, vtkCommand *, float priority=0.0f)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object.
 
unsigned long AddObserver (const char *event, vtkCommand *, float priority=0.0f)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object.
 
vtkCommandGetCommand (unsigned long tag)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object.
 
void RemoveObserver (vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object.
 
void RemoveObservers (unsigned long event, vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object.
 
void RemoveObservers (const char *event, vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object.
 
vtkTypeBool HasObserver (unsigned long event, vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object.
 
vtkTypeBool HasObserver (const char *event, vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object.
 
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.
 
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.
 
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.
 
int InvokeEvent (unsigned long event, void *callData)
 This method invokes an event and return whether the event was aborted or not.
 
int InvokeEvent (const char *event, void *callData)
 This method invokes an event and return whether the event was aborted or not.
 
- Public Member Functions inherited from vtkObjectBase
const char * GetClassName () const
 Return the class name as a string.
 
virtual vtkTypeBool IsA (const char *name)
 Return 1 if this class is the same type of (or a subclass of) the named class.
 
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).
 
virtual void Delete ()
 Delete a VTK object.
 
virtual void FastDelete ()
 Delete a reference to this object.
 
void InitializeObjectBase ()
 
void Print (ostream &os)
 Print an object to an ostream.
 
virtual void Register (vtkObjectBase *o)
 Increase the reference count (mark as used by another object).
 
virtual void UnRegister (vtkObjectBase *o)
 Decrease the reference count (release by another object).
 
int GetReferenceCount ()
 Return the current reference count of this object.
 
void SetReferenceCount (int)
 Sets the reference count.
 
bool GetIsInMemkind () const
 A local state flag that remembers whether this object lives in the normal or extended memory space.
 
virtual void PrintHeader (ostream &os, vtkIndent indent)
 Methods invoked by print to print information about the object including superclasses.
 
virtual void PrintTrailer (ostream &os, vtkIndent indent)
 Methods invoked by print to print information about the object including superclasses.
 

Protected Member Functions

 vtkCellIterator ()
 
 ~vtkCellIterator () override
 
virtual void ResetToFirstCell ()=0
 Update internal state to point to the first cell.
 
virtual void IncrementToNextCell ()=0
 Update internal state to point to the next cell.
 
virtual void FetchCellType ()=0
 Lookup the cell type in the data set and store it in this->CellType.
 
virtual void FetchPointIds ()=0
 Lookup the cell point ids in the data set and store them in this->PointIds.
 
virtual void FetchPoints ()=0
 Lookup the cell points in the data set and store them in this->Points.
 
virtual void FetchFaces ()
 Lookup the cell faces in the data set and store them in this->Faces.
 
- 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.
 
void InternalReleaseFocus ()
 These methods allow a command to exclusively grab all events.
 
- Protected Member Functions inherited from vtkObjectBase
 vtkObjectBase ()
 
virtual ~vtkObjectBase ()
 
virtual void RegisterInternal (vtkObjectBase *, vtkTypeBool check)
 
virtual void UnRegisterInternal (vtkObjectBase *, vtkTypeBool check)
 
virtual void ReportReferences (vtkGarbageCollector *)
 
 vtkObjectBase (const vtkObjectBase &)
 
void operator= (const vtkObjectBase &)
 

Protected Attributes

int CellType
 
vtkPointsPoints
 
vtkIdListPointIds
 
vtkIdListFaces
 
- Protected Attributes inherited from vtkObject
bool Debug
 
vtkTimeStamp MTime
 
vtkSubjectHelper * SubjectHelper
 
- Protected Attributes inherited from vtkObjectBase
std::atomic< int32_t > ReferenceCount
 
vtkWeakPointerBase ** WeakPointers
 

Additional Inherited Members

- 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.
 
static void BreakOnError ()
 This method is called when vtkErrorMacro executes.
 
static void SetGlobalWarningDisplay (int val)
 This is a global flag that controls whether any debug, warning or error messages are displayed.
 
static void GlobalWarningDisplayOn ()
 This is a global flag that controls whether any debug, warning or error messages are displayed.
 
static void GlobalWarningDisplayOff ()
 This is a global flag that controls whether any debug, warning or error messages are displayed.
 
static int GetGlobalWarningDisplay ()
 This is a global flag that controls whether any debug, warning or error messages are displayed.
 
- 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.
 
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).
 
static vtkObjectBaseNew ()
 Create an object with Debug turned off, modified time initialized to zero, and reference counting on.
 
static void SetMemkindDirectory (const char *directoryname)
 The name of a directory, ideally mounted -o dax, to memory map an extended memory space within.
 
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.
 
- Static Protected Member Functions inherited from vtkObjectBase
static vtkMallocingFunction GetCurrentMallocFunction ()
 
static vtkReallocingFunction GetCurrentReallocFunction ()
 
static vtkFreeingFunction GetCurrentFreeFunction ()
 
static vtkFreeingFunction GetAlternateFreeFunction ()
 

Detailed Description

Efficient cell iterator for vtkDataSet topologies.

vtkCellIterator provides a method for traversing cells in a data set. Call the vtkDataSet::NewCellIterator() method to use this class.

The cell is represented as a set of three pieces of information: The cell type, the ids of the points constituting the cell, and the points themselves. This iterator fetches these as needed. If only the cell type is used, the type is not looked up until GetCellType is called, and the point information is left uninitialized. This allows efficient screening of cells, since expensive point lookups may be skipped depending on the cell type/etc.

An example usage of this class:

void myWorkerFunction(vtkDataSet *ds)
{
for (it->InitTraversal(); !it->IsDoneWithTraversal(); it->GoToNextCell())
{
if (it->GetCellType() != VTK_TETRA)
{
continue; // Skip non-tetrahedral cells
}
vtkIdList *pointIds = it->GetPointIds();
// Do screening on the point ids, maybe figure out scalar range and skip
cells that do not lie in a certain range?
vtkPoints *points = it->GetPoints();
// Do work using the cell points, or ...
vtkGenericCell *cell = ...;
it->GetCell(cell);
// ... do work with a vtkCell.
}
it->Delete();
}
Efficient cell iterator for vtkDataSet topologies.
void GetCell(vtkGenericCell *cell)
Write the current full cell information into the argument.
void InitTraversal()
Reset to the first cell.
vtkPoints * GetPoints()
Get the points in the current cell.
vtkIdList * GetPointIds()
Get the ids of the points in the current cell.
int GetCellType()
Get the current cell type (e.g.
void GoToNextCell()
Increment to next cell.
virtual bool IsDoneWithTraversal()=0
Returns false while the iterator is valid.
abstract class to specify dataset behavior
Definition: vtkDataSet.h:166
virtual vtkCellIterator * NewCellIterator()
Return an iterator that traverses the cells in this data set.
provides thread-safe access to cells
list of point or cell ids
Definition: vtkIdList.h:140
virtual void Delete()
Delete a VTK object.
represent and manipulate 3D points
Definition: vtkPoints.h:143
@ VTK_TETRA
Definition: vtkCellType.h:95

The example above pulls in bits of information as needed to filter out cells that aren't relevant. The least expensive lookups are performed first (cell type, then point ids, then points/full cell) to prevent wasted cycles fetching unnecessary data. Also note that at the end of the loop, the iterator must be deleted as these iterators are vtkObject subclasses.

Online Examples:

Definition at line 176 of file vtkCellIterator.h.

Constructor & Destructor Documentation

◆ vtkCellIterator()

vtkCellIterator::vtkCellIterator ( )
protected

◆ ~vtkCellIterator()

vtkCellIterator::~vtkCellIterator ( )
overrideprotected

Member Function Documentation

◆ PrintSelf()

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

Methods invoked by print to print information about the object including superclasses.

Typically not called by the user (use Print() instead) but used in the hierarchical print process to combine the output of several classes.

Reimplemented from vtkObject.

Reimplemented in vtkCPExodusIIElementBlockCellIterator, vtkDataSetCellIterator, vtkMappedUnstructuredGridCellIterator< Implementation >, vtkPointSetCellIterator, and vtkUnstructuredGridCellIterator.

◆ vtkAbstractTypeMacro()

vtkCellIterator::vtkAbstractTypeMacro ( vtkCellIterator  ,
vtkObject   
)

◆ InitTraversal()

void vtkCellIterator::InitTraversal ( )
inline

Reset to the first cell.

Definition at line 324 of file vtkCellIterator.h.

◆ GoToNextCell()

void vtkCellIterator::GoToNextCell ( )
inline

Increment to next cell.

Always safe to call.

Definition at line 331 of file vtkCellIterator.h.

◆ IsDoneWithTraversal()

virtual bool vtkCellIterator::IsDoneWithTraversal ( )
pure virtual

Returns false while the iterator is valid.

Always safe to call.

Implemented in vtkDataSetCellIterator, vtkMappedUnstructuredGridCellIterator< Implementation >, vtkPointSetCellIterator, and vtkUnstructuredGridCellIterator.

◆ GetCellType()

int vtkCellIterator::GetCellType ( )
inline

Get the current cell type (e.g.

VTK_LINE, VTK_VERTEX, VTK_TETRA, etc). This should only be called when IsDoneWithTraversal() returns false.

Definition at line 338 of file vtkCellIterator.h.

◆ GetCellDimension()

int vtkCellIterator::GetCellDimension ( )

Get the current cell dimension (0, 1, 2, or 3).

This should only be called when IsDoneWithTraversal() returns false.

◆ GetCellId()

virtual vtkIdType vtkCellIterator::GetCellId ( )
pure virtual

◆ GetPointIds()

vtkIdList * vtkCellIterator::GetPointIds ( )
inline

Get the ids of the points in the current cell.

This should only be called when IsDoneWithTraversal() returns false.

Definition at line 349 of file vtkCellIterator.h.

◆ GetPoints()

vtkPoints * vtkCellIterator::GetPoints ( )
inline

Get the points in the current cell.

This is usually a very expensive call, and should be avoided when possible. This should only be called when IsDoneWithTraversal() returns false.

Definition at line 360 of file vtkCellIterator.h.

◆ GetFaces()

vtkIdList * vtkCellIterator::GetFaces ( )
inline

Get the faces for a polyhedral cell.

This is only valid when CellType is VTK_POLYHEDRON.

Definition at line 371 of file vtkCellIterator.h.

◆ GetCell()

void vtkCellIterator::GetCell ( vtkGenericCell cell)

Write the current full cell information into the argument.

This is usually a very expensive call, and should be avoided when possible. This should only be called when IsDoneWithTraversal() returns false.

◆ GetNumberOfPoints()

vtkIdType vtkCellIterator::GetNumberOfPoints ( )
inline

Return the number of points in the current cell.

This should only be called when IsDoneWithTraversal() returns false.

Definition at line 382 of file vtkCellIterator.h.

◆ GetNumberOfFaces()

vtkIdType vtkCellIterator::GetNumberOfFaces ( )
inline

Return the number of faces in the current cell.

This should only be called when IsDoneWithTraversal() returns false.

Definition at line 393 of file vtkCellIterator.h.

◆ ResetToFirstCell()

virtual void vtkCellIterator::ResetToFirstCell ( )
protectedpure virtual

◆ IncrementToNextCell()

virtual void vtkCellIterator::IncrementToNextCell ( )
protectedpure virtual

◆ FetchCellType()

virtual void vtkCellIterator::FetchCellType ( )
protectedpure virtual

◆ FetchPointIds()

virtual void vtkCellIterator::FetchPointIds ( )
protectedpure virtual

◆ FetchPoints()

virtual void vtkCellIterator::FetchPoints ( )
protectedpure virtual

◆ FetchFaces()

virtual void vtkCellIterator::FetchFaces ( )
inlineprotectedvirtual

Lookup the cell faces in the data set and store them in this->Faces.

Few data sets support faces, so this method has a no-op default implementation. See vtkUnstructuredGrid::GetFaceStream for a description of the layout that Faces should have.

Reimplemented in vtkUnstructuredGridCellIterator.

Definition at line 287 of file vtkCellIterator.h.

Member Data Documentation

◆ CellType

int vtkCellIterator::CellType
protected

Definition at line 289 of file vtkCellIterator.h.

◆ Points

vtkPoints* vtkCellIterator::Points
protected

Definition at line 290 of file vtkCellIterator.h.

◆ PointIds

vtkIdList* vtkCellIterator::PointIds
protected

Definition at line 291 of file vtkCellIterator.h.

◆ Faces

vtkIdList* vtkCellIterator::Faces
protected

Definition at line 292 of file vtkCellIterator.h.


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