VTK  9.3.20240328
Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
vtkDelaunay2D Class Reference

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

#include <vtkDelaunay2D.h>

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

Public Types

typedef vtkPolyDataAlgorithm Superclass
 
- Public Types inherited from vtkPolyDataAlgorithm
typedef vtkAlgorithm Superclass
 
- Public Types inherited from vtkAlgorithm
enum  DesiredOutputPrecision { SINGLE_PRECISION , DOUBLE_PRECISION , DEFAULT_PRECISION }
 Values used for setting the desired output precision for various algorithms. More...
 
typedef vtkObject Superclass
 

Public Member Functions

virtual vtkTypeBool IsA (const char *type)
 Return 1 if this class is the same type of (or a subclass of) the named class. More...
 
vtkDelaunay2DNewInstance () const
 
void PrintSelf (ostream &os, vtkIndent indent) override
 Methods invoked by print to print information about the object including superclasses. More...
 
void SetSourceData (vtkPolyData *)
 Specify the source object used to specify constrained edges and loops. More...
 
void SetSourceConnection (vtkAlgorithmOutput *algOutput)
 Specify the source object used to specify constrained edges and loops. More...
 
vtkPolyDataGetSource ()
 Get a pointer to the source object. More...
 
virtual void SetAlpha (double)
 Specify alpha (or distance) value to control output of this filter. More...
 
virtual double GetAlpha ()
 Specify alpha (or distance) value to control output of this filter. More...
 
virtual void SetTolerance (double)
 Specify a tolerance to control discarding of closely spaced points. More...
 
virtual double GetTolerance ()
 Specify a tolerance to control discarding of closely spaced points. More...
 
virtual void SetOffset (double)
 Specify a multiplier to control the size of the initial, bounding Delaunay triangulation. More...
 
virtual double GetOffset ()
 Specify a multiplier to control the size of the initial, bounding Delaunay triangulation. More...
 
virtual void SetBoundingTriangulation (vtkTypeBool)
 Boolean controls whether bounding triangulation points (and associated triangles) are included in the output. More...
 
virtual vtkTypeBool GetBoundingTriangulation ()
 Boolean controls whether bounding triangulation points (and associated triangles) are included in the output. More...
 
virtual void BoundingTriangulationOn ()
 Boolean controls whether bounding triangulation points (and associated triangles) are included in the output. More...
 
virtual void BoundingTriangulationOff ()
 Boolean controls whether bounding triangulation points (and associated triangles) are included in the output. More...
 
 vtkSetSmartPointerMacro (Transform, vtkAbstractTransform)
 Set / get the transform which is applied to points to generate a 2D problem. More...
 
 vtkGetSmartPointerMacro (Transform, vtkAbstractTransform)
 Set / get the transform which is applied to points to generate a 2D problem. More...
 
virtual void SetProjectionPlaneMode (int)
 Define the method to project the input 3D points into a 2D plane for triangulation. More...
 
virtual int GetProjectionPlaneMode ()
 Define the method to project the input 3D points into a 2D plane for triangulation. More...
 
virtual void SetRandomPointInsertion (vtkTypeBool)
 Indicate whether to insert the points in given order, or pseudo-random order. More...
 
virtual vtkTypeBool GetRandomPointInsertion ()
 Indicate whether to insert the points in given order, or pseudo-random order. More...
 
virtual void RandomPointInsertionOn ()
 Indicate whether to insert the points in given order, or pseudo-random order. More...
 
virtual void RandomPointInsertionOff ()
 Indicate whether to insert the points in given order, or pseudo-random order. More...
 
- Public Member Functions inherited from vtkPolyDataAlgorithm
vtkPolyDataAlgorithmNewInstance () const
 
vtkTypeBool ProcessRequest (vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
 see vtkAlgorithm for details More...
 
vtkDataObjectGetInput ()
 
vtkDataObjectGetInput (int port)
 
vtkPolyDataGetPolyDataInput (int port)
 
vtkPolyDataGetOutput ()
 Get the output data object for a port on this algorithm. More...
 
vtkPolyDataGetOutput (int)
 Get the output data object for a port on this algorithm. More...
 
virtual void SetOutput (vtkDataObject *d)
 Get the output data object for a port on this algorithm. More...
 
void SetInputData (vtkDataObject *)
 Assign a data object as input. More...
 
void SetInputData (int, vtkDataObject *)
 Assign a data object as input. More...
 
void AddInputData (vtkDataObject *)
 Assign a data object as input. More...
 
void AddInputData (int, vtkDataObject *)
 Assign a data object as input. More...
 
- Public Member Functions inherited from vtkAlgorithm
vtkAlgorithmNewInstance () const
 
vtkTypeBool HasExecutive ()
 Check whether this algorithm has an assigned executive. More...
 
vtkExecutiveGetExecutive ()
 Get this algorithm's executive. More...
 
virtual void SetExecutive (vtkExecutive *executive)
 Set this algorithm's executive. More...
 
vtkTypeBool ProcessRequest (vtkInformation *request, vtkCollection *inInfo, vtkInformationVector *outInfo)
 Version of ProcessRequest() that is wrapped. More...
 
virtual int ComputePipelineMTime (vtkInformation *request, vtkInformationVector **inInfoVec, vtkInformationVector *outInfoVec, int requestFromOutputPort, vtkMTimeType *mtime)
 A special version of ProcessRequest meant specifically for the pipeline modified time request. More...
 
virtual int ModifyRequest (vtkInformation *request, int when)
 This method gives the algorithm a chance to modify the contents of a request before or after (specified in the when argument) it is forwarded. More...
 
vtkInformationGetInputPortInformation (int port)
 Get the information object associated with an input port. More...
 
vtkInformationGetOutputPortInformation (int port)
 Get the information object associated with an output port. More...
 
int GetNumberOfInputPorts ()
 Get the number of input ports used by the algorithm. More...
 
int GetNumberOfOutputPorts ()
 Get the number of output ports provided by the algorithm. More...
 
void SetAbortExecuteAndUpdateTime ()
 Set AbortExecute Flag and update LastAbortTime. More...
 
void UpdateProgress (double amount)
 Update the progress of the process object. More...
 
bool CheckAbort ()
 Checks to see if this filter should abort. More...
 
virtual void SetInputArrayToProcess (int idx, int port, int connection, const char *fieldAssociation, const char *attributeTypeorName)
 String based versions of SetInputArrayToProcess(). More...
 
vtkInformationGetInputArrayInformation (int idx)
 Get the info object for the specified input array to this algorithm. More...
 
void RemoveAllInputs ()
 Remove all the input data. More...
 
vtkDataObjectGetOutputDataObject (int port)
 Get the data object that will contain the algorithm output for the given port. More...
 
vtkDataObjectGetInputDataObject (int port, int connection)
 Get the data object that will contain the algorithm input for the given port and given connection. More...
 
virtual void RemoveInputConnection (int port, vtkAlgorithmOutput *input)
 Remove a connection from the given input port index. More...
 
virtual void RemoveInputConnection (int port, int idx)
 Remove a connection given by index idx. More...
 
virtual void RemoveAllInputConnections (int port)
 Removes all input connections. More...
 
virtual void SetInputDataObject (int port, vtkDataObject *data)
 Sets the data-object as an input on the given port index. More...
 
virtual void SetInputDataObject (vtkDataObject *data)
 
virtual void AddInputDataObject (int port, vtkDataObject *data)
 Add the data-object as an input to this given port. More...
 
virtual void AddInputDataObject (vtkDataObject *data)
 
vtkAlgorithmOutputGetOutputPort (int index)
 Get a proxy object corresponding to the given output port of this algorithm. More...
 
vtkAlgorithmOutputGetOutputPort ()
 
int GetNumberOfInputConnections (int port)
 Get the number of inputs currently connected to a port. More...
 
int GetTotalNumberOfInputConnections ()
 Get the total number of inputs for this algorithm. More...
 
vtkAlgorithmOutputGetInputConnection (int port, int index)
 Get the algorithm output port connected to an input port. More...
 
vtkAlgorithmGetInputAlgorithm (int port, int index, int &algPort)
 Returns the algorithm and the output port index of that algorithm connected to a port-index pair. More...
 
vtkAlgorithmGetInputAlgorithm (int port, int index)
 Returns the algorithm connected to a port-index pair. More...
 
vtkAlgorithmGetInputAlgorithm ()
 Equivalent to GetInputAlgorithm(0, 0). More...
 
vtkExecutiveGetInputExecutive (int port, int index)
 Returns the executive associated with a particular input connection. More...
 
vtkExecutiveGetInputExecutive ()
 Equivalent to GetInputExecutive(0, 0) More...
 
vtkInformationGetInputInformation (int port, int index)
 Return the information object that is associated with a particular input connection. More...
 
vtkInformationGetInputInformation ()
 Equivalent to GetInputInformation(0, 0) More...
 
vtkInformationGetOutputInformation (int port)
 Return the information object that is associated with a particular output port. More...
 
virtual vtkTypeBool Update (int port, vtkInformationVector *requests)
 This method enables the passing of data requests to the algorithm to be used during execution (in addition to bringing a particular port up-to-date). More...
 
virtual vtkTypeBool Update (vtkInformation *requests)
 Convenience method to update an algorithm after passing requests to its first output port. More...
 
virtual int UpdatePiece (int piece, int numPieces, int ghostLevels, const int extents[6]=nullptr)
 Convenience method to update an algorithm after passing requests to its first output port. More...
 
virtual int UpdateExtent (const int extents[6])
 Convenience method to update an algorithm after passing requests to its first output port. More...
 
virtual int UpdateTimeStep (double time, int piece=-1, int numPieces=1, int ghostLevels=0, const int extents[6]=nullptr)
 Convenience method to update an algorithm after passing requests to its first output port. More...
 
virtual void UpdateInformation ()
 Bring the algorithm's information up-to-date. More...
 
virtual void UpdateDataObject ()
 Create output object(s). More...
 
virtual void PropagateUpdateExtent ()
 Propagate meta-data upstream. More...
 
virtual void UpdateWholeExtent ()
 Bring this algorithm's outputs up-to-date. More...
 
void ConvertTotalInputToPortConnection (int ind, int &port, int &conn)
 Convenience routine to convert from a linear ordering of input connections to a port/connection pair. More...
 
void RemoveNoPriorTemporalAccessInformationKey ()
 Removes any information key vtkStreamingDemandDrivenPipeline::NO_PRIOR_TEMPORAL_ACCESS() to all output ports of this vtkAlgorithm. More...
 
virtual vtkInformationGetInformation ()
 Set/Get the information object associated with this algorithm. More...
 
virtual void SetInformation (vtkInformation *)
 Set/Get the information object associated with this algorithm. More...
 
bool UsesGarbageCollector () const override
 Participate in garbage collection. More...
 
virtual void SetAbortExecute (vtkTypeBool)
 Set/Get the AbortExecute flag for the process object. More...
 
virtual vtkTypeBool GetAbortExecute ()
 Set/Get the AbortExecute flag for the process object. More...
 
virtual void AbortExecuteOn ()
 Set/Get the AbortExecute flag for the process object. More...
 
virtual void AbortExecuteOff ()
 Set/Get the AbortExecute flag for the process object. More...
 
virtual double GetProgress ()
 Get the execution progress of a process object. More...
 
void SetContainerAlgorithm (vtkAlgorithm *containerAlg)
 Set/get a Container algorithm for this algorithm. More...
 
vtkAlgorithmGetContainerAlgorithm ()
 Set/get a Container algorithm for this algorithm. More...
 
virtual void SetAbortOutput (bool)
 Set/Get an internal variable used to communicate between the algorithm and executive. More...
 
virtual bool GetAbortOutput ()
 Set/Get an internal variable used to communicate between the algorithm and executive. More...
 
void SetProgressShiftScale (double shift, double scale)
 Specify the shift and scale values to use to apply to the progress amount when UpdateProgress is called. More...
 
virtual double GetProgressShift ()
 Specify the shift and scale values to use to apply to the progress amount when UpdateProgress is called. More...
 
virtual double GetProgressScale ()
 Specify the shift and scale values to use to apply to the progress amount when UpdateProgress is called. More...
 
void SetProgressText (const char *ptext)
 Set the current text message associated with the progress state. More...
 
virtual char * GetProgressText ()
 Set the current text message associated with the progress state. More...
 
virtual unsigned long GetErrorCode ()
 The error code contains a possible error that occurred while reading or writing the file. More...
 
virtual void SetInputArrayToProcess (int idx, int port, int connection, int fieldAssociation, const char *name)
 Set the input data arrays that this algorithm will process. More...
 
virtual void SetInputArrayToProcess (int idx, int port, int connection, int fieldAssociation, int fieldAttributeType)
 Set the input data arrays that this algorithm will process. More...
 
virtual void SetInputArrayToProcess (int idx, vtkInformation *info)
 Set the input data arrays that this algorithm will process. More...
 
virtual void SetInputConnection (int port, vtkAlgorithmOutput *input)
 Set the connection for the given input port index. More...
 
virtual void SetInputConnection (vtkAlgorithmOutput *input)
 Set the connection for the given input port index. More...
 
virtual void AddInputConnection (int port, vtkAlgorithmOutput *input)
 Add a connection to the given input port index. More...
 
virtual void AddInputConnection (vtkAlgorithmOutput *input)
 Add a connection to the given input port index. More...
 
virtual void Update (int port)
 Bring this algorithm's outputs up-to-date. More...
 
virtual void Update ()
 Bring this algorithm's outputs up-to-date. More...
 
virtual void SetReleaseDataFlag (vtkTypeBool)
 Turn release data flag on or off for all output ports. More...
 
virtual vtkTypeBool GetReleaseDataFlag ()
 Turn release data flag on or off for all output ports. More...
 
void ReleaseDataFlagOn ()
 Turn release data flag on or off for all output ports. More...
 
void ReleaseDataFlagOff ()
 Turn release data flag on or off for all output ports. More...
 
int UpdateExtentIsEmpty (vtkInformation *pinfo, vtkDataObject *output)
 This detects when the UpdateExtent will generate no data This condition is satisfied when the UpdateExtent has zero volume (0,-1,...) or the UpdateNumberOfPieces is 0. More...
 
int UpdateExtentIsEmpty (vtkInformation *pinfo, int extentType)
 This detects when the UpdateExtent will generate no data This condition is satisfied when the UpdateExtent has zero volume (0,-1,...) or the UpdateNumberOfPieces is 0. More...
 
int * GetUpdateExtent ()
 These functions return the update extent for output ports that use 3D extents. More...
 
int * GetUpdateExtent (int port)
 These functions return the update extent for output ports that use 3D extents. More...
 
void GetUpdateExtent (int &x0, int &x1, int &y0, int &y1, int &z0, int &z1)
 These functions return the update extent for output ports that use 3D extents. More...
 
void GetUpdateExtent (int port, int &x0, int &x1, int &y0, int &y1, int &z0, int &z1)
 These functions return the update extent for output ports that use 3D extents. More...
 
void GetUpdateExtent (int extent[6])
 These functions return the update extent for output ports that use 3D extents. More...
 
void GetUpdateExtent (int port, int extent[6])
 These functions return the update extent for output ports that use 3D extents. More...
 
int GetUpdatePiece ()
 These functions return the update extent for output ports that use piece extents. More...
 
int GetUpdatePiece (int port)
 These functions return the update extent for output ports that use piece extents. More...
 
int GetUpdateNumberOfPieces ()
 These functions return the update extent for output ports that use piece extents. More...
 
int GetUpdateNumberOfPieces (int port)
 These functions return the update extent for output ports that use piece extents. More...
 
int GetUpdateGhostLevel ()
 These functions return the update extent for output ports that use piece extents. More...
 
int GetUpdateGhostLevel (int port)
 These functions return the update extent for output ports that use piece extents. More...
 
void SetProgressObserver (vtkProgressObserver *)
 If an ProgressObserver is set, the algorithm will report progress through it rather than directly. More...
 
virtual vtkProgressObserverGetProgressObserver ()
 If an ProgressObserver is set, the algorithm will report progress through it rather than directly. More...
 
void SetNoPriorTemporalAccessInformationKey (int key)
 Set to all output ports of this algorithm the information key vtkStreamingDemandDrivenPipeline::NO_PRIOR_TEMPORAL_ACCESS(). More...
 
void SetNoPriorTemporalAccessInformationKey ()
 Set to all output ports of this algorithm the information key vtkStreamingDemandDrivenPipeline::NO_PRIOR_TEMPORAL_ACCESS(). More...
 
- Public Member Functions inherited from vtkObject
 vtkBaseTypeMacro (vtkObject, vtkObjectBase)
 
virtual void DebugOn ()
 Turn debugging output on. More...
 
virtual void DebugOff ()
 Turn debugging output off. More...
 
bool GetDebug ()
 Get the value of the debug flag. More...
 
void SetDebug (bool debugFlag)
 Set the value of the debug flag. More...
 
virtual void Modified ()
 Update the modification time for this object. More...
 
virtual vtkMTimeType GetMTime ()
 Return this object's modified time. More...
 
void RemoveObserver (unsigned long tag)
 
void RemoveObservers (unsigned long event)
 
void RemoveObservers (const char *event)
 
void RemoveAllObservers ()
 
vtkTypeBool HasObserver (unsigned long event)
 
vtkTypeBool HasObserver (const char *event)
 
vtkTypeBool InvokeEvent (unsigned long event)
 
vtkTypeBool InvokeEvent (const char *event)
 
std::string GetObjectDescription () const override
 The object description printed in messages and PrintSelf output. More...
 
unsigned long AddObserver (unsigned long event, vtkCommand *, float priority=0.0f)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
unsigned long AddObserver (const char *event, vtkCommand *, float priority=0.0f)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
vtkCommandGetCommand (unsigned long tag)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
void RemoveObserver (vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
void RemoveObservers (unsigned long event, vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
void RemoveObservers (const char *event, vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
vtkTypeBool HasObserver (unsigned long event, vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
vtkTypeBool HasObserver (const char *event, vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object. More...
 
template<class U , class T >
unsigned long AddObserver (unsigned long event, U observer, void(T::*callback)(), float priority=0.0f)
 Overloads to AddObserver that allow developers to add class member functions as callbacks for events. More...
 
template<class U , class T >
unsigned long AddObserver (unsigned long event, U observer, void(T::*callback)(vtkObject *, unsigned long, void *), float priority=0.0f)
 Overloads to AddObserver that allow developers to add class member functions as callbacks for events. More...
 
template<class U , class T >
unsigned long AddObserver (unsigned long event, U observer, bool(T::*callback)(vtkObject *, unsigned long, void *), float priority=0.0f)
 Allow user to set the AbortFlagOn() with the return value of the callback method. More...
 
vtkTypeBool InvokeEvent (unsigned long event, void *callData)
 This method invokes an event and return whether the event was aborted or not. More...
 
vtkTypeBool InvokeEvent (const char *event, void *callData)
 This method invokes an event and return whether the event was aborted or not. More...
 
virtual void SetObjectName (const std::string &objectName)
 Set/get the name of this object for reporting purposes. More...
 
virtual std::string GetObjectName () const
 Set/get the name of this object for reporting purposes. More...
 
- Public Member Functions inherited from vtkObjectBase
const char * GetClassName () const
 Return the class name as a string. More...
 
virtual vtkIdType GetNumberOfGenerationsFromBase (const char *name)
 Given the name of a base class of this class type, return the distance of inheritance between this class type and the named class (how many generations of inheritance are there between this class and the named class). More...
 
virtual void Delete ()
 Delete a VTK object. More...
 
virtual void FastDelete ()
 Delete a reference to this object. More...
 
void InitializeObjectBase ()
 
void Print (ostream &os)
 Print an object to an ostream. More...
 
void Register (vtkObjectBase *o)
 Increase the reference count (mark as used by another object). More...
 
virtual void UnRegister (vtkObjectBase *o)
 Decrease the reference count (release by another object). More...
 
int GetReferenceCount ()
 Return the current reference count of this object. More...
 
void SetReferenceCount (int)
 Sets the reference count. More...
 
bool GetIsInMemkind () const
 A local state flag that remembers whether this object lives in the normal or extended memory space. More...
 
virtual void PrintHeader (ostream &os, vtkIndent indent)
 Methods invoked by print to print information about the object including superclasses. More...
 
virtual void PrintTrailer (ostream &os, vtkIndent indent)
 Methods invoked by print to print information about the object including superclasses. More...
 

Static Public Member Functions

static vtkTypeBool IsTypeOf (const char *type)
 
static vtkDelaunay2DSafeDownCast (vtkObjectBase *o)
 
static vtkDelaunay2DNew ()
 Construct object with Alpha = 0.0; Tolerance = 0.001; Offset = 1.25; BoundingTriangulation turned off. More...
 
static vtkAbstractTransformComputeBestFittingPlane (vtkPointSet *input)
 This method computes the best fit plane to a set of points represented by a vtkPointSet. More...
 
- Static Public Member Functions inherited from vtkPolyDataAlgorithm
static vtkPolyDataAlgorithmNew ()
 
static vtkTypeBool IsTypeOf (const char *type)
 
static vtkPolyDataAlgorithmSafeDownCast (vtkObjectBase *o)
 
- Static Public Member Functions inherited from vtkAlgorithm
static vtkAlgorithmNew ()
 
static vtkTypeBool IsTypeOf (const char *type)
 
static vtkAlgorithmSafeDownCast (vtkObjectBase *o)
 
static vtkInformationIntegerKeyINPUT_IS_OPTIONAL ()
 Keys used to specify input port requirements. More...
 
static vtkInformationIntegerKeyINPUT_IS_REPEATABLE ()
 
static vtkInformationInformationVectorKeyINPUT_REQUIRED_FIELDS ()
 
static vtkInformationStringVectorKeyINPUT_REQUIRED_DATA_TYPE ()
 
static vtkInformationInformationVectorKeyINPUT_ARRAYS_TO_PROCESS ()
 
static vtkInformationIntegerKeyINPUT_PORT ()
 
static vtkInformationIntegerKeyINPUT_CONNECTION ()
 
static vtkInformationIntegerKeyCAN_PRODUCE_SUB_EXTENT ()
 This key tells the executive that a particular output port is capable of producing an arbitrary subextent of the whole extent. More...
 
static vtkInformationIntegerKeyCAN_HANDLE_PIECE_REQUEST ()
 Key that tells the pipeline that a particular algorithm can or cannot handle piece request. More...
 
static vtkInformationIntegerKeyABORTED ()
 
static void SetDefaultExecutivePrototype (vtkExecutive *proto)
 If the DefaultExecutivePrototype is set, a copy of it is created in CreateDefaultExecutive() using NewInstance(). More...
 
- Static Public Member Functions inherited from vtkObject
static vtkObjectNew ()
 Create an object with Debug turned off, modified time initialized to zero, and reference counting on. More...
 
static void BreakOnError ()
 This method is called when vtkErrorMacro executes. More...
 
static void SetGlobalWarningDisplay (vtkTypeBool val)
 This is a global flag that controls whether any debug, warning or error messages are displayed. More...
 
static void GlobalWarningDisplayOn ()
 This is a global flag that controls whether any debug, warning or error messages are displayed. More...
 
static void GlobalWarningDisplayOff ()
 This is a global flag that controls whether any debug, warning or error messages are displayed. More...
 
static vtkTypeBool GetGlobalWarningDisplay ()
 This is a global flag that controls whether any debug, warning or error messages are displayed. More...
 
- Static Public Member Functions inherited from vtkObjectBase
static vtkTypeBool IsTypeOf (const char *name)
 Return 1 if this class type is the same type of (or a subclass of) the named class. More...
 
static vtkIdType GetNumberOfGenerationsFromBaseType (const char *name)
 Given a the name of a base class of this class type, return the distance of inheritance between this class type and the named class (how many generations of inheritance are there between this class and the named class). More...
 
static vtkObjectBaseNew ()
 Create an object with Debug turned off, modified time initialized to zero, and reference counting on. More...
 
static void SetMemkindDirectory (const char *directoryname)
 The name of a directory, ideally mounted -o dax, to memory map an extended memory space within. More...
 
static bool GetUsingMemkind ()
 A global state flag that controls whether vtkObjects are constructed in the usual way (the default) or within the extended memory space. More...
 

Protected Member Functions

virtual vtkObjectBaseNewInstanceInternal () const
 
 vtkDelaunay2D ()
 
int RequestData (vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
 This is called by the superclass. More...
 
- Protected Member Functions inherited from vtkPolyDataAlgorithm
 vtkPolyDataAlgorithm ()
 
 ~vtkPolyDataAlgorithm () override
 
virtual int RequestInformation (vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
 
virtual int RequestUpdateExtent (vtkInformation *, vtkInformationVector **, vtkInformationVector *)
 This is called by the superclass. More...
 
virtual int RequestUpdateTime (vtkInformation *, vtkInformationVector **, vtkInformationVector *)
 This is called by the superclass. More...
 
int FillOutputPortInformation (int port, vtkInformation *info) override
 Fill the output port information objects for this algorithm. More...
 
- Protected Member Functions inherited from vtkAlgorithm
 vtkAlgorithm ()
 
 ~vtkAlgorithm () override
 
bool CheckUpstreamAbort ()
 Checks to see if an upstream filter has been aborted. More...
 
virtual void SetNumberOfInputPorts (int n)
 Set the number of input ports used by the algorithm. More...
 
virtual void SetNumberOfOutputPorts (int n)
 Set the number of output ports provided by the algorithm. More...
 
int InputPortIndexInRange (int index, const char *action)
 
int OutputPortIndexInRange (int index, const char *action)
 
int GetInputArrayAssociation (int idx, vtkInformationVector **inputVector)
 Get the association of the actual data array for the input array specified by idx, this is only reasonable during the REQUEST_DATA pass. More...
 
vtkInformationGetInputArrayFieldInformation (int idx, vtkInformationVector **inputVector)
 This method takes in an index (as specified in SetInputArrayToProcess) and a pipeline information vector. More...
 
virtual vtkExecutiveCreateDefaultExecutive ()
 Create a default executive. More...
 
void ReportReferences (vtkGarbageCollector *) override
 
virtual void SetNthInputConnection (int port, int index, vtkAlgorithmOutput *input)
 Replace the Nth connection on the given input port. More...
 
virtual void SetNumberOfInputConnections (int port, int n)
 Set the number of input connections on the given input port. More...
 
void SetInputDataInternal (int port, vtkDataObject *input)
 These methods are used by subclasses to implement methods to set data objects directly as input. More...
 
void AddInputDataInternal (int port, vtkDataObject *input)
 
int GetInputArrayAssociation (int idx, int connection, vtkInformationVector **inputVector)
 Filters that have multiple connections on one port can use this signature. More...
 
int GetInputArrayAssociation (int idx, vtkDataObject *input)
 Filters that have multiple connections on one port can use this signature. More...
 
vtkDataArrayGetInputArrayToProcess (int idx, vtkInformationVector **inputVector)
 Get the actual data array for the input array specified by idx, this is only reasonable during the REQUEST_DATA pass. More...
 
vtkDataArrayGetInputArrayToProcess (int idx, vtkInformationVector **inputVector, int &association)
 Get the actual data array for the input array specified by idx, this is only reasonable during the REQUEST_DATA pass. More...
 
vtkDataArrayGetInputArrayToProcess (int idx, int connection, vtkInformationVector **inputVector)
 Filters that have multiple connections on one port can use this signature. More...
 
vtkDataArrayGetInputArrayToProcess (int idx, int connection, vtkInformationVector **inputVector, int &association)
 Filters that have multiple connections on one port can use this signature. More...
 
vtkDataArrayGetInputArrayToProcess (int idx, vtkDataObject *input)
 Filters that have multiple connections on one port can use this signature. More...
 
vtkDataArrayGetInputArrayToProcess (int idx, vtkDataObject *input, int &association)
 Filters that have multiple connections on one port can use this signature. More...
 
vtkAbstractArrayGetInputAbstractArrayToProcess (int idx, vtkInformationVector **inputVector)
 Get the actual data array for the input array specified by idx, this is only reasonable during the REQUEST_DATA pass. More...
 
vtkAbstractArrayGetInputAbstractArrayToProcess (int idx, vtkInformationVector **inputVector, int &association)
 Get the actual data array for the input array specified by idx, this is only reasonable during the REQUEST_DATA pass. More...
 
vtkAbstractArrayGetInputAbstractArrayToProcess (int idx, int connection, vtkInformationVector **inputVector)
 Filters that have multiple connections on one port can use this signature. More...
 
vtkAbstractArrayGetInputAbstractArrayToProcess (int idx, int connection, vtkInformationVector **inputVector, int &association)
 Filters that have multiple connections on one port can use this signature. More...
 
vtkAbstractArrayGetInputAbstractArrayToProcess (int idx, vtkDataObject *input)
 Filters that have multiple connections on one port can use this signature. More...
 
vtkAbstractArrayGetInputAbstractArrayToProcess (int idx, vtkDataObject *input, int &association)
 Filters that have multiple connections on one port can use this signature. More...
 
virtual void SetErrorCode (unsigned long)
 The error code contains a possible error that occurred while reading or writing the file. More...
 
- Protected Member Functions inherited from vtkObject
 vtkObject ()
 
 ~vtkObject () override
 
void RegisterInternal (vtkObjectBase *, vtkTypeBool check) override
 
void UnRegisterInternal (vtkObjectBase *, vtkTypeBool check) override
 
void InternalGrabFocus (vtkCommand *mouseEvents, vtkCommand *keypressEvents=nullptr)
 These methods allow a command to exclusively grab all events. More...
 
void InternalReleaseFocus ()
 These methods allow a command to exclusively grab all events. More...
 
- Protected Member Functions inherited from vtkObjectBase
 vtkObjectBase ()
 
virtual ~vtkObjectBase ()
 
 vtkObjectBase (const vtkObjectBase &)
 
void operator= (const vtkObjectBase &)
 

Protected Attributes

double Alpha
 
double Tolerance
 
vtkTypeBool BoundingTriangulation
 
double Offset
 
vtkTypeBool RandomPointInsertion
 
vtkSmartPointer< vtkAbstractTransformTransform
 
int ProjectionPlaneMode
 
- Protected Attributes inherited from vtkAlgorithm
vtkTimeStamp LastAbortCheckTime
 
vtkInformationInformation
 
double Progress
 
char * ProgressText
 
vtkProgressObserverProgressObserver
 
unsigned long ErrorCode
 The error code contains a possible error that occurred while reading or writing the file. More...
 
- Protected Attributes inherited from vtkObject
bool Debug
 
vtkTimeStamp MTime
 
vtkSubjectHelper * SubjectHelper
 
std::string ObjectName
 
- Protected Attributes inherited from vtkObjectBase
std::atomic< int32_t > ReferenceCount
 
vtkWeakPointerBase ** WeakPointers
 

Additional Inherited Members

- Public Attributes inherited from vtkAlgorithm
std::atomic< vtkTypeBoolAbortExecute
 
- Static Protected Member Functions inherited from vtkAlgorithm
static vtkInformationIntegerKeyPORT_REQUIREMENTS_FILLED ()
 
- Static Protected Member Functions inherited from vtkObjectBase
static vtkMallocingFunction GetCurrentMallocFunction ()
 
static vtkReallocingFunction GetCurrentReallocFunction ()
 
static vtkFreeingFunction GetCurrentFreeFunction ()
 
static vtkFreeingFunction GetAlternateFreeFunction ()
 
- Static Protected Attributes inherited from vtkAlgorithm
static vtkTimeStamp LastAbortTime
 
static vtkExecutiveDefaultExecutivePrototype
 

Detailed Description

create 2D Delaunay triangulation of input points

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 laying 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 can use the vtkTransformFilter to transform the points into and out of the x-y plane or you can specify a transform to the Delaunay2D directly. In the latter case, the input points are transformed, the transformed points are triangulated, and the output will use the triangulated topology for the original (non-transformed) points. This avoids transforming the data back as would be required when using the vtkTransformFilter method. Specifying a transform directly also allows any transform to be used: rigid, non-rigid, non-invertible, etc.

If an input transform is used, then alpha values are applied (for the most part) in the original data space. The exception is when BoundingTriangulation is on. In this case, alpha values are applied in the original data space unless a cell uses a bounding vertex.

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. Note also that the RandomPointInsertion mode can be set which will insert the points in pseudo-random order.

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. Only the topology of (lines and polygons) from this second input are used. The topology is assumed to reference points in the input point set (the one to be triangulated). In other words, the lines and polygons use point ids from the first input point set. 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.)

If an input transform is used, constraints are defined in the "transformed" space. So when the right hand rule is used for a polygon constraint, that operation is applied using the transformed points. Since the input transform can be any transformation (rigid or non-rigid), care must be taken in constructing constraints when an input transform is used.

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.
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.
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
Online Examples:

Tests:
vtkDelaunay2D (Tests)

Definition at line 243 of file vtkDelaunay2D.h.

Member Typedef Documentation

◆ Superclass

Definition at line 246 of file vtkDelaunay2D.h.

Constructor & Destructor Documentation

◆ vtkDelaunay2D()

vtkDelaunay2D::vtkDelaunay2D ( )
protected

Member Function Documentation

◆ IsTypeOf()

static vtkTypeBool vtkDelaunay2D::IsTypeOf ( const char *  type)
static

◆ IsA()

virtual vtkTypeBool vtkDelaunay2D::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 vtkPolyDataAlgorithm.

◆ SafeDownCast()

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

◆ NewInstanceInternal()

virtual vtkObjectBase* vtkDelaunay2D::NewInstanceInternal ( ) const
protectedvirtual

Reimplemented from vtkPolyDataAlgorithm.

◆ NewInstance()

vtkDelaunay2D* vtkDelaunay2D::NewInstance ( ) const

◆ PrintSelf()

void vtkDelaunay2D::PrintSelf ( ostream &  os,
vtkIndent  indent 
)
overridevirtual

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 vtkPolyDataAlgorithm.

◆ New()

static vtkDelaunay2D* vtkDelaunay2D::New ( )
static

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

◆ SetSourceData()

void vtkDelaunay2D::SetSourceData ( 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. The lines/polygons are assumed to reference points in the input point set (i.e. point ids are identical in the input and source). Note that this method does not connect the pipeline. See SetSourceConnection for connecting the pipeline.

◆ SetSourceConnection()

void vtkDelaunay2D::SetSourceConnection ( vtkAlgorithmOutput algOutput)

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. The lines/polygons are assumed to reference points in the input point set (i.e. point ids are identical in the input and source). New style. This method is equivalent to SetInputConnection(1, algOutput).

◆ GetSource()

vtkPolyData* vtkDelaunay2D::GetSource ( )

Get a pointer to the source object.

◆ SetAlpha()

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.

◆ GetAlpha()

virtual double vtkDelaunay2D::GetAlpha ( )
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.

◆ SetTolerance()

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.

◆ GetTolerance()

virtual double vtkDelaunay2D::GetTolerance ( )
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.

◆ SetOffset()

virtual void vtkDelaunay2D::SetOffset ( double  )
virtual

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

◆ GetOffset()

virtual double vtkDelaunay2D::GetOffset ( )
virtual

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

◆ SetBoundingTriangulation()

virtual void vtkDelaunay2D::SetBoundingTriangulation ( vtkTypeBool  )
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.)

◆ GetBoundingTriangulation()

virtual vtkTypeBool vtkDelaunay2D::GetBoundingTriangulation ( )
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.)

◆ BoundingTriangulationOn()

virtual void vtkDelaunay2D::BoundingTriangulationOn ( )
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.)

◆ BoundingTriangulationOff()

virtual void vtkDelaunay2D::BoundingTriangulationOff ( )
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.)

◆ vtkSetSmartPointerMacro()

vtkDelaunay2D::vtkSetSmartPointerMacro ( Transform  ,
vtkAbstractTransform   
)

Set / get the transform which is applied to points to generate a 2D problem.

This maps a 3D dataset into a 2D dataset where triangulation can be done on the XY plane. The points are transformed and triangulated. The topology of triangulated points is used as the output topology. The output points are the original (untransformed) points. The transform can be any subclass of vtkAbstractTransform (thus it does not need to be a linear or invertible transform).

◆ vtkGetSmartPointerMacro()

vtkDelaunay2D::vtkGetSmartPointerMacro ( Transform  ,
vtkAbstractTransform   
)

Set / get the transform which is applied to points to generate a 2D problem.

This maps a 3D dataset into a 2D dataset where triangulation can be done on the XY plane. The points are transformed and triangulated. The topology of triangulated points is used as the output topology. The output points are the original (untransformed) points. The transform can be any subclass of vtkAbstractTransform (thus it does not need to be a linear or invertible transform).

◆ SetProjectionPlaneMode()

virtual void vtkDelaunay2D::SetProjectionPlaneMode ( int  )
virtual

Define the method to project the input 3D points into a 2D plane for triangulation.

When the VTK_DELAUNAY_XY_PLANE is set, the z-coordinate is simply ignored. When VTK_SET_TRANSFORM_PLANE is set, then a transform must be supplied and the points are transformed using it. Finally, if VTK_BEST_FITTING_PLANE is set, then the filter computes a best fitting plane and projects the points onto it.

◆ GetProjectionPlaneMode()

virtual int vtkDelaunay2D::GetProjectionPlaneMode ( )
virtual

Define the method to project the input 3D points into a 2D plane for triangulation.

When the VTK_DELAUNAY_XY_PLANE is set, the z-coordinate is simply ignored. When VTK_SET_TRANSFORM_PLANE is set, then a transform must be supplied and the points are transformed using it. Finally, if VTK_BEST_FITTING_PLANE is set, then the filter computes a best fitting plane and projects the points onto it.

◆ ComputeBestFittingPlane()

static vtkAbstractTransform* vtkDelaunay2D::ComputeBestFittingPlane ( vtkPointSet input)
static

This method computes the best fit plane to a set of points represented by a vtkPointSet.

The method constructs a transform and returns it on successful completion (null otherwise). The user is responsible for deleting the transform instance.

◆ SetRandomPointInsertion()

virtual void vtkDelaunay2D::SetRandomPointInsertion ( vtkTypeBool  )
virtual

Indicate whether to insert the points in given order, or pseudo-random order.

Inserting in random order can improve performance and numerics in many circumstances.

◆ GetRandomPointInsertion()

virtual vtkTypeBool vtkDelaunay2D::GetRandomPointInsertion ( )
virtual

Indicate whether to insert the points in given order, or pseudo-random order.

Inserting in random order can improve performance and numerics in many circumstances.

◆ RandomPointInsertionOn()

virtual void vtkDelaunay2D::RandomPointInsertionOn ( )
virtual

Indicate whether to insert the points in given order, or pseudo-random order.

Inserting in random order can improve performance and numerics in many circumstances.

◆ RandomPointInsertionOff()

virtual void vtkDelaunay2D::RandomPointInsertionOff ( )
virtual

Indicate whether to insert the points in given order, or pseudo-random order.

Inserting in random order can improve performance and numerics in many circumstances.

◆ RequestData()

int vtkDelaunay2D::RequestData ( vtkInformation request,
vtkInformationVector **  inputVector,
vtkInformationVector outputVector 
)
overrideprotectedvirtual

This is called by the superclass.

This is the method you should override.

Reimplemented from vtkPolyDataAlgorithm.

Member Data Documentation

◆ Alpha

double vtkDelaunay2D::Alpha
protected

Definition at line 375 of file vtkDelaunay2D.h.

◆ Tolerance

double vtkDelaunay2D::Tolerance
protected

Definition at line 376 of file vtkDelaunay2D.h.

◆ BoundingTriangulation

vtkTypeBool vtkDelaunay2D::BoundingTriangulation
protected

Definition at line 377 of file vtkDelaunay2D.h.

◆ Offset

double vtkDelaunay2D::Offset
protected

Definition at line 378 of file vtkDelaunay2D.h.

◆ RandomPointInsertion

vtkTypeBool vtkDelaunay2D::RandomPointInsertion
protected

Definition at line 379 of file vtkDelaunay2D.h.

◆ Transform

vtkSmartPointer<vtkAbstractTransform> vtkDelaunay2D::Transform
protected

Definition at line 382 of file vtkDelaunay2D.h.

◆ ProjectionPlaneMode

int vtkDelaunay2D::ProjectionPlaneMode
protected

Definition at line 384 of file vtkDelaunay2D.h.


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