VTK  9.3.20240328
vtkPolyDataEdgeConnectivityFilter.h
Go to the documentation of this file.
1 // SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2 // SPDX-License-Identifier: BSD-3-Clause
73 #ifndef vtkPolyDataEdgeConnectivityFilter_h
74 #define vtkPolyDataEdgeConnectivityFilter_h
75 
76 #include "vtkFiltersCoreModule.h" // For export macro
77 #include "vtkIdTypeArray.h" //
78 #include "vtkPolyDataAlgorithm.h"
79 
80 #define VTK_EXTRACT_POINT_SEEDED_REGIONS 1
81 #define VTK_EXTRACT_CELL_SEEDED_REGIONS 2
82 #define VTK_EXTRACT_SPECIFIED_REGIONS 3
83 #define VTK_EXTRACT_LARGEST_REGION 4
84 #define VTK_EXTRACT_ALL_REGIONS 5
85 #define VTK_EXTRACT_CLOSEST_POINT_REGION 6
86 #define VTK_EXTRACT_LARGE_REGIONS 7
87 
88 VTK_ABI_NAMESPACE_BEGIN
89 class vtkDataArray;
90 class vtkCharArray;
91 class vtkIdList;
92 class vtkIdTypeArray;
93 class vtkEdgeTable;
94 
95 class VTKFILTERSCORE_EXPORT vtkPolyDataEdgeConnectivityFilter : public vtkPolyDataAlgorithm
96 {
97 public:
99 
104  void PrintSelf(ostream& os, vtkIndent indent) override;
106 
108 
111  vtkSetClampMacro(
113  vtkGetMacro(ExtractionMode, int);
115  {
116  this->SetExtractionMode(VTK_EXTRACT_POINT_SEEDED_REGIONS);
117  }
119  {
120  this->SetExtractionMode(VTK_EXTRACT_CELL_SEEDED_REGIONS);
121  }
124  {
125  this->SetExtractionMode(VTK_EXTRACT_SPECIFIED_REGIONS);
126  }
128  {
129  this->SetExtractionMode(VTK_EXTRACT_CLOSEST_POINT_REGION);
130  }
131  void SetExtractionModeToLargeRegions() { this->SetExtractionMode(VTK_EXTRACT_LARGE_REGIONS); }
132  void SetExtractionModeToAllRegions() { this->SetExtractionMode(VTK_EXTRACT_ALL_REGIONS); }
133  const char* GetExtractionModeAsString();
135 
137 
142  vtkSetMacro(BarrierEdges, vtkTypeBool);
143  vtkGetMacro(BarrierEdges, vtkTypeBool);
144  vtkBooleanMacro(BarrierEdges, vtkTypeBool);
146 
148 
159 
161 
169  vtkSetVector2Macro(BarrierEdgeLength, double);
170  vtkGetVector2Macro(BarrierEdgeLength, double);
172 
174 
179  vtkSetMacro(ScalarConnectivity, vtkTypeBool);
180  vtkGetMacro(ScalarConnectivity, vtkTypeBool);
181  vtkBooleanMacro(ScalarConnectivity, vtkTypeBool);
183 
185 
188  vtkSetVector2Macro(ScalarRange, double);
189  vtkGetVector2Macro(ScalarRange, double);
191 
193 
197  vtkGetObjectMacro(RegionSizes, vtkIdTypeArray);
199 
204 
208  void AddSeed(int id);
209 
213  void DeleteSeed(int id);
214 
219 
223  void AddSpecifiedRegion(int id);
224 
228  void DeleteSpecifiedRegion(int id);
229 
234 
236 
240  vtkSetVector3Macro(ClosestPoint, double);
241  vtkGetVectorMacro(ClosestPoint, double, 3);
243 
244  // Control the region growing process.
246  {
247  RegionGrowingOff = 0,
248  LargeRegions = 1,
249  SmallRegions = 2
250  };
251 
253 
263  vtkSetClampMacro(RegionGrowing, int, RegionGrowingOff, SmallRegions);
264  vtkGetMacro(RegionGrowing, int);
265  void SetRegionGrowingOff() { this->SetRegionGrowing(RegionGrowingOff); }
266  void GrowLargeRegionsOff() { this->SetRegionGrowing(RegionGrowingOff); }
267  void GrowSmallRegionsOff() { this->SetRegionGrowing(RegionGrowingOff); }
268  void SetRegionGrowingToLargeRegions() { this->SetRegionGrowing(LargeRegions); }
269  void GrowLargeRegionsOn() { this->SetRegionGrowing(LargeRegions); }
270  void SetRegionGrowingToSmallRegions() { this->SetRegionGrowing(SmallRegions); }
271  void GrowSmallRegionsOn() { this->SetRegionGrowing(SmallRegions); }
273 
275 
283  vtkSetClampMacro(LargeRegionThreshold, double, 0.0, 1.0);
284  vtkGetMacro(LargeRegionThreshold, double);
286 
291  int GetNumberOfExtractedRegions() { return this->NumberOfExtractedRegions; }
292 
296  double GetTotalArea() { return this->TotalArea; }
297 
299 
304  vtkSetMacro(ColorRegions, vtkTypeBool);
305  vtkGetMacro(ColorRegions, vtkTypeBool);
306  vtkBooleanMacro(ColorRegions, vtkTypeBool);
308 
310 
315  vtkSetMacro(CellRegionAreas, vtkTypeBool);
316  vtkGetMacro(CellRegionAreas, vtkTypeBool);
317  vtkBooleanMacro(CellRegionAreas, vtkTypeBool);
319 
321 
326  vtkSetMacro(OutputPointsPrecision, int);
327  vtkGetMacro(OutputPointsPrecision, int);
329 
330 protected:
333 
334  // Usual data generation method
336 
337  // Optional second input
339 
340  // Filter data members
341  vtkTypeBool ColorRegions; // boolean turns on/off scalar generation for separate regions
342  vtkTypeBool CellRegionAreas; // for each cell, the area of the region the cell is associated with
343 
344  int ExtractionMode; // how to extract regions
345  vtkTypeBool BarrierEdges; // enable barrier edges
346  double BarrierEdgeLength[2]; // edges of length within this range are barrier edges
348  double ScalarRange[2];
349  std::vector<vtkIdType> Seeds; // id's of points or cells used to seed regions
350  std::vector<vtkIdType> SpecifiedRegionIds; // regions specified for extraction
351  vtkSmartPointer<vtkIdTypeArray> RegionSizes; // size (in cells) of each region extracted
352  double ClosestPoint[3];
354 
355  // Methods for iterative traversal and marking cells
358  vtkIdType cellId, vtkIdType npts, const vtkIdType* pts, vtkIdList* neis);
361 
362  // Methods implementing iterative region growing
365  int CurrentGrowPass; // region growing is a multiple-pass process
367  void ExchangeRegions(vtkIdType currentRegionId, vtkIdType neiId, vtkIdType neiRegId);
370  int AssimilateCell(vtkIdType cellId, vtkIdType npts, const vtkIdType* pts);
373 
374  double TotalArea; // the total area of the input mesh
375  std::vector<double> CellAreas; // the area of each polygonal cell
376  std::vector<double> RegionAreas; // the total area of each region
377  std::vector<char> RegionClassification; // indicate whether the region is large or small
378 
379  // used to support algorithm execution
380  std::vector<vtkIdType> RegionIds;
381  std::vector<vtkIdType> PointMap;
389  std::vector<vtkIdType> Wave;
390  std::vector<vtkIdType> Wave2;
395  double BRange2[2]; // BarrierEdgeLength[0,1]**2 of edge lengths defining barriers
396 
397 private:
399  void operator=(const vtkPolyDataEdgeConnectivityFilter&) = delete;
400 };
401 
406 {
408  {
409  return "ExtractPointSeededRegions";
410  }
412  {
413  return "ExtractCellSeededRegions";
414  }
416  {
417  return "ExtractSpecifiedRegions";
418  }
419  else if (this->ExtractionMode == VTK_EXTRACT_ALL_REGIONS)
420  {
421  return "ExtractAllRegions";
422  }
424  {
425  return "ExtractClosestPointRegion";
426  }
427  else if (this->ExtractionMode == VTK_EXTRACT_LARGE_REGIONS)
428  {
429  return "ExtractLargeRegions";
430  }
431  else
432  {
433  return "ExtractLargestRegion";
434  }
435 }
436 
437 VTK_ABI_NAMESPACE_END
438 #endif
Proxy object to connect input/output ports.
dynamic, self-adjusting array of char
Definition: vtkCharArray.h:60
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:154
keep track of edges (edge is pair of integer id's)
Definition: vtkEdgeTable.h:30
list of point or cell ids
Definition: vtkIdList.h:132
dynamic, self-adjusting array of vtkIdType
a simple class to control print indentation
Definition: vtkIndent.h:108
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
Superclass for algorithms that produce only polydata as output.
segment polygonal mesh based on shared edge connectivity
void GrowLargeRegionsOff()
Specify a strategy for region growing.
void SetExtractionModeToLargeRegions()
Control the extraction of connected surfaces.
void GrowSmallRegionsOn()
Specify a strategy for region growing.
void SetRegionGrowingToLargeRegions()
Specify a strategy for region growing.
void GrowLargeRegionsOn()
Specify a strategy for region growing.
int GetNumberOfSpecifiedRegions()
Get number of specified regions.
vtkSmartPointer< vtkIdTypeArray > RegionSizes
int AssimilateCell(vtkIdType cellId, vtkIdType npts, const vtkIdType *pts)
int IsScalarConnected(vtkIdType cellId, vtkIdType neiId)
void GrowSmallRegionsOff()
Specify a strategy for region growing.
const char * GetExtractionModeAsString()
Return the method of extraction as a string.
int FillInputPortInformation(int, vtkInformation *) override
Fill the input port information objects for this algorithm.
void SetExtractionModeToLargestRegion()
Control the extraction of connected surfaces.
void SetSourceData(vtkPolyData *)
Specify the source vtkPolyData object used to specify barrier edges (this is an optional connection....
void SetExtractionModeToAllRegions()
Control the extraction of connected surfaces.
void SetRegionGrowingToSmallRegions()
Specify a strategy for region growing.
void ExchangeRegions(vtkIdType currentRegionId, vtkIdType neiId, vtkIdType neiRegId)
static vtkPolyDataEdgeConnectivityFilter * New()
Standard methods to instantiate, get type information, and print the object.
void SetExtractionModeToSpecifiedRegions()
Control the extraction of connected surfaces.
vtkPolyData * GetSource()
Specify the source vtkPolyData object used to specify barrier edges (this is an optional connection....
void SetExtractionModeToPointSeededRegions()
Control the extraction of connected surfaces.
void SetSourceConnection(vtkAlgorithmOutput *algOutput)
Specify the source vtkPolyData object used to specify barrier edges (this is an optional connection....
void InitializeSeedList()
Initialize list of point ids/cell ids used to seed regions.
void SetExtractionModeToCellSeededRegions()
Control the extraction of connected surfaces.
void SetRegionGrowingOff()
Specify a strategy for region growing.
void DeleteSpecifiedRegion(int id)
Delete a region id to extract.
bool IsBarrierEdge(vtkIdType p0, vtkIdType p1)
int GetNumberOfExtractedRegions()
Obtain the number of connected regions found.
void SetExtractionModeToClosestPointRegion()
Control the extraction of connected surfaces.
void DeleteSeed(int id)
Delete a seed id (point or cell id).
void GetConnectedNeighbors(vtkIdType cellId, vtkIdType npts, const vtkIdType *pts, vtkIdList *neis)
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
void AddSpecifiedRegion(int id)
Add a region id to extract.
void AddSeed(int id)
Add a seed id (point or cell id).
void InitializeSpecifiedRegionList()
Initialize list of region ids to extract.
void PrintSelf(ostream &os, vtkIndent indent) override
Standard methods to instantiate, get type information, and print the object.
double GetTotalArea()
Obtain the total area of all regions combined.
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:180
int vtkTypeBool
Definition: vtkABI.h:64
#define VTK_EXTRACT_CLOSEST_POINT_REGION
#define VTK_EXTRACT_POINT_SEEDED_REGIONS
#define VTK_EXTRACT_ALL_REGIONS
#define VTK_EXTRACT_CELL_SEEDED_REGIONS
#define VTK_EXTRACT_LARGE_REGIONS
#define VTK_EXTRACT_SPECIFIED_REGIONS
#define VTK_EXTRACT_LARGEST_REGION
int vtkIdType
Definition: vtkType.h:315