VTK  9.3.0
vtkGeometryFilter.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
81#ifndef vtkGeometryFilter_h
82#define vtkGeometryFilter_h
83
84#include "vtkFiltersGeometryModule.h" // For export macro
86
87#include <array> // For std::array
88
89VTK_ABI_NAMESPACE_BEGIN
96
97// Used to coordinate delegation to vtkDataSetSurfaceFilter
98struct VTKFILTERSGEOMETRY_EXPORT vtkGeometryFilterHelper
99{
101 {
102 VERTS = 0,
103 LINES = 1,
104 POLYS = 2,
105 STRIPS = 3,
106 OTHER_LINEAR_CELLS = 4,
107 NON_LINEAR_CELLS = 5,
108 NUM_CELL_TYPES
109 };
110 using CellTypesInformation = std::array<bool, NUM_CELL_TYPES>;
112 unsigned char IsLinear;
117 {
118 return this->CellTypesInfo[VERTS] && !this->CellTypesInfo[LINES] &&
119 !this->CellTypesInfo[POLYS] && !this->CellTypesInfo[STRIPS] &&
120 !this->CellTypesInfo[OTHER_LINEAR_CELLS] && !this->CellTypesInfo[NON_LINEAR_CELLS];
121 }
123 {
124 return !this->CellTypesInfo[VERTS] && this->CellTypesInfo[LINES] &&
125 !this->CellTypesInfo[POLYS] && !this->CellTypesInfo[STRIPS] &&
126 !this->CellTypesInfo[OTHER_LINEAR_CELLS] && !this->CellTypesInfo[NON_LINEAR_CELLS];
127 }
129 {
130 return !this->CellTypesInfo[VERTS] && !this->CellTypesInfo[LINES] &&
131 this->CellTypesInfo[POLYS] && !this->CellTypesInfo[STRIPS] &&
132 !this->CellTypesInfo[OTHER_LINEAR_CELLS] && !this->CellTypesInfo[NON_LINEAR_CELLS];
133 }
135 {
136 return !this->CellTypesInfo[VERTS] && !this->CellTypesInfo[LINES] &&
137 !this->CellTypesInfo[POLYS] && this->CellTypesInfo[STRIPS] &&
138 !this->CellTypesInfo[OTHER_LINEAR_CELLS] && !this->CellTypesInfo[NON_LINEAR_CELLS];
139 }
140};
141
142class VTKFILTERSGEOMETRY_EXPORT vtkGeometryFilter : public vtkPolyDataAlgorithm
143{
144public:
146
151 void PrintSelf(ostream& os, vtkIndent indent) override;
153
155
158 vtkSetMacro(PointClipping, bool);
159 vtkGetMacro(PointClipping, bool);
160 vtkBooleanMacro(PointClipping, bool);
162
164
167 vtkSetMacro(CellClipping, bool);
168 vtkGetMacro(CellClipping, bool);
169 vtkBooleanMacro(CellClipping, bool);
171
173
176 vtkSetMacro(ExtentClipping, bool);
177 vtkGetMacro(ExtentClipping, bool);
178 vtkBooleanMacro(ExtentClipping, bool);
180
182
185 vtkSetClampMacro(PointMinimum, vtkIdType, 0, VTK_ID_MAX);
186 vtkGetMacro(PointMinimum, vtkIdType);
188
190
193 vtkSetClampMacro(PointMaximum, vtkIdType, 0, VTK_ID_MAX);
194 vtkGetMacro(PointMaximum, vtkIdType);
196
198
201 vtkSetClampMacro(CellMinimum, vtkIdType, 0, VTK_ID_MAX);
202 vtkGetMacro(CellMinimum, vtkIdType);
204
206
209 vtkSetClampMacro(CellMaximum, vtkIdType, 0, VTK_ID_MAX);
210 vtkGetMacro(CellMaximum, vtkIdType);
212
216 void SetExtent(double xMin, double xMax, double yMin, double yMax, double zMin, double zMax);
217
219
222 void SetExtent(double extent[6]);
223 double* GetExtent() VTK_SIZEHINT(6) { return this->Extent; }
225
227
235 vtkSetMacro(Merging, bool);
236 vtkGetMacro(Merging, bool);
237 vtkBooleanMacro(Merging, bool);
239
241
248 void SetOutputPointsPrecision(int precision);
251
253
260 vtkSetMacro(FastMode, bool);
261 vtkGetMacro(FastMode, bool);
262 vtkBooleanMacro(FastMode, bool);
264
266
273 VTK_DEPRECATED_IN_9_2_0("This method is no longer used and has no effect.")
274 virtual void SetDegree(unsigned int vtkNotUsed(arg)) {}
275 VTK_DEPRECATED_IN_9_2_0("This method is no longer used and has no effect.")
276 virtual unsigned int GetDegreeMinValue() { return 1; }
277 VTK_DEPRECATED_IN_9_2_0("This method is no longer used and has no effect.")
278 virtual unsigned int GetDegreeMaxValue() { return static_cast<int>(~0u >> 1); }
279 VTK_DEPRECATED_IN_9_2_0("This method is no longer used and has no effect.")
280 virtual unsigned int GetDegree() { return 4; }
282
284
289 VTK_DEPRECATED_IN_9_2_0("This method is no longer used and has no effect.")
290 void SetLocator(vtkIncrementalPointLocator* locator);
291 VTK_DEPRECATED_IN_9_2_0("This method is no longer used and has no effect.")
292 vtkGetObjectMacro(Locator, vtkIncrementalPointLocator);
294
299 VTK_DEPRECATED_IN_9_2_0("This method is no longer used and has no effect.")
300 void CreateDefaultLocator();
301
302 // The following are methods compatible with vtkDataSetSurfaceFilter.
303
305
310 vtkSetMacro(PieceInvariant, int);
311 vtkGetMacro(PieceInvariant, int);
313
315
322 vtkSetMacro(PassThroughCellIds, vtkTypeBool);
323 vtkGetMacro(PassThroughCellIds, vtkTypeBool);
324 vtkBooleanMacro(PassThroughCellIds, vtkTypeBool);
326
328
335 vtkSetMacro(PassThroughPointIds, vtkTypeBool);
336 vtkGetMacro(PassThroughPointIds, vtkTypeBool);
337 vtkBooleanMacro(PassThroughPointIds, vtkTypeBool);
339
341
347 vtkSetStringMacro(OriginalCellIdsName);
348 virtual const char* GetOriginalCellIdsName()
349 {
350 return (this->OriginalCellIdsName ? this->OriginalCellIdsName : "vtkOriginalCellIds");
351 }
352 vtkSetStringMacro(OriginalPointIdsName);
353 virtual const char* GetOriginalPointIdsName()
354 {
355 return (this->OriginalPointIdsName ? this->OriginalPointIdsName : "vtkOriginalPointIds");
356 }
358
360
375
377
388 vtkSetMacro(NonlinearSubdivisionLevel, int);
389 vtkGetMacro(NonlinearSubdivisionLevel, int);
391
393
396 vtkSetMacro(Delegation, vtkTypeBool);
397 vtkGetMacro(Delegation, vtkTypeBool);
398 vtkBooleanMacro(Delegation, vtkTypeBool);
400
402
413 vtkSetMacro(RemoveGhostInterfaces, bool);
414 vtkBooleanMacro(RemoveGhostInterfaces, bool);
415 vtkGetMacro(RemoveGhostInterfaces, bool);
417
419
426
428 vtkDataSet* input, vtkPolyData* output, vtkGeometryFilterHelper* info, vtkPolyData* exc);
429 virtual int UnstructuredGridExecute(vtkDataSet* input, vtkPolyData* output);
430
431 VTK_DEPRECATED_IN_9_3_0("Use the new version that has int* instead of vtkInformation*")
432 int StructuredExecute(vtkDataSet* input, vtkPolyData* output, vtkInformation* inInfo,
433 vtkPolyData* exc, bool* extractFace = nullptr);
434 int StructuredExecute(vtkDataSet* input, vtkPolyData* output, int* wholeExtent, vtkPolyData* exc,
435 bool* extractFace = nullptr);
436 VTK_DEPRECATED_IN_9_3_0("Use the new version that has int* instead of vtkInformation*")
437 virtual int StructuredExecute(
438 vtkDataSet* input, vtkPolyData* output, vtkInformation* inInfo, bool* extractFace = nullptr);
439 virtual int StructuredExecute(
440 vtkDataSet* input, vtkPolyData* output, int* wholeExt, bool* extractFace = nullptr);
441
442 int DataSetExecute(vtkDataSet* input, vtkPolyData* output, vtkPolyData* exc);
443 virtual int DataSetExecute(vtkDataSet* input, vtkPolyData* output);
445
446protected:
448 ~vtkGeometryFilter() override;
449
451 int FillInputPortInformation(int port, vtkInformation* info) override;
452
453 // special cases for performance
454 int RequestUpdateExtent(vtkInformation*, vtkInformationVector**, vtkInformationVector*) override;
455
456 vtkIdType PointMaximum;
457 vtkIdType PointMinimum;
458 vtkIdType CellMinimum;
459 vtkIdType CellMaximum;
460 double Extent[6];
461 bool PointClipping;
462 bool CellClipping;
463 bool ExtentClipping;
464 int OutputPointsPrecision;
465 bool RemoveGhostInterfaces;
466
467 bool Merging;
469
470 bool FastMode;
471
472 // These methods support compatibility with vtkDataSetSurfaceFilter
473 int PieceInvariant;
474 vtkTypeBool PassThroughCellIds;
475 char* OriginalCellIdsName;
476
477 vtkTypeBool PassThroughPointIds;
478 char* OriginalPointIdsName;
479
480 int NonlinearSubdivisionLevel;
481
482 vtkTypeBool Delegation;
483
484private:
485 vtkGeometryFilter(const vtkGeometryFilter&) = delete;
486 void operator=(const vtkGeometryFilter&) = delete;
487};
488
489VTK_ABI_NAMESPACE_END
490#endif
Proxy object to connect input/output ports.
Extracts outer surface (as vtkPolyData) of any dataset.
abstract class to specify dataset behavior
Definition vtkDataSet.h:53
extract boundary geometry from dataset (or convert data to polygonal type)
virtual const char * GetOriginalPointIdsName()
If PassThroughCellIds or PassThroughPointIds is on, then these ivars control the name given to the fi...
virtual int PolyDataExecute(vtkDataSet *, vtkPolyData *)
Direct access methods so that this class can be used as an algorithm without using it as a filter (i....
vtkPolyData * GetExcludedFaces()
If a second, vtkPolyData input is provided, this second input specifies a list of faces to be exclude...
int GetOutputPointsPrecision() const
Set/get the desired precision for the output types.
void SetOutputPointsPrecision(int precision)
Set/get the desired precision for the output types.
static vtkGeometryFilter * New()
Standard methods for instantiation, type information, and printing.
void PrintSelf(ostream &os, vtkIndent indent) override
Standard methods for instantiation, type information, and printing.
int UnstructuredGridExecute(vtkDataSet *input, vtkPolyData *output, vtkGeometryFilterHelper *info, vtkPolyData *exc)
Direct access methods so that this class can be used as an algorithm without using it as a filter (i....
virtual int UnstructuredGridExecute(vtkDataSet *input, vtkPolyData *output)
Direct access methods so that this class can be used as an algorithm without using it as a filter (i....
void SetExcludedFacesConnection(vtkAlgorithmOutput *algOutput)
If a second, vtkPolyData input is provided, this second input specifies a list of faces to be exclude...
int PolyDataExecute(vtkDataSet *input, vtkPolyData *output, vtkPolyData *exc)
Direct access methods so that this class can be used as an algorithm without using it as a filter (i....
void SetExtent(double extent[6])
Set / get a (xmin,xmax, ymin,ymax, zmin,zmax) bounding box to clip data.
void SetExcludedFacesData(vtkPolyData *)
If a second, vtkPolyData input is provided, this second input specifies a list of faces to be exclude...
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.
double * GetExtent()
Set / get a (xmin,xmax, ymin,ymax, zmin,zmax) bounding box to clip data.
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition vtkIndent.h:29
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
Superclass for algorithms that produce only polydata as output.
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition vtkPolyData.h:80
topologically regular array of data
dataset represents arbitrary combinations of all possible cell types.
static vtkGeometryFilterHelper * CharacterizeUnstructuredGrid(vtkUnstructuredGridBase *)
CellTypesInformation CellTypesInfo
static void CopyFilterParams(vtkDataSetSurfaceFilter *dssf, vtkGeometryFilter *gf)
std::array< bool, NUM_CELL_TYPES > CellTypesInformation
static void CopyFilterParams(vtkGeometryFilter *gf, vtkDataSetSurfaceFilter *dssf)
int vtkTypeBool
Definition vtkABI.h:64
#define VTK_DEPRECATED_IN_9_3_0(reason)
#define VTK_DEPRECATED_IN_9_2_0(reason)
int vtkIdType
Definition vtkType.h:315
#define VTK_ID_MAX
Definition vtkType.h:319
#define VTK_SIZEHINT(...)