/*========================================================================= Program: Visualization Toolkit Module: vtkSurfaceNets2D.h Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen All rights reserved. See Copyright.txt or http://www.kitware.com/Copyright.htm for details. This software is distributed WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the above copyright notice for more information. =========================================================================*/ /** * @class vtkSurfaceNets2D * @brief generate smoothed constours from segmented 2D image data (i.e., "label maps") * * vtkSurfaceNets2D creates boundary/isocontour lines from a label map (e.g., * a segmented image) using a threaded, 2D version of the multiple * regions/labels surface nets algorithm. The input is a 2D image where each * pixel is labeled (integer labels are preferred to real values), and the * output data is polygonal data separating labeled regions via line * segments. (Note that on output each region [corresponding to a different * segmented object] will share points/edges on a common boundary, i.e., two * objects next to each other will share the boundary that separates them.) * * While this filter is similar to a contouring operation, classic contouring * methods assume a continuous scalar field. In comparision, label maps are * not continuous in scalar function value, meaning that usual data * interpolation (e.g., along edges) is not possible. Instead, when the edge * endpoint pixels are labeled in differing regions, the edge is split and * transected by a line segment that connects the center points of the squares * on either side of the edge. Later, using a energy minimization smoothing * process, these split edges will be 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 published literature describing the surface nets * algorithm.) * * The SurfaceNets 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.) * * The filter can optionally output a two-component, cell data array * indicating the labels/regions on either side of the line segments * 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". * * 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 lines. The smoothing process is also threaded using a * double-buffering approach. * * @warning * This filter is specialized to 2D images. * * @warning * Subtle differences in the output may result when the number of objects / * labels 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 regions / labels, then the smoothing operation * can produce different results. * * @warning * The filters vtkDiscreteMarchingCubes, vtkDiscreteFlyingEdges2D, * vtkDiscreteFlyingEdges3D, and vtkDiscreteFlyingEdgesClipper2D also * perform isocontour extraction. However these filters produce output * that may not share common boundary cells, and may produce "gaps" * between segmented regions. For example, vtkDiscreteMarchingCubes will * share points between adjacent regions, but not triangle cells (which * will be conincident). Also, no center point is inserted into voxels, * meaning that intermittent gaps may form between regions. * * @warning * 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. * * @sa * vtkSurfaceNets3D vtkDiscreteFlyingEdges2D vtkDiscreteFlyingEdgesClipper2D * vtkConstrainedSmoothingFilter vtkFlyingEdges2D vtkFlyingEdges3D * vtkWindowedSincPolyDataFilter */ #ifndef vtkSurfaceNets2D_h #define vtkSurfaceNets2D_h #include "vtkConstrainedSmoothingFilter.h" // Perform mesh smoothing #include "vtkContourValues.h" // Needed for direct access to ContourValues #include "vtkFiltersCoreModule.h" // For export macro #include "vtkPolyData.h" // To support data caching #include "vtkPolyDataAlgorithm.h" class vtkImageData; class VTKFILTERSCORE_EXPORT vtkSurfaceNets2D : public vtkPolyDataAlgorithm { public: ///@{ /** * Standard methods for instantiation, printing, and type information. */ static vtkSurfaceNets2D* New(); vtkTypeMacro(vtkSurfaceNets2D, vtkPolyDataAlgorithm); void PrintSelf(ostream& os, vtkIndent indent) override; ///@} /** * The modified time is also a function of the label values and * the smoothing filter. */ vtkMTimeType GetMTime() override; ///@{ /** * Set a particular label value at label number i. The index i ranges * between 0 <= i Labels->SetValue(i, value); } void SetLabel(int i, double value) { this->Labels->SetValue(i, value); } ///@} ///@{ /** * Get the ith label value. */ double GetValue(int i) { return this->Labels->GetValue(i); } double GetLabel(int i) { return this->Labels->GetValue(i); } ///@} ///@{ /** * Get a pointer to an array of labels. There will be * GetNumberOfLabels() values in the list. */ double* GetValues() { return this->Labels->GetValues(); } double* GetLabels() { return this->Labels->GetValues(); } ///@} ///@{ /** * 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. */ void GetValues(double* contourValues) { this->Labels->GetValues(contourValues); } void GetLabels(double* contourValues) { this->Labels->GetValues(contourValues); } ///@} ///@{ /** * 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. */ void SetNumberOfLabels(int number) { this->Labels->SetNumberOfContours(number); } void SetNumberOfContours(int number) { this->Labels->SetNumberOfContours(number); } ///@} ///@{ /** * Get the number of labels in the list of label values. */ vtkIdType GetNumberOfLabels() { return this->Labels->GetNumberOfContours(); } vtkIdType GetNumberOfContours() { return this->Labels->GetNumberOfContours(); } ///@} ///@{ /** * Generate numLabels equally spaced labels between the specified * range. The labels will include the min/max range values. */ void GenerateLabels(int numLabels, double range[2]) { this->Labels->GenerateValues(numLabels, range); } void GenerateValues(int numContours, double range[2]) { this->Labels->GenerateValues(numContours, range); } void GenerateLabels(int numLabels, double rangeStart, double rangeEnd) { this->Labels->GenerateValues(numLabels, rangeStart, rangeEnd); } void GenerateValues(int numContours, double rangeStart, double rangeEnd) { this->Labels->GenerateValues(numContours, rangeStart, rangeEnd); } ///@} ///@{ /** * Enable/disable an option to generate cell scalars in the output. The * cell scalars are a two-tuple that indicates which labels (i.e., * segmented regions) are on either side of each (line) cell (l0,l1) with * l0 Labels; bool ComputeScalars; double BackgroundLabel; int ArrayComponent; bool Smoothing; vtkSmartPointer Smoother; // Support data caching of the extracted surface nets. This is used to // avoid repeated surface extraction when only smoothing filter // parameters are modified. bool DataCaching; vtkSmartPointer GeometryCache; vtkSmartPointer StencilsCache; vtkTimeStamp SmoothingTime; bool IsCacheEmpty(); void CacheData(vtkPolyData* pd, vtkCellArray* ca); private: vtkSurfaceNets2D(const vtkSurfaceNets2D&) = delete; void operator=(const vtkSurfaceNets2D&) = delete; }; #endif