VTK  9.4.20241218
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
vtkGeometryFilter Class Reference

extract boundary geometry from dataset (or convert data to polygonal type) More...

#include <vtkGeometryFilter.h>

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

Public Member Functions

void SetExtent (double xMin, double xMax, double yMin, double yMax, double zMin, double zMax)
 Specify a (xmin,xmax, ymin,ymax, zmin,zmax) bounding box to clip data.
 
virtual void SetPointClipping (bool)
 Turn on/off selection of geometry by point id.
 
virtual bool GetPointClipping ()
 Turn on/off selection of geometry by point id.
 
virtual void PointClippingOn ()
 Turn on/off selection of geometry by point id.
 
virtual void PointClippingOff ()
 Turn on/off selection of geometry by point id.
 
virtual void SetCellClipping (bool)
 Turn on/off selection of geometry by cell id.
 
virtual bool GetCellClipping ()
 Turn on/off selection of geometry by cell id.
 
virtual void CellClippingOn ()
 Turn on/off selection of geometry by cell id.
 
virtual void CellClippingOff ()
 Turn on/off selection of geometry by cell id.
 
virtual void SetExtentClipping (bool)
 Turn on/off selection of geometry via bounding box.
 
virtual bool GetExtentClipping ()
 Turn on/off selection of geometry via bounding box.
 
virtual void ExtentClippingOn ()
 Turn on/off selection of geometry via bounding box.
 
virtual void ExtentClippingOff ()
 Turn on/off selection of geometry via bounding box.
 
virtual void SetPointMinimum (vtkIdType)
 Specify the minimum point id for point id selection.
 
virtual vtkIdType GetPointMinimum ()
 Specify the minimum point id for point id selection.
 
virtual void SetPointMaximum (vtkIdType)
 Specify the maximum point id for point id selection.
 
virtual vtkIdType GetPointMaximum ()
 Specify the maximum point id for point id selection.
 
virtual void SetCellMinimum (vtkIdType)
 Specify the minimum cell id for point id selection.
 
virtual vtkIdType GetCellMinimum ()
 Specify the minimum cell id for point id selection.
 
virtual void SetCellMaximum (vtkIdType)
 Specify the maximum cell id for point id selection.
 
virtual vtkIdType GetCellMaximum ()
 Specify the maximum cell id for point id selection.
 
void SetExtent (double extent[6])
 Set / get a (xmin,xmax, ymin,ymax, zmin,zmax) bounding box to clip data.
 
double * GetExtent ()
 Set / get a (xmin,xmax, ymin,ymax, zmin,zmax) bounding box to clip data.
 
virtual void SetMerging (bool)
 Turn on/off merging of points.
 
virtual bool GetMerging ()
 Turn on/off merging of points.
 
virtual void MergingOn ()
 Turn on/off merging of points.
 
virtual void MergingOff ()
 Turn on/off merging of points.
 
void SetOutputPointsPrecision (int precision)
 Set/get the desired precision for the output types.
 
int GetOutputPointsPrecision () const
 Set/get the desired precision for the output types.
 
virtual void SetFastMode (bool)
 Turn on/off fast mode execution.
 
virtual bool GetFastMode ()
 Turn on/off fast mode execution.
 
virtual void FastModeOn ()
 Turn on/off fast mode execution.
 
virtual void FastModeOff ()
 Turn on/off fast mode execution.
 
virtual void SetPieceInvariant (int)
 If PieceInvariant is true, vtkGeometryFilter requests 1 ghost level from input in order to remove internal surface that are between processes.
 
virtual int GetPieceInvariant ()
 If PieceInvariant is true, vtkGeometryFilter requests 1 ghost level from input in order to remove internal surface that are between processes.
 
virtual void SetPassThroughCellIds (vtkTypeBool)
 This parameter drives the generation or not of a CellData array for the output polygonal dataset that holds the cell index of the original 3D cell that produced each output cell.
 
virtual vtkTypeBool GetPassThroughCellIds ()
 This parameter drives the generation or not of a CellData array for the output polygonal dataset that holds the cell index of the original 3D cell that produced each output cell.
 
virtual void PassThroughCellIdsOn ()
 This parameter drives the generation or not of a CellData array for the output polygonal dataset that holds the cell index of the original 3D cell that produced each output cell.
 
virtual void PassThroughCellIdsOff ()
 This parameter drives the generation or not of a CellData array for the output polygonal dataset that holds the cell index of the original 3D cell that produced each output cell.
 
virtual void SetPassThroughPointIds (vtkTypeBool)
 This parameter drives the generation or not of a PointData array for the output polygonal dataset that holds the cell/point index of the original point that produced each output point.
 
virtual vtkTypeBool GetPassThroughPointIds ()
 This parameter drives the generation or not of a PointData array for the output polygonal dataset that holds the cell/point index of the original point that produced each output point.
 
virtual void PassThroughPointIdsOn ()
 This parameter drives the generation or not of a PointData array for the output polygonal dataset that holds the cell/point index of the original point that produced each output point.
 
virtual void PassThroughPointIdsOff ()
 This parameter drives the generation or not of a PointData array for the output polygonal dataset that holds the cell/point index of the original point that produced each output point.
 
virtual void SetOriginalCellIdsName (const char *)
 If PassThroughCellIds or PassThroughPointIds is on, then these ivars control the name given to the field in which the ids are written into.
 
virtual const char * GetOriginalCellIdsName ()
 If PassThroughCellIds or PassThroughPointIds is on, then these ivars control the name given to the field in which the ids are written into.
 
virtual void SetOriginalPointIdsName (const char *)
 If PassThroughCellIds or PassThroughPointIds is on, then these ivars control the name given to the field in which the ids are written into.
 
virtual const char * GetOriginalPointIdsName ()
 If PassThroughCellIds or PassThroughPointIds is on, then these ivars control the name given to the field in which the ids are written into.
 
void SetExcludedFacesData (vtkPolyData *)
 If a second, vtkPolyData input is provided, this second input specifies a list of faces to be excluded from the output (in the vtkPolyData::Polys attribute).
 
void SetExcludedFacesConnection (vtkAlgorithmOutput *algOutput)
 If a second, vtkPolyData input is provided, this second input specifies a list of faces to be excluded from the output (in the vtkPolyData::Polys attribute).
 
vtkPolyDataGetExcludedFaces ()
 If a second, vtkPolyData input is provided, this second input specifies a list of faces to be excluded from the output (in the vtkPolyData::Polys attribute).
 
virtual void SetNonlinearSubdivisionLevel (int)
 If the input is an unstructured grid with nonlinear faces, this parameter determines how many times the face is subdivided into linear faces.
 
virtual int GetNonlinearSubdivisionLevel ()
 If the input is an unstructured grid with nonlinear faces, this parameter determines how many times the face is subdivided into linear faces.
 
virtual void SetMatchBoundariesIgnoringCellOrder (int)
 When two volumetric cells of different order are connected by their corners (for instance, a quadratic hexahedron next to a linear hexahedron ), the internal face is rendered and is not considered as a ghost cell.
 
virtual int GetMatchBoundariesIgnoringCellOrder ()
 When two volumetric cells of different order are connected by their corners (for instance, a quadratic hexahedron next to a linear hexahedron ), the internal face is rendered and is not considered as a ghost cell.
 
virtual void SetDelegation (vtkTypeBool)
 Disable delegation to an internal vtkDataSetSurfaceFilter.
 
virtual vtkTypeBool GetDelegation ()
 Disable delegation to an internal vtkDataSetSurfaceFilter.
 
virtual void DelegationOn ()
 Disable delegation to an internal vtkDataSetSurfaceFilter.
 
virtual void DelegationOff ()
 Disable delegation to an internal vtkDataSetSurfaceFilter.
 
virtual void SetRemoveGhostInterfaces (bool)
 Set/Get if Ghost interfaces will be removed.
 
virtual void RemoveGhostInterfacesOn ()
 Set/Get if Ghost interfaces will be removed.
 
virtual void RemoveGhostInterfacesOff ()
 Set/Get if Ghost interfaces will be removed.
 
virtual bool GetRemoveGhostInterfaces ()
 Set/Get if Ghost interfaces will be removed.
 
int PolyDataExecute (vtkDataSet *input, vtkPolyData *output, vtkPolyData *exc)
 Direct access methods so that this class can be used as an algorithm without using it as a filter (i.e., no pipeline updates).
 
virtual int PolyDataExecute (vtkDataSet *, vtkPolyData *)
 Direct access methods so that this class can be used as an algorithm without using it as a filter (i.e., no pipeline updates).
 
int UnstructuredGridExecute (vtkDataSet *input, vtkPolyData *output, vtkGeometryFilterHelper *info, vtkPolyData *exc)
 Direct access methods so that this class can be used as an algorithm without using it as a filter (i.e., no pipeline updates).
 
virtual int UnstructuredGridExecute (vtkDataSet *input, vtkPolyData *output)
 Direct access methods so that this class can be used as an algorithm without using it as a filter (i.e., no pipeline updates).
 
int StructuredExecute (vtkDataSet *input, vtkPolyData *output, int *wholeExtent, vtkPolyData *exc, bool *extractFace=nullptr)
 Direct access methods so that this class can be used as an algorithm without using it as a filter (i.e., no pipeline updates).
 
virtual int StructuredExecute (vtkDataSet *input, vtkPolyData *output, int *wholeExt, bool *extractFace=nullptr)
 Direct access methods so that this class can be used as an algorithm without using it as a filter (i.e., no pipeline updates).
 
int DataSetExecute (vtkDataSet *input, vtkPolyData *output, vtkPolyData *exc)
 Direct access methods so that this class can be used as an algorithm without using it as a filter (i.e., no pipeline updates).
 
virtual int DataSetExecute (vtkDataSet *input, vtkPolyData *output)
 Direct access methods so that this class can be used as an algorithm without using it as a filter (i.e., no pipeline updates).
 
- Public Member Functions inherited from vtkPolyDataAlgorithm
virtual vtkTypeBool IsA (const char *type)
 Return 1 if this class is the same type of (or a subclass of) the named class.
 
vtkPolyDataAlgorithmNewInstance () const
 
void PrintSelf (ostream &os, vtkIndent indent) override
 Methods invoked by print to print information about the object including superclasses.
 
vtkTypeBool ProcessRequest (vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
 see vtkAlgorithm for details
 
vtkDataObjectGetInput ()
 
vtkDataObjectGetInput (int port)
 
vtkPolyDataGetPolyDataInput (int port)
 
vtkPolyDataGetOutput ()
 Get the output data object for a port on this algorithm.
 
vtkPolyDataGetOutput (int)
 Get the output data object for a port on this algorithm.
 
virtual void SetOutput (vtkDataObject *d)
 Get the output data object for a port on this algorithm.
 
void SetInputData (vtkDataObject *)
 Assign a data object as input.
 
void SetInputData (int, vtkDataObject *)
 Assign a data object as input.
 
void AddInputData (vtkDataObject *)
 Assign a data object as input.
 
void AddInputData (int, vtkDataObject *)
 Assign a data object as input.
 
- Public Member Functions inherited from vtkAlgorithm
virtual vtkTypeBool IsA (const char *type)
 Return 1 if this class is the same type of (or a subclass of) the named class.
 
vtkAlgorithmNewInstance () const
 
void PrintSelf (ostream &os, vtkIndent indent) override
 Methods invoked by print to print information about the object including superclasses.
 
vtkTypeBool HasExecutive ()
 Check whether this algorithm has an assigned executive.
 
vtkExecutiveGetExecutive ()
 Get this algorithm's executive.
 
virtual void SetExecutive (vtkExecutive *executive)
 Set this algorithm's executive.
 
virtual vtkTypeBool ProcessRequest (vtkInformation *request, vtkInformationVector **inInfo, vtkInformationVector *outInfo)
 Upstream/Downstream requests form the generalized interface through which executives invoke a algorithm's functionality.
 
vtkTypeBool ProcessRequest (vtkInformation *request, vtkCollection *inInfo, vtkInformationVector *outInfo)
 Version of ProcessRequest() that is wrapped.
 
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.
 
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.
 
vtkInformationGetInputPortInformation (int port)
 Get the information object associated with an input port.
 
vtkInformationGetOutputPortInformation (int port)
 Get the information object associated with an output port.
 
int GetNumberOfInputPorts ()
 Get the number of input ports used by the algorithm.
 
int GetNumberOfOutputPorts ()
 Get the number of output ports provided by the algorithm.
 
void SetAbortExecuteAndUpdateTime ()
 Set AbortExecute Flag and update LastAbortTime.
 
void UpdateProgress (double amount)
 Update the progress of the process object.
 
bool CheckAbort ()
 Checks to see if this filter should abort.
 
vtkInformationGetInputArrayInformation (int idx)
 Get the info object for the specified input array to this algorithm.
 
void RemoveAllInputs ()
 Remove all the input data.
 
vtkDataObjectGetOutputDataObject (int port)
 Get the data object that will contain the algorithm output for the given port.
 
vtkDataObjectGetInputDataObject (int port, int connection)
 Get the data object that will contain the algorithm input for the given port and given connection.
 
virtual void RemoveInputConnection (int port, vtkAlgorithmOutput *input)
 Remove a connection from the given input port index.
 
virtual void RemoveInputConnection (int port, int idx)
 Remove a connection given by index idx.
 
virtual void RemoveAllInputConnections (int port)
 Removes all input connections.
 
virtual void SetInputDataObject (int port, vtkDataObject *data)
 Sets the data-object as an input on the given port index.
 
virtual void SetInputDataObject (vtkDataObject *data)
 
virtual void AddInputDataObject (int port, vtkDataObject *data)
 Add the data-object as an input to this given port.
 
virtual void AddInputDataObject (vtkDataObject *data)
 
vtkAlgorithmOutputGetOutputPort (int index)
 Get a proxy object corresponding to the given output port of this algorithm.
 
vtkAlgorithmOutputGetOutputPort ()
 
int GetNumberOfInputConnections (int port)
 Get the number of inputs currently connected to a port.
 
int GetTotalNumberOfInputConnections ()
 Get the total number of inputs for this algorithm.
 
vtkAlgorithmOutputGetInputConnection (int port, int index)
 Get the algorithm output port connected to an input port.
 
vtkAlgorithmGetInputAlgorithm (int port, int index, int &algPort)
 Returns the algorithm and the output port index of that algorithm connected to a port-index pair.
 
vtkAlgorithmGetInputAlgorithm (int port, int index)
 Returns the algorithm connected to a port-index pair.
 
vtkAlgorithmGetInputAlgorithm ()
 Equivalent to GetInputAlgorithm(0, 0).
 
vtkExecutiveGetInputExecutive (int port, int index)
 Returns the executive associated with a particular input connection.
 
vtkExecutiveGetInputExecutive ()
 Equivalent to GetInputExecutive(0, 0)
 
vtkInformationGetInputInformation (int port, int index)
 Return the information object that is associated with a particular input connection.
 
vtkInformationGetInputInformation ()
 Equivalent to GetInputInformation(0, 0)
 
vtkInformationGetOutputInformation (int port)
 Return the information object that is associated with a particular output port.
 
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).
 
virtual vtkTypeBool Update (vtkInformation *requests)
 Convenience method to update an algorithm after passing requests to its first output port.
 
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.
 
virtual VTK_UNBLOCKTHREADS int UpdateExtent (const int extents[6])
 Convenience method to update an algorithm after passing requests to its first output port.
 
virtual VTK_UNBLOCKTHREADS 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.
 
virtual VTK_UNBLOCKTHREADS void UpdateInformation ()
 Bring the algorithm's information up-to-date.
 
virtual void UpdateDataObject ()
 Create output object(s).
 
virtual void PropagateUpdateExtent ()
 Propagate meta-data upstream.
 
virtual VTK_UNBLOCKTHREADS void UpdateWholeExtent ()
 Bring this algorithm's outputs up-to-date.
 
void ConvertTotalInputToPortConnection (int ind, int &port, int &conn)
 Convenience routine to convert from a linear ordering of input connections to a port/connection pair.
 
void RemoveNoPriorTemporalAccessInformationKey ()
 Removes any information key vtkStreamingDemandDrivenPipeline::NO_PRIOR_TEMPORAL_ACCESS() to all output ports of this vtkAlgorithm.
 
virtual vtkInformationGetInformation ()
 Set/Get the information object associated with this algorithm.
 
virtual void SetInformation (vtkInformation *)
 Set/Get the information object associated with this algorithm.
 
bool UsesGarbageCollector () const override
 Participate in garbage collection.
 
virtual void SetAbortExecute (vtkTypeBool)
 Set/Get the AbortExecute flag for the process object.
 
virtual vtkTypeBool GetAbortExecute ()
 Set/Get the AbortExecute flag for the process object.
 
virtual void AbortExecuteOn ()
 Set/Get the AbortExecute flag for the process object.
 
virtual void AbortExecuteOff ()
 Set/Get the AbortExecute flag for the process object.
 
virtual double GetProgress ()
 Get the execution progress of a process object.
 
void SetContainerAlgorithm (vtkAlgorithm *containerAlg)
 Set/get a Container algorithm for this algorithm.
 
vtkAlgorithmGetContainerAlgorithm ()
 Set/get a Container algorithm for this algorithm.
 
virtual void SetAbortOutput (bool)
 Set/Get an internal variable used to communicate between the algorithm and executive.
 
virtual bool GetAbortOutput ()
 Set/Get an internal variable used to communicate between the algorithm and executive.
 
void SetProgressShiftScale (double shift, double scale)
 Specify the shift and scale values to use to apply to the progress amount when UpdateProgress is called.
 
virtual double GetProgressShift ()
 Specify the shift and scale values to use to apply to the progress amount when UpdateProgress is called.
 
virtual double GetProgressScale ()
 Specify the shift and scale values to use to apply to the progress amount when UpdateProgress is called.
 
void SetProgressText (const char *ptext)
 Set the current text message associated with the progress state.
 
virtual char * GetProgressText ()
 Set the current text message associated with the progress state.
 
virtual unsigned long GetErrorCode ()
 The error code contains a possible error that occurred while reading or writing the file.
 
void SetInputArrayToProcess (const char *name, int fieldAssociation)
 Set the input data arrays that this algorithm will process.
 
virtual void SetInputArrayToProcess (int idx, int port, int connection, int fieldAssociation, const char *name)
 Set the input data arrays that this algorithm will process.
 
virtual void SetInputArrayToProcess (int idx, int port, int connection, int fieldAssociation, int fieldAttributeType)
 Set the input data arrays that this algorithm will process.
 
virtual void SetInputArrayToProcess (int idx, vtkInformation *info)
 Set the input data arrays that this algorithm will process.
 
virtual void SetInputArrayToProcess (int idx, int port, int connection, const char *fieldAssociation, const char *attributeTypeorName)
 Set the input data arrays that this algorithm will process.
 
virtual void SetInputConnection (int port, vtkAlgorithmOutput *input)
 Set the connection for the given input port index.
 
virtual void SetInputConnection (vtkAlgorithmOutput *input)
 Set the connection for the given input port index.
 
virtual void AddInputConnection (int port, vtkAlgorithmOutput *input)
 Add a connection to the given input port index.
 
virtual void AddInputConnection (vtkAlgorithmOutput *input)
 Add a connection to the given input port index.
 
virtual VTK_UNBLOCKTHREADS void Update (int port)
 Bring this algorithm's outputs up-to-date.
 
virtual VTK_UNBLOCKTHREADS void Update ()
 Bring this algorithm's outputs up-to-date.
 
virtual void SetReleaseDataFlag (vtkTypeBool)
 Turn release data flag on or off for all output ports.
 
virtual vtkTypeBool GetReleaseDataFlag ()
 Turn release data flag on or off for all output ports.
 
void ReleaseDataFlagOn ()
 Turn release data flag on or off for all output ports.
 
void ReleaseDataFlagOff ()
 Turn release data flag on or off for all output ports.
 
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.
 
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.
 
int * GetUpdateExtent ()
 These functions return the update extent for output ports that use 3D extents.
 
int * GetUpdateExtent (int port)
 These functions return the update extent for output ports that use 3D extents.
 
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.
 
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.
 
void GetUpdateExtent (int extent[6])
 These functions return the update extent for output ports that use 3D extents.
 
void GetUpdateExtent (int port, int extent[6])
 These functions return the update extent for output ports that use 3D extents.
 
int GetUpdatePiece ()
 These functions return the update extent for output ports that use piece extents.
 
int GetUpdatePiece (int port)
 These functions return the update extent for output ports that use piece extents.
 
int GetUpdateNumberOfPieces ()
 These functions return the update extent for output ports that use piece extents.
 
int GetUpdateNumberOfPieces (int port)
 These functions return the update extent for output ports that use piece extents.
 
int GetUpdateGhostLevel ()
 These functions return the update extent for output ports that use piece extents.
 
int GetUpdateGhostLevel (int port)
 These functions return the update extent for output ports that use piece extents.
 
void SetProgressObserver (vtkProgressObserver *)
 If an ProgressObserver is set, the algorithm will report progress through it rather than directly.
 
virtual vtkProgressObserverGetProgressObserver ()
 If an ProgressObserver is set, the algorithm will report progress through it rather than directly.
 
void SetNoPriorTemporalAccessInformationKey (int key)
 Set to all output ports of this algorithm the information key vtkStreamingDemandDrivenPipeline::NO_PRIOR_TEMPORAL_ACCESS().
 
void SetNoPriorTemporalAccessInformationKey ()
 Set to all output ports of this algorithm the information key vtkStreamingDemandDrivenPipeline::NO_PRIOR_TEMPORAL_ACCESS().
 
- Public Member Functions inherited from vtkObject
 vtkBaseTypeMacro (vtkObject, vtkObjectBase)
 
virtual void DebugOn ()
 Turn debugging output on.
 
virtual void DebugOff ()
 Turn debugging output off.
 
bool GetDebug ()
 Get the value of the debug flag.
 
void SetDebug (bool debugFlag)
 Set the value of the debug flag.
 
virtual void Modified ()
 Update the modification time for this object.
 
virtual vtkMTimeType GetMTime ()
 Return this object's modified time.
 
void PrintSelf (ostream &os, vtkIndent indent) override
 Methods invoked by print to print information about the object including superclasses.
 
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.
 
unsigned long AddObserver (unsigned long event, vtkCommand *, float priority=0.0f)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object.
 
unsigned long AddObserver (const char *event, vtkCommand *, float priority=0.0f)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object.
 
vtkCommandGetCommand (unsigned long tag)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object.
 
void RemoveObserver (vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object.
 
void RemoveObservers (unsigned long event, vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object.
 
void RemoveObservers (const char *event, vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object.
 
vtkTypeBool HasObserver (unsigned long event, vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object.
 
vtkTypeBool HasObserver (const char *event, vtkCommand *)
 Allow people to add/remove/invoke observers (callbacks) to any VTK object.
 
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.
 
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.
 
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.
 
vtkTypeBool InvokeEvent (unsigned long event, void *callData)
 This method invokes an event and return whether the event was aborted or not.
 
vtkTypeBool InvokeEvent (const char *event, void *callData)
 This method invokes an event and return whether the event was aborted or not.
 
virtual void SetObjectName (const std::string &objectName)
 Set/get the name of this object for reporting purposes.
 
virtual std::string GetObjectName () const
 Set/get the name of this object for reporting purposes.
 
- Public Member Functions inherited from vtkObjectBase
const char * GetClassName () const
 Return the class name as a string.
 
virtual std::string GetObjectDescription () const
 The object description printed in messages and PrintSelf output.
 
virtual vtkTypeBool IsA (const char *name)
 Return 1 if this class is the same type of (or a subclass of) the named class.
 
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).
 
virtual void Delete ()
 Delete a VTK object.
 
virtual void FastDelete ()
 Delete a reference to this object.
 
void InitializeObjectBase ()
 
void Print (ostream &os)
 Print an object to an ostream.
 
void Register (vtkObjectBase *o)
 Increase the reference count (mark as used by another object).
 
virtual void UnRegister (vtkObjectBase *o)
 Decrease the reference count (release by another object).
 
int GetReferenceCount ()
 Return the current reference count of this object.
 
void SetReferenceCount (int)
 Sets the reference count.
 
bool GetIsInMemkind () const
 A local state flag that remembers whether this object lives in the normal or extended memory space.
 
virtual void PrintHeader (ostream &os, vtkIndent indent)
 Methods invoked by print to print information about the object including superclasses.
 
virtual void PrintTrailer (ostream &os, vtkIndent indent)
 Methods invoked by print to print information about the object including superclasses.
 

Protected Member Functions

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

Protected Attributes

vtkIdType PointMaximum
 
vtkIdType PointMinimum
 
vtkIdType CellMinimum
 
vtkIdType CellMaximum
 
double Extent [6]
 
bool PointClipping
 
bool CellClipping
 
bool ExtentClipping
 
int OutputPointsPrecision
 
bool RemoveGhostInterfaces
 
bool Merging
 
vtkIncrementalPointLocatorLocator
 
bool FastMode
 
int PieceInvariant
 
vtkTypeBool PassThroughCellIds
 
char * OriginalCellIdsName
 
vtkTypeBool PassThroughPointIds
 
char * OriginalPointIdsName
 
int NonlinearSubdivisionLevel
 
int MatchBoundariesIgnoringCellOrder
 
vtkTypeBool Delegation
 
- 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.
 
- 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
 
typedef vtkPolyDataAlgorithm Superclass
 Standard methods for instantiation, type information, and printing.
 
static vtkGeometryFilterNew ()
 Standard methods for instantiation, type information, and printing.
 
static vtkTypeBool IsTypeOf (const char *type)
 Standard methods for instantiation, type information, and printing.
 
static vtkGeometryFilterSafeDownCast (vtkObjectBase *o)
 Standard methods for instantiation, type information, and printing.
 
virtual vtkTypeBool IsA (const char *type)
 Standard methods for instantiation, type information, and printing.
 
vtkGeometryFilterNewInstance () const
 Standard methods for instantiation, type information, and printing.
 
void PrintSelf (ostream &os, vtkIndent indent) override
 Standard methods for instantiation, type information, and printing.
 
virtual vtkObjectBaseNewInstanceInternal () const
 Standard methods for instantiation, type information, and printing.
 

Additional Inherited Members

- 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
 
- 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.
 
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.
 
static vtkInformationIntegerKeyCAN_HANDLE_PIECE_REQUEST ()
 Key that tells the pipeline that a particular algorithm can or cannot handle piece request.
 
static vtkInformationIntegerKeyABORTED ()
 
static void SetDefaultExecutivePrototype (vtkExecutive *proto)
 If the DefaultExecutivePrototype is set, a copy of it is created in CreateDefaultExecutive() using NewInstance().
 
- 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.
 
static void BreakOnError ()
 This method is called when vtkErrorMacro executes.
 
static void SetGlobalWarningDisplay (vtkTypeBool val)
 This is a global flag that controls whether any debug, warning or error messages are displayed.
 
static void GlobalWarningDisplayOn ()
 This is a global flag that controls whether any debug, warning or error messages are displayed.
 
static void GlobalWarningDisplayOff ()
 This is a global flag that controls whether any debug, warning or error messages are displayed.
 
static vtkTypeBool GetGlobalWarningDisplay ()
 This is a global flag that controls whether any debug, warning or error messages are displayed.
 
- 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.
 
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).
 
static vtkObjectBaseNew ()
 Create an object with Debug turned off, modified time initialized to zero, and reference counting on.
 
static void SetMemkindDirectory (const char *directoryname)
 The name of a directory, ideally mounted -o dax, to memory map an extended memory space within.
 
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.
 
- 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

extract boundary geometry from dataset (or convert data to polygonal type)

vtkGeometryFilter is a general-purpose filter to extract dataset boundary geometry, topology, and associated attribute data from any type of dataset. Geometry is obtained as follows: all 0D, 1D, and 2D cells are extracted. All 2D faces that are used by only one 3D cell (i.e., boundary faces) are extracted. It also is possible to specify conditions on point ids, cell ids, and on a bounding box (referred to as "Extent") to control the extraction process. This point and cell id- and extent-based clipping is a powerful way to "see inside" datasets; however it may impact performance significantly.

This filter may also be used to convert any type of data to polygonal type. This is particularly useful for surface rendering. The conversion process may be less than satisfactory for some 3D datasets. For example, this filter will extract the outer surface of a volume or structured grid dataset (if point, cell, and extent clipping is disabled). (For structured data you may want to use vtkImageDataGeometryFilter, vtkStructuredGridGeometryFilter, vtkExtractUnstructuredGrid, vtkRectilinearGridGeometryFilter, or vtkExtractVOI.)

Another important feature of vtkGeometryFilter is that it preserves topological connectivity. This enables filters that depend on correct connectivity (e.g., vtkQuadricDecimation, vtkFeatureEdges, etc.) to operate properly . It is possible to label the output polydata with an originating cell (PassThroughCellIds) or point id (PassThroughPointIds). The output precision of created points (if they need to be created) can also be specified.

Finally, this filter takes an optional second, vtkPolyData input. This input represents a list of faces that are to be excluded from the output of vtkGeometryFilter.

Warning
While vtkGeometryFilter and vtkDataSetSurfaceFilter perform similar operations, there are important differences as follows:
  1. vtkGeometryFilter can preserve (using RemoveGhostInterfaces) topological connectivity. vtkDataSetSurfaceFilter produces output primitives which may be disconnected from one another.
  2. vtkGeometryFilter can generate output based on cell ids, point ids, and/or extent (bounding box) clipping. vtkDataSetSurfaceFilter strictly extracts the boundary surface of a dataset.
  3. vtkGeometryFilter is much faster than vtkDataSetSurfaceFilter, because it's multi-threaded. As a result, vtkDataSetSurfaceFilter will delegate the processing of linear unstructured grids to vtkGeometryFilter.
  4. vtkGeometryFilter can (currently) only handle linear cells. The filter will delegate to vtkDataSetSurfaceFilter for higher-order cells. (This is a historical artifact and may be rectified in the future.)
If point merging (MergingOff) is disabled, the filter will (if possible) use the input points and point attributes. This can result in a lot of unused points in the output, at some gain in filter performance. If enabled, point merging will generate only new points that are used by the output polydata cells.
This class is templated. It may run slower than serial execution if the code is not optimized during compilation. Build in Release or ReleaseWithDebugInfo.
This class has been threaded with vtkSMPTools. Using TBB or other non-sequential type (set in the CMake variable VTK_SMP_IMPLEMENTATION_TYPE) may improve performance significantly.
See also
vtkDataSetSurfaceFilter vtkImageDataGeometryFilter vtkStructuredGridGeometryFilter vtkExtractGeometry vtkExtractVOI vtkMarkBoundaryFilter vtkRemovePolyData
Examples:
vtkGeometryFilter (Examples)
Online Examples:

Tests:
vtkGeometryFilter (Tests)

Definition at line 211 of file vtkGeometryFilter.h.

Member Typedef Documentation

◆ Superclass

Standard methods for instantiation, type information, and printing.

Definition at line 219 of file vtkGeometryFilter.h.

Constructor & Destructor Documentation

◆ vtkGeometryFilter()

vtkGeometryFilter::vtkGeometryFilter ( )
protected

◆ ~vtkGeometryFilter()

vtkGeometryFilter::~vtkGeometryFilter ( )
overrideprotected

Member Function Documentation

◆ New()

static vtkGeometryFilter * vtkGeometryFilter::New ( )
static

Standard methods for instantiation, type information, and printing.

◆ IsTypeOf()

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

Standard methods for instantiation, type information, and printing.

◆ IsA()

virtual vtkTypeBool vtkGeometryFilter::IsA ( const char *  type)
virtual

Standard methods for instantiation, type information, and printing.

Reimplemented from vtkPolyDataAlgorithm.

Reimplemented in vtkAdaptiveDataSetSurfaceFilter.

◆ SafeDownCast()

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

Standard methods for instantiation, type information, and printing.

◆ NewInstanceInternal()

virtual vtkObjectBase * vtkGeometryFilter::NewInstanceInternal ( ) const
protectedvirtual

Standard methods for instantiation, type information, and printing.

Reimplemented from vtkPolyDataAlgorithm.

Reimplemented in vtkAdaptiveDataSetSurfaceFilter.

◆ NewInstance()

vtkGeometryFilter * vtkGeometryFilter::NewInstance ( ) const

Standard methods for instantiation, type information, and printing.

◆ PrintSelf()

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

Standard methods for instantiation, type information, and printing.

Reimplemented from vtkAlgorithm.

◆ SetPointClipping()

virtual void vtkGeometryFilter::SetPointClipping ( bool  )
virtual

Turn on/off selection of geometry by point id.

◆ GetPointClipping()

virtual bool vtkGeometryFilter::GetPointClipping ( )
virtual

Turn on/off selection of geometry by point id.

◆ PointClippingOn()

virtual void vtkGeometryFilter::PointClippingOn ( )
virtual

Turn on/off selection of geometry by point id.

◆ PointClippingOff()

virtual void vtkGeometryFilter::PointClippingOff ( )
virtual

Turn on/off selection of geometry by point id.

◆ SetCellClipping()

virtual void vtkGeometryFilter::SetCellClipping ( bool  )
virtual

Turn on/off selection of geometry by cell id.

◆ GetCellClipping()

virtual bool vtkGeometryFilter::GetCellClipping ( )
virtual

Turn on/off selection of geometry by cell id.

◆ CellClippingOn()

virtual void vtkGeometryFilter::CellClippingOn ( )
virtual

Turn on/off selection of geometry by cell id.

◆ CellClippingOff()

virtual void vtkGeometryFilter::CellClippingOff ( )
virtual

Turn on/off selection of geometry by cell id.

◆ SetExtentClipping()

virtual void vtkGeometryFilter::SetExtentClipping ( bool  )
virtual

Turn on/off selection of geometry via bounding box.

◆ GetExtentClipping()

virtual bool vtkGeometryFilter::GetExtentClipping ( )
virtual

Turn on/off selection of geometry via bounding box.

◆ ExtentClippingOn()

virtual void vtkGeometryFilter::ExtentClippingOn ( )
virtual

Turn on/off selection of geometry via bounding box.

◆ ExtentClippingOff()

virtual void vtkGeometryFilter::ExtentClippingOff ( )
virtual

Turn on/off selection of geometry via bounding box.

◆ SetPointMinimum()

virtual void vtkGeometryFilter::SetPointMinimum ( vtkIdType  )
virtual

Specify the minimum point id for point id selection.

◆ GetPointMinimum()

virtual vtkIdType vtkGeometryFilter::GetPointMinimum ( )
virtual

Specify the minimum point id for point id selection.

◆ SetPointMaximum()

virtual void vtkGeometryFilter::SetPointMaximum ( vtkIdType  )
virtual

Specify the maximum point id for point id selection.

◆ GetPointMaximum()

virtual vtkIdType vtkGeometryFilter::GetPointMaximum ( )
virtual

Specify the maximum point id for point id selection.

◆ SetCellMinimum()

virtual void vtkGeometryFilter::SetCellMinimum ( vtkIdType  )
virtual

Specify the minimum cell id for point id selection.

◆ GetCellMinimum()

virtual vtkIdType vtkGeometryFilter::GetCellMinimum ( )
virtual

Specify the minimum cell id for point id selection.

◆ SetCellMaximum()

virtual void vtkGeometryFilter::SetCellMaximum ( vtkIdType  )
virtual

Specify the maximum cell id for point id selection.

◆ GetCellMaximum()

virtual vtkIdType vtkGeometryFilter::GetCellMaximum ( )
virtual

Specify the maximum cell id for point id selection.

◆ SetExtent() [1/2]

void vtkGeometryFilter::SetExtent ( double  xMin,
double  xMax,
double  yMin,
double  yMax,
double  zMin,
double  zMax 
)

Specify a (xmin,xmax, ymin,ymax, zmin,zmax) bounding box to clip data.

◆ SetExtent() [2/2]

void vtkGeometryFilter::SetExtent ( double  extent[6])

Set / get a (xmin,xmax, ymin,ymax, zmin,zmax) bounding box to clip data.

◆ GetExtent()

double * vtkGeometryFilter::GetExtent ( )
inline

Set / get a (xmin,xmax, ymin,ymax, zmin,zmax) bounding box to clip data.

Definition at line 292 of file vtkGeometryFilter.h.

◆ SetMerging()

virtual void vtkGeometryFilter::SetMerging ( bool  )
virtual

Turn on/off merging of points.

This will reduce the number of output points, at some cost to performance. If Merging is off, then if possible (i.e., if the point representation is explicit), the filter will reuse the input points to create the output polydata. Certain input dataset types (with implicit point representations) will always create new points (effectively performing a merge operation).

◆ GetMerging()

virtual bool vtkGeometryFilter::GetMerging ( )
virtual

Turn on/off merging of points.

This will reduce the number of output points, at some cost to performance. If Merging is off, then if possible (i.e., if the point representation is explicit), the filter will reuse the input points to create the output polydata. Certain input dataset types (with implicit point representations) will always create new points (effectively performing a merge operation).

◆ MergingOn()

virtual void vtkGeometryFilter::MergingOn ( )
virtual

Turn on/off merging of points.

This will reduce the number of output points, at some cost to performance. If Merging is off, then if possible (i.e., if the point representation is explicit), the filter will reuse the input points to create the output polydata. Certain input dataset types (with implicit point representations) will always create new points (effectively performing a merge operation).

◆ MergingOff()

virtual void vtkGeometryFilter::MergingOff ( )
virtual

Turn on/off merging of points.

This will reduce the number of output points, at some cost to performance. If Merging is off, then if possible (i.e., if the point representation is explicit), the filter will reuse the input points to create the output polydata. Certain input dataset types (with implicit point representations) will always create new points (effectively performing a merge operation).

◆ SetOutputPointsPrecision()

void vtkGeometryFilter::SetOutputPointsPrecision ( int  precision)

Set/get the desired precision for the output types.

See the documentation for the vtkAlgorithm::DesiredOutputPrecision enum for an explanation of the available precision settings. This only applies for data types where we create points (merging) as opposed to passing them from input to output, such as unstructured grids.

◆ GetOutputPointsPrecision()

int vtkGeometryFilter::GetOutputPointsPrecision ( ) const

Set/get the desired precision for the output types.

See the documentation for the vtkAlgorithm::DesiredOutputPrecision enum for an explanation of the available precision settings. This only applies for data types where we create points (merging) as opposed to passing them from input to output, such as unstructured grids.

◆ SetFastMode()

virtual void vtkGeometryFilter::SetFastMode ( bool  )
virtual

Turn on/off fast mode execution.

If enabled, fast mode typically runs much faster (2-3x) than the standard algorithm, however the output is an approximation to the correct result. FastMode is only meaningful when the input is vtkImageData/vtkRectilinearGrid/vtkStructuredGrid and there are blank cells.

◆ GetFastMode()

virtual bool vtkGeometryFilter::GetFastMode ( )
virtual

Turn on/off fast mode execution.

If enabled, fast mode typically runs much faster (2-3x) than the standard algorithm, however the output is an approximation to the correct result. FastMode is only meaningful when the input is vtkImageData/vtkRectilinearGrid/vtkStructuredGrid and there are blank cells.

◆ FastModeOn()

virtual void vtkGeometryFilter::FastModeOn ( )
virtual

Turn on/off fast mode execution.

If enabled, fast mode typically runs much faster (2-3x) than the standard algorithm, however the output is an approximation to the correct result. FastMode is only meaningful when the input is vtkImageData/vtkRectilinearGrid/vtkStructuredGrid and there are blank cells.

◆ FastModeOff()

virtual void vtkGeometryFilter::FastModeOff ( )
virtual

Turn on/off fast mode execution.

If enabled, fast mode typically runs much faster (2-3x) than the standard algorithm, however the output is an approximation to the correct result. FastMode is only meaningful when the input is vtkImageData/vtkRectilinearGrid/vtkStructuredGrid and there are blank cells.

◆ SetPieceInvariant()

virtual void vtkGeometryFilter::SetPieceInvariant ( int  )
virtual

If PieceInvariant is true, vtkGeometryFilter requests 1 ghost level from input in order to remove internal surface that are between processes.

False by default.

◆ GetPieceInvariant()

virtual int vtkGeometryFilter::GetPieceInvariant ( )
virtual

If PieceInvariant is true, vtkGeometryFilter requests 1 ghost level from input in order to remove internal surface that are between processes.

False by default.

◆ SetPassThroughCellIds()

virtual void vtkGeometryFilter::SetPassThroughCellIds ( vtkTypeBool  )
virtual

This parameter drives the generation or not of a CellData array for the output polygonal dataset that holds the cell index of the original 3D cell that produced each output cell.

This is useful for cell picking. The default is off to conserve memory.

Note: Use SetOriginalCellIdsName() to set the name of the CellData array.

◆ GetPassThroughCellIds()

virtual vtkTypeBool vtkGeometryFilter::GetPassThroughCellIds ( )
virtual

This parameter drives the generation or not of a CellData array for the output polygonal dataset that holds the cell index of the original 3D cell that produced each output cell.

This is useful for cell picking. The default is off to conserve memory.

Note: Use SetOriginalCellIdsName() to set the name of the CellData array.

◆ PassThroughCellIdsOn()

virtual void vtkGeometryFilter::PassThroughCellIdsOn ( )
virtual

This parameter drives the generation or not of a CellData array for the output polygonal dataset that holds the cell index of the original 3D cell that produced each output cell.

This is useful for cell picking. The default is off to conserve memory.

Note: Use SetOriginalCellIdsName() to set the name of the CellData array.

◆ PassThroughCellIdsOff()

virtual void vtkGeometryFilter::PassThroughCellIdsOff ( )
virtual

This parameter drives the generation or not of a CellData array for the output polygonal dataset that holds the cell index of the original 3D cell that produced each output cell.

This is useful for cell picking. The default is off to conserve memory.

Note: Use SetOriginalCellIdsName() to set the name of the CellData array.

◆ SetPassThroughPointIds()

virtual void vtkGeometryFilter::SetPassThroughPointIds ( vtkTypeBool  )
virtual

This parameter drives the generation or not of a PointData array for the output polygonal dataset that holds the cell/point index of the original point that produced each output point.

This is useful for point picking. The default is off to conserve memory.

Note: Use SetOriginalPointIdsName() to set the name of the PointData array.

◆ GetPassThroughPointIds()

virtual vtkTypeBool vtkGeometryFilter::GetPassThroughPointIds ( )
virtual

This parameter drives the generation or not of a PointData array for the output polygonal dataset that holds the cell/point index of the original point that produced each output point.

This is useful for point picking. The default is off to conserve memory.

Note: Use SetOriginalPointIdsName() to set the name of the PointData array.

◆ PassThroughPointIdsOn()

virtual void vtkGeometryFilter::PassThroughPointIdsOn ( )
virtual

This parameter drives the generation or not of a PointData array for the output polygonal dataset that holds the cell/point index of the original point that produced each output point.

This is useful for point picking. The default is off to conserve memory.

Note: Use SetOriginalPointIdsName() to set the name of the PointData array.

◆ PassThroughPointIdsOff()

virtual void vtkGeometryFilter::PassThroughPointIdsOff ( )
virtual

This parameter drives the generation or not of a PointData array for the output polygonal dataset that holds the cell/point index of the original point that produced each output point.

This is useful for point picking. The default is off to conserve memory.

Note: Use SetOriginalPointIdsName() to set the name of the PointData array.

◆ SetOriginalCellIdsName()

virtual void vtkGeometryFilter::SetOriginalCellIdsName ( const char *  )
virtual

If PassThroughCellIds or PassThroughPointIds is on, then these ivars control the name given to the field in which the ids are written into.

If set to nullptr, then vtkOriginalCellIds or vtkOriginalPointIds (the default) is used, respectively.

◆ GetOriginalCellIdsName()

virtual const char * vtkGeometryFilter::GetOriginalCellIdsName ( )
inlinevirtual

If PassThroughCellIds or PassThroughPointIds is on, then these ivars control the name given to the field in which the ids are written into.

If set to nullptr, then vtkOriginalCellIds or vtkOriginalPointIds (the default) is used, respectively.

Definition at line 380 of file vtkGeometryFilter.h.

◆ SetOriginalPointIdsName()

virtual void vtkGeometryFilter::SetOriginalPointIdsName ( const char *  )
virtual

If PassThroughCellIds or PassThroughPointIds is on, then these ivars control the name given to the field in which the ids are written into.

If set to nullptr, then vtkOriginalCellIds or vtkOriginalPointIds (the default) is used, respectively.

◆ GetOriginalPointIdsName()

virtual const char * vtkGeometryFilter::GetOriginalPointIdsName ( )
inlinevirtual

If PassThroughCellIds or PassThroughPointIds is on, then these ivars control the name given to the field in which the ids are written into.

If set to nullptr, then vtkOriginalCellIds or vtkOriginalPointIds (the default) is used, respectively.

Definition at line 385 of file vtkGeometryFilter.h.

◆ SetExcludedFacesData()

void vtkGeometryFilter::SetExcludedFacesData ( vtkPolyData )

If a second, vtkPolyData input is provided, this second input specifies a list of faces to be excluded from the output (in the vtkPolyData::Polys attribute).

This is useful to prevent the same face to be output multiple times in complex pipelines. (A candidate output boundary face is the same as a face in the excluded face list if it uses the same point ids as one of the polygons defined in the second input.) For example, a face may be extracted separately via a threshold filter; thus this face should not be also extracted via the vtkGeometryFilter. (This functionality is related to vtkRemovePolyData.)

◆ SetExcludedFacesConnection()

void vtkGeometryFilter::SetExcludedFacesConnection ( vtkAlgorithmOutput algOutput)

If a second, vtkPolyData input is provided, this second input specifies a list of faces to be excluded from the output (in the vtkPolyData::Polys attribute).

This is useful to prevent the same face to be output multiple times in complex pipelines. (A candidate output boundary face is the same as a face in the excluded face list if it uses the same point ids as one of the polygons defined in the second input.) For example, a face may be extracted separately via a threshold filter; thus this face should not be also extracted via the vtkGeometryFilter. (This functionality is related to vtkRemovePolyData.)

◆ GetExcludedFaces()

vtkPolyData * vtkGeometryFilter::GetExcludedFaces ( )

If a second, vtkPolyData input is provided, this second input specifies a list of faces to be excluded from the output (in the vtkPolyData::Polys attribute).

This is useful to prevent the same face to be output multiple times in complex pipelines. (A candidate output boundary face is the same as a face in the excluded face list if it uses the same point ids as one of the polygons defined in the second input.) For example, a face may be extracted separately via a threshold filter; thus this face should not be also extracted via the vtkGeometryFilter. (This functionality is related to vtkRemovePolyData.)

◆ SetNonlinearSubdivisionLevel()

virtual void vtkGeometryFilter::SetNonlinearSubdivisionLevel ( int  )
virtual

If the input is an unstructured grid with nonlinear faces, this parameter determines how many times the face is subdivided into linear faces.

If 0, the output is the equivalent of its linear counterpart (and the midpoints determining the nonlinear interpolation are discarded). If 1 (the default), the nonlinear face is triangulated based on the midpoints. If greater than 1, the triangulated pieces are recursively subdivided to reach the desired subdivision. Setting the value to greater than 1 may cause some point data to not be passed even if no nonlinear faces exist. This option has no effect if the input is not an unstructured grid.

◆ GetNonlinearSubdivisionLevel()

virtual int vtkGeometryFilter::GetNonlinearSubdivisionLevel ( )
virtual

If the input is an unstructured grid with nonlinear faces, this parameter determines how many times the face is subdivided into linear faces.

If 0, the output is the equivalent of its linear counterpart (and the midpoints determining the nonlinear interpolation are discarded). If 1 (the default), the nonlinear face is triangulated based on the midpoints. If greater than 1, the triangulated pieces are recursively subdivided to reach the desired subdivision. Setting the value to greater than 1 may cause some point data to not be passed even if no nonlinear faces exist. This option has no effect if the input is not an unstructured grid.

◆ SetMatchBoundariesIgnoringCellOrder()

virtual void vtkGeometryFilter::SetMatchBoundariesIgnoringCellOrder ( int  )
virtual

When two volumetric cells of different order are connected by their corners (for instance, a quadratic hexahedron next to a linear hexahedron ), the internal face is rendered and is not considered as a ghost cell.

To remove these faces, switch MatchBoundariesIgnoringCellOrder to 1 (default is 0).

◆ GetMatchBoundariesIgnoringCellOrder()

virtual int vtkGeometryFilter::GetMatchBoundariesIgnoringCellOrder ( )
virtual

When two volumetric cells of different order are connected by their corners (for instance, a quadratic hexahedron next to a linear hexahedron ), the internal face is rendered and is not considered as a ghost cell.

To remove these faces, switch MatchBoundariesIgnoringCellOrder to 1 (default is 0).

◆ SetDelegation()

virtual void vtkGeometryFilter::SetDelegation ( vtkTypeBool  )
virtual

Disable delegation to an internal vtkDataSetSurfaceFilter.

◆ GetDelegation()

virtual vtkTypeBool vtkGeometryFilter::GetDelegation ( )
virtual

Disable delegation to an internal vtkDataSetSurfaceFilter.

◆ DelegationOn()

virtual void vtkGeometryFilter::DelegationOn ( )
virtual

Disable delegation to an internal vtkDataSetSurfaceFilter.

◆ DelegationOff()

virtual void vtkGeometryFilter::DelegationOff ( )
virtual

Disable delegation to an internal vtkDataSetSurfaceFilter.

◆ SetRemoveGhostInterfaces()

virtual void vtkGeometryFilter::SetRemoveGhostInterfaces ( bool  )
virtual

Set/Get if Ghost interfaces will be removed.

When you are rendering you want to remove ghost interfaces that originate from duplicate cells.

There are certain algorithms though that need the ghost interfaces, such as GhostCellGenerator and FeatureEdges.

Since Rendering is the most common case, the Default is on.

Note: DON'T change it if there are no ghost cells.

◆ RemoveGhostInterfacesOn()

virtual void vtkGeometryFilter::RemoveGhostInterfacesOn ( )
virtual

Set/Get if Ghost interfaces will be removed.

When you are rendering you want to remove ghost interfaces that originate from duplicate cells.

There are certain algorithms though that need the ghost interfaces, such as GhostCellGenerator and FeatureEdges.

Since Rendering is the most common case, the Default is on.

Note: DON'T change it if there are no ghost cells.

◆ RemoveGhostInterfacesOff()

virtual void vtkGeometryFilter::RemoveGhostInterfacesOff ( )
virtual

Set/Get if Ghost interfaces will be removed.

When you are rendering you want to remove ghost interfaces that originate from duplicate cells.

There are certain algorithms though that need the ghost interfaces, such as GhostCellGenerator and FeatureEdges.

Since Rendering is the most common case, the Default is on.

Note: DON'T change it if there are no ghost cells.

◆ GetRemoveGhostInterfaces()

virtual bool vtkGeometryFilter::GetRemoveGhostInterfaces ( )
virtual

Set/Get if Ghost interfaces will be removed.

When you are rendering you want to remove ghost interfaces that originate from duplicate cells.

There are certain algorithms though that need the ghost interfaces, such as GhostCellGenerator and FeatureEdges.

Since Rendering is the most common case, the Default is on.

Note: DON'T change it if there are no ghost cells.

◆ PolyDataExecute() [1/2]

int vtkGeometryFilter::PolyDataExecute ( vtkDataSet input,
vtkPolyData output,
vtkPolyData exc 
)

Direct access methods so that this class can be used as an algorithm without using it as a filter (i.e., no pipeline updates).

Also some internal methods with additional options.

◆ PolyDataExecute() [2/2]

virtual int vtkGeometryFilter::PolyDataExecute ( vtkDataSet ,
vtkPolyData  
)
virtual

Direct access methods so that this class can be used as an algorithm without using it as a filter (i.e., no pipeline updates).

Also some internal methods with additional options.

◆ UnstructuredGridExecute() [1/2]

int vtkGeometryFilter::UnstructuredGridExecute ( vtkDataSet input,
vtkPolyData output,
vtkGeometryFilterHelper info,
vtkPolyData exc 
)

Direct access methods so that this class can be used as an algorithm without using it as a filter (i.e., no pipeline updates).

Also some internal methods with additional options.

◆ UnstructuredGridExecute() [2/2]

virtual int vtkGeometryFilter::UnstructuredGridExecute ( vtkDataSet input,
vtkPolyData output 
)
virtual

Direct access methods so that this class can be used as an algorithm without using it as a filter (i.e., no pipeline updates).

Also some internal methods with additional options.

◆ StructuredExecute() [1/2]

int vtkGeometryFilter::StructuredExecute ( vtkDataSet input,
vtkPolyData output,
int *  wholeExtent,
vtkPolyData exc,
bool *  extractFace = nullptr 
)

Direct access methods so that this class can be used as an algorithm without using it as a filter (i.e., no pipeline updates).

Also some internal methods with additional options.

◆ StructuredExecute() [2/2]

virtual int vtkGeometryFilter::StructuredExecute ( vtkDataSet input,
vtkPolyData output,
int *  wholeExt,
bool *  extractFace = nullptr 
)
virtual

Direct access methods so that this class can be used as an algorithm without using it as a filter (i.e., no pipeline updates).

Also some internal methods with additional options.

◆ DataSetExecute() [1/2]

int vtkGeometryFilter::DataSetExecute ( vtkDataSet input,
vtkPolyData output,
vtkPolyData exc 
)

Direct access methods so that this class can be used as an algorithm without using it as a filter (i.e., no pipeline updates).

Also some internal methods with additional options.

◆ DataSetExecute() [2/2]

virtual int vtkGeometryFilter::DataSetExecute ( vtkDataSet input,
vtkPolyData output 
)
virtual

Direct access methods so that this class can be used as an algorithm without using it as a filter (i.e., no pipeline updates).

Also some internal methods with additional options.

◆ RequestData()

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

This is called by the superclass.

This is the method you should override.

Reimplemented from vtkPolyDataAlgorithm.

◆ FillInputPortInformation()

int vtkGeometryFilter::FillInputPortInformation ( int  port,
vtkInformation info 
)
overrideprotectedvirtual

Fill the input port information objects for this algorithm.

This is invoked by the first call to GetInputPortInformation for each port so subclasses can specify what they can handle.

Reimplemented from vtkAlgorithm.

◆ RequestUpdateExtent()

int vtkGeometryFilter::RequestUpdateExtent ( vtkInformation ,
vtkInformationVector **  ,
vtkInformationVector  
)
overrideprotectedvirtual

This is called by the superclass.

This is the method you should override.

Reimplemented from vtkPolyDataAlgorithm.

Member Data Documentation

◆ PointMaximum

vtkIdType vtkGeometryFilter::PointMaximum
protected

Definition at line 493 of file vtkGeometryFilter.h.

◆ PointMinimum

vtkIdType vtkGeometryFilter::PointMinimum
protected

Definition at line 494 of file vtkGeometryFilter.h.

◆ CellMinimum

vtkIdType vtkGeometryFilter::CellMinimum
protected

Definition at line 495 of file vtkGeometryFilter.h.

◆ CellMaximum

vtkIdType vtkGeometryFilter::CellMaximum
protected

Definition at line 496 of file vtkGeometryFilter.h.

◆ Extent

double vtkGeometryFilter::Extent[6]
protected

Definition at line 497 of file vtkGeometryFilter.h.

◆ PointClipping

bool vtkGeometryFilter::PointClipping
protected

Definition at line 498 of file vtkGeometryFilter.h.

◆ CellClipping

bool vtkGeometryFilter::CellClipping
protected

Definition at line 499 of file vtkGeometryFilter.h.

◆ ExtentClipping

bool vtkGeometryFilter::ExtentClipping
protected

Definition at line 500 of file vtkGeometryFilter.h.

◆ OutputPointsPrecision

int vtkGeometryFilter::OutputPointsPrecision
protected

Definition at line 501 of file vtkGeometryFilter.h.

◆ RemoveGhostInterfaces

bool vtkGeometryFilter::RemoveGhostInterfaces
protected

Definition at line 502 of file vtkGeometryFilter.h.

◆ Merging

bool vtkGeometryFilter::Merging
protected

Definition at line 504 of file vtkGeometryFilter.h.

◆ Locator

vtkIncrementalPointLocator* vtkGeometryFilter::Locator
protected

Definition at line 505 of file vtkGeometryFilter.h.

◆ FastMode

bool vtkGeometryFilter::FastMode
protected

Definition at line 507 of file vtkGeometryFilter.h.

◆ PieceInvariant

int vtkGeometryFilter::PieceInvariant
protected

Definition at line 510 of file vtkGeometryFilter.h.

◆ PassThroughCellIds

vtkTypeBool vtkGeometryFilter::PassThroughCellIds
protected

Definition at line 511 of file vtkGeometryFilter.h.

◆ OriginalCellIdsName

char* vtkGeometryFilter::OriginalCellIdsName
protected

Definition at line 512 of file vtkGeometryFilter.h.

◆ PassThroughPointIds

vtkTypeBool vtkGeometryFilter::PassThroughPointIds
protected

Definition at line 514 of file vtkGeometryFilter.h.

◆ OriginalPointIdsName

char* vtkGeometryFilter::OriginalPointIdsName
protected

Definition at line 515 of file vtkGeometryFilter.h.

◆ NonlinearSubdivisionLevel

int vtkGeometryFilter::NonlinearSubdivisionLevel
protected

Definition at line 517 of file vtkGeometryFilter.h.

◆ MatchBoundariesIgnoringCellOrder

int vtkGeometryFilter::MatchBoundariesIgnoringCellOrder
protected

Definition at line 518 of file vtkGeometryFilter.h.

◆ Delegation

vtkTypeBool vtkGeometryFilter::Delegation
protected

Definition at line 520 of file vtkGeometryFilter.h.


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