Main Page | Directories | Class Hierarchy | Alphabetical List | Class List | File List | Class Members | File Members | Related Pages

vtkPolyData Class Reference

#include <vtkPolyData.h>

Inheritance diagram for vtkPolyData:

Inheritance graph
[legend]
Collaboration diagram for vtkPolyData:

Collaboration graph
[legend]
List of all members.

Detailed Description

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

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.

Created by:
  • Bertel, Francois
CVS contributions (if > 5%):
  • Bertel, Francois (95%)
CVS logs (CVSweb):
  • .cxx (/Filtering/vtkPolyData.cxx)
  • .h (/Filtering/vtkPolyData.h)
Examples:
vtkPolyData (Examples)
Tests:
vtkPolyData (Tests)

Definition at line 70 of file vtkPolyData.h.

Public Types

typedef vtkPointSet Superclass

Public Member Functions

virtual const char * GetClassName ()
virtual int IsA (const char *type)
void PrintSelf (ostream &os, vtkIndent indent)
int GetDataObjectType ()
void CopyStructure (vtkDataSet *ds)
void GetCellPoints (vtkIdType cellId, vtkIdList *ptIds)
void GetPointCells (vtkIdType 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 ()
void Allocate (vtkIdType numCells=1000, int extSize=1000)
int InsertNextCell (int type, int npts, vtkIdType *pts)
int InsertNextCell (int type, vtkIdList *pts)
void Reset ()
void BuildCells ()
void BuildLinks (int initialSize=0)
void DeleteCells ()
void DeleteLinks ()
void GetCellPoints (vtkIdType cellId, vtkIdType &npts, vtkIdType *&pts)
int IsTriangle (int v1, int v2, int v3)
int IsEdge (vtkIdType p1, vtkIdType p2)
int IsPointUsedByCell (vtkIdType ptId, vtkIdType cellId)
void ReplaceCell (vtkIdType cellId, int npts, vtkIdType *pts)
void ReverseCell (vtkIdType cellId)
int InsertNextLinkedCell (int type, int npts, vtkIdType *pts)
void ReplaceLinkedCell (vtkIdType cellId, int npts, vtkIdType *pts)
void RemoveCellReference (vtkIdType cellId)
void AddCellReference (vtkIdType cellId)
void RemoveReferenceToCell (vtkIdType ptId, vtkIdType cellId)
void AddReferenceToCell (vtkIdType ptId, vtkIdType cellId)
void ResizeCellList (vtkIdType ptId, int size)
virtual void Initialize ()
virtual int GetGhostLevel ()
unsigned long GetActualMemorySize ()
void RemoveGhostCells (int level)
vtkIdType GetNumberOfCells ()
vtkCellGetCell (vtkIdType cellId)
void GetCell (vtkIdType cellId, vtkGenericCell *cell)
int GetCellType (vtkIdType cellId)
void GetCellBounds (vtkIdType cellId, double bounds[6])
void GetCellNeighbors (vtkIdType cellId, vtkIdList *ptIds, vtkIdList *cellIds)
void CopyCells (vtkPolyData *pd, vtkIdList *idList, vtkPointLocator *locator=NULL)
vtkIdType GetNumberOfVerts ()
vtkIdType GetNumberOfLines ()
vtkIdType GetNumberOfPolys ()
vtkIdType GetNumberOfStrips ()
void Allocate (vtkPolyData *inPolyData, vtkIdType numCells=1000, int extSize=1000)
void GetPointCells (vtkIdType ptId, unsigned short &ncells, vtkIdType *&cells)
void GetCellEdgeNeighbors (vtkIdType cellId, vtkIdType p1, vtkIdType p2, vtkIdList *cellIds)
void ReplaceCellPoint (vtkIdType cellId, vtkIdType oldPtId, vtkIdType newPtId)
void DeletePoint (vtkIdType ptId)
void DeleteCell (vtkIdType cellId)
int InsertNextLinkedPoint (int numLinks)
int InsertNextLinkedPoint (double x[3], int numLinks)
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 &x0, int &x1, int &y0, int &y1, int &z0, int &z1)
virtual void GetUpdateExtent (int extent[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 ()
void ShallowCopy (vtkDataObject *src)
void DeepCopy (vtkDataObject *src)

Static Public Member Functions

vtkPolyDataNew ()
int IsTypeOf (const char *type)
vtkPolyDataSafeDownCast (vtkObject *o)

Protected Member Functions

 vtkPolyData ()
 ~vtkPolyData ()
virtual void Crop ()

Protected Attributes

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

Static Protected Attributes

vtkCellArrayDummy


Member Typedef Documentation

typedef vtkPointSet vtkPolyData::Superclass
 

Reimplemented from vtkPointSet.

Definition at line 75 of file vtkPolyData.h.


Constructor & Destructor Documentation

vtkPolyData::vtkPolyData  )  [protected]
 

vtkPolyData::~vtkPolyData  )  [protected]
 


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.

virtual const char* vtkPolyData::GetClassName  )  [virtual]
 

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 vtkTypeRevisionMacro 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 vtkTypeRevisionMacro found in vtkSetGet.h.

Reimplemented from vtkPointSet.

vtkPolyData* vtkPolyData::SafeDownCast vtkObject o  )  [static]
 

Reimplemented from vtkPointSet.

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.

int vtkPolyData::GetDataObjectType  )  [inline, virtual]
 

Return what type of dataset this is.

Reimplemented from vtkDataSet.

Definition at line 79 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.

vtkIdType vtkPolyData::GetNumberOfCells  )  [virtual]
 

Standard vtkDataSet interface.

Implements vtkDataSet.

vtkCell* vtkPolyData::GetCell vtkIdType  cellId  )  [virtual]
 

Standard vtkDataSet interface.

Implements vtkDataSet.

void vtkPolyData::GetCell vtkIdType  cellId,
vtkGenericCell cell
[virtual]
 

Standard vtkDataSet interface.

Implements vtkDataSet.

int vtkPolyData::GetCellType vtkIdType  cellId  )  [virtual]
 

Standard vtkDataSet interface.

Implements vtkDataSet.

void vtkPolyData::GetCellBounds vtkIdType  cellId,
double  bounds[6]
[virtual]
 

Standard vtkDataSet interface.

Reimplemented from vtkDataSet.

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

Standard vtkDataSet interface.

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 vtkIdType  cellId,
vtkIdList ptIds
[virtual]
 

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

Implements vtkDataSet.

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

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

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

Implements vtkDataSet.

Referenced by 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.

Implements 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).

vtkIdType vtkPolyData::GetNumberOfVerts  ) 
 

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

vtkIdType vtkPolyData::GetNumberOfLines  ) 
 

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

vtkIdType vtkPolyData::GetNumberOfPolys  ) 
 

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

vtkIdType vtkPolyData::GetNumberOfStrips  ) 
 

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

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

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

Similar to the method above, this method allocates initial storage for vertex, line, polygon, and triangle strip arrays. It does this more intelligently, examining the supplied inPolyData to determine whether to allocate the verts, lines, polys, and strips arrays. (These arrays are allocated only if there is data in the corresponding arrays in the inPolyData.) Caution: if the inPolyData has no verts, and after allocating with this method an PolyData::InsertNextCell() is invoked where a vertex is inserted, bad things will happen.

int vtkPolyData::InsertNextCell int  type,
int  npts,
vtkIdType 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.

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 int  initialSize = 0  ) 
 

Create upward links from points to cells that use each point. Enables topologically complex queries. Normally the links array is allocated based on the number of points in the vtkPolyData. The optional initialSize parameter can be used to allocate a larger size initially.

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 vtkIdType  ptId,
unsigned short &  ncells,
vtkIdType *&  cells
[inline]
 

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

Definition at line 438 of file vtkPolyData.h.

References vtkCellLinks::GetCells(), vtkCellLinks::GetNcells(), Links, and vtkIdType.

void vtkPolyData::GetCellEdgeNeighbors vtkIdType  cellId,
vtkIdType  p1,
vtkIdType  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 vtkIdType  cellId,
vtkIdType npts,
vtkIdType *&  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 445 of file vtkPolyData.h.

References GetCellPoints(), GetPointCells(), and vtkIdType.

int vtkPolyData::IsEdge vtkIdType  p1,
vtkIdType  p2
 

Determine whether two points form an edge. If they do, return non-zero. By definition PolyVertex and PolyLine have no edges since 1-dimensional edges are only found on cells 2D and higher. Edges are defined as 1-D boundary entities to cells. Make sure BuildLinks() has been called first.

int vtkPolyData::IsPointUsedByCell vtkIdType  ptId,
vtkIdType  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 475 of file vtkPolyData.h.

References GetCellPoints(), and vtkIdType.

void vtkPolyData::ReplaceCell vtkIdType  cellId,
int  npts,
vtkIdType 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 vtkIdType  cellId,
vtkIdType  oldPtId,
vtkIdType  newPtId
[inline]
 

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

Definition at line 528 of file vtkPolyData.h.

References GetCellPoints(), and vtkIdType.

void vtkPolyData::ReverseCell vtkIdType  cellId  ) 
 

Reverse the order of point ids defining the cell.

void vtkPolyData::DeletePoint vtkIdType  ptId  )  [inline]
 

Mark a point/cell as deleted from this vtkPolyData.

Definition at line 491 of file vtkPolyData.h.

References vtkCellLinks::DeletePoint(), Links, and vtkIdType.

void vtkPolyData::DeleteCell vtkIdType  cellId  )  [inline]
 

Mark a point/cell as deleted from this vtkPolyData.

Definition at line 496 of file vtkPolyData.h.

References Cells, vtkCellTypes::DeleteCell(), and vtkIdType.

int vtkPolyData::InsertNextLinkedPoint 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.) Of the two methods below, one inserts a point coordinate and the other just makes room for cell links.

int vtkPolyData::InsertNextLinkedPoint double  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.) Of the two methods below, one inserts a point coordinate and the other just makes room for cell links.

int vtkPolyData::InsertNextLinkedCell int  type,
int  npts,
vtkIdType 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 vtkIdType  cellId,
int  npts,
vtkIdType 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 vtkIdType  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 501 of file vtkPolyData.h.

References GetCellPoints(), Links, vtkCellLinks::RemoveCellReference(), and vtkIdType.

void vtkPolyData::AddCellReference vtkIdType  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 512 of file vtkPolyData.h.

References vtkCellLinks::AddCellReference(), GetCellPoints(), Links, and vtkIdType.

void vtkPolyData::RemoveReferenceToCell vtkIdType  ptId,
vtkIdType  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 vtkIdType  ptId,
vtkIdType  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 vtkIdType  ptId,
int  size
[inline]
 

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

Definition at line 523 of file vtkPolyData.h.

References Links, vtkCellLinks::ResizeCellList(), and vtkIdType.

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]
 

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

Reimplemented from vtkDataObject.

Definition at line 337 of file vtkPolyData.h.

References vtkDataObject::SetUpdateExtent().

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

virtual int* vtkPolyData::GetUpdateExtent  )  [virtual]
 

We need this here to avoid hiding superclass method

Reimplemented from vtkDataObject.

virtual void vtkPolyData::GetUpdateExtent int &  x0,
int &  x1,
int &  y0,
int &  y1,
int &  z0,
int &  z1
[virtual]
 

We need this here to avoid hiding superclass method

Reimplemented from vtkDataObject.

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

We need this here to avoid hiding superclass method

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 355 of file vtkPolyData.h.

References vtkDataObject::SetUpdateExtent().

void vtkPolyData::SetUpdateExtent int  ext[6]  )  [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 357 of file vtkPolyData.h.

References vtkDataObject::SetUpdateExtent().

virtual int vtkPolyData::GetPiece  )  [virtual]
 

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

virtual int vtkPolyData::GetNumberOfPieces  )  [virtual]
 

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

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]
 

Shallow and Deep copy.

Reimplemented from vtkPointSet.

void vtkPolyData::RemoveGhostCells int  level  ) 
 

This method will remove any cell that has a ghost level array value greater or equal to level. It does not remove unused points (yet).

virtual void vtkPolyData::Crop  )  [protected, virtual]
 

This method crops the data object (if necesary) so that the extent matches the update extent.

Reimplemented from vtkDataObject.


Member Data Documentation

vtkVertex* vtkPolyData::Vertex [protected]
 

Definition at line 392 of file vtkPolyData.h.

vtkPolyVertex* vtkPolyData::PolyVertex [protected]
 

Definition at line 393 of file vtkPolyData.h.

vtkLine* vtkPolyData::Line [protected]
 

Definition at line 394 of file vtkPolyData.h.

vtkPolyLine* vtkPolyData::PolyLine [protected]
 

Definition at line 395 of file vtkPolyData.h.

vtkTriangle* vtkPolyData::Triangle [protected]
 

Definition at line 396 of file vtkPolyData.h.

vtkQuad* vtkPolyData::Quad [protected]
 

Definition at line 397 of file vtkPolyData.h.

vtkPolygon* vtkPolyData::Polygon [protected]
 

Definition at line 398 of file vtkPolyData.h.

vtkTriangleStrip* vtkPolyData::TriangleStrip [protected]
 

Definition at line 399 of file vtkPolyData.h.

vtkEmptyCell* vtkPolyData::EmptyCell [protected]
 

Definition at line 400 of file vtkPolyData.h.

vtkCellArray* vtkPolyData::Verts [protected]
 

Definition at line 404 of file vtkPolyData.h.

vtkCellArray* vtkPolyData::Lines [protected]
 

Definition at line 405 of file vtkPolyData.h.

vtkCellArray* vtkPolyData::Polys [protected]
 

Definition at line 406 of file vtkPolyData.h.

vtkCellArray* vtkPolyData::Strips [protected]
 

Definition at line 407 of file vtkPolyData.h.

vtkCellArray* vtkPolyData::Dummy [static, protected]
 

Definition at line 410 of file vtkPolyData.h.

vtkCellTypes* vtkPolyData::Cells [protected]
 

Definition at line 414 of file vtkPolyData.h.

Referenced by DeleteCell().

vtkCellLinks* vtkPolyData::Links [protected]
 

Definition at line 415 of file vtkPolyData.h.

Referenced by AddCellReference(), DeletePoint(), GetPointCells(), RemoveCellReference(), and ResizeCellList().


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