VTK  9.4.20241212
Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
vtkSurfaceNets3D Class Reference

generate smoothed isocontours from segmented 3D image data (i.e., "label maps") More...

#include <vtkSurfaceNets3D.h>

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

Public Types

enum  MeshType { MESH_TYPE_DEFAULT = 0 , MESH_TYPE_TRIANGLES , MESH_TYPE_QUADS }
 This enum is used to control the type of the output polygonal mesh. More...
 
enum  OutputType { OUTPUT_STYLE_DEFAULT = 0 , OUTPUT_STYLE_BOUNDARY , OUTPUT_STYLE_SELECTED }
 This enum is used to control the production of the filter output. More...
 
enum  TriangulationType { TRIANGULATION_GREEDY = 0 , TRIANGULATION_MIN_EDGE , TRIANGULATION_MIN_AREA }
 This enum is used to control how quadrilaterals are triangulated. More...
 
- Public Types inherited from vtkPolyDataAlgorithm
typedef vtkAlgorithm Superclass
 
- Public Types inherited from vtkAlgorithm
enum  DesiredOutputPrecision { SINGLE_PRECISION , DOUBLE_PRECISION , DEFAULT_PRECISION }
 Values used for setting the desired output precision for various algorithms. More...
 
typedef vtkObject Superclass
 

Public Member Functions

vtkMTimeType GetMTime () override
 The modified time is also a function of the label values and the smoothing filter.
 
void SetValue (int i, double value)
 Set a particular label value at label number i.
 
void SetLabel (int i, double value)
 Set a particular label value at label number i.
 
double GetValue (int i)
 Get the ith label value.
 
double GetLabel (int i)
 Get the ith label value.
 
double * GetValues ()
 Get a pointer to an array of labels.
 
double * GetLabels ()
 Get a pointer to an array of labels.
 
void GetValues (double *contourValues)
 Fill a supplied list with label values.
 
void GetLabels (double *contourValues)
 Fill a supplied list with label values.
 
void SetNumberOfLabels (int number)
 Set the number of labels to place into the list.
 
void SetNumberOfContours (int number)
 Set the number of labels to place into the list.
 
vtkIdType GetNumberOfLabels ()
 Get the number of labels in the list of label values.
 
vtkIdType GetNumberOfContours ()
 Get the number of labels in the list of label values.
 
void GenerateLabels (int numLabels, double range[2])
 Generate numLabels equally spaced labels between the specified range.
 
void GenerateValues (int numContours, double range[2])
 Generate numLabels equally spaced labels between the specified range.
 
void GenerateLabels (int numLabels, double rangeStart, double rangeEnd)
 Generate numLabels equally spaced labels between the specified range.
 
void GenerateValues (int numContours, double rangeStart, double rangeEnd)
 Generate numLabels equally spaced labels between the specified range.
 
virtual void SetBackgroundLabel (double)
 This value specifies the label value to use when referencing the background region outside of any of the specified regions.
 
virtual double GetBackgroundLabel ()
 This value specifies the label value to use when referencing the background region outside of any of the specified regions.
 
virtual void SetArrayComponent (int)
 Set/get which component of a input multi-component scalar array to contour with; defaults to component 0.
 
virtual int GetArrayComponent ()
 Set/get which component of a input multi-component scalar array to contour with; defaults to component 0.
 
virtual void SetOutputMeshType (int)
 Control the type of output mesh.
 
virtual int GetOutputMeshType ()
 Control the type of output mesh.
 
void SetOutputMeshTypeToDefault ()
 Control the type of output mesh.
 
void SetOutputMeshTypeToTriangles ()
 Control the type of output mesh.
 
void SetOutputMeshTypeToQuads ()
 Control the type of output mesh.
 
virtual void SetSmoothing (bool)
 Indicate whether smoothing should be enabled.
 
virtual bool GetSmoothing ()
 Indicate whether smoothing should be enabled.
 
virtual void SmoothingOn ()
 Indicate whether smoothing should be enabled.
 
virtual void SmoothingOff ()
 Indicate whether smoothing should be enabled.
 
void SetNumberOfIterations (int n)
 Convenience methods that delegate to the internal smoothing filter follow below.
 
int GetNumberOfIterations ()
 Convenience methods that delegate to the internal smoothing filter follow below.
 
void SetRelaxationFactor (double f)
 Convenience methods that delegate to the internal smoothing filter follow below.
 
double GetRelaxationFactor ()
 Convenience methods that delegate to the internal smoothing filter follow below.
 
void SetConstraintDistance (double d)
 Convenience methods that delegate to the internal smoothing filter follow below.
 
double GetConstraintDistance ()
 Convenience methods that delegate to the internal smoothing filter follow below.
 
void SetConstraintBox (double sx, double sy, double sz)
 Convenience methods that delegate to the internal smoothing filter follow below.
 
void SetConstraintBox (double s[3])
 Convenience methods that delegate to the internal smoothing filter follow below.
 
double * GetConstraintBox ()
 Convenience methods that delegate to the internal smoothing filter follow below.
 
void GetConstraintBox (double s[3])
 Convenience methods that delegate to the internal smoothing filter follow below.
 
void SetConstraintStrategyToConstraintDistance ()
 Convenience methods that delegate to the internal smoothing filter follow below.
 
void SetConstraintStrategyToConstraintBox ()
 Convenience methods that delegate to the internal smoothing filter follow below.
 
int GetConstraintStrategy ()
 Convenience methods that delegate to the internal smoothing filter follow below.
 
virtual void SetAutomaticSmoothingConstraints (bool)
 Specify whether to set the smoothing constraints automatically.
 
virtual bool GetAutomaticSmoothingConstraints ()
 Specify whether to set the smoothing constraints automatically.
 
virtual void AutomaticSmoothingConstraintsOn ()
 Specify whether to set the smoothing constraints automatically.
 
virtual void AutomaticSmoothingConstraintsOff ()
 Specify whether to set the smoothing constraints automatically.
 
virtual void SetConstraintScale (double)
 Specify whether to set the smoothing constraints automatically.
 
virtual double GetConstraintScale ()
 Specify whether to set the smoothing constraints automatically.
 
virtual void SetOptimizedSmoothingStencils (bool)
 Indicate whether to use optimized smoothing stencils.
 
virtual bool GetOptimizedSmoothingStencils ()
 Indicate whether to use optimized smoothing stencils.
 
virtual void OptimizedSmoothingStencilsOn ()
 Indicate whether to use optimized smoothing stencils.
 
virtual void OptimizedSmoothingStencilsOff ()
 Indicate whether to use optimized smoothing stencils.
 
 vtkGetSmartPointerMacro (Smoother, vtkConstrainedSmoothingFilter)
 Get the instance of vtkConstrainedSmoothingFilter used to smooth the extracted surface net.
 
virtual void SetOutputStyle (int)
 Specify the form (i.e., the style) of the output.
 
virtual int GetOutputStyle ()
 Specify the form (i.e., the style) of the output.
 
void SetOutputStyleToDefault ()
 Specify the form (i.e., the style) of the output.
 
void SetOutputStyleToBoundary ()
 Specify the form (i.e., the style) of the output.
 
void SetOutputStyleToSelected ()
 Specify the form (i.e., the style) of the output.
 
void InitializeSelectedLabelsList ()
 When the OutputStyle is set to OUTPUT_STYLE_SELECTED, these methods are used to specify the labeled regions to output.
 
void AddSelectedLabel (double label)
 When the OutputStyle is set to OUTPUT_STYLE_SELECTED, these methods are used to specify the labeled regions to output.
 
void DeleteSelectedLabel (double label)
 When the OutputStyle is set to OUTPUT_STYLE_SELECTED, these methods are used to specify the labeled regions to output.
 
vtkIdType GetNumberOfSelectedLabels ()
 When the OutputStyle is set to OUTPUT_STYLE_SELECTED, these methods are used to specify the labeled regions to output.
 
double GetSelectedLabel (vtkIdType ithLabel)
 When the OutputStyle is set to OUTPUT_STYLE_SELECTED, these methods are used to specify the labeled regions to output.
 
virtual void SetTriangulationStrategy (int)
 Specify the strategy to triangulate the quads (not applicable if the output mesh type is set to MESH_TYPE_QUADS).
 
virtual int GetTriangulationStrategy ()
 Specify the strategy to triangulate the quads (not applicable if the output mesh type is set to MESH_TYPE_QUADS).
 
void SetTriangulationStrategyToGreedy ()
 Specify the strategy to triangulate the quads (not applicable if the output mesh type is set to MESH_TYPE_QUADS).
 
void SetTriangulationStrategyToMinEdge ()
 Specify the strategy to triangulate the quads (not applicable if the output mesh type is set to MESH_TYPE_QUADS).
 
void SetTriangulationStrategyToMinArea ()
 Specify the strategy to triangulate the quads (not applicable if the output mesh type is set to MESH_TYPE_QUADS).
 
virtual void SetDataCaching (bool)
 Enable caching of intermediate data.
 
virtual bool GetDataCaching ()
 Enable caching of intermediate data.
 
virtual void DataCachingOn ()
 Enable caching of intermediate data.
 
virtual void DataCachingOff ()
 Enable caching of intermediate data.
 
- 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.
 
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 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.
 
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 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

 vtkSurfaceNets3D ()
 
 ~vtkSurfaceNets3D () override=default
 
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.
 
bool IsCacheEmpty ()
 
void CacheData (vtkPolyData *pd, vtkCellArray *ca)
 
- 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

vtkSmartPointer< vtkContourValuesLabels
 
bool ComputeScalars
 
double BackgroundLabel
 
int ArrayComponent
 
int OutputMeshType
 
bool Smoothing
 
bool OptimizedSmoothingStencils
 
vtkSmartPointer< vtkConstrainedSmoothingFilterSmoother
 
bool AutomaticSmoothingConstraints
 
double ConstraintScale
 
bool DataCaching
 
vtkSmartPointer< vtkPolyDataGeometryCache
 
vtkSmartPointer< vtkCellArrayStencilsCache
 
vtkTimeStamp SmoothingTime
 
int OutputStyle
 
std::vector< double > SelectedLabels
 
vtkTimeStamp SelectedLabelsTime
 
int TriangulationStrategy
 
- 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, printing, and obtaining type information.
 
static vtkSurfaceNets3DNew ()
 Standard methods for instantiation, printing, and obtaining type information.
 
static vtkTypeBool IsTypeOf (const char *type)
 Standard methods for instantiation, printing, and obtaining type information.
 
static vtkSurfaceNets3DSafeDownCast (vtkObjectBase *o)
 Standard methods for instantiation, printing, and obtaining type information.
 
virtual vtkTypeBool IsA (const char *type)
 Standard methods for instantiation, printing, and obtaining type information.
 
vtkSurfaceNets3DNewInstance () const
 Standard methods for instantiation, printing, and obtaining type information.
 
void PrintSelf (ostream &os, vtkIndent indent) override
 Standard methods for instantiation, printing, and obtaining type information.
 
virtual vtkObjectBaseNewInstanceInternal () const
 Standard methods for instantiation, printing, and obtaining type information.
 

Additional Inherited Members

- 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

generate smoothed isocontours from segmented 3D image data (i.e., "label maps")

vtkSurfaceNets3D creates boundary/isocontour surfaces from a label map (e.g., a segmented image) using a threaded, 3D version of the multiple objects/labels Surface Nets algorithm. The input is a 3D image (i.e., volume) where each voxel is labeled (integer labels are preferred to real values), and the output data is a polygonal mesh separating labeled regions / objects. (Note that on output each region [corresponding to a different segmented object] will share points/edges on a common boundary, i.e., two neighboring objects will share the boundary that separates them.) This threaded implementation uses concepts from Flying Edges to achieve high performance and scalability.

The filter implements a contouring operation over a non-continuous scalar field. In comparison, classic contouring methods (like Flying Edges or Marching Cubes) presume a continuous scalar field. In comparison, this method processes non-continuous label maps, which corresponds to discrete regions in an input 3D image (i.e., volume). With a non-continuous scalar function, the usual data interpolation across a continuous function (e.g., interpolation along cell edges) is not possible. Instead, when the edge endpoint voxels are labeled in differing regions, the edge is split and transected by a quad polygon that connects the center points of the voxels on either side of the edge. Later, using a energy minimization smoothing process, the resulting polygonal mesh is adjusted to produce a smoother result. (Constraints on smoothing displacements may be specified to prevent excessive shrinkage and/or object distortion.)

The smoothing process is controlled by setting a convergence measure, the number of smoothing iterations, the step size, and the allowed (constraint) distance that points may move. These can be adjusted to provide the desired result. This class provides a method to access an internal instance of vtkConstrainedSmoothingFilter, through which these smoothing parameters may be specified, and which actually performs the smoothing operation. (Note: it is possible to skip the smoothing process altogether by disabling smoothing [e.g., invoking SmoothingOff()] or setting the number of smoothing iterations to zero. This can be useful when using a different smoothing filter like vtkWindowedSincPolyDataFilter; or if an unsmoothed, aliased output is desired. The reason the smoothing is built in to this filter is to remain faithful to the original published literature describing the Surface Nets algorithm, and for performance reasons since smoothing stencils can be generated on the fly.)

See the following reference for more details about the implementation: W. Schroeder, S. Tsalikis, M. Halle, S. Frisken. A High-Performance SurfaceNets Discrete Isocontouring Algorithm. arXiv:2401.14906. 2024. (http://arxiv.org/abs/2401.14906).

The Surface Nets algorithm was first proposed by Sarah Frisken. Two important papers include the description of surface nets for binary objects (i.e., extracting just one segmented object from a volume) and multi-label (multiple object extraction).

S. Frisken (Gibson), “Constrained Elastic SurfaceNets: Generating Smooth Surfaces from Binary Segmented Data”, Proc. MICCAI, 1998, pp. 888-898.

S. Frisken, “SurfaceNets for Multi-Label Segmentations with Preservation of Sharp Boundaries”, J. Computer Graphics Techniques, 2022.

Note that one nice feature of this filter is that algorithm execution occurs only once no matter the number of object labels / contour values. In many contouring-like algorithms, each separate contour value requires an additional algorithm execution with a new contour value. So in this filter large numbers of contour values do not significantly affect overall speed. The user can specify which objects (i.e., labels) are to be output to the filter. (Unspecified labels are treated as background and not output.)

Besides output geometry defining the surface net, the filter outputs a two-component, cell data array indicating the labels/regions on either side of the polygons composing the output vtkPolyData. (This can be used for advanced operations like extracting shared/contacting boundaries between two objects. The name of this cell data array is "BoundaryLabels".)

Note also that the content of the filter's output can be controlled by specifying the OutputStyle. This produces different output which may better serve a particular workflow. For example, it is possible to produce just exterior boundary faces, or extract selected objects/ labeled regions from the surface net.

Implementation note: For performance reasons, this filter is internally implemented quite differently than described in the literature. The main difference is that concepts from the Flying Edges parallel isocontouring algorithm are used. Namely, parallel, edge-by-edge processing is used to define cell cases, generate smoothing stencils, and produce points and output polygons. Plus the constrained smoothing process is also threaded using a double-buffering approach. For more information on Flying Edges see the paper:

"Flying Edges: A High-Performance Scalable Isocontouring Algorithm" by Schroeder, Maynard, Geveci. Proc. of LDAV 2015. Chicago, IL.

or visit VTK's FE implementation vtkFlyingEdges3D.

Warning
This filter is specialized to 3D images.
The output of this filter is a polygonal mesh. By default when smoothing is disabled, the output is quad polygons. However, once smoothing is enabled, the quads are typically decomposed into triangles since the quads are typically no longer planar. A filter option is available to force the type of output polygonal mesh (quads, or triangles).
Subtle differences in the output may result when the number of objects / labels extracted changes. This is because the smoothing operation operates on all of the boundaries simultaneously. If the boundaries change due to a difference in the number of extracted regions / labels, then the smoothing operation can produce slightly different results.
The filters vtkDiscreteMarchingCubes and vtkDiscreteFlyingEdges3D also perform contouring of label maps. However these filters produce output that may not share coincident points and/or cells, or may produce "gaps" between segmented regions. For example, vtkDiscreteMarchingCubes will share points between adjacent regions, but not triangle cells (which will be coincident). Also, no center point is inserted into voxels, meaning that intermittent gaps may form between regions. This Surface Nets implementation fully shares the boundary (points and cells) between adjacent objects; and no gaps between objects are formed (if the objects are neighbors to one another).
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. Note that for "small" volumes, serial execution may be faster due to the cost of managing threads. To force serial execution set VTK_IMPLEMENTATION_TYPE to "Sequential".
See also vtkPackLabels which is a utility class for renumbering the labels found in the input segmentation mask to contiguous forms of smaller type.
See also
vtkSurfaceNets2D vtkDiscreteMarchingCubes vtkDiscreteFlyingEdges3D vtkConstrainedSmoothingFilter vtkFlyingEdges3D vtkWindowedSincPolyDataFilter vtkPackLabels
Tests:
vtkSurfaceNets3D (Tests)

Definition at line 165 of file vtkSurfaceNets3D.h.

Member Typedef Documentation

◆ Superclass

Standard methods for instantiation, printing, and obtaining type information.

Definition at line 174 of file vtkSurfaceNets3D.h.

Member Enumeration Documentation

◆ MeshType

This enum is used to control the type of the output polygonal mesh.

Enumerator
MESH_TYPE_DEFAULT 
MESH_TYPE_TRIANGLES 
MESH_TYPE_QUADS 

Definition at line 299 of file vtkSurfaceNets3D.h.

◆ OutputType

This enum is used to control the production of the filter output.

Different output styles are used to transform the data so they can be used in different workflows, providing tradeoffs between speed, memory, and auxiliary information. By default (OUTPUT_STYLE_DEFAULT) the filter produces a mesh with shared points (i.e., points are not duplicated), and all mesh polygons, both interior and exterior, are produced. OUTPUT_STYLE_BOUNDARY is similar to OUTPUT_STYLE_DEFAULT except that only mesh polygons that are on the boundary are produced (i.e., only polygons that border the background region) - thus no interior polygons are produced. OUTPUT_STYLE_SELECTED is used to extract faces bounding selected regions.

Enumerator
OUTPUT_STYLE_DEFAULT 
OUTPUT_STYLE_BOUNDARY 
OUTPUT_STYLE_SELECTED 

Definition at line 429 of file vtkSurfaceNets3D.h.

◆ TriangulationType

This enum is used to control how quadrilaterals are triangulated.

Enumerator
TRIANGULATION_GREEDY 
TRIANGULATION_MIN_EDGE 
TRIANGULATION_MIN_AREA 

Definition at line 471 of file vtkSurfaceNets3D.h.

Constructor & Destructor Documentation

◆ vtkSurfaceNets3D()

vtkSurfaceNets3D::vtkSurfaceNets3D ( )
protected

◆ ~vtkSurfaceNets3D()

vtkSurfaceNets3D::~vtkSurfaceNets3D ( )
overrideprotecteddefault

Member Function Documentation

◆ New()

static vtkSurfaceNets3D * vtkSurfaceNets3D::New ( )
static

Standard methods for instantiation, printing, and obtaining type information.

◆ IsTypeOf()

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

Standard methods for instantiation, printing, and obtaining type information.

◆ IsA()

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

Standard methods for instantiation, printing, and obtaining type information.

Reimplemented from vtkPolyDataAlgorithm.

◆ SafeDownCast()

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

Standard methods for instantiation, printing, and obtaining type information.

◆ NewInstanceInternal()

virtual vtkObjectBase * vtkSurfaceNets3D::NewInstanceInternal ( ) const
protectedvirtual

Standard methods for instantiation, printing, and obtaining type information.

Reimplemented from vtkPolyDataAlgorithm.

◆ NewInstance()

vtkSurfaceNets3D * vtkSurfaceNets3D::NewInstance ( ) const

Standard methods for instantiation, printing, and obtaining type information.

◆ PrintSelf()

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

Standard methods for instantiation, printing, and obtaining type information.

Reimplemented from vtkPolyDataAlgorithm.

◆ GetMTime()

vtkMTimeType vtkSurfaceNets3D::GetMTime ( )
overridevirtual

The modified time is also a function of the label values and the smoothing filter.

Reimplemented from vtkObject.

◆ SetValue()

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

Set a particular label value at label number i.

The index i ranges between (0 <= i < NumberOfLabels). (Note: while labels values are expressed as doubles, the underlying scalar data may be a different type. During execution the label values are cast to the type of the scalar data.) Note the use of "Value" and "Label" when specifying regions to extract. The use of "Value" is consistent with other VTK continuous-scalar field isocontouring algorithms; however the term "Label" is more consistent with label maps. Warning: make sure that the value of the background label (see definition below) is different than any of the specified labels, otherwise the generated cell scalars may be incorrect.

Definition at line 198 of file vtkSurfaceNets3D.h.

◆ SetLabel()

void vtkSurfaceNets3D::SetLabel ( int  i,
double  value 
)
inline

Set a particular label value at label number i.

The index i ranges between (0 <= i < NumberOfLabels). (Note: while labels values are expressed as doubles, the underlying scalar data may be a different type. During execution the label values are cast to the type of the scalar data.) Note the use of "Value" and "Label" when specifying regions to extract. The use of "Value" is consistent with other VTK continuous-scalar field isocontouring algorithms; however the term "Label" is more consistent with label maps. Warning: make sure that the value of the background label (see definition below) is different than any of the specified labels, otherwise the generated cell scalars may be incorrect.

Definition at line 199 of file vtkSurfaceNets3D.h.

◆ GetValue()

double vtkSurfaceNets3D::GetValue ( int  i)
inline

Get the ith label value.

Definition at line 206 of file vtkSurfaceNets3D.h.

◆ GetLabel()

double vtkSurfaceNets3D::GetLabel ( int  i)
inline

Get the ith label value.

Definition at line 207 of file vtkSurfaceNets3D.h.

◆ GetValues() [1/2]

double * vtkSurfaceNets3D::GetValues ( )
inline

Get a pointer to an array of labels.

There will be GetNumberOfLabels() values in the list.

Definition at line 215 of file vtkSurfaceNets3D.h.

◆ GetLabels() [1/2]

double * vtkSurfaceNets3D::GetLabels ( )
inline

Get a pointer to an array of labels.

There will be GetNumberOfLabels() values in the list.

Definition at line 216 of file vtkSurfaceNets3D.h.

◆ GetValues() [2/2]

void vtkSurfaceNets3D::GetValues ( double *  contourValues)
inline

Fill a supplied list with label values.

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

Definition at line 225 of file vtkSurfaceNets3D.h.

◆ GetLabels() [2/2]

void vtkSurfaceNets3D::GetLabels ( double *  contourValues)
inline

Fill a supplied list with label values.

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

Definition at line 226 of file vtkSurfaceNets3D.h.

◆ SetNumberOfLabels()

void vtkSurfaceNets3D::SetNumberOfLabels ( int  number)
inline

Set the number of labels 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. Note that for consistency with other isocountoring-related algorithms, some methods use "Labels" and "Contours" interchangeably.

Definition at line 237 of file vtkSurfaceNets3D.h.

◆ SetNumberOfContours()

void vtkSurfaceNets3D::SetNumberOfContours ( int  number)
inline

Set the number of labels 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. Note that for consistency with other isocountoring-related algorithms, some methods use "Labels" and "Contours" interchangeably.

Definition at line 238 of file vtkSurfaceNets3D.h.

◆ GetNumberOfLabels()

vtkIdType vtkSurfaceNets3D::GetNumberOfLabels ( )
inline

Get the number of labels in the list of label values.

Definition at line 245 of file vtkSurfaceNets3D.h.

◆ GetNumberOfContours()

vtkIdType vtkSurfaceNets3D::GetNumberOfContours ( )
inline

Get the number of labels in the list of label values.

Definition at line 246 of file vtkSurfaceNets3D.h.

◆ GenerateLabels() [1/2]

void vtkSurfaceNets3D::GenerateLabels ( int  numLabels,
double  range[2] 
)
inline

Generate numLabels equally spaced labels between the specified range.

The labels will include the min/max range values.

Definition at line 254 of file vtkSurfaceNets3D.h.

◆ GenerateValues() [1/2]

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

Generate numLabels equally spaced labels between the specified range.

The labels will include the min/max range values.

Definition at line 258 of file vtkSurfaceNets3D.h.

◆ GenerateLabels() [2/2]

void vtkSurfaceNets3D::GenerateLabels ( int  numLabels,
double  rangeStart,
double  rangeEnd 
)
inline

Generate numLabels equally spaced labels between the specified range.

The labels will include the min/max range values.

Definition at line 262 of file vtkSurfaceNets3D.h.

◆ GenerateValues() [2/2]

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

Generate numLabels equally spaced labels between the specified range.

The labels will include the min/max range values.

Definition at line 266 of file vtkSurfaceNets3D.h.

◆ SetBackgroundLabel()

virtual void vtkSurfaceNets3D::SetBackgroundLabel ( double  )
virtual

This value specifies the label value to use when referencing the background region outside of any of the specified regions.

(This value is used when producing cell scalars.) By default this value is zero. Be very careful of the value being used here, it should not overlap an extracted label value, and because it is the same type as the input image scalars, make sure the value can be properly represented (e.g., if the input scalars are an unsigned type, then BackgroundLabel should not be negative).

◆ GetBackgroundLabel()

virtual double vtkSurfaceNets3D::GetBackgroundLabel ( )
virtual

This value specifies the label value to use when referencing the background region outside of any of the specified regions.

(This value is used when producing cell scalars.) By default this value is zero. Be very careful of the value being used here, it should not overlap an extracted label value, and because it is the same type as the input image scalars, make sure the value can be properly represented (e.g., if the input scalars are an unsigned type, then BackgroundLabel should not be negative).

◆ SetArrayComponent()

virtual void vtkSurfaceNets3D::SetArrayComponent ( int  )
virtual

Set/get which component of a input multi-component scalar array to contour with; defaults to component 0.

◆ GetArrayComponent()

virtual int vtkSurfaceNets3D::GetArrayComponent ( )
virtual

Set/get which component of a input multi-component scalar array to contour with; defaults to component 0.

◆ SetOutputMeshType()

virtual void vtkSurfaceNets3D::SetOutputMeshType ( int  )
virtual

Control the type of output mesh.

By default, if smoothing is off, the output mesh is a polygonal mesh consisting of quadrilaterals (quads). However, if smoothing is enabled, then the output mesh type is a polygonal mesh consisting of triangles. It is possible to force the output mesh type to be of a certain type (triangles, or quads) regardless whether smoothing is enabled or not. Note that if an output mesh is forced to be quads, and smoothing is enabled, the resulting quads may not be planar.

◆ GetOutputMeshType()

virtual int vtkSurfaceNets3D::GetOutputMeshType ( )
virtual

Control the type of output mesh.

By default, if smoothing is off, the output mesh is a polygonal mesh consisting of quadrilaterals (quads). However, if smoothing is enabled, then the output mesh type is a polygonal mesh consisting of triangles. It is possible to force the output mesh type to be of a certain type (triangles, or quads) regardless whether smoothing is enabled or not. Note that if an output mesh is forced to be quads, and smoothing is enabled, the resulting quads may not be planar.

◆ SetOutputMeshTypeToDefault()

void vtkSurfaceNets3D::SetOutputMeshTypeToDefault ( )
inline

Control the type of output mesh.

By default, if smoothing is off, the output mesh is a polygonal mesh consisting of quadrilaterals (quads). However, if smoothing is enabled, then the output mesh type is a polygonal mesh consisting of triangles. It is possible to force the output mesh type to be of a certain type (triangles, or quads) regardless whether smoothing is enabled or not. Note that if an output mesh is forced to be quads, and smoothing is enabled, the resulting quads may not be planar.

Definition at line 319 of file vtkSurfaceNets3D.h.

◆ SetOutputMeshTypeToTriangles()

void vtkSurfaceNets3D::SetOutputMeshTypeToTriangles ( )
inline

Control the type of output mesh.

By default, if smoothing is off, the output mesh is a polygonal mesh consisting of quadrilaterals (quads). However, if smoothing is enabled, then the output mesh type is a polygonal mesh consisting of triangles. It is possible to force the output mesh type to be of a certain type (triangles, or quads) regardless whether smoothing is enabled or not. Note that if an output mesh is forced to be quads, and smoothing is enabled, the resulting quads may not be planar.

Definition at line 320 of file vtkSurfaceNets3D.h.

◆ SetOutputMeshTypeToQuads()

void vtkSurfaceNets3D::SetOutputMeshTypeToQuads ( )
inline

Control the type of output mesh.

By default, if smoothing is off, the output mesh is a polygonal mesh consisting of quadrilaterals (quads). However, if smoothing is enabled, then the output mesh type is a polygonal mesh consisting of triangles. It is possible to force the output mesh type to be of a certain type (triangles, or quads) regardless whether smoothing is enabled or not. Note that if an output mesh is forced to be quads, and smoothing is enabled, the resulting quads may not be planar.

Definition at line 321 of file vtkSurfaceNets3D.h.

◆ SetSmoothing()

virtual void vtkSurfaceNets3D::SetSmoothing ( bool  )
virtual

Indicate whether smoothing should be enabled.

By default, after the surface net is extracted, smoothing occurs using the built-in smoother. To disable smoothing, invoke SmoothingOff().

◆ GetSmoothing()

virtual bool vtkSurfaceNets3D::GetSmoothing ( )
virtual

Indicate whether smoothing should be enabled.

By default, after the surface net is extracted, smoothing occurs using the built-in smoother. To disable smoothing, invoke SmoothingOff().

◆ SmoothingOn()

virtual void vtkSurfaceNets3D::SmoothingOn ( )
virtual

Indicate whether smoothing should be enabled.

By default, after the surface net is extracted, smoothing occurs using the built-in smoother. To disable smoothing, invoke SmoothingOff().

◆ SmoothingOff()

virtual void vtkSurfaceNets3D::SmoothingOff ( )
virtual

Indicate whether smoothing should be enabled.

By default, after the surface net is extracted, smoothing occurs using the built-in smoother. To disable smoothing, invoke SmoothingOff().

◆ SetNumberOfIterations()

void vtkSurfaceNets3D::SetNumberOfIterations ( int  n)
inline

Convenience methods that delegate to the internal smoothing filter follow below.

See the documentation for vtkConstrainedSmoothingAlgorithm for more information.

Definition at line 346 of file vtkSurfaceNets3D.h.

◆ GetNumberOfIterations()

int vtkSurfaceNets3D::GetNumberOfIterations ( )
inline

Convenience methods that delegate to the internal smoothing filter follow below.

See the documentation for vtkConstrainedSmoothingAlgorithm for more information.

Definition at line 347 of file vtkSurfaceNets3D.h.

◆ SetRelaxationFactor()

void vtkSurfaceNets3D::SetRelaxationFactor ( double  f)
inline

Convenience methods that delegate to the internal smoothing filter follow below.

See the documentation for vtkConstrainedSmoothingAlgorithm for more information.

Definition at line 348 of file vtkSurfaceNets3D.h.

◆ GetRelaxationFactor()

double vtkSurfaceNets3D::GetRelaxationFactor ( )
inline

Convenience methods that delegate to the internal smoothing filter follow below.

See the documentation for vtkConstrainedSmoothingAlgorithm for more information.

Definition at line 349 of file vtkSurfaceNets3D.h.

◆ SetConstraintDistance()

void vtkSurfaceNets3D::SetConstraintDistance ( double  d)
inline

Convenience methods that delegate to the internal smoothing filter follow below.

See the documentation for vtkConstrainedSmoothingAlgorithm for more information.

Definition at line 350 of file vtkSurfaceNets3D.h.

◆ GetConstraintDistance()

double vtkSurfaceNets3D::GetConstraintDistance ( )
inline

Convenience methods that delegate to the internal smoothing filter follow below.

See the documentation for vtkConstrainedSmoothingAlgorithm for more information.

Definition at line 351 of file vtkSurfaceNets3D.h.

◆ SetConstraintBox() [1/2]

void vtkSurfaceNets3D::SetConstraintBox ( double  sx,
double  sy,
double  sz 
)
inline

Convenience methods that delegate to the internal smoothing filter follow below.

See the documentation for vtkConstrainedSmoothingAlgorithm for more information.

Definition at line 352 of file vtkSurfaceNets3D.h.

◆ SetConstraintBox() [2/2]

void vtkSurfaceNets3D::SetConstraintBox ( double  s[3])
inline

Convenience methods that delegate to the internal smoothing filter follow below.

See the documentation for vtkConstrainedSmoothingAlgorithm for more information.

Definition at line 356 of file vtkSurfaceNets3D.h.

◆ GetConstraintBox() [1/2]

double * vtkSurfaceNets3D::GetConstraintBox ( )
inline

Convenience methods that delegate to the internal smoothing filter follow below.

See the documentation for vtkConstrainedSmoothingAlgorithm for more information.

Definition at line 357 of file vtkSurfaceNets3D.h.

◆ GetConstraintBox() [2/2]

void vtkSurfaceNets3D::GetConstraintBox ( double  s[3])
inline

Convenience methods that delegate to the internal smoothing filter follow below.

See the documentation for vtkConstrainedSmoothingAlgorithm for more information.

Definition at line 358 of file vtkSurfaceNets3D.h.

◆ SetConstraintStrategyToConstraintDistance()

void vtkSurfaceNets3D::SetConstraintStrategyToConstraintDistance ( )
inline

Convenience methods that delegate to the internal smoothing filter follow below.

See the documentation for vtkConstrainedSmoothingAlgorithm for more information.

Definition at line 359 of file vtkSurfaceNets3D.h.

◆ SetConstraintStrategyToConstraintBox()

void vtkSurfaceNets3D::SetConstraintStrategyToConstraintBox ( )
inline

Convenience methods that delegate to the internal smoothing filter follow below.

See the documentation for vtkConstrainedSmoothingAlgorithm for more information.

Definition at line 363 of file vtkSurfaceNets3D.h.

◆ GetConstraintStrategy()

int vtkSurfaceNets3D::GetConstraintStrategy ( )
inline

Convenience methods that delegate to the internal smoothing filter follow below.

See the documentation for vtkConstrainedSmoothingAlgorithm for more information.

Definition at line 367 of file vtkSurfaceNets3D.h.

◆ SetAutomaticSmoothingConstraints()

virtual void vtkSurfaceNets3D::SetAutomaticSmoothingConstraints ( bool  )
virtual

Specify whether to set the smoothing constraints automatically.

If automatic is on, the constraint distance and constraint box will calculated and set (based on the input size of the volume voxel). Note that the ConstraintScale is used to adjust the size of the constraint distance or box when set automatically. (Typically the constraint distance defines a circumscribing sphere around a voxel, and the constraint box is a box with voxel spacing.) If constraints are not set automatically, then the constraint distance and/or constraint box should be set manually.) By default, automatic smoothing constraints are enabled.

◆ GetAutomaticSmoothingConstraints()

virtual bool vtkSurfaceNets3D::GetAutomaticSmoothingConstraints ( )
virtual

Specify whether to set the smoothing constraints automatically.

If automatic is on, the constraint distance and constraint box will calculated and set (based on the input size of the volume voxel). Note that the ConstraintScale is used to adjust the size of the constraint distance or box when set automatically. (Typically the constraint distance defines a circumscribing sphere around a voxel, and the constraint box is a box with voxel spacing.) If constraints are not set automatically, then the constraint distance and/or constraint box should be set manually.) By default, automatic smoothing constraints are enabled.

◆ AutomaticSmoothingConstraintsOn()

virtual void vtkSurfaceNets3D::AutomaticSmoothingConstraintsOn ( )
virtual

Specify whether to set the smoothing constraints automatically.

If automatic is on, the constraint distance and constraint box will calculated and set (based on the input size of the volume voxel). Note that the ConstraintScale is used to adjust the size of the constraint distance or box when set automatically. (Typically the constraint distance defines a circumscribing sphere around a voxel, and the constraint box is a box with voxel spacing.) If constraints are not set automatically, then the constraint distance and/or constraint box should be set manually.) By default, automatic smoothing constraints are enabled.

◆ AutomaticSmoothingConstraintsOff()

virtual void vtkSurfaceNets3D::AutomaticSmoothingConstraintsOff ( )
virtual

Specify whether to set the smoothing constraints automatically.

If automatic is on, the constraint distance and constraint box will calculated and set (based on the input size of the volume voxel). Note that the ConstraintScale is used to adjust the size of the constraint distance or box when set automatically. (Typically the constraint distance defines a circumscribing sphere around a voxel, and the constraint box is a box with voxel spacing.) If constraints are not set automatically, then the constraint distance and/or constraint box should be set manually.) By default, automatic smoothing constraints are enabled.

◆ SetConstraintScale()

virtual void vtkSurfaceNets3D::SetConstraintScale ( double  )
virtual

Specify whether to set the smoothing constraints automatically.

If automatic is on, the constraint distance and constraint box will calculated and set (based on the input size of the volume voxel). Note that the ConstraintScale is used to adjust the size of the constraint distance or box when set automatically. (Typically the constraint distance defines a circumscribing sphere around a voxel, and the constraint box is a box with voxel spacing.) If constraints are not set automatically, then the constraint distance and/or constraint box should be set manually.) By default, automatic smoothing constraints are enabled.

◆ GetConstraintScale()

virtual double vtkSurfaceNets3D::GetConstraintScale ( )
virtual

Specify whether to set the smoothing constraints automatically.

If automatic is on, the constraint distance and constraint box will calculated and set (based on the input size of the volume voxel). Note that the ConstraintScale is used to adjust the size of the constraint distance or box when set automatically. (Typically the constraint distance defines a circumscribing sphere around a voxel, and the constraint box is a box with voxel spacing.) If constraints are not set automatically, then the constraint distance and/or constraint box should be set manually.) By default, automatic smoothing constraints are enabled.

◆ SetOptimizedSmoothingStencils()

virtual void vtkSurfaceNets3D::SetOptimizedSmoothingStencils ( bool  )
virtual

Indicate whether to use optimized smoothing stencils.

Optimized stencils (which are on by default) are designed to better smooth sharp edges across the surface net. In some cases it may be desired to disable the use of optimized smoothing stencils.

◆ GetOptimizedSmoothingStencils()

virtual bool vtkSurfaceNets3D::GetOptimizedSmoothingStencils ( )
virtual

Indicate whether to use optimized smoothing stencils.

Optimized stencils (which are on by default) are designed to better smooth sharp edges across the surface net. In some cases it may be desired to disable the use of optimized smoothing stencils.

◆ OptimizedSmoothingStencilsOn()

virtual void vtkSurfaceNets3D::OptimizedSmoothingStencilsOn ( )
virtual

Indicate whether to use optimized smoothing stencils.

Optimized stencils (which are on by default) are designed to better smooth sharp edges across the surface net. In some cases it may be desired to disable the use of optimized smoothing stencils.

◆ OptimizedSmoothingStencilsOff()

virtual void vtkSurfaceNets3D::OptimizedSmoothingStencilsOff ( )
virtual

Indicate whether to use optimized smoothing stencils.

Optimized stencils (which are on by default) are designed to better smooth sharp edges across the surface net. In some cases it may be desired to disable the use of optimized smoothing stencils.

◆ vtkGetSmartPointerMacro()

vtkSurfaceNets3D::vtkGetSmartPointerMacro ( Smoother  ,
vtkConstrainedSmoothingFilter   
)

Get the instance of vtkConstrainedSmoothingFilter used to smooth the extracted surface net.

To control smoothing, access this instance and specify its parameters such as number of smoothing iterations and constraint distance. If you wish to disable smoothing, set SmoothingOff().

◆ SetOutputStyle()

virtual void vtkSurfaceNets3D::SetOutputStyle ( int  )
virtual

Specify the form (i.e., the style) of the output.

Different styles are meant to support different workflows. OUTPUT_STYLE_DEFAULT provides the basic information defining the output surface net. OUTPUT_STYLE_BOUNDARY produces much smaller output since the interior polygon faces are not produced. Finally, OUTPUT_STYLE_SELECTED enables the user to extract a subset of the labeled regions. This is useful because the smoothing operation will occur across all the specified input regions, meaning that the selected regions do not change shape due to changes in the specified input regions. You must specify the selected regions (i.e., labels) to output.

◆ GetOutputStyle()

virtual int vtkSurfaceNets3D::GetOutputStyle ( )
virtual

Specify the form (i.e., the style) of the output.

Different styles are meant to support different workflows. OUTPUT_STYLE_DEFAULT provides the basic information defining the output surface net. OUTPUT_STYLE_BOUNDARY produces much smaller output since the interior polygon faces are not produced. Finally, OUTPUT_STYLE_SELECTED enables the user to extract a subset of the labeled regions. This is useful because the smoothing operation will occur across all the specified input regions, meaning that the selected regions do not change shape due to changes in the specified input regions. You must specify the selected regions (i.e., labels) to output.

◆ SetOutputStyleToDefault()

void vtkSurfaceNets3D::SetOutputStyleToDefault ( )
inline

Specify the form (i.e., the style) of the output.

Different styles are meant to support different workflows. OUTPUT_STYLE_DEFAULT provides the basic information defining the output surface net. OUTPUT_STYLE_BOUNDARY produces much smaller output since the interior polygon faces are not produced. Finally, OUTPUT_STYLE_SELECTED enables the user to extract a subset of the labeled regions. This is useful because the smoothing operation will occur across all the specified input regions, meaning that the selected regions do not change shape due to changes in the specified input regions. You must specify the selected regions (i.e., labels) to output.

Definition at line 451 of file vtkSurfaceNets3D.h.

◆ SetOutputStyleToBoundary()

void vtkSurfaceNets3D::SetOutputStyleToBoundary ( )
inline

Specify the form (i.e., the style) of the output.

Different styles are meant to support different workflows. OUTPUT_STYLE_DEFAULT provides the basic information defining the output surface net. OUTPUT_STYLE_BOUNDARY produces much smaller output since the interior polygon faces are not produced. Finally, OUTPUT_STYLE_SELECTED enables the user to extract a subset of the labeled regions. This is useful because the smoothing operation will occur across all the specified input regions, meaning that the selected regions do not change shape due to changes in the specified input regions. You must specify the selected regions (i.e., labels) to output.

Definition at line 452 of file vtkSurfaceNets3D.h.

◆ SetOutputStyleToSelected()

void vtkSurfaceNets3D::SetOutputStyleToSelected ( )
inline

Specify the form (i.e., the style) of the output.

Different styles are meant to support different workflows. OUTPUT_STYLE_DEFAULT provides the basic information defining the output surface net. OUTPUT_STYLE_BOUNDARY produces much smaller output since the interior polygon faces are not produced. Finally, OUTPUT_STYLE_SELECTED enables the user to extract a subset of the labeled regions. This is useful because the smoothing operation will occur across all the specified input regions, meaning that the selected regions do not change shape due to changes in the specified input regions. You must specify the selected regions (i.e., labels) to output.

Definition at line 453 of file vtkSurfaceNets3D.h.

◆ InitializeSelectedLabelsList()

void vtkSurfaceNets3D::InitializeSelectedLabelsList ( )

When the OutputStyle is set to OUTPUT_STYLE_SELECTED, these methods are used to specify the labeled regions to output.

◆ AddSelectedLabel()

void vtkSurfaceNets3D::AddSelectedLabel ( double  label)

When the OutputStyle is set to OUTPUT_STYLE_SELECTED, these methods are used to specify the labeled regions to output.

◆ DeleteSelectedLabel()

void vtkSurfaceNets3D::DeleteSelectedLabel ( double  label)

When the OutputStyle is set to OUTPUT_STYLE_SELECTED, these methods are used to specify the labeled regions to output.

◆ GetNumberOfSelectedLabels()

vtkIdType vtkSurfaceNets3D::GetNumberOfSelectedLabels ( )

When the OutputStyle is set to OUTPUT_STYLE_SELECTED, these methods are used to specify the labeled regions to output.

◆ GetSelectedLabel()

double vtkSurfaceNets3D::GetSelectedLabel ( vtkIdType  ithLabel)

When the OutputStyle is set to OUTPUT_STYLE_SELECTED, these methods are used to specify the labeled regions to output.

◆ SetTriangulationStrategy()

virtual void vtkSurfaceNets3D::SetTriangulationStrategy ( int  )
virtual

Specify the strategy to triangulate the quads (not applicable if the output mesh type is set to MESH_TYPE_QUADS).

If TRIANGULATE_GREEDY is specified, then quads are triangulated in no particular order. If TRIANGULATED_MIN_EDGE is specified, then trianglate the quad using a minimum-edge-length diagonal. If TRIANGULATED_MIN_AREA is specified, then trianglate the quad to produce a minimum surface area. By default, TRIANGULATE_MIN_EDGE is used. (Slight performance affects may occur, with TRIANGULATION_GREEDY generally the fastest.)

◆ GetTriangulationStrategy()

virtual int vtkSurfaceNets3D::GetTriangulationStrategy ( )
virtual

Specify the strategy to triangulate the quads (not applicable if the output mesh type is set to MESH_TYPE_QUADS).

If TRIANGULATE_GREEDY is specified, then quads are triangulated in no particular order. If TRIANGULATED_MIN_EDGE is specified, then trianglate the quad using a minimum-edge-length diagonal. If TRIANGULATED_MIN_AREA is specified, then trianglate the quad to produce a minimum surface area. By default, TRIANGULATE_MIN_EDGE is used. (Slight performance affects may occur, with TRIANGULATION_GREEDY generally the fastest.)

◆ SetTriangulationStrategyToGreedy()

void vtkSurfaceNets3D::SetTriangulationStrategyToGreedy ( )
inline

Specify the strategy to triangulate the quads (not applicable if the output mesh type is set to MESH_TYPE_QUADS).

If TRIANGULATE_GREEDY is specified, then quads are triangulated in no particular order. If TRIANGULATED_MIN_EDGE is specified, then trianglate the quad using a minimum-edge-length diagonal. If TRIANGULATED_MIN_AREA is specified, then trianglate the quad to produce a minimum surface area. By default, TRIANGULATE_MIN_EDGE is used. (Slight performance affects may occur, with TRIANGULATION_GREEDY generally the fastest.)

Definition at line 491 of file vtkSurfaceNets3D.h.

◆ SetTriangulationStrategyToMinEdge()

void vtkSurfaceNets3D::SetTriangulationStrategyToMinEdge ( )
inline

Specify the strategy to triangulate the quads (not applicable if the output mesh type is set to MESH_TYPE_QUADS).

If TRIANGULATE_GREEDY is specified, then quads are triangulated in no particular order. If TRIANGULATED_MIN_EDGE is specified, then trianglate the quad using a minimum-edge-length diagonal. If TRIANGULATED_MIN_AREA is specified, then trianglate the quad to produce a minimum surface area. By default, TRIANGULATE_MIN_EDGE is used. (Slight performance affects may occur, with TRIANGULATION_GREEDY generally the fastest.)

Definition at line 492 of file vtkSurfaceNets3D.h.

◆ SetTriangulationStrategyToMinArea()

void vtkSurfaceNets3D::SetTriangulationStrategyToMinArea ( )
inline

Specify the strategy to triangulate the quads (not applicable if the output mesh type is set to MESH_TYPE_QUADS).

If TRIANGULATE_GREEDY is specified, then quads are triangulated in no particular order. If TRIANGULATED_MIN_EDGE is specified, then trianglate the quad using a minimum-edge-length diagonal. If TRIANGULATED_MIN_AREA is specified, then trianglate the quad to produce a minimum surface area. By default, TRIANGULATE_MIN_EDGE is used. (Slight performance affects may occur, with TRIANGULATION_GREEDY generally the fastest.)

Definition at line 496 of file vtkSurfaceNets3D.h.

◆ SetDataCaching()

virtual void vtkSurfaceNets3D::SetDataCaching ( bool  )
virtual

Enable caching of intermediate data.

A common workflow using this filter requires extracting object boundaries (i.e., the isocontour), and then repeatedly rerunning the smoothing process with different parameters. To improve performance by avoiding repeated extraction of the boundary, the filter can cache intermediate data prior to the smoothing process. In this way, the boundary is only extracted once, and as long as only the internal constrained smoothing filter is modified, then boundary extraction will not be reexecuted. By default this is enabled.

◆ GetDataCaching()

virtual bool vtkSurfaceNets3D::GetDataCaching ( )
virtual

Enable caching of intermediate data.

A common workflow using this filter requires extracting object boundaries (i.e., the isocontour), and then repeatedly rerunning the smoothing process with different parameters. To improve performance by avoiding repeated extraction of the boundary, the filter can cache intermediate data prior to the smoothing process. In this way, the boundary is only extracted once, and as long as only the internal constrained smoothing filter is modified, then boundary extraction will not be reexecuted. By default this is enabled.

◆ DataCachingOn()

virtual void vtkSurfaceNets3D::DataCachingOn ( )
virtual

Enable caching of intermediate data.

A common workflow using this filter requires extracting object boundaries (i.e., the isocontour), and then repeatedly rerunning the smoothing process with different parameters. To improve performance by avoiding repeated extraction of the boundary, the filter can cache intermediate data prior to the smoothing process. In this way, the boundary is only extracted once, and as long as only the internal constrained smoothing filter is modified, then boundary extraction will not be reexecuted. By default this is enabled.

◆ DataCachingOff()

virtual void vtkSurfaceNets3D::DataCachingOff ( )
virtual

Enable caching of intermediate data.

A common workflow using this filter requires extracting object boundaries (i.e., the isocontour), and then repeatedly rerunning the smoothing process with different parameters. To improve performance by avoiding repeated extraction of the boundary, the filter can cache intermediate data prior to the smoothing process. In this way, the boundary is only extracted once, and as long as only the internal constrained smoothing filter is modified, then boundary extraction will not be reexecuted. By default this is enabled.

◆ RequestData()

int vtkSurfaceNets3D::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 vtkSurfaceNets3D::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 vtkPolyDataAlgorithm.

◆ IsCacheEmpty()

bool vtkSurfaceNets3D::IsCacheEmpty ( )
protected

◆ CacheData()

void vtkSurfaceNets3D::CacheData ( vtkPolyData pd,
vtkCellArray ca 
)
protected

Member Data Documentation

◆ Labels

vtkSmartPointer<vtkContourValues> vtkSurfaceNets3D::Labels
protected

Definition at line 527 of file vtkSurfaceNets3D.h.

◆ ComputeScalars

bool vtkSurfaceNets3D::ComputeScalars
protected

Definition at line 528 of file vtkSurfaceNets3D.h.

◆ BackgroundLabel

double vtkSurfaceNets3D::BackgroundLabel
protected

Definition at line 529 of file vtkSurfaceNets3D.h.

◆ ArrayComponent

int vtkSurfaceNets3D::ArrayComponent
protected

Definition at line 530 of file vtkSurfaceNets3D.h.

◆ OutputMeshType

int vtkSurfaceNets3D::OutputMeshType
protected

Definition at line 531 of file vtkSurfaceNets3D.h.

◆ Smoothing

bool vtkSurfaceNets3D::Smoothing
protected

Definition at line 534 of file vtkSurfaceNets3D.h.

◆ OptimizedSmoothingStencils

bool vtkSurfaceNets3D::OptimizedSmoothingStencils
protected

Definition at line 535 of file vtkSurfaceNets3D.h.

◆ Smoother

vtkSmartPointer<vtkConstrainedSmoothingFilter> vtkSurfaceNets3D::Smoother
protected

Definition at line 536 of file vtkSurfaceNets3D.h.

◆ AutomaticSmoothingConstraints

bool vtkSurfaceNets3D::AutomaticSmoothingConstraints
protected

Definition at line 537 of file vtkSurfaceNets3D.h.

◆ ConstraintScale

double vtkSurfaceNets3D::ConstraintScale
protected

Definition at line 538 of file vtkSurfaceNets3D.h.

◆ DataCaching

bool vtkSurfaceNets3D::DataCaching
protected

Definition at line 543 of file vtkSurfaceNets3D.h.

◆ GeometryCache

vtkSmartPointer<vtkPolyData> vtkSurfaceNets3D::GeometryCache
protected

Definition at line 544 of file vtkSurfaceNets3D.h.

◆ StencilsCache

vtkSmartPointer<vtkCellArray> vtkSurfaceNets3D::StencilsCache
protected

Definition at line 545 of file vtkSurfaceNets3D.h.

◆ SmoothingTime

vtkTimeStamp vtkSurfaceNets3D::SmoothingTime
protected

Definition at line 546 of file vtkSurfaceNets3D.h.

◆ OutputStyle

int vtkSurfaceNets3D::OutputStyle
protected

Definition at line 551 of file vtkSurfaceNets3D.h.

◆ SelectedLabels

std::vector<double> vtkSurfaceNets3D::SelectedLabels
protected

Definition at line 552 of file vtkSurfaceNets3D.h.

◆ SelectedLabelsTime

vtkTimeStamp vtkSurfaceNets3D::SelectedLabelsTime
protected

Definition at line 553 of file vtkSurfaceNets3D.h.

◆ TriangulationStrategy

int vtkSurfaceNets3D::TriangulationStrategy
protected

Definition at line 556 of file vtkSurfaceNets3D.h.


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