VTK
|
Allows datasets with arbitrary storage layouts to be used with VTK. More...
#include <vtkMappedUnstructuredGrid.h>
Allows datasets with arbitrary storage layouts to be used with VTK.
This class fulfills the vtkUnstructuredGridBase API while delegating to a arbitrary implementation of the dataset topology. The purpose of vtkMappedUnstructuredGrid is to allow external data structures to be used directly in a VTK pipeline, e.g. for in-situ analysis of a running simulation.
When introducing an external data structure into VTK, there are 3 principle components of the dataset to consider:
Points and attributes can be handled by subclassing vtkMappedDataArray and implementing that interface to adapt the external data structures into VTK. The vtkMappedDataArray subclasses can then be used as the vtkPoints's Data member (for points/nodes) or added directly to vtkPointData, vtkCellData, or vtkFieldData for attribute information. Filters used in the pipeline will need to be modified to remove calls to vtkDataArray::GetVoidPointer and use vtkTypedDataArrayIterators instead. See vtkDataArrayIteratorMacro.h for a (relatively) simple way to write processing algorithms that will use efficient raw memory accesses for standard VTK data arrays and safe iterators for non-standard data arrays in a single templated implementation.
Introducing an arbitrary topology implementation into VTK requires the use of the vtkMappedUnstructuredGrid class. Unlike the data array counterpart, the mapped unstructured grid is not subclassed, but rather takes an adaptor class as a template argument. This is to allow cheap shallow copies of the data by passing the reference-counted implementation object to new instances of vtkMappedUnstructuredGrid.
The implementation class should derive from vtkObject (for reference counting) and implement the usual vtkObject API requirements, such as a static New() method and PrintSelf function. The following methods must also be implemented:
These methods should provide the same functionality as defined in vtkUnstructuredGrid. See that class's documentation for more information.
Note that since the implementation class is used as a compile-time template parameter in vtkMappedUnstructuredGrid, the above methods do not need be virtuals. The compiler will statically bind the calls, making dynamic vtable lookups unneccessary and giving a slight performance boost.
Adapting a filter or algorithm to safely traverse the vtkMappedUnstructuredGrid's topology requires removing calls the following implementation-dependent vtkUnstructuredGrid methods:
A custom vtkCellIterator implementation may be specified for a particular vtkMappedUnstructuredGrid as the second template parameter. By default, vtkMappedUnstructuredGridCellIterator will be used, which increments an internal cell id counter and performs random-access lookup as-needed. More efficient implementations may be used with data structures better suited for sequential access, see vtkUnstructuredGridCellIterator for an example.
A set of four macros are provided to generate a concrete subclass of vtkMappedUnstructuredGrid with a specified implementation, cell iterator, and export declaration. They are:
To instantiate a vtkMappedUnstructuredGrid subclass created by the above macro, the follow pattern is encouraged:
MyGrid.h: ---------------------------------------------------------------------- class MyGridImplementation : public vtkObject { public: ... (vtkObject required API) ... ... (vtkMappedUnstructuredGrid Implementation required API) ... void SetImplementationDetails(...raw data from external source...); }; vtkMakeMappedUnstructuredGrid(MyGrid, MyGridImplementation) SomeSource.cxx ---------------------------------------------------------------------- vtkNew<MyGrid> grid; grid->GetImplementation()->SetImplementationDetails(...); /* grid is now ready to use. */
The vtkCPExodusIIElementBlock class provides an example of vtkMappedUnstructuredGrid usage, adapting the Exodus II data structures for the VTK pipeline.
Definition at line 154 of file vtkMappedUnstructuredGrid.h.
typedef vtkTypeTemplate<vtkMappedUnstructuredGrid<Implementation, CellIterator>, vtkUnstructuredGridBase> vtkMappedUnstructuredGrid< Implementation, CellIterator >::Superclass |
Reimplemented from vtkTypeTemplate< vtkMappedUnstructuredGrid< Implementation, CellIterator >, vtkUnstructuredGridBase >.
Definition at line 161 of file vtkMappedUnstructuredGrid.h.
typedef Implementation vtkMappedUnstructuredGrid< Implementation, CellIterator >::ImplementationType |
Definition at line 162 of file vtkMappedUnstructuredGrid.h.
typedef CellIterator vtkMappedUnstructuredGrid< Implementation, CellIterator >::CellIteratorType |
Definition at line 163 of file vtkMappedUnstructuredGrid.h.
typedef vtkMappedUnstructuredGrid<Implementation, CellIterator> vtkMappedUnstructuredGrid< Implementation, CellIterator >::ThisType [protected] |
Definition at line 195 of file vtkMappedUnstructuredGrid.h.
vtkMappedUnstructuredGrid< Implementation, CellIterator >::vtkMappedUnstructuredGrid | ( | ) | [protected] |
vtkMappedUnstructuredGrid< Implementation, CellIterator >::~vtkMappedUnstructuredGrid | ( | ) | [protected] |
void vtkMappedUnstructuredGrid< Implementation, CellIterator >::PrintSelf | ( | ostream & | os, |
vtkIndent | indent | ||
) | [virtual] |
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 vtkPointSet.
void vtkMappedUnstructuredGrid< Implementation, CellIterator >::CopyStructure | ( | vtkDataSet * | pd | ) | [virtual] |
Copy the geometric structure of an input point set object.
Reimplemented from vtkPointSet.
void vtkMappedUnstructuredGrid< Implementation, CellIterator >::ShallowCopy | ( | vtkDataObject * | src | ) | [virtual] |
Shallow and Deep copy.
Reimplemented from vtkPointSet.
vtkIdType vtkMappedUnstructuredGrid< Implementation, CellIterator >::GetNumberOfCells | ( | ) | [virtual] |
Determine the number of cells composing the dataset. THIS METHOD IS THREAD SAFE
Implements vtkDataSet.
vtkCell* vtkMappedUnstructuredGrid< Implementation, CellIterator >::GetCell | ( | vtkIdType | cellId | ) | [virtual] |
Get cell with cellId such that: 0 <= cellId < NumberOfCells. THIS METHOD IS NOT THREAD SAFE.
Implements vtkDataSet.
void vtkMappedUnstructuredGrid< Implementation, CellIterator >::GetCell | ( | vtkIdType | cellId, |
vtkGenericCell * | cell | ||
) | [virtual] |
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
Implements vtkDataSet.
int vtkMappedUnstructuredGrid< Implementation, CellIterator >::GetCellType | ( | vtkIdType | cellId | ) | [virtual] |
Get type of cell with cellId such that: 0 <= cellId < NumberOfCells. THIS METHOD IS THREAD SAFE IF FIRST CALLED FROM A SINGLE THREAD AND THE DATASET IS NOT MODIFIED
Implements vtkDataSet.
void vtkMappedUnstructuredGrid< Implementation, CellIterator >::GetCellPoints | ( | vtkIdType | cellId, |
vtkIdList * | ptIds | ||
) | [virtual] |
Topological inquiry to get points defining cell. THIS METHOD IS THREAD SAFE IF FIRST CALLED FROM A SINGLE THREAD AND THE DATASET IS NOT MODIFIED
Implements vtkDataSet.
vtkCellIterator* vtkMappedUnstructuredGrid< Implementation, CellIterator >::NewCellIterator | ( | ) | [virtual] |
Return an iterator that traverses the cells in this data set.
Reimplemented from vtkPointSet.
void vtkMappedUnstructuredGrid< Implementation, CellIterator >::GetPointCells | ( | vtkIdType | ptId, |
vtkIdList * | cellIds | ||
) | [virtual] |
Topological inquiry to get cells using point. THIS METHOD IS THREAD SAFE IF FIRST CALLED FROM A SINGLE THREAD AND THE DATASET IS NOT MODIFIED
Implements vtkDataSet.
int vtkMappedUnstructuredGrid< Implementation, CellIterator >::GetMaxCellSize | ( | ) | [virtual] |
Convenience method returns largest cell size in dataset. This is generally used to allocate memory for supporting data structures. THIS METHOD IS THREAD SAFE
Implements vtkDataSet.
void vtkMappedUnstructuredGrid< Implementation, CellIterator >::GetIdsOfCellsOfType | ( | int | type, |
vtkIdTypeArray * | array | ||
) | [virtual] |
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.
int vtkMappedUnstructuredGrid< Implementation, CellIterator >::IsHomogeneous | ( | ) | [virtual] |
Traverse cells and determine if cells are all of the same type.
Implements vtkUnstructuredGridBase.
void vtkMappedUnstructuredGrid< Implementation, CellIterator >::Allocate | ( | vtkIdType | numCells, |
int | extSize = 1000 |
||
) |
vtkIdType vtkMappedUnstructuredGrid< Implementation, CellIterator >::InsertNextCell | ( | int | type, |
vtkIdList * | ptIds | ||
) | [virtual] |
Insert/create cell in object by a list of point ids defining cell topology. Most cells require just a type which implicitly defines a set of points and their ordering. For non-polyhedron cell type, ptIds is the list of global Ids of unique cell points. For polyhedron cell, a special ptIds input format is required: (numCellFaces, numFace0Pts, id1, id2, id3, numFace1Pts,id1, id2, id3, ...)
Implements vtkUnstructuredGridBase.
vtkIdType vtkMappedUnstructuredGrid< Implementation, CellIterator >::InsertNextCell | ( | int | type, |
vtkIdType | npts, | ||
vtkIdType * | ptIds | ||
) | [virtual] |
Insert/create cell in object by type and list of point ids defining cell topology. Most cells require just a type which implicitly defines a set of points and their ordering. For non-polyhedron cell type, npts is the number of unique points in the cell. pts are the list of global point Ids. For polyhedron cell, a special input format is required. npts is the number of faces in the cell. ptIds is the list of face stream: (numFace0Pts, id1, id2, id3, numFace1Pts,id1, id2, id3, ...)
Implements vtkUnstructuredGridBase.
vtkIdType vtkMappedUnstructuredGrid< Implementation, CellIterator >::InsertNextCell | ( | int | type, |
vtkIdType | npts, | ||
vtkIdType * | ptIds, | ||
vtkIdType | nfaces, | ||
vtkIdType * | faces | ||
) | [virtual] |
Implements vtkUnstructuredGridBase.
void vtkMappedUnstructuredGrid< Implementation, CellIterator >::ReplaceCell | ( | vtkIdType | cellId, |
int | npts, | ||
vtkIdType * | pts | ||
) | [virtual] |
Replace the points defining cell "cellId" with a new set of points. This operator is (typically) used when links from points to cells have not been built (i.e., BuildLinks() has not been executed). Use the operator ReplaceLinkedCell() to replace a cell when cell structure has been built.
Implements vtkUnstructuredGridBase.
unsigned long vtkMappedUnstructuredGrid< Implementation, CellIterator >::GetMTime | ( | ) | [virtual] |
Get MTime which also considers its vtkPoints MTime.
Reimplemented from vtkPointSet.
void vtkMappedUnstructuredGrid< Implementation, CellIterator >::SetImplementation | ( | ImplementationType * | impl | ) |
ImplementationType* vtkMappedUnstructuredGrid< Implementation, CellIterator >::GetImplementation | ( | ) |
vtkSmartPointer<ImplementationType> vtkMappedUnstructuredGrid< Implementation, CellIterator >::Impl [protected] |
Definition at line 197 of file vtkMappedUnstructuredGrid.h.