VTK  9.1.0
Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
vtkContour3DLinearGrid Class Reference

fast generation of isosurface from 3D linear cells More...

#include <vtkContour3DLinearGrid.h>

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

Public Member Functions

virtual void SetComputeScalars (vtkTypeBool)
 Set/Get flag to compute scalars.
 
virtual vtkTypeBool GetComputeScalars ()
 
virtual void ComputeScalarsOn ()
 
virtual void ComputeScalarsOff ()
 
vtkMTimeType GetMTime () override
 Overloaded GetMTime() because of delegation to the internal vtkContourValues class.
 
int GetNumberOfThreadsUsed ()
 Return the number of threads actually used during execution.
 
bool GetLargeIds ()
 Inform the user as to whether large ids were used during filter execution.
 
void SetValue (int i, double value)
 Methods to set / get contour values.
 
double GetValue (int i)
 Get the ith contour value.
 
double * GetValues ()
 Get a pointer to an array of contour values.
 
void GetValues (double *contourValues)
 Fill a supplied list with contour values.
 
void SetNumberOfContours (int number)
 Set the number of contours to place into the list.
 
vtkIdType GetNumberOfContours ()
 Get the number of contours in the list of contour values.
 
void GenerateValues (int numContours, double range[2])
 Generate numContours equally spaced contour values between specified range.
 
void GenerateValues (int numContours, double rangeStart, double rangeEnd)
 Generate numContours equally spaced contour values between specified range.
 
virtual void SetMergePoints (vtkTypeBool)
 Indicate whether to merge coincident points.
 
virtual vtkTypeBool GetMergePoints ()
 Indicate whether to merge coincident points.
 
virtual void MergePointsOn ()
 Indicate whether to merge coincident points.
 
virtual void MergePointsOff ()
 Indicate whether to merge coincident points.
 
virtual void SetInterpolateAttributes (vtkTypeBool)
 Indicate whether to interpolate input attributes onto the isosurface.
 
virtual vtkTypeBool GetInterpolateAttributes ()
 Indicate whether to interpolate input attributes onto the isosurface.
 
virtual void InterpolateAttributesOn ()
 Indicate whether to interpolate input attributes onto the isosurface.
 
virtual void InterpolateAttributesOff ()
 Indicate whether to interpolate input attributes onto the isosurface.
 
virtual void SetComputeNormals (vtkTypeBool)
 Indicate whether to compute output point normals.
 
virtual vtkTypeBool GetComputeNormals ()
 Indicate whether to compute output point normals.
 
virtual void ComputeNormalsOn ()
 Indicate whether to compute output point normals.
 
virtual void ComputeNormalsOff ()
 Indicate whether to compute output point normals.
 
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 SetUseScalarTree (vtkTypeBool)
 Enable the use of a scalar tree to accelerate contour extraction.
 
virtual vtkTypeBool GetUseScalarTree ()
 Enable the use of a scalar tree to accelerate contour extraction.
 
virtual void UseScalarTreeOn ()
 Enable the use of a scalar tree to accelerate contour extraction.
 
virtual void UseScalarTreeOff ()
 Enable the use of a scalar tree to accelerate contour extraction.
 
virtual void SetScalarTree (vtkScalarTree *)
 Specify the scalar tree to use.
 
virtual vtkScalarTreeGetScalarTree ()
 Specify the scalar tree to use.
 
virtual void SetSequentialProcessing (vtkTypeBool)
 Force sequential processing (i.e.
 
virtual vtkTypeBool GetSequentialProcessing ()
 Force sequential processing (i.e.
 
virtual void SequentialProcessingOn ()
 Force sequential processing (i.e.
 
virtual void SequentialProcessingOff ()
 Force sequential processing (i.e.
 
- Public Member Functions inherited from vtkDataObjectAlgorithm
virtual vtkTypeBool IsA (const char *type)
 Return 1 if this class is the same type of (or a subclass of) the named class.
 
vtkDataObjectAlgorithmNewInstance () 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)
 
vtkDataObjectGetOutput ()
 Get the output data object for a port on this algorithm.
 
vtkDataObjectGetOutput (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.
 
int 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 SetProgress (double)
 SetProgress is deprecated.
 
void UpdateProgress (double amount)
 Update the progress of the process object.
 
virtual void SetInputArrayToProcess (int idx, int port, int connection, const char *fieldAssociation, const char *attributeTypeorName)
 String based versions of SetInputArrayToProcess().
 
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 int UpdateExtent (const int extents[6])
 Convenience method to update an algorithm after passing requests to its first output port.
 
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.
 
virtual 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 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.
 
virtual vtkInformationGetInformation ()
 Set/Get the information object associated with this algorithm.
 
virtual void SetInformation (vtkInformation *)
 Set/Get the information object associated with this algorithm.
 
void Register (vtkObjectBase *o) override
 Participate in garbage collection.
 
void UnRegister (vtkObjectBase *o) 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 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.
 
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 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 void Update (int port)
 Bring this algorithm's outputs up-to-date.
 
virtual void Update ()
 Bring this algorithm's outputs up-to-date.
 
virtual void SetReleaseDataFlag (int)
 Turn release data flag on or off for all output ports.
 
virtual int 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.
 
- 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)
 
int InvokeEvent (unsigned long event)
 
int InvokeEvent (const char *event)
 
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.
 
int InvokeEvent (unsigned long event, void *callData)
 This method invokes an event and return whether the event was aborted or not.
 
int InvokeEvent (const char *event, void *callData)
 This method invokes an event and return whether the event was aborted or not.
 
- Public Member Functions inherited from vtkObjectBase
const char * GetClassName () const
 Return the class name as a string.
 
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.
 
virtual 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.
 

Static Public Member Functions

static bool CanFullyProcessDataObject (vtkDataObject *object, const char *scalarArrayName)
 Returns true if the data object passed in is fully supported by this filter, i.e., all cell types are linear.
 
- 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.
 
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 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 (int 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 int 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.
 

Protected Member Functions

 vtkContour3DLinearGrid ()
 
 ~vtkContour3DLinearGrid () override
 
void ProcessPiece (vtkUnstructuredGrid *input, vtkDataArray *inScalars, vtkPolyData *output)
 
int RequestDataObject (vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override
 This is called by the superclass.
 
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.
 
- Protected Member Functions inherited from vtkDataObjectAlgorithm
virtual vtkObjectBaseNewInstanceInternal () const
 
 vtkDataObjectAlgorithm ()
 
 ~vtkDataObjectAlgorithm () override
 
virtual int RequestInformation (vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
 
virtual int RequestDataObject (vtkInformation *, vtkInformationVector **, vtkInformationVector *)
 This is called by the superclass.
 
virtual int RequestUpdateExtent (vtkInformation *, vtkInformationVector **, vtkInformationVector *)
 This is called by the superclass.
 
virtual int RequestData (vtkInformation *, vtkInformationVector **, vtkInformationVector *)
 
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
 
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 assocition 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 *)
 
 vtkObjectBase (const vtkObjectBase &)
 
void operator= (const vtkObjectBase &)
 

Protected Attributes

vtkContourValuesContourValues
 
int OutputPointsPrecision
 
vtkTypeBool MergePoints
 
vtkTypeBool InterpolateAttributes
 
vtkTypeBool ComputeNormals
 
vtkTypeBool ComputeScalars
 
vtkTypeBool SequentialProcessing
 
int NumberOfThreadsUsed
 
bool LargeIds
 
vtkTypeBool UseScalarTree
 
vtkScalarTreeScalarTree
 
struct vtkScalarTreeMap * ScalarTreeMap
 
- Protected Attributes inherited from vtkAlgorithm
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
 
- Protected Attributes inherited from vtkObjectBase
std::atomic< int32_t > ReferenceCount
 
vtkWeakPointerBase ** WeakPointers
 
typedef vtkDataObjectAlgorithm Superclass
 Standard methods for construction, type info, and printing.
 
static vtkContour3DLinearGridNew ()
 Standard methods for construction, type info, and printing.
 
static vtkTypeBool IsTypeOf (const char *type)
 Standard methods for construction, type info, and printing.
 
static vtkContour3DLinearGridSafeDownCast (vtkObjectBase *o)
 Standard methods for construction, type info, and printing.
 
virtual vtkTypeBool IsA (const char *type)
 Standard methods for construction, type info, and printing.
 
vtkContour3DLinearGridNewInstance () const
 Standard methods for construction, type info, and printing.
 
void PrintSelf (ostream &os, vtkIndent indent) override
 Standard methods for construction, type info, and printing.
 
virtual vtkObjectBaseNewInstanceInternal () const
 Standard methods for construction, type info, and printing.
 

Additional Inherited Members

- 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 Attributes inherited from vtkAlgorithm
vtkTypeBool AbortExecute
 
- Static Protected Member Functions inherited from vtkDataObjectAlgorithm
static bool SetOutputDataObject (int dataType, vtkInformation *outputInformation, bool exact=false)
 A helper method that can be used by subclasses in RequestDataObject to create an output data object of the given type if not already present.
 
- 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 vtkExecutiveDefaultExecutivePrototype
 

Detailed Description

fast generation of isosurface from 3D linear cells

vtkContour3DLinearGrid is a specialized filter that generates isocontours from 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. (Note: the filter will also process input vtkCompositeDataSets containing vtkUnstructuredGrids.)

To use this filter you must specify an input unstructured grid or vtkCompositeDataSet, and one or more contour values. You can either use the method SetValue() to specify each contour value, or use GenerateValues() to generate a series of evenly spaced contours.

The filter performance varies depending on optional output information. Basically if point merging is required (when PointMerging, InterpolateAttributes, and/or ComputeNormals is enabled), a sorting process is required to eliminate duplicate output points in the isosurface. Otherwise when point merging is not required, a fast path process produces independent triangles representing the isosurface. In many situations the results of the fast path are quite good and do not require additional processing.

Note that another performance option exists, using a vtkScalarTree, which is an object that accelerates isosurface extraction, at the initial cost of building the scalar tree. (This feature is useful for exploratory isosurface extraction when the isovalue is frequently changed.) In some cases this can improve performance, however this algorithm is so highly tuned that random memory jumps (due to random access of cells provided by the scalar tree) can actually negatively impact performance, especially if the input dataset type consists of homogeneous cell types.

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.
The fast path simply produces output points and triangles (the fast path executes when MergePoints if off; InterpolateAttributes is off; and ComputeNormals is off). Since the fast path does not merge points, it produces many more output points, typically on the order of 5-6x more than when MergePoints is enabled. Adding in the other options point merging, field interpolation, and normal generation results in additional performance impacts. By default the fast path is enabled.
When a vtkCompositeDataSet is provided as input, and UseScalarTree is enabled and a ScalarTree specified, then the specified scalar tree is cloned to create new ones for each dataset in the composite dataset. Otherwise (i.e., when vtkUnstructuredGrid input) the specified scalar tree is directly used (no cloning required).
Internal to this filter, a caching iterator is used to traverse the cells that compose the vtkUnstructuredGrid. Maximum performance is obtained if the cells are all of one type (i.e., input grid of homogeneous cell types); repeated switching from different types may have detrimental effects on performance.
For unstructured data, gradients are not computed. Normals are computed if requested; they are "pseudo-normals" in that the normals of output triangles that use a common point are averaged at the point. Alternatively use vtkPolyDataNormals to compute the surface normals.
The output of this filter is subtly different than the more general filter vtkContourGrid. vtkContourGrid eliminates small, degenerate triangles with concident points which are consequently not sent to the output. In practice this makes little impact on visual appearance but may have repercussions if the output is used for modelling and/or analysis.
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, and a limited subset of scalar types.
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
vtkContourGrid vtkContourFilter vtkFlyingEdges3D vtkMarchingCubes vtkPolyDataNormals vtkStaticEdgeLocatorTemplate.h vtkScalarTree vtkSpanSpace

Definition at line 125 of file vtkContour3DLinearGrid.h.

Member Typedef Documentation

◆ Superclass

Standard methods for construction, type info, and printing.

Definition at line 133 of file vtkContour3DLinearGrid.h.

Constructor & Destructor Documentation

◆ vtkContour3DLinearGrid()

vtkContour3DLinearGrid::vtkContour3DLinearGrid ( )
protected

◆ ~vtkContour3DLinearGrid()

vtkContour3DLinearGrid::~vtkContour3DLinearGrid ( )
overrideprotected

Member Function Documentation

◆ New()

static vtkContour3DLinearGrid * vtkContour3DLinearGrid::New ( )
static

Standard methods for construction, type info, and printing.

◆ IsTypeOf()

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

Standard methods for construction, type info, and printing.

◆ IsA()

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

Standard methods for construction, type info, and printing.

Reimplemented from vtkDataObjectAlgorithm.

◆ SafeDownCast()

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

Standard methods for construction, type info, and printing.

◆ NewInstanceInternal()

virtual vtkObjectBase * vtkContour3DLinearGrid::NewInstanceInternal ( ) const
protectedvirtual

Standard methods for construction, type info, and printing.

Reimplemented from vtkDataObjectAlgorithm.

◆ NewInstance()

vtkContour3DLinearGrid * vtkContour3DLinearGrid::NewInstance ( ) const

Standard methods for construction, type info, and printing.

◆ PrintSelf()

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

Standard methods for construction, type info, and printing.

Reimplemented from vtkDataObjectAlgorithm.

◆ SetValue()

void vtkContour3DLinearGrid::SetValue ( int  i,
double  value 
)
inline

Methods to set / get contour values.

Set a particular contour value at contour number i.

The index i ranges between 0<=i<NumberOfContours.

Definition at line 306 of file vtkContour3DLinearGrid.h.

◆ GetValue()

double vtkContour3DLinearGrid::GetValue ( int  i)
inline

Get the ith contour value.

Definition at line 314 of file vtkContour3DLinearGrid.h.

◆ GetValues() [1/2]

double * vtkContour3DLinearGrid::GetValues ( )
inline

Get a pointer to an array of contour values.

There will be GetNumberOfContours() values in the list.

Definition at line 323 of file vtkContour3DLinearGrid.h.

◆ GetValues() [2/2]

void vtkContour3DLinearGrid::GetValues ( double *  contourValues)
inline

Fill a supplied list with contour values.

There will be GetNumberOfContours() values in the list. Make sure you allocate enough memory to hold the list.

Definition at line 333 of file vtkContour3DLinearGrid.h.

◆ SetNumberOfContours()

void vtkContour3DLinearGrid::SetNumberOfContours ( int  number)
inline

Set the number of contours to place into the list.

You only really need to use this method to reduce list size. The method SetValue() will automatically increase list size as needed.

Definition at line 343 of file vtkContour3DLinearGrid.h.

◆ GetNumberOfContours()

vtkIdType vtkContour3DLinearGrid::GetNumberOfContours ( )
inline

Get the number of contours in the list of contour values.

Definition at line 351 of file vtkContour3DLinearGrid.h.

◆ GenerateValues() [1/2]

void vtkContour3DLinearGrid::GenerateValues ( int  numContours,
double  range[2] 
)
inline

Generate numContours equally spaced contour values between specified range.

Contour values will include min/max range values.

Definition at line 360 of file vtkContour3DLinearGrid.h.

◆ GenerateValues() [2/2]

void vtkContour3DLinearGrid::GenerateValues ( int  numContours,
double  rangeStart,
double  rangeEnd 
)
inline

Generate numContours equally spaced contour values between specified range.

Contour values will include min/max range values.

Definition at line 369 of file vtkContour3DLinearGrid.h.

◆ SetMergePoints()

virtual void vtkContour3DLinearGrid::SetMergePoints ( vtkTypeBool  )
virtual

Indicate whether to merge coincident points.

This takes extra time and produces fewer output points, creating a "watertight" contour surface. By default this is off.

◆ GetMergePoints()

virtual vtkTypeBool vtkContour3DLinearGrid::GetMergePoints ( )
virtual

Indicate whether to merge coincident points.

This takes extra time and produces fewer output points, creating a "watertight" contour surface. By default this is off.

◆ MergePointsOn()

virtual void vtkContour3DLinearGrid::MergePointsOn ( )
virtual

Indicate whether to merge coincident points.

This takes extra time and produces fewer output points, creating a "watertight" contour surface. By default this is off.

◆ MergePointsOff()

virtual void vtkContour3DLinearGrid::MergePointsOff ( )
virtual

Indicate whether to merge coincident points.

This takes extra time and produces fewer output points, creating a "watertight" contour surface. By default this is off.

◆ SetInterpolateAttributes()

virtual void vtkContour3DLinearGrid::SetInterpolateAttributes ( vtkTypeBool  )
virtual

Indicate whether to interpolate input attributes onto the isosurface.

By default this option is off.

◆ GetInterpolateAttributes()

virtual vtkTypeBool vtkContour3DLinearGrid::GetInterpolateAttributes ( )
virtual

Indicate whether to interpolate input attributes onto the isosurface.

By default this option is off.

◆ InterpolateAttributesOn()

virtual void vtkContour3DLinearGrid::InterpolateAttributesOn ( )
virtual

Indicate whether to interpolate input attributes onto the isosurface.

By default this option is off.

◆ InterpolateAttributesOff()

virtual void vtkContour3DLinearGrid::InterpolateAttributesOff ( )
virtual

Indicate whether to interpolate input attributes onto the isosurface.

By default this option is off.

◆ SetComputeNormals()

virtual void vtkContour3DLinearGrid::SetComputeNormals ( vtkTypeBool  )
virtual

Indicate whether to compute output point normals.

An averaging method is used to average shared triangle normals. By default this if off. This is a relatively expensive option so use judiciously.

◆ GetComputeNormals()

virtual vtkTypeBool vtkContour3DLinearGrid::GetComputeNormals ( )
virtual

Indicate whether to compute output point normals.

An averaging method is used to average shared triangle normals. By default this if off. This is a relatively expensive option so use judiciously.

◆ ComputeNormalsOn()

virtual void vtkContour3DLinearGrid::ComputeNormalsOn ( )
virtual

Indicate whether to compute output point normals.

An averaging method is used to average shared triangle normals. By default this if off. This is a relatively expensive option so use judiciously.

◆ ComputeNormalsOff()

virtual void vtkContour3DLinearGrid::ComputeNormalsOff ( )
virtual

Indicate whether to compute output point normals.

An averaging method is used to average shared triangle normals. By default this if off. This is a relatively expensive option so use judiciously.

◆ SetComputeScalars()

virtual void vtkContour3DLinearGrid::SetComputeScalars ( vtkTypeBool  )
virtual

Set/Get flag to compute scalars.

When enabled, and when the InterpolateAttributes option is on, vtkContour3DLinearGrid will add an array corresponding to the array used to compute the contour and populate it with values.

◆ GetComputeScalars()

virtual vtkTypeBool vtkContour3DLinearGrid::GetComputeScalars ( )
virtual

◆ ComputeScalarsOn()

virtual void vtkContour3DLinearGrid::ComputeScalarsOn ( )
virtual

◆ ComputeScalarsOff()

virtual void vtkContour3DLinearGrid::ComputeScalarsOff ( )
virtual

◆ SetOutputPointsPrecision()

void vtkContour3DLinearGrid::SetOutputPointsPrecision ( int  precision)

Set/get the desired precision for the output types.

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

◆ GetOutputPointsPrecision()

int vtkContour3DLinearGrid::GetOutputPointsPrecision ( ) const

Set/get the desired precision for the output types.

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

◆ GetMTime()

vtkMTimeType vtkContour3DLinearGrid::GetMTime ( )
overridevirtual

Overloaded GetMTime() because of delegation to the internal vtkContourValues class.

Reimplemented from vtkObject.

◆ SetUseScalarTree()

virtual void vtkContour3DLinearGrid::SetUseScalarTree ( vtkTypeBool  )
virtual

Enable the use of a scalar tree to accelerate contour extraction.

By default this is off. If enabled, and a scalar tree is not specified, then a vtkSpanSpace instance will be constructed and used.

◆ GetUseScalarTree()

virtual vtkTypeBool vtkContour3DLinearGrid::GetUseScalarTree ( )
virtual

Enable the use of a scalar tree to accelerate contour extraction.

By default this is off. If enabled, and a scalar tree is not specified, then a vtkSpanSpace instance will be constructed and used.

◆ UseScalarTreeOn()

virtual void vtkContour3DLinearGrid::UseScalarTreeOn ( )
virtual

Enable the use of a scalar tree to accelerate contour extraction.

By default this is off. If enabled, and a scalar tree is not specified, then a vtkSpanSpace instance will be constructed and used.

◆ UseScalarTreeOff()

virtual void vtkContour3DLinearGrid::UseScalarTreeOff ( )
virtual

Enable the use of a scalar tree to accelerate contour extraction.

By default this is off. If enabled, and a scalar tree is not specified, then a vtkSpanSpace instance will be constructed and used.

◆ SetScalarTree()

virtual void vtkContour3DLinearGrid::SetScalarTree ( vtkScalarTree )
virtual

Specify the scalar tree to use.

By default a vtkSpanSpace scalar tree is used.

◆ GetScalarTree()

virtual vtkScalarTree * vtkContour3DLinearGrid::GetScalarTree ( )
virtual

Specify the scalar tree to use.

By default a vtkSpanSpace scalar tree is used.

◆ SetSequentialProcessing()

virtual void vtkContour3DLinearGrid::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 vtkContour3DLinearGrid::GetSequentialProcessing ( )
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.

◆ SequentialProcessingOn()

virtual void vtkContour3DLinearGrid::SequentialProcessingOn ( )
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.

◆ SequentialProcessingOff()

virtual void vtkContour3DLinearGrid::SequentialProcessingOff ( )
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.

◆ GetNumberOfThreadsUsed()

int vtkContour3DLinearGrid::GetNumberOfThreadsUsed ( )
inline

Return the number of threads actually used during execution.

This is valid only after algorithm execution.

Definition at line 249 of file vtkContour3DLinearGrid.h.

◆ GetLargeIds()

bool vtkContour3DLinearGrid::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 259 of file vtkContour3DLinearGrid.h.

◆ CanFullyProcessDataObject()

static bool vtkContour3DLinearGrid::CanFullyProcessDataObject ( vtkDataObject object,
const char *  scalarArrayName 
)
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. The second array is the name of the array to process.

◆ ProcessPiece()

void vtkContour3DLinearGrid::ProcessPiece ( vtkUnstructuredGrid input,
vtkDataArray inScalars,
vtkPolyData output 
)
protected

◆ RequestDataObject()

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

This is called by the superclass.

This is the method you should override.

Reimplemented from vtkDataObjectAlgorithm.

◆ RequestData()

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

Reimplemented from vtkDataObjectAlgorithm.

◆ FillInputPortInformation()

int vtkContour3DLinearGrid::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

◆ ContourValues

vtkContourValues* vtkContour3DLinearGrid::ContourValues
protected

Definition at line 273 of file vtkContour3DLinearGrid.h.

◆ OutputPointsPrecision

int vtkContour3DLinearGrid::OutputPointsPrecision
protected

Definition at line 274 of file vtkContour3DLinearGrid.h.

◆ MergePoints

vtkTypeBool vtkContour3DLinearGrid::MergePoints
protected

Definition at line 275 of file vtkContour3DLinearGrid.h.

◆ InterpolateAttributes

vtkTypeBool vtkContour3DLinearGrid::InterpolateAttributes
protected

Definition at line 276 of file vtkContour3DLinearGrid.h.

◆ ComputeNormals

vtkTypeBool vtkContour3DLinearGrid::ComputeNormals
protected

Definition at line 277 of file vtkContour3DLinearGrid.h.

◆ ComputeScalars

vtkTypeBool vtkContour3DLinearGrid::ComputeScalars
protected

Definition at line 278 of file vtkContour3DLinearGrid.h.

◆ SequentialProcessing

vtkTypeBool vtkContour3DLinearGrid::SequentialProcessing
protected

Definition at line 279 of file vtkContour3DLinearGrid.h.

◆ NumberOfThreadsUsed

int vtkContour3DLinearGrid::NumberOfThreadsUsed
protected

Definition at line 280 of file vtkContour3DLinearGrid.h.

◆ LargeIds

bool vtkContour3DLinearGrid::LargeIds
protected

Definition at line 281 of file vtkContour3DLinearGrid.h.

◆ UseScalarTree

vtkTypeBool vtkContour3DLinearGrid::UseScalarTree
protected

Definition at line 284 of file vtkContour3DLinearGrid.h.

◆ ScalarTree

vtkScalarTree* vtkContour3DLinearGrid::ScalarTree
protected

Definition at line 285 of file vtkContour3DLinearGrid.h.

◆ ScalarTreeMap

struct vtkScalarTreeMap* vtkContour3DLinearGrid::ScalarTreeMap
protected

Definition at line 286 of file vtkContour3DLinearGrid.h.


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