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

defines dataset interface More...

#include <vtkGenericDataSet.h>

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

Public Member Functions

virtual vtkIdType GetNumberOfPoints ()=0
 
virtual vtkIdType GetNumberOfCells (int dim=-1)=0
 
virtual int GetCellDimension ()=0
 
virtual void GetCellTypes (vtkCellTypes *types)
 
virtual vtkGenericCellIteratorNewCellIterator (int dim=-1)=0
 
virtual vtkGenericPointIteratorNewPointIterator ()=0
 
virtual unsigned long int GetMTime ()
 
virtual void ComputeBounds ()=0
 
virtual doubleGetBounds ()
 
virtual void GetBounds (double bounds[6])
 
virtual doubleGetCenter ()
 
virtual void GetCenter (double center[3])
 
virtual double GetLength ()
 
virtual unsigned long GetActualMemorySize ()
 
int GetDataObjectType ()
 
virtual vtkIdType GetEstimatedSize ()=0
 
virtual vtkGenericCellIteratorNewBoundaryIterator (int dim=-1, int exteriorOnly=0)=0
 
virtual int FindCell (double x[3], vtkGenericCellIterator *&cell, double tol2, int &subId, double pcoords[3])=0
 
virtual void FindPoint (double x[3], vtkGenericPointIterator *p)=0
 
virtual vtkGenericAttributeCollectionGetAttributes ()
 
virtual vtkDataSetAttributesGetAttributes (int type)
 
virtual void SetTessellator (vtkGenericCellTessellator *tessellator)
 
virtual vtkGenericCellTessellatorGetTessellator ()
 
- Public Member Functions inherited from vtkDataObject
vtkDataObjectNewInstance () const
 
virtual void Initialize ()
 
void ReleaseData ()
 
unsigned long GetUpdateTime ()
 
virtual void CopyInformationToPipeline (vtkInformation *vtkNotUsed(info))
 
void DataHasBeenGenerated ()
 
virtual void PrepareForNewData ()
 
virtual int GetExtentType ()
 
virtual void Crop (const int *updateExtent)
 
virtual vtkFieldDataGetAttributesAsFieldData (int type)
 
virtual int GetAttributeTypeForArray (vtkAbstractArray *arr)
 
virtual vtkIdType GetNumberOfElements (int type)
 
virtual vtkInformationGetInformation ()
 
virtual void SetInformation (vtkInformation *)
 
virtual int GetDataReleased ()
 
virtual void SetFieldData (vtkFieldData *)
 
virtual vtkFieldDataGetFieldData ()
 
virtual void CopyInformationFromPipeline (vtkInformation *vtkNotUsed(info))
 
virtual void ShallowCopy (vtkDataObject *src)
 
virtual void DeepCopy (vtkDataObject *src)
 
void GlobalReleaseDataFlagOn ()
 
void GlobalReleaseDataFlagOff ()
 
- Public Member Functions inherited from vtkObject
vtkObjectNewInstance () const
 
virtual void DebugOn ()
 
virtual void DebugOff ()
 
bool GetDebug ()
 
void SetDebug (bool debugFlag)
 
virtual void Modified ()
 
unsigned long AddObserver (unsigned long event, vtkCommand *, float priority=0.0f)
 
unsigned long AddObserver (const char *event, vtkCommand *, float priority=0.0f)
 
vtkCommandGetCommand (unsigned long tag)
 
void RemoveObserver (vtkCommand *)
 
void RemoveObservers (unsigned long event, vtkCommand *)
 
void RemoveObservers (const char *event, vtkCommand *)
 
int HasObserver (unsigned long event, vtkCommand *)
 
int HasObserver (const char *event, vtkCommand *)
 
void RemoveObserver (unsigned long tag)
 
void RemoveObservers (unsigned long event)
 
void RemoveObservers (const char *event)
 
void RemoveAllObservers ()
 
int HasObserver (unsigned long event)
 
int HasObserver (const char *event)
 
template<class U , class T >
unsigned long AddObserver (unsigned long event, U observer, void(T::*callback)(), float priority=0.0f)
 
template<class U , class T >
unsigned long AddObserver (unsigned long event, U observer, void(T::*callback)(vtkObject *, unsigned long, void *), float priority=0.0f)
 
template<class U , class T >
unsigned long AddObserver (unsigned long event, U observer, bool(T::*callback)(vtkObject *, unsigned long, void *), float priority=0.0f)
 
int InvokeEvent (unsigned long event, void *callData)
 
int InvokeEvent (const char *event, void *callData)
 
int InvokeEvent (unsigned long event)
 
int InvokeEvent (const char *event)
 
- Public Member Functions inherited from vtkObjectBase
const char * GetClassName () const
 
virtual void Delete ()
 
virtual void FastDelete ()
 
void Print (ostream &os)
 
virtual void Register (vtkObjectBase *o)
 
virtual void UnRegister (vtkObjectBase *o)
 
void SetReferenceCount (int)
 
void PrintRevisions (ostream &)
 
virtual void PrintHeader (ostream &os, vtkIndent indent)
 
virtual void PrintTrailer (ostream &os, vtkIndent indent)
 
int GetReferenceCount ()
 

Static Public Member Functions

static vtkGenericDataSetGetData (vtkInformation *info)
 
static vtkGenericDataSetGetData (vtkInformationVector *v, int i=0)
 
- Static Public Member Functions inherited from vtkDataObject
static vtkDataObjectNew ()
 
static int IsTypeOf (const char *type)
 
static vtkDataObjectSafeDownCast (vtkObjectBase *o)
 
static const char * GetAssociationTypeAsString (int associationType)
 
static int GetAssociationTypeFromString (const char *associationType)
 
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 vtkInformationDoubleVectorKeyBOUNDING_BOX ()
 
static vtkInformationDataObjectKeySIL ()
 
static vtkInformationGetActiveFieldInformation (vtkInformation *info, int fieldAssociation, int attributeType)
 
static vtkInformationGetNamedFieldInformation (vtkInformation *info, int fieldAssociation, const char *name)
 
static void RemoveNamedFieldInformation (vtkInformation *info, int fieldAssociation, const char *name)
 
static vtkInformationSetActiveAttribute (vtkInformation *info, int fieldAssociation, const char *attributeName, int attributeType)
 
static void SetActiveAttributeInfo (vtkInformation *info, int fieldAssociation, int attributeType, const char *name, int arrayType, int numComponents, int numTuples)
 
static void SetPointDataActiveScalarInfo (vtkInformation *info, int arrayType, int numComponents)
 
static vtkDataObjectGetData (vtkInformation *info)
 
static vtkDataObjectGetData (vtkInformationVector *v, int i=0)
 
static void SetGlobalReleaseDataFlag (int val)
 
static int GetGlobalReleaseDataFlag ()
 
- Static Public Member Functions inherited from vtkObject
static int IsTypeOf (const char *type)
 
static vtkObjectSafeDownCast (vtkObjectBase *o)
 
static vtkObjectNew ()
 
static void BreakOnError ()
 
static void SetGlobalWarningDisplay (int val)
 
static void GlobalWarningDisplayOn ()
 
static void GlobalWarningDisplayOff ()
 
static int GetGlobalWarningDisplay ()
 
- Static Public Member Functions inherited from vtkObjectBase
static int IsTypeOf (const char *name)
 
static vtkObjectBaseNew ()
 

Protected Member Functions

 vtkGenericDataSet ()
 
virtual ~vtkGenericDataSet ()
 
- Protected Member Functions inherited from vtkDataObject
 vtkDataObject ()
 
 ~vtkDataObject ()
 
- Protected Member Functions inherited from vtkObject
 vtkObject ()
 
virtual ~vtkObject ()
 
virtual void RegisterInternal (vtkObjectBase *, int check)
 
virtual void UnRegisterInternal (vtkObjectBase *, int check)
 
void InternalGrabFocus (vtkCommand *mouseEvents, vtkCommand *keypressEvents=NULL)
 
void InternalReleaseFocus ()
 
- Protected Member Functions inherited from vtkObjectBase
 vtkObjectBase ()
 
virtual ~vtkObjectBase ()
 
virtual void CollectRevisions (ostream &)
 
virtual void ReportReferences (vtkGarbageCollector *)
 
 vtkObjectBase (const vtkObjectBase &)
 
void operator= (const vtkObjectBase &)
 

Protected Attributes

vtkGenericAttributeCollectionAttributes
 
vtkGenericCellTessellatorTessellator
 
double Bounds [6]
 
double Center [3]
 
vtkTimeStamp ComputeTime
 
- Protected Attributes inherited from vtkDataObject
vtkFieldDataFieldData
 
int DataReleased
 
vtkTimeStamp UpdateTime
 
vtkInformationInformation
 
- Protected Attributes inherited from vtkObject
bool Debug
 
vtkTimeStamp MTime
 
vtkSubjectHelper * SubjectHelper
 
- Protected Attributes inherited from vtkObjectBase
vtkAtomicInt32 ReferenceCount
 
vtkWeakPointerBase ** WeakPointers
 
typedef vtkDataObject Superclass
 
static int IsTypeOf (const char *type)
 
static vtkGenericDataSetSafeDownCast (vtkObjectBase *o)
 
virtual int IsA (const char *type)
 
vtkGenericDataSetNewInstance () const
 
void PrintSelf (ostream &os, vtkIndent indent)
 
virtual vtkObjectBaseNewInstanceInternal () const
 

Additional Inherited Members

- Public Types inherited from vtkDataObject
typedef vtkObject Superclass
 
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
}
 
enum  AttributeTypes {
  POINT, CELL, FIELD, POINT_THEN_CELL,
  VERTEX, EDGE, ROW, NUMBER_OF_ATTRIBUTE_TYPES
}
 
enum  FieldOperations { FIELD_OPERATION_PRESERVED, FIELD_OPERATION_REINTERPOLATED, FIELD_OPERATION_MODIFIED, FIELD_OPERATION_REMOVED }
 
- Public Types inherited from vtkObject
typedef vtkObjectBase Superclass
 

Detailed Description

defines dataset interface

In VTK, spatial-temporal data is defined in terms of a dataset. The dataset consists of geometry (e.g., points), topology (e.g., cells), and attributes (e.g., scalars, vectors, etc.) vtkGenericDataSet is an abstract class defining this abstraction.

Since vtkGenericDataSet provides a general interface to manipulate data, algorithms that process it tend to be slower than those specialized for a particular data type. For this reason, there are concrete, non-abstract subclasses that represent and provide access to data more efficiently. Note that filters to process this dataset type are currently found in the VTK/GenericFiltering/ subdirectory.

Unlike the vtkDataSet class, vtkGenericDataSet provides a more flexible interface including support for iterators. vtkGenericDataSet is also designed to interface VTK to external simulation packages without the penalty of copying memory (see VTK/GenericFiltering/README.html) for more information. Thus vtkGenericDataSet plays a central role in the adaptor framework.

Please note that this class introduces the concepts of "boundary cells". This refers to the boundaries of a cell (e.g., face of a tetrahedron) which may in turn be represented as a cell. Boundary cells are derivative topological features of cells, and are therefore never explicitly represented in the dataset. Often in visualization algorithms, looping over boundaries (edges or faces) is employed, while the actual dataset cells may not traversed. Thus there are methods to loop over these boundary cells.

Finally, as a point of clarification, points are not the same as vertices. Vertices refer to points, and points specify a position is space. Vertices are a type of 0-D cell. Also, the concept of a DOFNode, which is where coefficients for higher-order cells are kept, is a new concept introduced by the adaptor framework (see vtkGenericAdaptorCell for more information).

See also
vtkGenericAdaptorCell vtkDataSet

Definition at line 68 of file vtkGenericDataSet.h.

Member Typedef Documentation

Standard VTK type and print macros.

Definition at line 73 of file vtkGenericDataSet.h.

Constructor & Destructor Documentation

vtkGenericDataSet::vtkGenericDataSet ( )
protected

Constructor with uninitialized bounds (1,-1, 1,-1, 1,-1), empty attribute collection and default tessellator.

virtual vtkGenericDataSet::~vtkGenericDataSet ( )
protectedvirtual

Member Function Documentation

static int vtkGenericDataSet::IsTypeOf ( const char *  type)
static

Standard VTK type and print macros.

virtual int vtkGenericDataSet::IsA ( const char *  type)
virtual

Standard VTK type and print macros.

Reimplemented from vtkDataObject.

Reimplemented in vtkBridgeDataSet.

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

Standard VTK type and print macros.

virtual vtkObjectBase* vtkGenericDataSet::NewInstanceInternal ( ) const
protectedvirtual

Standard VTK type and print macros.

Reimplemented from vtkDataObject.

Reimplemented in vtkBridgeDataSet.

vtkGenericDataSet* vtkGenericDataSet::NewInstance ( ) const

Standard VTK type and print macros.

void vtkGenericDataSet::PrintSelf ( ostream &  os,
vtkIndent  indent 
)
virtual

Standard VTK type and print macros.

Reimplemented from vtkDataObject.

Reimplemented in vtkBridgeDataSet.

virtual vtkIdType vtkGenericDataSet::GetNumberOfPoints ( )
pure virtual

Return the number of points composing the dataset. See NewPointIterator() for more details.

Postcondition
positive_result: result>=0

Implemented in vtkBridgeDataSet.

virtual vtkIdType vtkGenericDataSet::GetNumberOfCells ( int  dim = -1)
pure virtual

Return the number of cells that explicitly define the dataset. See NewCellIterator() for more details.

Precondition
valid_dim_range: (dim>=-1) && (dim<=3)
Postcondition
positive_result: result>=0

Implemented in vtkBridgeDataSet.

virtual int vtkGenericDataSet::GetCellDimension ( )
pure virtual

Return -1 if the dataset is explicitly defined by cells of varying dimensions or if there are no cells. If the dataset is explicitly defined by cells of a unique dimension, return this dimension.

Postcondition
valid_range: (result>=-1) && (result<=3)

Implemented in vtkBridgeDataSet.

virtual void vtkGenericDataSet::GetCellTypes ( vtkCellTypes types)
virtual

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 5 triangles, 3 lines, and 100 hexahedra would result a list of three entries, corresponding to the types VTK_TRIANGLE, VTK_LINE, and VTK_HEXAHEDRON. THIS METHOD IS THREAD SAFE IF FIRST CALLED FROM A SINGLE THREAD AND THE DATASET IS NOT MODIFIED

Precondition
types_exist: types!=0

Reimplemented in vtkBridgeDataSet.

virtual vtkGenericCellIterator* vtkGenericDataSet::NewCellIterator ( int  dim = -1)
pure virtual

Return an iterator to traverse cells of dimension `dim' (or all dimensions if -1) that explicitly define the dataset. For instance, it will return only tetrahedra if the mesh is defined by tetrahedra. If the mesh is composed of two parts, one with tetrahedra and another part with triangles, it will return both, but will not return the boundary edges and vertices of these cells. The user is responsible for deleting the iterator.

Precondition
valid_dim_range: (dim>=-1) && (dim<=3)
Postcondition
result_exists: result!=0

Implemented in vtkBridgeDataSet.

virtual vtkGenericCellIterator* vtkGenericDataSet::NewBoundaryIterator ( int  dim = -1,
int  exteriorOnly = 0 
)
pure virtual

Return an iterator to traverse cell boundaries of dimension `dim' (or all dimensions if -1) of the dataset. If `exteriorOnly' is true, only the exterior cell boundaries of the dataset will be returned, otherwise it will return exterior and interior cell boundaries. The user is responsible for deleting the iterator.

Precondition
valid_dim_range: (dim>=-1) && (dim<=2)
Postcondition
result_exists: result!=0

Implemented in vtkBridgeDataSet.

virtual vtkGenericPointIterator* vtkGenericDataSet::NewPointIterator ( )
pure virtual

Return an iterator to traverse the points composing the dataset; they can be points that define a cell (corner points) or isolated points. The user is responsible for deleting the iterator.

Postcondition
result_exists: result!=0

Implemented in vtkBridgeDataSet.

virtual int vtkGenericDataSet::FindCell ( double  x[3],
vtkGenericCellIterator *&  cell,
double  tol2,
int subId,
double  pcoords[3] 
)
pure virtual

Locate the closest cell to position `x' (global coordinates) with respect to a tolerance squared `tol2' and an initial guess `cell' (if valid). The result consists in the `cell', the `subId' of the sub-cell (0 if primary cell), the parametric coordinates `pcoord' of the position. It returns whether the position is inside the cell or not (boolean). Tolerance is used to control how close the point is to be considered "in" the cell. THIS METHOD IS NOT THREAD SAFE.

Precondition
not_empty: GetNumberOfCells()>0
cell_exists: cell!=0
positive_tolerance: tol2>0

Implemented in vtkBridgeDataSet.

virtual void vtkGenericDataSet::FindPoint ( double  x[3],
vtkGenericPointIterator p 
)
pure virtual

Locate the closest point `p' to position `x' (global coordinates).

Precondition
not_empty: GetNumberOfPoints()>0
p_exists: p!=0

Implemented in vtkBridgeDataSet.

virtual unsigned long int vtkGenericDataSet::GetMTime ( )
virtual

Datasets are composite objects and need to check each part for their modified time.

Reimplemented from vtkDataObject.

Reimplemented in vtkBridgeDataSet.

virtual void vtkGenericDataSet::ComputeBounds ( )
pure virtual

Compute the geometry bounding box.

Implemented in vtkBridgeDataSet.

virtual double* vtkGenericDataSet::GetBounds ( )
virtual

Return a pointer to the geometry bounding box in the form (xmin,xmax, ymin,ymax, zmin,zmax). The return value is VOLATILE.

Postcondition
result_exists: result!=0
virtual void vtkGenericDataSet::GetBounds ( double  bounds[6])
virtual

Return the geometry bounding box in global coordinates in the form (xmin,xmax, ymin,ymax, zmin,zmax) in the `bounds' array.

virtual double* vtkGenericDataSet::GetCenter ( )
virtual

Get the center of the bounding box in global coordinates. The return value is VOLATILE.

Postcondition
result_exists: result!=0
virtual void vtkGenericDataSet::GetCenter ( double  center[3])
virtual

Get the center of the bounding box in global coordinates.

virtual double vtkGenericDataSet::GetLength ( )
virtual

Return the length of the diagonal of the bounding box.

Postcondition
positive_result: result>=0
virtual vtkGenericAttributeCollection* vtkGenericDataSet::GetAttributes ( )
virtual

Get the collection of attributes associated with this dataset.

virtual vtkDataSetAttributes* vtkGenericDataSet::GetAttributes ( int  type)
inlinevirtual

Returns the attributes of the data object of the specified attribute type. The type may be:

The other attribute type, FIELD, will return NULL since field data is stored as a vtkFieldData instance, not a vtkDataSetAttributes instance. To retrieve field data, use GetAttributesAsFieldData.

Reimplemented from vtkDataObject.

Definition at line 193 of file vtkGenericDataSet.h.

virtual void vtkGenericDataSet::SetTessellator ( vtkGenericCellTessellator tessellator)
virtual

Set/Get a cell tessellator if cells must be tessellated during processing.

Precondition
tessellator_exists: tessellator!=0
virtual vtkGenericCellTessellator* vtkGenericDataSet::GetTessellator ( )
virtual

Set/Get a cell tessellator if cells must be tessellated during processing.

Precondition
tessellator_exists: tessellator!=0
virtual unsigned long vtkGenericDataSet::GetActualMemorySize ( )
virtual

Actual size of the data in kibibytes (1024 bytes); only valid after the pipeline has updated. It is guaranteed to be greater than or equal to the memory required to represent the data.

Reimplemented from vtkDataObject.

int vtkGenericDataSet::GetDataObjectType ( )
virtual

Return the type of data object.

Reimplemented from vtkDataObject.

virtual vtkIdType vtkGenericDataSet::GetEstimatedSize ( )
pure virtual

Estimated size needed after tessellation (or special operation)

Implemented in vtkBridgeDataSet.

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

Retrieve an instance of this class from an information object.

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

Retrieve an instance of this class from an information object.

Member Data Documentation

vtkGenericAttributeCollection* vtkGenericDataSet::Attributes
protected

Definition at line 230 of file vtkGenericDataSet.h.

vtkGenericCellTessellator* vtkGenericDataSet::Tessellator
protected

Definition at line 233 of file vtkGenericDataSet.h.

double vtkGenericDataSet::Bounds[6]
protected

Definition at line 235 of file vtkGenericDataSet.h.

double vtkGenericDataSet::Center[3]
protected

Definition at line 236 of file vtkGenericDataSet.h.

vtkTimeStamp vtkGenericDataSet::ComputeTime
protected

Definition at line 237 of file vtkGenericDataSet.h.


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