VTK
Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes
vtkAbstractCellLocator Class Reference

an abstract base class for locators which find cells More...

#include <vtkAbstractCellLocator.h>

Inheritance diagram for vtkAbstractCellLocator:
Inheritance graph
[legend]
Collaboration diagram for vtkAbstractCellLocator:
Collaboration graph
[legend]

List of all members.

Public Types

typedef vtkLocator Superclass

Public Member Functions

virtual const char * GetClassName ()
virtual int IsA (const char *type)
void PrintSelf (ostream &os, vtkIndent indent)
virtual void FindCellsWithinBounds (double *bbox, vtkIdList *cells)
virtual vtkIdType FindCell (double x[3])
virtual bool InsideCellBounds (double x[3], vtkIdType cell_ID)
virtual void SetNumberOfCellsPerNode (int)
virtual int GetNumberOfCellsPerNode ()
virtual void SetCacheCellBounds (int)
virtual int GetCacheCellBounds ()
virtual void CacheCellBoundsOn ()
virtual void CacheCellBoundsOff ()
virtual void SetRetainCellLists (int)
virtual int GetRetainCellLists ()
virtual void RetainCellListsOn ()
virtual void RetainCellListsOff ()
virtual void SetLazyEvaluation (int)
virtual int GetLazyEvaluation ()
virtual void LazyEvaluationOn ()
virtual void LazyEvaluationOff ()
virtual void SetUseExistingSearchStructure (int)
virtual int GetUseExistingSearchStructure ()
virtual void UseExistingSearchStructureOn ()
virtual void UseExistingSearchStructureOff ()
virtual int IntersectWithLine (double p1[3], double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId)
virtual int IntersectWithLine (double p1[3], double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId, vtkIdType &cellId)
virtual int IntersectWithLine (double p1[3], double p2[3], double tol, double &t, double x[3], double pcoords[3], int &subId, vtkIdType &cellId, vtkGenericCell *cell)
virtual int IntersectWithLine (const double p1[3], const double p2[3], vtkPoints *points, vtkIdList *cellIds)
virtual void FindClosestPoint (double x[3], double closestPoint[3], vtkIdType &cellId, int &subId, double &dist2)
virtual void FindClosestPoint (double x[3], double closestPoint[3], vtkGenericCell *cell, vtkIdType &cellId, int &subId, double &dist2)
virtual vtkIdType FindClosestPointWithinRadius (double x[3], double radius, double closestPoint[3], vtkIdType &cellId, int &subId, double &dist2)
virtual vtkIdType FindClosestPointWithinRadius (double x[3], double radius, double closestPoint[3], vtkGenericCell *cell, vtkIdType &cellId, int &subId, double &dist2)
virtual vtkIdType FindClosestPointWithinRadius (double x[3], double radius, double closestPoint[3], vtkGenericCell *cell, vtkIdType &cellId, int &subId, double &dist2, int &inside)
virtual void FindCellsAlongLine (double p1[3], double p2[3], double tolerance, vtkIdList *cells)
virtual vtkIdType FindCell (double x[3], double tol2, vtkGenericCell *GenCell, double pcoords[3], double *weights)

Static Public Member Functions

static int IsTypeOf (const char *type)
static vtkAbstractCellLocatorSafeDownCast (vtkObject *o)

Protected Member Functions

 vtkAbstractCellLocator ()
 ~vtkAbstractCellLocator ()
virtual bool StoreCellBounds ()
virtual void FreeCellBounds ()

Protected Attributes

int NumberOfCellsPerNode
int RetainCellLists
int CacheCellBounds
int LazyEvaluation
int UseExistingSearchStructure
vtkGenericCellGenericCell
double(* CellBounds )[6]

Detailed Description

an abstract base class for locators which find cells

vtkAbstractCellLocator is a spatial search object to quickly locate cells in 3D. vtkAbstractCellLocator supplies a basic interface which concrete subclasses should implement.

Warning:
When deriving a class from vtkAbstractCellLocator, one should include the 'hidden' member functions by the following construct in the derived class
 //BTX
  using vtkAbstractCellLocator::IntersectWithLine;
  using vtkAbstractCellLocator::FindClosestPoint;
  using vtkAbstractCellLocator::FindClosestPointWithinRadius;
 //ETX
 
See also:
vtkLocator vtkPointLocator vtkOBBTree vtkCellLocator

Definition at line 48 of file vtkAbstractCellLocator.h.


Member Typedef Documentation

Reimplemented from vtkLocator.

Reimplemented in vtkModifiedBSPTree, vtkOBBTree, vtkCellTreeLocator, and vtkCellLocator.

Definition at line 51 of file vtkAbstractCellLocator.h.


Constructor & Destructor Documentation

vtkAbstractCellLocator::vtkAbstractCellLocator ( ) [protected]
vtkAbstractCellLocator::~vtkAbstractCellLocator ( ) [protected]

Member Function Documentation

virtual const char* vtkAbstractCellLocator::GetClassName ( ) [virtual]

Reimplemented from vtkLocator.

Reimplemented in vtkModifiedBSPTree, vtkOBBTree, vtkCellTreeLocator, and vtkCellLocator.

static int vtkAbstractCellLocator::IsTypeOf ( const char *  name) [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 vtkLocator.

Reimplemented in vtkModifiedBSPTree, vtkOBBTree, vtkCellTreeLocator, and vtkCellLocator.

virtual int vtkAbstractCellLocator::IsA ( const char *  name) [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 vtkLocator.

Reimplemented in vtkModifiedBSPTree, vtkOBBTree, vtkCellTreeLocator, and vtkCellLocator.

static vtkAbstractCellLocator* vtkAbstractCellLocator::SafeDownCast ( vtkObject o) [static]

Reimplemented from vtkLocator.

Reimplemented in vtkModifiedBSPTree, vtkOBBTree, vtkCellTreeLocator, and vtkCellLocator.

void vtkAbstractCellLocator::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 vtkLocator.

Reimplemented in vtkModifiedBSPTree, vtkOBBTree, vtkCellTreeLocator, and vtkCellLocator.

virtual void vtkAbstractCellLocator::SetNumberOfCellsPerNode ( int  ) [virtual]

Specify the preferred/maximum number of cells in each node/bucket. Default 32. Locators generally operate by subdividing space into smaller regions until the number of cells in each region (or node) reaches the desired level.

virtual int vtkAbstractCellLocator::GetNumberOfCellsPerNode ( ) [virtual]

Specify the preferred/maximum number of cells in each node/bucket. Default 32. Locators generally operate by subdividing space into smaller regions until the number of cells in each region (or node) reaches the desired level.

virtual void vtkAbstractCellLocator::SetCacheCellBounds ( int  ) [virtual]

Boolean controls whether the bounds of each cell are computed only once and then saved. Should be 10 to 20% faster if repeatedly calling any of the Intersect/Find routines and the extra memory won't cause disk caching (24 extra bytes per cell are required to save the bounds).

virtual int vtkAbstractCellLocator::GetCacheCellBounds ( ) [virtual]

Boolean controls whether the bounds of each cell are computed only once and then saved. Should be 10 to 20% faster if repeatedly calling any of the Intersect/Find routines and the extra memory won't cause disk caching (24 extra bytes per cell are required to save the bounds).

virtual void vtkAbstractCellLocator::CacheCellBoundsOn ( ) [virtual]

Boolean controls whether the bounds of each cell are computed only once and then saved. Should be 10 to 20% faster if repeatedly calling any of the Intersect/Find routines and the extra memory won't cause disk caching (24 extra bytes per cell are required to save the bounds).

virtual void vtkAbstractCellLocator::CacheCellBoundsOff ( ) [virtual]

Boolean controls whether the bounds of each cell are computed only once and then saved. Should be 10 to 20% faster if repeatedly calling any of the Intersect/Find routines and the extra memory won't cause disk caching (24 extra bytes per cell are required to save the bounds).

virtual void vtkAbstractCellLocator::SetRetainCellLists ( int  ) [virtual]

Boolean controls whether to maintain list of cells in each node. not applicable to all implementations, but if the locator is being used as a geometry simplification technique, there is no need to keep them.

virtual int vtkAbstractCellLocator::GetRetainCellLists ( ) [virtual]

Boolean controls whether to maintain list of cells in each node. not applicable to all implementations, but if the locator is being used as a geometry simplification technique, there is no need to keep them.

virtual void vtkAbstractCellLocator::RetainCellListsOn ( ) [virtual]

Boolean controls whether to maintain list of cells in each node. not applicable to all implementations, but if the locator is being used as a geometry simplification technique, there is no need to keep them.

virtual void vtkAbstractCellLocator::RetainCellListsOff ( ) [virtual]

Boolean controls whether to maintain list of cells in each node. not applicable to all implementations, but if the locator is being used as a geometry simplification technique, there is no need to keep them.

virtual void vtkAbstractCellLocator::SetLazyEvaluation ( int  ) [virtual]

Most Locators build their search structures during BuildLocator but some may delay construction until it is actually needed. If LazyEvaluation is supported, this turns on/off the feature. if not supported, it is ignored.

virtual int vtkAbstractCellLocator::GetLazyEvaluation ( ) [virtual]

Most Locators build their search structures during BuildLocator but some may delay construction until it is actually needed. If LazyEvaluation is supported, this turns on/off the feature. if not supported, it is ignored.

virtual void vtkAbstractCellLocator::LazyEvaluationOn ( ) [virtual]

Most Locators build their search structures during BuildLocator but some may delay construction until it is actually needed. If LazyEvaluation is supported, this turns on/off the feature. if not supported, it is ignored.

virtual void vtkAbstractCellLocator::LazyEvaluationOff ( ) [virtual]

Most Locators build their search structures during BuildLocator but some may delay construction until it is actually needed. If LazyEvaluation is supported, this turns on/off the feature. if not supported, it is ignored.

virtual void vtkAbstractCellLocator::SetUseExistingSearchStructure ( int  ) [virtual]

Some locators support querying a new dataset without rebuilding the search structure (typically this may occur when a dataset changes due to a time update, but is actually the same topology) Turning on this flag enables some locators to skip the rebuilding phase

virtual int vtkAbstractCellLocator::GetUseExistingSearchStructure ( ) [virtual]

Some locators support querying a new dataset without rebuilding the search structure (typically this may occur when a dataset changes due to a time update, but is actually the same topology) Turning on this flag enables some locators to skip the rebuilding phase

virtual void vtkAbstractCellLocator::UseExistingSearchStructureOn ( ) [virtual]

Some locators support querying a new dataset without rebuilding the search structure (typically this may occur when a dataset changes due to a time update, but is actually the same topology) Turning on this flag enables some locators to skip the rebuilding phase

virtual void vtkAbstractCellLocator::UseExistingSearchStructureOff ( ) [virtual]

Some locators support querying a new dataset without rebuilding the search structure (typically this may occur when a dataset changes due to a time update, but is actually the same topology) Turning on this flag enables some locators to skip the rebuilding phase

virtual int vtkAbstractCellLocator::IntersectWithLine ( double  p1[3],
double  p2[3],
double  tol,
double t,
double  x[3],
double  pcoords[3],
int subId 
) [virtual]

Return intersection point (if any) of finite line with cells contained in cell locator.

Reimplemented in vtkModifiedBSPTree, vtkOBBTree, vtkCellTreeLocator, and vtkCellLocator.

virtual int vtkAbstractCellLocator::IntersectWithLine ( double  p1[3],
double  p2[3],
double  tol,
double t,
double  x[3],
double  pcoords[3],
int subId,
vtkIdType cellId 
) [virtual]

Return intersection point (if any) AND the cell which was intersected by the finite line.

Reimplemented in vtkModifiedBSPTree, vtkOBBTree, vtkCellTreeLocator, and vtkCellLocator.

virtual int vtkAbstractCellLocator::IntersectWithLine ( double  p1[3],
double  p2[3],
double  tol,
double t,
double  x[3],
double  pcoords[3],
int subId,
vtkIdType cellId,
vtkGenericCell cell 
) [virtual]

Return intersection point (if any) AND the cell which was intersected by the finite line. The cell is returned as a cell id and as a generic cell.

Reimplemented in vtkModifiedBSPTree, vtkOBBTree, vtkCellLocator, and vtkCellTreeLocator.

virtual int vtkAbstractCellLocator::IntersectWithLine ( const double  p1[3],
const double  p2[3],
vtkPoints points,
vtkIdList cellIds 
) [virtual]

Take the passed line segment and intersect it with the data set. This method assumes that the data set is a vtkPolyData that describes a closed surface, and the intersection points that are returned in 'points' alternate between entrance points and exit points. The return value of the function is 0 if no intersections were found, -1 if point 'a0' lies inside the closed surface, or +1 if point 'a0' lies outside the closed surface. Either 'points' or 'cellIds' can be set to NULL if you don't want to receive that information. This method is currently only implemented in vtkOBBTree

Reimplemented in vtkModifiedBSPTree, vtkOBBTree, vtkCellTreeLocator, and vtkCellLocator.

virtual void vtkAbstractCellLocator::FindClosestPoint ( double  x[3],
double  closestPoint[3],
vtkIdType cellId,
int subId,
double dist2 
) [virtual]

Return the closest point and the cell which is closest to the point x. The closest point is somewhere on a cell, it need not be one of the vertices of the cell.

Reimplemented in vtkOBBTree, and vtkCellLocator.

virtual void vtkAbstractCellLocator::FindClosestPoint ( double  x[3],
double  closestPoint[3],
vtkGenericCell cell,
vtkIdType cellId,
int subId,
double dist2 
) [virtual]

Return the closest point and the cell which is closest to the point x. The closest point is somewhere on a cell, it need not be one of the vertices of the cell. This version takes in a vtkGenericCell to avoid allocating and deallocating the cell. This is much faster than the version which does not take a *cell, especially when this function is called many times in a row such as by a for loop, where the allocation and deallocation can be done only once outside the for loop. If a cell is found, "cell" contains the points and ptIds for the cell "cellId" upon exit.

Reimplemented in vtkOBBTree, and vtkCellLocator.

virtual vtkIdType vtkAbstractCellLocator::FindClosestPointWithinRadius ( double  x[3],
double  radius,
double  closestPoint[3],
vtkIdType cellId,
int subId,
double dist2 
) [virtual]

Return the closest point within a specified radius and the cell which is closest to the point x. The closest point is somewhere on a cell, it need not be one of the vertices of the cell. This method returns 1 if a point is found within the specified radius. If there are no cells within the specified radius, the method returns 0 and the values of closestPoint, cellId, subId, and dist2 are undefined.

Reimplemented in vtkOBBTree, and vtkCellLocator.

virtual vtkIdType vtkAbstractCellLocator::FindClosestPointWithinRadius ( double  x[3],
double  radius,
double  closestPoint[3],
vtkGenericCell cell,
vtkIdType cellId,
int subId,
double dist2 
) [virtual]

Return the closest point within a specified radius and the cell which is closest to the point x. The closest point is somewhere on a cell, it need not be one of the vertices of the cell. This method returns 1 if a point is found within the specified radius. If there are no cells within the specified radius, the method returns 0 and the values of closestPoint, cellId, subId, and dist2 are undefined. This version takes in a vtkGenericCell to avoid allocating and deallocating the cell. This is much faster than the version which does not take a *cell, especially when this function is called many times in a row such as by a for loop, where the allocation and deallocation can be done only once outside the for loop. If a closest point is found, "cell" contains the points and ptIds for the cell "cellId" upon exit.

Reimplemented in vtkOBBTree, and vtkCellLocator.

virtual vtkIdType vtkAbstractCellLocator::FindClosestPointWithinRadius ( double  x[3],
double  radius,
double  closestPoint[3],
vtkGenericCell cell,
vtkIdType cellId,
int subId,
double dist2,
int inside 
) [virtual]

Return the closest point within a specified radius and the cell which is closest to the point x. The closest point is somewhere on a cell, it need not be one of the vertices of the cell. This method returns 1 if a point is found within the specified radius. If there are no cells within the specified radius, the method returns 0 and the values of closestPoint, cellId, subId, and dist2 are undefined. This version takes in a vtkGenericCell to avoid allocating and deallocating the cell. This is much faster than the version which does not take a *cell, especially when this function is called many times in a row such as by a for loop, where the allocation and dealloction can be done only once outside the for loop. If a closest point is found, "cell" contains the points and ptIds for the cell "cellId" upon exit. If a closest point is found, inside returns the return value of the EvaluatePosition call to the closest cell; inside(=1) or outside(=0).

Reimplemented in vtkOBBTree, and vtkCellLocator.

virtual void vtkAbstractCellLocator::FindCellsWithinBounds ( double bbox,
vtkIdList cells 
) [virtual]

Return a list of unique cell ids inside of a given bounding box. The user must provide the vtkIdList to populate. This method returns data only after the locator has been built.

Reimplemented in vtkCellLocator, and vtkCellTreeLocator.

virtual void vtkAbstractCellLocator::FindCellsAlongLine ( double  p1[3],
double  p2[3],
double  tolerance,
vtkIdList cells 
) [virtual]

Given a finite line defined by the two points (p1,p2), return the list of unique cell ids in the buckets containing the line. It is possible that an empty cell list is returned. The user must provide the vtkIdList to populate. This method returns data only after the locator has been built.

Reimplemented in vtkCellLocator.

virtual vtkIdType vtkAbstractCellLocator::FindCell ( double  x[3]) [virtual]

Returns the Id of the cell containing the point, returns -1 if no cell found. This interface uses a tolerance of zero

Reimplemented in vtkModifiedBSPTree, vtkCellLocator, and vtkCellTreeLocator.

virtual vtkIdType vtkAbstractCellLocator::FindCell ( double  x[3],
double  tol2,
vtkGenericCell GenCell,
double  pcoords[3],
double weights 
) [virtual]

Find the cell containing a given point. returns -1 if no cell found the cell parameters are copied into the supplied variables, a cell must be provided to store the information.

Reimplemented in vtkModifiedBSPTree, vtkCellLocator, and vtkCellTreeLocator.

virtual bool vtkAbstractCellLocator::InsideCellBounds ( double  x[3],
vtkIdType  cell_ID 
) [virtual]

Quickly test if a point is inside the bounds of a particular cell. Some locators cache cell bounds and this function can make use of fast access to the data.

Reimplemented in vtkModifiedBSPTree.

virtual bool vtkAbstractCellLocator::StoreCellBounds ( ) [protected, virtual]

This command is used internally by the locator to copy all cell Bounds into the internal CellBounds array. Subsequent calls to InsideCellBounds(...) can make use of the data A valid dataset must be present for this to work. Returns true if bounds wre copied, false otherwise.

virtual void vtkAbstractCellLocator::FreeCellBounds ( ) [protected, virtual]

This command is used internally by the locator to copy all cell Bounds into the internal CellBounds array. Subsequent calls to InsideCellBounds(...) can make use of the data A valid dataset must be present for this to work. Returns true if bounds wre copied, false otherwise.


Member Data Documentation

Definition at line 273 of file vtkAbstractCellLocator.h.

Definition at line 274 of file vtkAbstractCellLocator.h.

Definition at line 275 of file vtkAbstractCellLocator.h.

Definition at line 276 of file vtkAbstractCellLocator.h.

Definition at line 277 of file vtkAbstractCellLocator.h.

Definition at line 278 of file vtkAbstractCellLocator.h.

Definition at line 280 of file vtkAbstractCellLocator.h.


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