VTK  9.1.0
vtkGeometryFilter.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkGeometryFilter.h
5
6 Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7 All rights reserved.
8 See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9
10 This software is distributed WITHOUT ANY WARRANTY; without even
11 the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12 PURPOSE. See the above copyright notice for more information.
13
14=========================================================================*/
96#ifndef vtkGeometryFilter_h
97#define vtkGeometryFilter_h
98
99#include "vtkFiltersGeometryModule.h" // For export macro
100#include "vtkPolyDataAlgorithm.h"
101
108struct vtkExcludedFaces;
109
110// Used to coordinate delegation to vtkDataSetSurfaceFilter
111struct VTKFILTERSGEOMETRY_EXPORT vtkGeometryFilterHelper
112{
113 unsigned char IsLinear;
117};
118
119class VTKFILTERSGEOMETRY_EXPORT vtkGeometryFilter : public vtkPolyDataAlgorithm
120{
121public:
123
128 void PrintSelf(ostream& os, vtkIndent indent) override;
130
132
135 vtkSetMacro(PointClipping, bool);
136 vtkGetMacro(PointClipping, bool);
137 vtkBooleanMacro(PointClipping, bool);
139
141
144 vtkSetMacro(CellClipping, bool);
145 vtkGetMacro(CellClipping, bool);
146 vtkBooleanMacro(CellClipping, bool);
148
150
153 vtkSetMacro(ExtentClipping, bool);
154 vtkGetMacro(ExtentClipping, bool);
155 vtkBooleanMacro(ExtentClipping, bool);
157
159
162 vtkSetClampMacro(PointMinimum, vtkIdType, 0, VTK_ID_MAX);
163 vtkGetMacro(PointMinimum, vtkIdType);
165
167
170 vtkSetClampMacro(PointMaximum, vtkIdType, 0, VTK_ID_MAX);
171 vtkGetMacro(PointMaximum, vtkIdType);
173
175
178 vtkSetClampMacro(CellMinimum, vtkIdType, 0, VTK_ID_MAX);
179 vtkGetMacro(CellMinimum, vtkIdType);
181
183
186 vtkSetClampMacro(CellMaximum, vtkIdType, 0, VTK_ID_MAX);
187 vtkGetMacro(CellMaximum, vtkIdType);
189
193 void SetExtent(double xMin, double xMax, double yMin, double yMax, double zMin, double zMax);
194
196
199 void SetExtent(double extent[6]);
200 double* GetExtent() VTK_SIZEHINT(6) { return this->Extent; }
202
204
212 vtkSetMacro(Merging, bool);
213 vtkGetMacro(Merging, bool);
214 vtkBooleanMacro(Merging, bool);
216
218
225 void SetOutputPointsPrecision(int precision);
228
230
236 vtkSetMacro(FastMode, bool);
237 vtkGetMacro(FastMode, bool);
238 vtkBooleanMacro(FastMode, bool);
240
242
249 vtkSetClampMacro(Degree, unsigned int, 1, VTK_INT_MAX);
250 vtkGetMacro(Degree, unsigned int);
252
254
260 vtkGetObjectMacro(Locator, vtkIncrementalPointLocator);
262
268
269 // The following are methods compatible with vtkDataSetSurfaceFilter.
270
272
277 vtkSetMacro(PieceInvariant, int);
278 vtkGetMacro(PieceInvariant, int);
280
282
290 vtkSetMacro(PassThroughCellIds, vtkTypeBool);
291 vtkGetMacro(PassThroughCellIds, vtkTypeBool);
292 vtkBooleanMacro(PassThroughCellIds, vtkTypeBool);
293 vtkSetMacro(PassThroughPointIds, vtkTypeBool);
294 vtkGetMacro(PassThroughPointIds, vtkTypeBool);
295 vtkBooleanMacro(PassThroughPointIds, vtkTypeBool);
297
299
305 vtkSetStringMacro(OriginalCellIdsName);
306 virtual const char* GetOriginalCellIdsName()
307 {
308 return (this->OriginalCellIdsName ? this->OriginalCellIdsName : "vtkOriginalCellIds");
309 }
310 vtkSetStringMacro(OriginalPointIdsName);
311 virtual const char* GetOriginalPointIdsName()
312 {
313 return (this->OriginalPointIdsName ? this->OriginalPointIdsName : "vtkOriginalPointIds");
314 }
316
318
333
335
346 vtkSetMacro(NonlinearSubdivisionLevel, int);
347 vtkGetMacro(NonlinearSubdivisionLevel, int);
349
351
354 vtkSetMacro(Delegation, vtkTypeBool);
355 vtkGetMacro(Delegation, vtkTypeBool);
356 vtkBooleanMacro(Delegation, vtkTypeBool);
358
360
365 int PolyDataExecute(vtkDataSet* input, vtkPolyData* output, vtkExcludedFaces* exc);
367
369 vtkDataSet* input, vtkPolyData* output, vtkGeometryFilterHelper* info, vtkExcludedFaces* exc);
370 virtual int UnstructuredGridExecute(vtkDataSet* input, vtkPolyData* output);
371
373 vtkDataSet* input, vtkPolyData* output, vtkInformation* inInfo, vtkExcludedFaces* exc);
374 virtual int StructuredExecute(vtkDataSet* input, vtkPolyData* output, vtkInformation* inInfo);
375
376 int DataSetExecute(vtkDataSet* input, vtkPolyData* output, vtkExcludedFaces* exc);
377 virtual int DataSetExecute(vtkDataSet* input, vtkPolyData* output);
379
380protected:
383
386
387 // special cases for performance
389
394 double Extent[6];
399
402
404 unsigned int Degree;
405
406 // This methods support compatability with vtkDataSetSurfaceFilter
410
413
415
417
418private:
419 vtkGeometryFilter(const vtkGeometryFilter&) = delete;
420 void operator=(const vtkGeometryFilter&) = delete;
421};
422
423#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:57
extract boundary geometry from dataset (or convert data to polygonal type)
int PolyDataExecute(vtkDataSet *input, vtkPolyData *output, vtkExcludedFaces *exc)
Direct access methods so that this class can be used as an algorithm without using it as a filter (i....
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
virtual const char * GetOriginalPointIdsName()
If PassThroughCellIds or PassThroughPointIds is on, then these ivars control the name given to the fi...
vtkTypeBool PassThroughCellIds
void CreateDefaultLocator()
Create default locator.
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....
virtual const char * GetOriginalCellIdsName()
If PassThroughCellIds or PassThroughPointIds is on, then these ivars control the name given to the fi...
vtkIncrementalPointLocator * Locator
vtkPolyData * GetExcludedFaces()
If a second, vtkPolyData input is provided, this second input specifies a list of faces to be exclude...
int RequestUpdateExtent(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
vtkTypeBool PassThroughPointIds
int GetOutputPointsPrecision() const
Set/get the desired precision for the output types.
int DataSetExecute(vtkDataSet *input, vtkPolyData *output, vtkExcludedFaces *exc)
Direct access methods so that this class can be used as an algorithm without using it as a filter (i....
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.
void SetLocator(vtkIncrementalPointLocator *locator)
Set / get a spatial locator for merging points.
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....
int UnstructuredGridExecute(vtkDataSet *input, vtkPolyData *output, vtkGeometryFilterHelper *info, vtkExcludedFaces *exc)
Direct access methods so that this class can be used as an algorithm without using it as a filter (i....
int StructuredExecute(vtkDataSet *input, vtkPolyData *output, vtkInformation *inInfo, vtkExcludedFaces *exc)
Direct access methods so that this class can be used as an algorithm without using it as a filter (i....
~vtkGeometryFilter() override
virtual int StructuredExecute(vtkDataSet *input, vtkPolyData *output, vtkInformation *inInfo)
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...
void SetExtent(double extent[6])
Set / get a (xmin,xmax, ymin,ymax, zmin,zmax) bounding box to clip data.
int FillInputPortInformation(int port, vtkInformation *info) override
Fill the input port information objects for this algorithm.
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.
virtual int DataSetExecute(vtkDataSet *input, vtkPolyData *output)
Direct access methods so that this class can be used as an algorithm without using it as a filter (i....
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition: vtkIndent.h:34
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:86
topologically regular array of data
dataset represents arbitrary combinations of all possible cell types
@ info
Definition: vtkX3D.h:382
@ port
Definition: vtkX3D.h:453
@ extent
Definition: vtkX3D.h:351
static vtkGeometryFilterHelper * CharacterizeUnstructuredGrid(vtkUnstructuredGrid *)
static void CopyFilterParams(vtkDataSetSurfaceFilter *dssf, vtkGeometryFilter *gf)
static void CopyFilterParams(vtkGeometryFilter *gf, vtkDataSetSurfaceFilter *dssf)
int vtkTypeBool
Definition: vtkABI.h:69
int vtkIdType
Definition: vtkType.h:332
#define VTK_ID_MAX
Definition: vtkType.h:336
#define VTK_INT_MAX
Definition: vtkType.h:155
#define VTK_SIZEHINT(...)