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

fast plane cutting of vtkUnstructuredGrid containing 3D linear cells More...

#include <vtk3DLinearGridPlaneCutter.h>

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

Public Types

typedef vtkDataObjectAlgorithm Superclass
 
- Public Types inherited from vtkDataObjectAlgorithm
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...
 
vtk3DLinearGridPlaneCutterNewInstance () const
 
void PrintSelf (ostream &os, vtkIndent indent) override
 Methods invoked by print to print information about the object including superclasses. More...
 
virtual void SetPlane (vtkPlane *)
 Specify the plane (an implicit function) to perform the cutting. More...
 
virtual vtkPlaneGetPlane ()
 
virtual void SetMergePoints (vtkTypeBool)
 Indicate whether to merge coincident points. More...
 
virtual vtkTypeBool GetMergePoints ()
 
virtual void MergePointsOn ()
 
virtual void MergePointsOff ()
 
virtual void SetInterpolateAttributes (vtkTypeBool)
 Indicate whether to interpolate input attributes onto the cut plane. More...
 
virtual vtkTypeBool GetInterpolateAttributes ()
 
virtual void InterpolateAttributesOn ()
 
virtual void InterpolateAttributesOff ()
 
virtual void SetComputeNormals (vtkTypeBool)
 Set/Get the computation of normals. More...
 
virtual vtkTypeBool GetComputeNormals ()
 
virtual void ComputeNormalsOn ()
 
virtual void ComputeNormalsOff ()
 
vtkMTimeType GetMTime () override
 Overloaded GetMTime() because of delegation to the helper vtkPlane. More...
 
void SetOutputPointsPrecision (int precision)
 Set/get the desired precision for the output points. More...
 
int GetOutputPointsPrecision () const
 
virtual void SetSequentialProcessing (vtkTypeBool)
 Force sequential processing (i.e. More...
 
virtual vtkTypeBool GetSequentialProcessing ()
 
virtual void SequentialProcessingOn ()
 
virtual void SequentialProcessingOff ()
 
int GetNumberOfThreadsUsed ()
 Return the number of threads actually used during execution. More...
 
bool GetLargeIds ()
 Inform the user as to whether large ids were used during filter execution. More...
 
- Public Member Functions inherited from vtkDataObjectAlgorithm
vtkDataObjectAlgorithmNewInstance () const
 
vtkDataObjectGetOutput ()
 Get the output data object for a port on this algorithm. More...
 
vtkDataObjectGetOutput (int)
 
virtual void SetOutput (vtkDataObject *d)
 
vtkTypeBool ProcessRequest (vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
 see vtkAlgorithm for details More...
 
vtkDataObjectGetInput ()
 
vtkDataObjectGetInput (int port)
 
void SetInputData (vtkDataObject *)
 Assign a data object as input. More...
 
void SetInputData (int, vtkDataObject *)
 
void AddInputData (vtkDataObject *)
 Assign a data object as input. More...
 
void AddInputData (int, vtkDataObject *)
 
- Public Member Functions inherited from vtkAlgorithm
vtkAlgorithmNewInstance () const
 
int 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...
 
virtual vtkInformationGetInformation ()
 Set/Get the information object associated with this algorithm. More...
 
virtual void SetInformation (vtkInformation *)
 
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 Register (vtkObjectBase *o) override
 Participate in garbage collection. More...
 
void UnRegister (vtkObjectBase *o) override
 Decrease the reference count (release by another object). More...
 
virtual void SetAbortExecute (vtkTypeBool)
 Set/Get the AbortExecute flag for the process object. More...
 
virtual vtkTypeBool GetAbortExecute ()
 
virtual void AbortExecuteOn ()
 
virtual void AbortExecuteOff ()
 
virtual double GetProgress ()
 Get the execution progress of a process object. More...
 
void SetProgress (double)
 SetProgress is deprecated. More...
 
void UpdateProgress (double amount)
 Update the progress of the process object. 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 ()
 
virtual double GetProgressScale ()
 
void SetProgressText (const char *ptext)
 Set the current text message associated with the progress state. More...
 
virtual char * GetProgressText ()
 
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)
 
virtual void SetInputArrayToProcess (int idx, vtkInformation *info)
 
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 SetInputConnection (int port, vtkAlgorithmOutput *input)
 Set the connection for the given input port index. More...
 
virtual void SetInputConnection (vtkAlgorithmOutput *input)
 
virtual void AddInputConnection (int port, vtkAlgorithmOutput *input)
 Add a connection to the given input port index. More...
 
virtual void AddInputConnection (vtkAlgorithmOutput *input)
 
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 void Update (int port)
 Bring this algorithm's outputs up-to-date. More...
 
virtual void Update ()
 
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...
 
virtual void SetReleaseDataFlag (int)
 Turn release data flag on or off for all output ports. More...
 
virtual int GetReleaseDataFlag ()
 
void ReleaseDataFlagOn ()
 
void ReleaseDataFlagOff ()
 
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)
 
intGetUpdateExtent ()
 These functions return the update extent for output ports that use 3D extents. More...
 
intGetUpdateExtent (int port)
 
void GetUpdateExtent (int &x0, int &x1, int &y0, int &y1, int &z0, int &z1)
 
void GetUpdateExtent (int port, int &x0, int &x1, int &y0, int &y1, int &z0, int &z1)
 
void GetUpdateExtent (int extent[6])
 
void GetUpdateExtent (int port, int extent[6])
 
int GetUpdatePiece ()
 These functions return the update extent for output ports that use piece extents. More...
 
int GetUpdatePiece (int port)
 
int GetUpdateNumberOfPieces ()
 
int GetUpdateNumberOfPieces (int port)
 
int GetUpdateGhostLevel ()
 
int GetUpdateGhostLevel (int port)
 
void SetProgressObserver (vtkProgressObserver *)
 If an ProgressObserver is set, the algorithm will report progress through it rather than directly. More...
 
virtual vtkProgressObserverGetProgressObserver ()
 
- 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...
 
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)
 
vtkCommandGetCommand (unsigned long tag)
 
void RemoveObserver (vtkCommand *)
 
void RemoveObservers (unsigned long event, vtkCommand *)
 
void RemoveObservers (const char *event, vtkCommand *)
 
vtkTypeBool HasObserver (unsigned long event, vtkCommand *)
 
vtkTypeBool HasObserver (const char *event, vtkCommand *)
 
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)
 
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)
 
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...
 
int InvokeEvent (unsigned long event, void *callData)
 This method invokes an event and return whether the event was aborted or not. More...
 
int InvokeEvent (const char *event, void *callData)
 
int InvokeEvent (unsigned long event)
 
int InvokeEvent (const char *event)
 
- Public Member Functions inherited from vtkObjectBase
const char * GetClassName () const
 Return the class name as a string. 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...
 
virtual void PrintHeader (ostream &os, vtkIndent indent)
 
virtual void PrintTrailer (ostream &os, vtkIndent indent)
 
int GetReferenceCount ()
 Return the current reference count of this object. More...
 
void SetReferenceCount (int)
 Sets the reference count. More...
 
void PrintRevisions (ostream &)
 Legacy. More...
 

Static Public Member Functions

static vtk3DLinearGridPlaneCutterNew ()
 Standard methods for construction, type info, and printing. More...
 
static vtkTypeBool IsTypeOf (const char *type)
 
static vtk3DLinearGridPlaneCutterSafeDownCast (vtkObjectBase *o)
 
static bool CanFullyProcessDataObject (vtkDataObject *object)
 Returns true if the data object passed in is fully supported by this filter, i.e., all cell types are linear. More...
 
- Static Public Member Functions inherited from vtkDataObjectAlgorithm
static vtkDataObjectAlgorithmNew ()
 
static vtkTypeBool IsTypeOf (const char *type)
 
static vtkDataObjectAlgorithmSafeDownCast (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 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 (int val)
 This is a global flag that controls whether any debug, warning or error messages are displayed. More...
 
static void GlobalWarningDisplayOn ()
 
static void GlobalWarningDisplayOff ()
 
static int GetGlobalWarningDisplay ()
 
- 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 vtkObjectBaseNew ()
 Create an object with Debug turned off, modified time initialized to zero, and reference counting on. More...
 

Protected Member Functions

virtual vtkObjectBaseNewInstanceInternal () const
 
 vtk3DLinearGridPlaneCutter ()
 
 ~vtk3DLinearGridPlaneCutter () override
 
int ProcessPiece (vtkUnstructuredGrid *input, vtkPlane *plane, vtkPolyData *output)
 
int RequestDataObject (vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
 This is called by the superclass. More...
 
int RequestData (vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
 
int FillInputPortInformation (int port, vtkInformation *info) override
 Fill the input port information objects for this algorithm. More...
 
- Protected Member Functions inherited from vtkDataObjectAlgorithm
 vtkDataObjectAlgorithm ()
 
 ~vtkDataObjectAlgorithm () override
 
virtual int RequestInformation (vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
 
virtual int RequestUpdateExtent (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
 
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 assocition of the actual data array for the input array specified by idx, this is only reasonable during the REQUEST_DATA pass. More...
 
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)
 
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)
 
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)
 
vtkDataArrayGetInputArrayToProcess (int idx, vtkDataObject *input)
 
vtkDataArrayGetInputArrayToProcess (int idx, vtkDataObject *input, int &association)
 
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)
 
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)
 
vtkAbstractArrayGetInputAbstractArrayToProcess (int idx, vtkDataObject *input)
 
vtkAbstractArrayGetInputAbstractArrayToProcess (int idx, vtkDataObject *input, int &association)
 
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...
 
virtual void SetErrorCode (unsigned long)
 The error code contains a possible error that occurred while reading or writing the file. 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)
 
- 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 ()
 
- Protected Member Functions inherited from vtkObjectBase
 vtkObjectBase ()
 
virtual ~vtkObjectBase ()
 
virtual void CollectRevisions (ostream &)
 
 vtkObjectBase (const vtkObjectBase &)
 
void operator= (const vtkObjectBase &)
 

Protected Attributes

vtkPlanePlane
 
vtkTypeBool MergePoints
 
vtkTypeBool InterpolateAttributes
 
vtkTypeBool ComputeNormals
 
int OutputPointsPrecision
 
vtkTypeBool SequentialProcessing
 
int NumberOfThreadsUsed
 
bool LargeIds
 
- Protected Attributes inherited from vtkAlgorithm
vtkInformationInformation
 
unsigned long ErrorCode
 
double Progress
 
char * ProgressText
 
vtkProgressObserverProgressObserver
 
- Protected Attributes inherited from vtkObject
bool Debug
 
vtkTimeStamp MTime
 
vtkSubjectHelper * SubjectHelper
 
- Protected Attributes inherited from vtkObjectBase
std::atomic< int32_t > ReferenceCount
 
vtkWeakPointerBase ** WeakPointers
 

Additional Inherited Members

- Public Attributes inherited from vtkAlgorithm
vtkTypeBool AbortExecute
 
- Static Protected Member Functions inherited from vtkAlgorithm
static vtkInformationIntegerKeyPORT_REQUIREMENTS_FILLED ()
 
- Static Protected Attributes inherited from vtkAlgorithm
static vtkExecutiveDefaultExecutivePrototype
 

Detailed Description

fast plane cutting of vtkUnstructuredGrid containing 3D linear cells

vtk3DLinearGridPlaneCutter is a specialized filter that cuts an input vtkUnstructuredGrid consisting of 3D linear cells: tetrahedra, hexahedra, voxels, pyramids, and/or wedges. (The cells are linear in the sense that each cell edge is a straight line.) The filter is designed for high-speed, specialized operation. All other cell types are skipped and produce no output.

To use this filter you must specify an input unstructured grid or vtkCompositeDataSet (containing unstructured grids) and a plane to cut with.

The filter performance varies depending on optional output information. Basically if point merging is required (when PointMerging is set) a sorting process is required to eliminate duplicate output points in the cut surface. Otherwise when point merging is not required, a fast path process produces independent triangles representing the cut surface.

This algorithm is fast because it is threaded, and may perform oeprations (in a preprocessing step) to accelerate the plane cutting.

Because this filter may build an initial data structure during a preprocessing step, the first execution of the filter may take longer than subsequent operations. Typically the first execution is still faster than vtkCutter (especially with threading enabled), but for certain types of data this may not be true. However if you are using the filter to cut a dataset multiple times (as in an exploratory or interactive workflow) this filter works well.

Warning
When the input is of type vtkCompositeDataSet the filter will process the unstructured grid(s) contained in the composite data set. As a result the output of this filter is then a vtkMultiBlockDataSet containing multiple vtkPolyData. When a vtkUnstructuredGrid is provided as input the output is a single vtkPolyData.
Input cells that are not of 3D linear type (tetrahedron, hexahedron, wedge, pyramid, and voxel) are simply skipped and not processed.
The filter is templated on types of input and output points, and input scalar type. To reduce object file bloat, only real points (float,double) are processed.
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
vtkCutter vtkFlyingEdgesPlaneCutter vtkPlaneCutter vtkPlane vtkSphereTree vtkContour3DLinearGrid

Definition at line 83 of file vtk3DLinearGridPlaneCutter.h.

Member Typedef Documentation

◆ Superclass

Definition at line 91 of file vtk3DLinearGridPlaneCutter.h.

Constructor & Destructor Documentation

◆ vtk3DLinearGridPlaneCutter()

vtk3DLinearGridPlaneCutter::vtk3DLinearGridPlaneCutter ( )
protected

◆ ~vtk3DLinearGridPlaneCutter()

vtk3DLinearGridPlaneCutter::~vtk3DLinearGridPlaneCutter ( )
overrideprotected

Member Function Documentation

◆ New()

static vtk3DLinearGridPlaneCutter* vtk3DLinearGridPlaneCutter::New ( )
static

Standard methods for construction, type info, and printing.

◆ IsTypeOf()

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

◆ IsA()

virtual vtkTypeBool vtk3DLinearGridPlaneCutter::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 vtkDataObjectAlgorithm.

◆ SafeDownCast()

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

◆ NewInstanceInternal()

virtual vtkObjectBase* vtk3DLinearGridPlaneCutter::NewInstanceInternal ( ) const
protectedvirtual

Reimplemented from vtkDataObjectAlgorithm.

◆ NewInstance()

vtk3DLinearGridPlaneCutter* vtk3DLinearGridPlaneCutter::NewInstance ( ) const

◆ PrintSelf()

void vtk3DLinearGridPlaneCutter::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 vtkDataObjectAlgorithm.

◆ SetPlane()

virtual void vtk3DLinearGridPlaneCutter::SetPlane ( vtkPlane )
virtual

Specify the plane (an implicit function) to perform the cutting.

The definition of the plane (its origin and normal) is controlled via this instance of vtkPlane.

◆ GetPlane()

virtual vtkPlane* vtk3DLinearGridPlaneCutter::GetPlane ( )
virtual

◆ SetMergePoints()

virtual void vtk3DLinearGridPlaneCutter::SetMergePoints ( vtkTypeBool  )
virtual

Indicate whether to merge coincident points.

Merging can take extra time and produces fewer output points, creating a "watertight" output surface. On the other hand, merging reduced output data size and may be just as fast especially for smaller data. By default this is off.

◆ GetMergePoints()

virtual vtkTypeBool vtk3DLinearGridPlaneCutter::GetMergePoints ( )
virtual

◆ MergePointsOn()

virtual void vtk3DLinearGridPlaneCutter::MergePointsOn ( )
virtual

◆ MergePointsOff()

virtual void vtk3DLinearGridPlaneCutter::MergePointsOff ( )
virtual

◆ SetInterpolateAttributes()

virtual void vtk3DLinearGridPlaneCutter::SetInterpolateAttributes ( vtkTypeBool  )
virtual

Indicate whether to interpolate input attributes onto the cut plane.

By default this option is on.

◆ GetInterpolateAttributes()

virtual vtkTypeBool vtk3DLinearGridPlaneCutter::GetInterpolateAttributes ( )
virtual

◆ InterpolateAttributesOn()

virtual void vtk3DLinearGridPlaneCutter::InterpolateAttributesOn ( )
virtual

◆ InterpolateAttributesOff()

virtual void vtk3DLinearGridPlaneCutter::InterpolateAttributesOff ( )
virtual

◆ SetComputeNormals()

virtual void vtk3DLinearGridPlaneCutter::SetComputeNormals ( vtkTypeBool  )
virtual

Set/Get the computation of normals.

The normal generated is simply the cut plane normal. The normal, if generated, is defined by cell data associated with the output polygons. By default computing of normals is off.

◆ GetComputeNormals()

virtual vtkTypeBool vtk3DLinearGridPlaneCutter::GetComputeNormals ( )
virtual

◆ ComputeNormalsOn()

virtual void vtk3DLinearGridPlaneCutter::ComputeNormalsOn ( )
virtual

◆ ComputeNormalsOff()

virtual void vtk3DLinearGridPlaneCutter::ComputeNormalsOff ( )
virtual

◆ GetMTime()

vtkMTimeType vtk3DLinearGridPlaneCutter::GetMTime ( )
overridevirtual

Overloaded GetMTime() because of delegation to the helper vtkPlane.

Reimplemented from vtkObject.

◆ SetOutputPointsPrecision()

void vtk3DLinearGridPlaneCutter::SetOutputPointsPrecision ( int  precision)

Set/get the desired precision for the output points.

See the documentation for the vtkAlgorithm::Precision enum for an explanation of the available precision settings.

◆ GetOutputPointsPrecision()

int vtk3DLinearGridPlaneCutter::GetOutputPointsPrecision ( ) const

◆ SetSequentialProcessing()

virtual void vtk3DLinearGridPlaneCutter::SetSequentialProcessing ( vtkTypeBool  )
virtual

Force sequential processing (i.e.

single thread) of the contouring process. By default, sequential processing is off. Note this flag only applies if the class has been compiled with VTK_SMP_IMPLEMENTATION_TYPE set to something other than Sequential. (If set to Sequential, then the filter always runs in serial mode.) This flag is typically used for benchmarking purposes.

◆ GetSequentialProcessing()

virtual vtkTypeBool vtk3DLinearGridPlaneCutter::GetSequentialProcessing ( )
virtual

◆ SequentialProcessingOn()

virtual void vtk3DLinearGridPlaneCutter::SequentialProcessingOn ( )
virtual

◆ SequentialProcessingOff()

virtual void vtk3DLinearGridPlaneCutter::SequentialProcessingOff ( )
virtual

◆ GetNumberOfThreadsUsed()

int vtk3DLinearGridPlaneCutter::GetNumberOfThreadsUsed ( )
inline

Return the number of threads actually used during execution.

This is valid only after algorithm execution.

Definition at line 173 of file vtk3DLinearGridPlaneCutter.h.

◆ GetLargeIds()

bool vtk3DLinearGridPlaneCutter::GetLargeIds ( )
inline

Inform the user as to whether large ids were used during filter execution.

This flag only has meaning after the filter has executed. Large ids are used when the id of the larges cell or point is greater than signed 32-bit precision. (Smaller ids reduce memory usage and speed computation. Note that LargeIds are only available on 64-bit architectures.)

Definition at line 184 of file vtk3DLinearGridPlaneCutter.h.

◆ CanFullyProcessDataObject()

static bool vtk3DLinearGridPlaneCutter::CanFullyProcessDataObject ( vtkDataObject object)
static

Returns true if the data object passed in is fully supported by this filter, i.e., all cell types are linear.

For composite datasets, this means all dataset leaves have only linear cell types that can be processed by this filter.

◆ ProcessPiece()

int vtk3DLinearGridPlaneCutter::ProcessPiece ( vtkUnstructuredGrid input,
vtkPlane plane,
vtkPolyData output 
)
protected

◆ RequestDataObject()

int vtk3DLinearGridPlaneCutter::RequestDataObject ( vtkInformation ,
vtkInformationVector **  ,
vtkInformationVector  
)
overrideprotectedvirtual

This is called by the superclass.

This is the method you should override.

Reimplemented from vtkDataObjectAlgorithm.

◆ RequestData()

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

Reimplemented from vtkDataObjectAlgorithm.

◆ FillInputPortInformation()

int vtk3DLinearGridPlaneCutter::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 vtkDataObjectAlgorithm.

Member Data Documentation

◆ Plane

vtkPlane* vtk3DLinearGridPlaneCutter::Plane
protected

Definition at line 199 of file vtk3DLinearGridPlaneCutter.h.

◆ MergePoints

vtkTypeBool vtk3DLinearGridPlaneCutter::MergePoints
protected

Definition at line 200 of file vtk3DLinearGridPlaneCutter.h.

◆ InterpolateAttributes

vtkTypeBool vtk3DLinearGridPlaneCutter::InterpolateAttributes
protected

Definition at line 201 of file vtk3DLinearGridPlaneCutter.h.

◆ ComputeNormals

vtkTypeBool vtk3DLinearGridPlaneCutter::ComputeNormals
protected

Definition at line 202 of file vtk3DLinearGridPlaneCutter.h.

◆ OutputPointsPrecision

int vtk3DLinearGridPlaneCutter::OutputPointsPrecision
protected

Definition at line 203 of file vtk3DLinearGridPlaneCutter.h.

◆ SequentialProcessing

vtkTypeBool vtk3DLinearGridPlaneCutter::SequentialProcessing
protected

Definition at line 204 of file vtk3DLinearGridPlaneCutter.h.

◆ NumberOfThreadsUsed

int vtk3DLinearGridPlaneCutter::NumberOfThreadsUsed
protected

Definition at line 205 of file vtk3DLinearGridPlaneCutter.h.

◆ LargeIds

bool vtk3DLinearGridPlaneCutter::LargeIds
protected

Definition at line 206 of file vtk3DLinearGridPlaneCutter.h.


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