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

vtkDelaunay2D Class Reference

create 2D Delaunay triangulation of input points. More...

#include <vtkDelaunay2D.h>

Inheritance diagram for vtkDelaunay2D:

Inheritance graph
[legend]
Collaboration diagram for vtkDelaunay2D:

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)
void SetSource (vtkPolyData *)
vtkPolyDataGetSource ()
virtual void SetAlpha (double)
virtual double GetAlpha ()
virtual void SetTolerance (double)
virtual double GetTolerance ()
virtual void SetOffset (double)
virtual double GetOffset ()
virtual void SetBoundingTriangulation (int)
virtual int GetBoundingTriangulation ()
virtual void BoundingTriangulationOn ()
virtual void BoundingTriangulationOff ()
virtual void SetInput (vtkPointSet *input)
vtkPointSetGetInput ()

Static Public Methods

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

Protected Methods

 vtkDelaunay2D ()
 ~vtkDelaunay2D ()
 vtkDelaunay2D (const vtkDelaunay2D &)
void operator= (const vtkDelaunay2D &)
void Execute ()

Protected Attributes

double Alpha
double Tolerance
int BoundingTriangulation
double Offset

Detailed Description

create 2D Delaunay triangulation of input points.

Date:
2000/12/10 20:08:35
Revision:
1.28

vtkDelaunay2D is a filter that constructs a 2D Delaunay triangulation from a list of input points. These points may be represented by any dataset of type vtkPointSet and subclasses. The output of the filter is a polygonal dataset. Usually the output is a triangle mesh, but if a non-zero alpha distance value is specified (called the "alpha" value), then only triangles, edges, and vertices lying within the alpha radius are output. In other words, non-zero alpha values may result in arbitrary combinations of triangles, lines, and vertices. (The notion of alpha value is derived from Edelsbrunner's work on "alpha shapes".) Also, it is possible to generate "constrained triangulations" using this filter. A constrained triangulation is one where edges and loops (i.e., polygons) can be defined and the triangulation will preserve them (read on for more information).

The 2D Delaunay triangulation is defined as the triangulation that satisfies the Delaunay criterion for n-dimensional simplexes (in this case n=2 and the simplexes are triangles). This criterion states that a circumsphere of each simplex in a triangulation contains only the n+1 defining points of the simplex. (See "The Visualization Toolkit" text for more information.) In two dimensions, this translates into an optimal triangulation. That is, the maximum interior angle of any triangle is less than or equal to that of any possible triangulation.

Delaunay triangulations are used to build topological structures from unorganized (or unstructured) points. The input to this filter is a list of points specified in 3D, even though the triangulation is 2D. Thus the triangulation is constructed in the x-y plane, and the z coordinate is ignored (although carried through to the output). (If you desire to triangulate in a different plane, you'll have to use the vtkTransformFilter to transform the points into and out of the x-y plane.)

The Delaunay triangulation can be numerically sensitive in some cases. To prevent problems, try to avoid injecting points that will result in triangles with bad aspect ratios (1000:1 or greater). In practice this means inserting points that are "widely dispersed", and enables smooth transition of triangle sizes throughout the mesh. (You may even want to add extra points to create a better point distribution.) If numerical problems are present, you will see a warning message to this effect at the end of the triangulation process.

To create constrained meshes, you must define an additional input. This input is an instance of vtkPolyData which contains lines, polylines, and/or polygons that define constrained edges and loops. Lines and polylines found in the input will be mesh edges in the output. Polygons define a loop with inside and outside regions. The inside of the polygon is determined by using the right-hand-rule, i.e., looking down the z-axis a polygon should be ordered counter-clockwise. Holes in a polygon should be ordered clockwise. If you choose to create a constrained triangulation, the final mesh may not satisfy the Delaunay criterion. (Noted: the lines/polygon edges must not intersect when projected onto the 2D plane. It may not be possible to recover all edges due to not enough points in the triangulation, or poorly defined edges (coincident or excessively long). The form of the lines or polygons is a list of point ids that correspond to the input point ids used to generate the triangulation.)

Warning:
Points arranged on a regular lattice (termed degenerate cases) can be triangulated in more than one way (at least according to the Delaunay criterion). The choice of triangulation (as implemented by this algorithm) depends on the order of the input points. The first three points will form a triangle; other degenerate points will not break this triangle.
Warning:
Points that are coincident (or nearly so) may be discarded by the algorithm. This is because the Delaunay triangulation requires unique input points. You can control the definition of coincidence with the "Tolerance" instance variable.
Warning:
The output of the Delaunay triangulation is supposedly a convex hull. In certain cases this implementation may not generate the convex hull. This behavior can be controlled by the Offset instance variable. Offset is a multiplier used to control the size of the initial triangulation. The larger the offset value, the more likely you will generate a convex hull; but the more likely you are to see numerical problems.
See also:
vtkDelaunay3D vtkTransformFilter vtkGaussianSplatter
Examples:
vtkDelaunay2D (examples)

Definition at line 137 of file vtkDelaunay2D.h.


Constructor & Destructor Documentation

vtkDelaunay2D::vtkDelaunay2D   [protected]
 

vtkDelaunay2D::~vtkDelaunay2D   [protected]
 

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

Definition at line 186 of file vtkDelaunay2D.h.


Member Function Documentation

virtual const char* vtkDelaunay2D::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 vtkPolyDataSource.

int vtkDelaunay2D::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 vtkPolyDataSource.

virtual int vtkDelaunay2D::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 vtkPolyDataSource.

vtkDelaunay2D* vtkDelaunay2D::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 vtkPolyDataSource.

void vtkDelaunay2D::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 vtkSource.

vtkDelaunay2D* vtkDelaunay2D::New   [static]
 

Construct object with Alpha = 0.0; Tolerance = 0.001; Offset = 1.25; BoundingTriangulation turned off.

Reimplemented from vtkPolyDataSource.

void vtkDelaunay2D::SetSource vtkPolyData  
 

Specify the source object used to specify constrained edges and loops. (This is optional.) If set, and lines/polygons are defined, a constrained triangulation is created.

vtkPolyData* vtkDelaunay2D::GetSource  
 

virtual void vtkDelaunay2D::SetAlpha double    [virtual]
 

Specify alpha (or distance) value to control output of this filter. For a non-zero alpha value, only edges or triangles contained within a sphere centered at mesh vertices will be output. Otherwise, only triangles will be output.

virtual double vtkDelaunay2D::GetAlpha   [virtual]
 

virtual void vtkDelaunay2D::SetTolerance double    [virtual]
 

Specify a tolerance to control discarding of closely spaced points. This tolerance is specified as a fraction of the diagonal length of the bounding box of the points.

virtual double vtkDelaunay2D::GetTolerance   [virtual]
 

virtual void vtkDelaunay2D::SetOffset double    [virtual]
 

Specify a multiplier to control the size of the initial, bounding Delaunay triangulation.

virtual double vtkDelaunay2D::GetOffset   [virtual]
 

virtual void vtkDelaunay2D::SetBoundingTriangulation int    [virtual]
 

Boolean controls whether bounding triangulation points (and associated triangles) are included in the output. (These are introduced as an initial triangulation to begin the triangulation process. This feature is nice for debugging output.)

virtual int vtkDelaunay2D::GetBoundingTriangulation   [virtual]
 

virtual void vtkDelaunay2D::BoundingTriangulationOn   [virtual]
 

virtual void vtkDelaunay2D::BoundingTriangulationOff   [virtual]
 

virtual void vtkDelaunay2D::SetInput vtkPointSet   input [virtual]
 

Set / get the input data or filter.

vtkPointSet* vtkDelaunay2D::GetInput  
 

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

Definition at line 187 of file vtkDelaunay2D.h.

void vtkDelaunay2D::Execute   [protected, virtual]
 

Reimplemented from vtkSource.


Member Data Documentation

double vtkDelaunay2D::Alpha [protected]
 

Definition at line 191 of file vtkDelaunay2D.h.

double vtkDelaunay2D::Tolerance [protected]
 

Definition at line 192 of file vtkDelaunay2D.h.

int vtkDelaunay2D::BoundingTriangulation [protected]
 

Definition at line 193 of file vtkDelaunay2D.h.

double vtkDelaunay2D::Offset [protected]
 

Definition at line 194 of file vtkDelaunay2D.h.


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