Main Page   Class Hierarchy   Alphabetical List   Compound List   File List   Compound Members   File Members   Related Pages  

vtkPolyData Class Reference

concrete dataset represents vertices, lines, polygons, and triangle strips. More...

#include <vtkPolyData.h>

Inheritance diagram for vtkPolyData:

Inheritance graph
[legend]
Collaboration diagram for vtkPolyData:

Collaboration graph
[legend]
List of all members.

Public Methods

virtual const char * GetClassName ()
virtual int IsA (const char *type)
void PrintSelf (ostream &os, vtkIndent indent)
vtkDataObjectMakeObject ()
int GetDataObjectType ()
void CopyStructure (vtkDataSet *ds)
int GetNumberOfCells ()
vtkCellGetCell (int cellId)
void GetCell (int cellId, vtkGenericCell *cell)
int GetCellType (int cellId)
void GetCellBounds (int cellId, float bounds[6])
void GetCellNeighbors (int cellId, vtkIdList *ptIds, vtkIdList *cellIds)
void CopyCells (vtkPolyData *pd, vtkIdList *idList, vtkPointLocator *locator=NULL)
void GetCellPoints (int cellId, vtkIdList *ptIds)
void GetPointCells (int ptId, vtkIdList *cellIds)
void ComputeBounds ()
void Squeeze ()
int GetMaxCellSize ()
void SetVerts (vtkCellArray *v)
vtkCellArrayGetVerts ()
void SetLines (vtkCellArray *l)
vtkCellArrayGetLines ()
void SetPolys (vtkCellArray *p)
vtkCellArrayGetPolys ()
void SetStrips (vtkCellArray *s)
vtkCellArrayGetStrips ()
int GetNumberOfVerts ()
int GetNumberOfLines ()
int GetNumberOfPolys ()
int GetNumberOfStrips ()
void Allocate (int numCells=1000, int extSize=1000)
int InsertNextCell (int type, int npts, int *pts)
int InsertNextCell (int type, vtkIdList *pts)
void Reset ()
void BuildCells ()
void BuildLinks ()
void DeleteCells ()
void DeleteLinks ()
void GetPointCells (int ptId, unsigned short &ncells, int *&cells)
void GetCellEdgeNeighbors (int cellId, int p1, int p2, vtkIdList *cellIds)
void GetCellPoints (int cellId, int &npts, int *&pts)
int IsTriangle (int v1, int v2, int v3)
int IsEdge (int v1, int v2)
int IsPointUsedByCell (int ptId, int cellId)
void ReplaceCell (int cellId, int npts, int *pts)
void ReplaceCellPoint (int cellId, int oldPtId, int newPtId)
void ReverseCell (int cellId)
void DeletePoint (int ptId)
void DeleteCell (int cellId)
int InsertNextLinkedPoint (float x[3], int numLinks)
int InsertNextLinkedCell (int type, int npts, int *pts)
void ReplaceLinkedCell (int cellId, int npts, int *pts)
void RemoveCellReference (int cellId)
void AddCellReference (int cellId)
void RemoveReferenceToCell (int ptId, int cellId)
void AddReferenceToCell (int ptId, int cellId)
void ResizeCellList (int ptId, int size)
virtual void Initialize ()
void SetUpdateExtent (int piece, int numPieces, int ghostLevel)
void SetUpdateExtent (int piece, int numPieces)
void GetUpdateExtent (int &piece, int &numPieces, int &ghostLevel)
virtual int * GetUpdateExtent ()
virtual void GetUpdateExtent (int &, int &, int &, int &, int &, int &)
virtual void GetUpdateExtent (int[6])
void SetUpdateExtent (int x1, int x2, int y1, int y2, int z1, int z2)
void SetUpdateExtent (int ext[6])
virtual int GetPiece ()
virtual int GetNumberOfPieces ()
virtual int GetGhostLevel ()
unsigned long GetActualMemorySize ()
void ShallowCopy (vtkDataObject *src)
void DeepCopy (vtkDataObject *src)
void GetCellPoints (int cellId, vtkIdList &ptIds)
void GetPointCells (int ptId, vtkIdList &cellIds)
int InsertNextCell (int type, vtkIdList &pts)
void GetCellEdgeNeighbors (int cellId, int p1, int p2, vtkIdList &cellIds)

Static Public Methods

vtkPolyData * New ()
int IsTypeOf (const char *type)
vtkPolyData * SafeDownCast (vtkObject *o)

Protected Methods

 vtkPolyData ()
 ~vtkPolyData ()
 vtkPolyData (const vtkPolyData &)
void operator= (const vtkPolyData &)

Protected Attributes

vtkVertexVertex
vtkPolyVertexPolyVertex
vtkLineLine
vtkPolyLinePolyLine
vtkTriangleTriangle
vtkQuadQuad
vtkPolygonPolygon
vtkTriangleStripTriangleStrip
vtkEmptyCellEmptyCell
vtkCellArrayVerts
vtkCellArrayLines
vtkCellArrayPolys
vtkCellArrayStrips
vtkCellTypesCells
vtkCellLinksLinks

Static Protected Attributes

vtkCellArrayDummy

Detailed Description

concrete dataset represents vertices, lines, polygons, and triangle strips.

Date:
2000/12/10 20:08:15
Revision:
1.106

vtkPolyData is a data object that is a concrete implementation of vtkDataSet. vtkPolyData represents a geometric structure consisting of vertices, lines, polygons, and triangle strips. Point attribute values (e.g., scalars, vectors, etc.) also are represented.

The actual cell types (CellType.h) supported by vtkPolyData are: vtkVertex, vtkPolyVertex, vtkLine, vtkPolyLine, vtkTriangle, vtkTriangleStrip, vtkPolygon, vtkPixel, and vtkQuad.

One important feature of vtkPolyData objects is that special traversal and data manipulation methods are available to process data. These methods are generally more efficient than vtkDataSet methods and should be used whenever possible. For example, traversing the cells in a dataset we would use GetCell(). To traverse cells with vtkPolyData we would retrieve the cell array object representing polygons (for example) and then use vtkCellArray's InitTraversal() and GetNextCell() methods.

Examples:
vtkPolyData (examples)

Definition at line 88 of file vtkPolyData.h.


Constructor & Destructor Documentation

vtkPolyData::vtkPolyData   [protected]
 

vtkPolyData::~vtkPolyData   [protected]
 

vtkPolyData::vtkPolyData const vtkPolyData &    [inline, protected]
 

Definition at line 364 of file vtkPolyData.h.


Member Function Documentation

vtkPolyData* vtkPolyData::New   [static]
 

Create an object with Debug turned off, modified time initialized to zero, and reference counting on.

Reimplemented from vtkDataObject.

Referenced by MakeObject().

virtual const char* vtkPolyData::GetClassName   [virtual]
 

Return the class name as a string. This method is defined in all subclasses of vtkObject with the vtkTypeMacro found in vtkSetGet.h.

Reimplemented from vtkPointSet.

int vtkPolyData::IsTypeOf const char *    type [static]
 

Return 1 if this class type is the same type of (or a subclass of) the named class. Returns 0 otherwise. This method works in combination with vtkTypeMacro found in vtkSetGet.h.

Reimplemented from vtkPointSet.

virtual int vtkPolyData::IsA const char *    type [virtual]
 

Return 1 if this class is the same type of (or a subclass of) the named class. Returns 0 otherwise. This method works in combination with vtkTypeMacro found in vtkSetGet.h.

Reimplemented from vtkPointSet.

vtkPolyData* vtkPolyData::SafeDownCast vtkObject   o [static]
 

Will cast the supplied object to vtkObject* is this is a safe operation (i.e., a safe downcast); otherwise NULL is returned. This method is defined in all subclasses of vtkObject with the vtkTypeMacro found in vtkSetGet.h.

Reimplemented from vtkPointSet.

Referenced by vtkPolyDataCollection::GetNextItem().

void vtkPolyData::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.

vtkDataObject* vtkPolyData::MakeObject   [inline, virtual]
 

Create a similar type object.

Reimplemented from vtkDataObject.

Definition at line 97 of file vtkPolyData.h.

int vtkPolyData::GetDataObjectType   [inline, virtual]
 

Return what type of dataset this is.

Reimplemented from vtkDataSet.

Definition at line 100 of file vtkPolyData.h.

void vtkPolyData::CopyStructure vtkDataSet   ds [virtual]
 

Copy the geometric and topological structure of an input poly data object.

Reimplemented from vtkPointSet.

int vtkPolyData::GetNumberOfCells   [virtual]
 

Standard vtkDataSet interface.

Reimplemented from vtkDataSet.

vtkCell* vtkPolyData::GetCell int    cellId [virtual]
 

Get cell with cellId such that: 0 <= cellId < NumberOfCells. THIS METHOD IS NOT THREAD SAFE.

Reimplemented from vtkDataSet.

void vtkPolyData::GetCell int    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

Reimplemented from vtkDataSet.

int vtkPolyData::GetCellType int    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

Reimplemented from vtkDataSet.

void vtkPolyData::GetCellBounds int    cellId,
float    bounds[6]
[virtual]
 

Get the bounds of the cell with cellId such that: 0 <= cellId < NumberOfCells. A subclass may be able to determine the bounds of cell without using an expensive GetCell() method. A default implementation is provided that actually uses a GetCell() call. This is to ensure the method is available to all datasets. Subclasses should override this method to provide an efficient implementation. THIS METHOD IS THREAD SAFE IF FIRST CALLED FROM A SINGLE THREAD AND THE DATASET IS NOT MODIFIED

Reimplemented from vtkDataSet.

void vtkPolyData::GetCellNeighbors int    cellId,
vtkIdList   ptIds,
vtkIdList   cellIds
[virtual]
 

Topological inquiry to get all cells using list of points exclusive of 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.

void vtkPolyData::CopyCells vtkPolyData *    pd,
vtkIdList   idList,
vtkPointLocator   locator = NULL
 

Copy cells listed in idList from pd, including points, point data, and cell data. This method assumes that point and cell data have been allocated. If you pass in a point locator, then the points won't be duplicated in the output.

void vtkPolyData::GetCellPoints int    cellId,
vtkIdList   ptIds
[virtual]
 

Copy a cells point ids into list provided. (Less efficient.)

Reimplemented from vtkDataSet.

Referenced by AddCellReference(), IsPointUsedByCell(), IsTriangle(), RemoveCellReference(), and ReplaceCellPoint().

void vtkPolyData::GetPointCells int    ptId,
vtkIdList   cellIds
[virtual]
 

Efficient method to obtain cells using a particular point. Make sure that routine BuildLinks() has been called.

Reimplemented from vtkDataSet.

Referenced by IsEdge(), and IsTriangle().

void vtkPolyData::ComputeBounds   [virtual]
 

Compute the (X, Y, Z) bounds of the data.

Reimplemented from vtkPointSet.

void vtkPolyData::Squeeze   [virtual]
 

Recover extra allocated memory when creating data whose initial size is unknown. Examples include using the InsertNextCell() method, or when using the CellArray::EstimateSize() method to create vertices, lines, polygons, or triangle strips.

Reimplemented from vtkPointSet.

int vtkPolyData::GetMaxCellSize   [virtual]
 

Return the maximum cell size in this poly data.

Reimplemented from vtkDataSet.

void vtkPolyData::SetVerts vtkCellArray   v
 

Set the cell array defining vertices.

vtkCellArray* vtkPolyData::GetVerts  
 

Get the cell array defining vertices. If there are no vertices, an empty array will be returned (convenience to simplify traversal).

void vtkPolyData::SetLines vtkCellArray   l
 

Set the cell array defining lines.

vtkCellArray* vtkPolyData::GetLines  
 

Get the cell array defining lines. If there are no lines, an empty array will be returned (convenience to simplify traversal).

void vtkPolyData::SetPolys vtkCellArray   p
 

Set the cell array defining polygons.

vtkCellArray* vtkPolyData::GetPolys  
 

Get the cell array defining polygons. If there are no polygons, an empty array will be returned (convenience to simplify traversal).

void vtkPolyData::SetStrips vtkCellArray   s
 

Set the cell array defining triangle strips.

vtkCellArray* vtkPolyData::GetStrips  
 

Get the cell array defining triangle strips. If there are no triangle strips, an empty array will be returned (convenience to simplify traversal).

int vtkPolyData::GetNumberOfVerts  
 

Return the number of primitives of a particular type held..

int vtkPolyData::GetNumberOfLines  
 

int vtkPolyData::GetNumberOfPolys  
 

int vtkPolyData::GetNumberOfStrips  
 

void vtkPolyData::Allocate int    numCells = 1000,
int    extSize = 1000
 

Method allocates initial storage for vertex, line, polygon, and triangle strip arrays. Use this method before the method PolyData::InsertNextCell(). (Or, provide vertex, line, polygon, and triangle strip cell arrays.)

int vtkPolyData::InsertNextCell int    type,
int    npts,
int *    pts
 

Insert a cell of type vtkVERTEX, vtkPOLY_VERTEX, vtkLINE, vtkPOLY_LINE, vtkTRIANGLE, vtkQUAD, vtkPOLYGON, or vtkTRIANGLE_STRIP. Make sure that the PolyData::Allocate() function has been called first or that vertex, line, polygon, and triangle strip arrays have been supplied. Note: will also insert vtkPIXEL, but converts it to vtkQUAD.

int vtkPolyData::InsertNextCell int    type,
vtkIdList   pts
 

Insert a cell of type VTK_VERTEX, VTK_POLY_VERTEX, VTK_LINE, VTK_POLY_LINE, VTK_TRIANGLE, VTK_QUAD, VTK_POLYGON, or VTK_TRIANGLE_STRIP. Make sure that the PolyData::Allocate() function has been called first or that vertex, line, polygon, and triangle strip arrays have been supplied. Note: will also insert VTK_PIXEL, but converts it to VTK_QUAD.

void vtkPolyData::Reset  
 

Begin inserting data all over again. Memory is not freed but otherwise objects are returned to their initial state.

void vtkPolyData::BuildCells  
 

Create data structure that allows random access of cells.

void vtkPolyData::BuildLinks  
 

Create upward links from points to cells that use each point. Enables topologically complex queries.

void vtkPolyData::DeleteCells  
 

Release data structure that allows random access of the cells. This must be done before a 2nd call to BuildLinks(). DeleteCells implicitly deletes the links as well since they are no longer valid.

void vtkPolyData::DeleteLinks  
 

Release the upward links from point to cells that use each point.

void vtkPolyData::GetPointCells int    ptId,
unsigned short &    ncells,
int *&    cells
[inline]
 

Special (efficient) operations on poly data. Use carefully.

Definition at line 402 of file vtkPolyData.h.

void vtkPolyData::GetCellEdgeNeighbors int    cellId,
int    p1,
int    p2,
vtkIdList   cellIds
 

Get the neighbors at an edge. More efficient than the general GetCellNeighbors(). Assumes links have been built (with BuildLinks()), and looks specifically for edge neighbors.

void vtkPolyData::GetCellPoints int    cellId,
int &    npts,
int *&    pts
 

Return a pointer to a list of point ids defining cell. (More efficient.) Assumes that cells have been built (with BuildCells()).

int vtkPolyData::IsTriangle int    v1,
int    v2,
int    v3
[inline]
 

Given three vertices, determine whether it's a triangle. Make sure BuildLinks() has been called first.

Definition at line 409 of file vtkPolyData.h.

int vtkPolyData::IsEdge int    p1,
int    p2
[inline]
 

Determine whether two points form an edge. If they do, return non-zero. Make sure BuildLinks() has been called first.

Definition at line 453 of file vtkPolyData.h.

int vtkPolyData::IsPointUsedByCell int    ptId,
int    cellId
[inline]
 

Determine whether a point is used by a particular cell. If it is, return non-zero. Make sure BuildCells() has been called first.

Definition at line 438 of file vtkPolyData.h.

Referenced by IsEdge().

void vtkPolyData::ReplaceCell int    cellId,
int    npts,
int *    pts
 

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.

void vtkPolyData::ReplaceCellPoint int    cellId,
int    oldPtId,
int    newPtId
[inline]
 

Replace a point in the cell connectivity list with a different point.

Definition at line 505 of file vtkPolyData.h.

void vtkPolyData::ReverseCell int    cellId
 

Reverse the order of point ids defining the cell.

void vtkPolyData::DeletePoint int    ptId [inline]
 

Mark a point/cell as deleted from this vtkPolyData.

Definition at line 470 of file vtkPolyData.h.

void vtkPolyData::DeleteCell int    cellId [inline]
 

Definition at line 475 of file vtkPolyData.h.

int vtkPolyData::InsertNextLinkedPoint float    x[3],
int    numLinks
 

Add a point to the cell data structure (after cell pointers have been built). This method adds the point and then allocates memory for the links to the cells. (To use this method, make sure points are available and BuildLinks() has been invoked.)

int vtkPolyData::InsertNextLinkedCell int    type,
int    npts,
int *    pts
 

Add a new cell to the cell data structure (after cell pointers have been built). This method adds the cell and then updates the links from the points to the cells. (Memory is allocated as necessary.)

void vtkPolyData::ReplaceLinkedCell int    cellId,
int    npts,
int *    pts
 

Replace one cell with another in cell structure. This operator updates the connectivity list and the point's link list. It does not delete references to the old cell in the point's link list. Use the operator RemoveCellReference() to delete all references from points to (old) cell. You may also want to consider using the operator ResizeCellList() if the link list is changing size.

void vtkPolyData::RemoveCellReference int    cellId [inline]
 

Remove all references to cell in cell structure. This means the links from the cell's points to the cell are deleted. Memory is not reclaimed. Use the method ResizeCellList() to resize the link list from a point to its using cells. (This operator assumes BuildLinks() has been called.)

Definition at line 480 of file vtkPolyData.h.

void vtkPolyData::AddCellReference int    cellId [inline]
 

Add references to cell in cell structure. This means the links from the cell's points to the cell are modified. Memory is not extended. Use the method ResizeCellList() to resize the link list from a point to its using cells. (This operator assumes BuildLinks() has been called.)

Definition at line 490 of file vtkPolyData.h.

void vtkPolyData::RemoveReferenceToCell int    ptId,
int    cellId
 

Remove a reference to a cell in a particular point's link list. You may also consider using RemoveCellReference() to remove the references from all the cell's points to the cell. This operator does not reallocate memory; use the operator ResizeCellList() to do this if necessary.

void vtkPolyData::AddReferenceToCell int    ptId,
int    cellId
 

Add a reference to a cell in a particular point's link list. (You may also consider using AddCellReference() to add the references from all the cell's points to the cell.) This operator does not realloc memory; use the operator ResizeCellList() to do this if necessary.

void vtkPolyData::ResizeCellList int    ptId,
int    size
[inline]
 

Resize the list of cells using a particular point. (This operator assumes that BuildLinks() has been called.)

Definition at line 500 of file vtkPolyData.h.

virtual void vtkPolyData::Initialize   [virtual]
 

Restore object to initial state. Release memory back to system.

Reimplemented from vtkPointSet.

void vtkPolyData::SetUpdateExtent int    piece,
int    numPieces,
int    ghostLevel
 

For streaming. User/next filter specifies which piece they want updated. The source of this poly data has to return exactly this piece.

void vtkPolyData::SetUpdateExtent int    piece,
int    numPieces
[inline]
 

Reimplemented from vtkDataObject.

Definition at line 314 of file vtkPolyData.h.

void vtkPolyData::GetUpdateExtent int &    piece,
int &    numPieces,
int &    ghostLevel
 

virtual int* vtkPolyData::GetUpdateExtent   [virtual]
 

We need this here to avoid hiding superclass method

Reimplemented from vtkDataObject.

virtual void vtkPolyData::GetUpdateExtent int &   ,
int &   ,
int &   ,
int &   ,
int &   ,
int &   
[virtual]
 

Reimplemented from vtkDataObject.

virtual void vtkPolyData::GetUpdateExtent int   [6] [virtual]
 

Reimplemented from vtkDataObject.

void vtkPolyData::SetUpdateExtent int    x1,
int    x2,
int    y1,
int    y2,
int    z1,
int    z2
[inline, virtual]
 

Call superclass method to avoid hiding Since this data type does not use 3D extents, this set method is useless but necessary since vtkDataSetToDataSetFilter does not know what type of data it is working on.

Reimplemented from vtkDataObject.

Definition at line 325 of file vtkPolyData.h.

void vtkPolyData::SetUpdateExtent int    ext[6] [inline, virtual]
 

Reimplemented from vtkDataObject.

Definition at line 327 of file vtkPolyData.h.

virtual int vtkPolyData::GetPiece   [virtual]
 

Set / Get the piece and the number of pieces. Similar to extent in 3D.

virtual int vtkPolyData::GetNumberOfPieces   [virtual]
 

virtual int vtkPolyData::GetGhostLevel   [virtual]
 

Get the ghost level.

unsigned long vtkPolyData::GetActualMemorySize   [virtual]
 

Return the actual size of the data in kilobytes. 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.

void vtkPolyData::ShallowCopy vtkDataObject   src [virtual]
 

Shallow and Deep copy.

Reimplemented from vtkPointSet.

void vtkPolyData::DeepCopy vtkDataObject   src [virtual]
 

Reimplemented from vtkPointSet.

void vtkPolyData::GetCellPoints int    cellId,
vtkIdList   ptIds
[inline]
 

For legacy compatibility. Do not use.

Reimplemented from vtkDataSet.

Definition at line 351 of file vtkPolyData.h.

void vtkPolyData::GetPointCells int    ptId,
vtkIdList   cellIds
[inline]
 

Reimplemented from vtkDataSet.

Definition at line 353 of file vtkPolyData.h.

int vtkPolyData::InsertNextCell int    type,
vtkIdList   pts
[inline]
 

Definition at line 355 of file vtkPolyData.h.

void vtkPolyData::GetCellEdgeNeighbors int    cellId,
int    p1,
int    p2,
vtkIdList   cellIds
[inline]
 

Definition at line 357 of file vtkPolyData.h.

void vtkPolyData::operator= const vtkPolyData &    [inline, protected]
 

Definition at line 365 of file vtkPolyData.h.


Member Data Documentation

vtkVertex* vtkPolyData::Vertex [protected]
 

Definition at line 368 of file vtkPolyData.h.

vtkPolyVertex* vtkPolyData::PolyVertex [protected]
 

Definition at line 369 of file vtkPolyData.h.

vtkLine* vtkPolyData::Line [protected]
 

Definition at line 370 of file vtkPolyData.h.

vtkPolyLine* vtkPolyData::PolyLine [protected]
 

Definition at line 371 of file vtkPolyData.h.

vtkTriangle* vtkPolyData::Triangle [protected]
 

Definition at line 372 of file vtkPolyData.h.

vtkQuad* vtkPolyData::Quad [protected]
 

Definition at line 373 of file vtkPolyData.h.

vtkPolygon* vtkPolyData::Polygon [protected]
 

Definition at line 374 of file vtkPolyData.h.

vtkTriangleStrip* vtkPolyData::TriangleStrip [protected]
 

Definition at line 375 of file vtkPolyData.h.

vtkEmptyCell* vtkPolyData::EmptyCell [protected]
 

Definition at line 376 of file vtkPolyData.h.

vtkCellArray* vtkPolyData::Verts [protected]
 

Definition at line 380 of file vtkPolyData.h.

vtkCellArray* vtkPolyData::Lines [protected]
 

Definition at line 381 of file vtkPolyData.h.

vtkCellArray* vtkPolyData::Polys [protected]
 

Definition at line 382 of file vtkPolyData.h.

vtkCellArray* vtkPolyData::Strips [protected]
 

Definition at line 383 of file vtkPolyData.h.

vtkCellArray* vtkPolyData::Dummy [static, protected]
 

Definition at line 386 of file vtkPolyData.h.

vtkCellTypes* vtkPolyData::Cells [protected]
 

Definition at line 390 of file vtkPolyData.h.

vtkCellLinks* vtkPolyData::Links [protected]
 

Definition at line 391 of file vtkPolyData.h.


The documentation for this class was generated from the following file:
Generated on Wed Nov 21 12:57:37 2001 for VTK by doxygen1.2.11.1 written by Dimitri van Heesch, © 1997-2001