VTK  9.3.20240417
vtkUnstructuredGridGeometryFilter.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
32 #ifndef vtkUnstructuredGridGeometryFilter_h
33 #define vtkUnstructuredGridGeometryFilter_h
34 
35 #include "vtkFiltersGeometryModule.h" // For export macro
37 
38 VTK_ABI_NAMESPACE_BEGIN
40 class vtkHashTableOfSurfels; // internal class
41 
42 class VTKFILTERSGEOMETRY_EXPORT vtkUnstructuredGridGeometryFilter
44 {
45 public:
48  void PrintSelf(ostream& os, vtkIndent indent) override;
49 
51 
54  vtkSetMacro(PointClipping, vtkTypeBool);
55  vtkGetMacro(PointClipping, vtkTypeBool);
56  vtkBooleanMacro(PointClipping, vtkTypeBool);
58 
60 
63  vtkSetMacro(CellClipping, vtkTypeBool);
64  vtkGetMacro(CellClipping, vtkTypeBool);
65  vtkBooleanMacro(CellClipping, vtkTypeBool);
67 
69 
72  vtkSetMacro(ExtentClipping, vtkTypeBool);
73  vtkGetMacro(ExtentClipping, vtkTypeBool);
74  vtkBooleanMacro(ExtentClipping, vtkTypeBool);
76 
78 
82  vtkSetMacro(DuplicateGhostCellClipping, vtkTypeBool);
83  vtkGetMacro(DuplicateGhostCellClipping, vtkTypeBool);
84  vtkBooleanMacro(DuplicateGhostCellClipping, vtkTypeBool);
86 
88 
91  vtkSetClampMacro(PointMinimum, vtkIdType, 0, VTK_ID_MAX);
92  vtkGetMacro(PointMinimum, vtkIdType);
94 
96 
99  vtkSetClampMacro(PointMaximum, vtkIdType, 0, VTK_ID_MAX);
100  vtkGetMacro(PointMaximum, vtkIdType);
102 
104 
107  vtkSetClampMacro(CellMinimum, vtkIdType, 0, VTK_ID_MAX);
108  vtkGetMacro(CellMinimum, vtkIdType);
110 
112 
115  vtkSetClampMacro(CellMaximum, vtkIdType, 0, VTK_ID_MAX);
116  vtkGetMacro(CellMaximum, vtkIdType);
118 
122  void SetExtent(double xMin, double xMax, double yMin, double yMax, double zMin, double zMax);
123 
125 
128  void SetExtent(double extent[6]);
129  double* GetExtent() { return this->Extent; }
131 
133 
138  vtkSetMacro(Merging, vtkTypeBool);
139  vtkGetMacro(Merging, vtkTypeBool);
140  vtkBooleanMacro(Merging, vtkTypeBool);
142 
144 
152  vtkSetMacro(PassThroughCellIds, vtkTypeBool);
153  vtkGetMacro(PassThroughCellIds, vtkTypeBool);
154  vtkBooleanMacro(PassThroughCellIds, vtkTypeBool);
155  vtkSetMacro(PassThroughPointIds, vtkTypeBool);
156  vtkGetMacro(PassThroughPointIds, vtkTypeBool);
157  vtkBooleanMacro(PassThroughPointIds, vtkTypeBool);
159 
161 
167  vtkSetMacro(MatchBoundariesIgnoringCellOrder, vtkTypeBool);
168  vtkGetMacro(MatchBoundariesIgnoringCellOrder, vtkTypeBool);
170 
172 
178  vtkSetStringMacro(OriginalCellIdsName);
179  virtual const char* GetOriginalCellIdsName()
180  {
181  return (this->OriginalCellIdsName ? this->OriginalCellIdsName : "vtkOriginalCellIds");
182  }
183  vtkSetStringMacro(OriginalPointIdsName);
184  virtual const char* GetOriginalPointIdsName()
185  {
186  return (this->OriginalPointIdsName ? this->OriginalPointIdsName : "vtkOriginalPointIds");
187  }
189 
191 
196  vtkGetObjectMacro(Locator, vtkIncrementalPointLocator);
198 
203 
207  vtkMTimeType GetMTime() override;
208 
209 protected:
212 
215 
217 
222  double Extent[6];
227 
233 
236 
237  vtkHashTableOfSurfels* HashTable;
238 
239 private:
241  void operator=(const vtkUnstructuredGridGeometryFilter&) = delete;
242 };
243 
244 VTK_ABI_NAMESPACE_END
245 #endif
Abstract class in support of both point location and point insertion.
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 vtkUnstructureGridBase subclasses as output.
extract geometry from an unstructured grid
void SetExtent(double xMin, double xMax, double yMin, double yMax, double zMin, double zMax)
Specify a (xmin,xmax, ymin,ymax, zmin,zmax) bounding box to clip data.
int RequestUpdateExtent(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
void SetExtent(double extent[6])
Set / get a (xmin,xmax, ymin,ymax, zmin,zmax) bounding box to clip data.
void SetLocator(vtkIncrementalPointLocator *locator)
Set / get a spatial locator for merging points.
static vtkUnstructuredGridGeometryFilter * New()
void CreateDefaultLocator()
Create default locator.
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
vtkMTimeType GetMTime() override
Return the MTime also considering the locator.
double * GetExtent()
Set / get a (xmin,xmax, ymin,ymax, zmin,zmax) bounding box to clip data.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
virtual const char * GetOriginalPointIdsName()
If PassThroughCellIds or PassThroughPointIds is on, then these ivars control the name given to the fi...
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
virtual const char * GetOriginalCellIdsName()
If PassThroughCellIds or PassThroughPointIds is on, then these ivars control the name given to the fi...
@ info
Definition: vtkX3D.h:376
@ port
Definition: vtkX3D.h:447
@ extent
Definition: vtkX3D.h:345
int vtkTypeBool
Definition: vtkABI.h:64
int vtkIdType
Definition: vtkType.h:315
#define VTK_ID_MAX
Definition: vtkType.h:319
vtkTypeUInt32 vtkMTimeType
Definition: vtkType.h:270