VTK  9.3.20240417
vtkStaticCleanUnstructuredGrid.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
61 #ifndef vtkStaticCleanUnstructuredGrid_h
62 #define vtkStaticCleanUnstructuredGrid_h
63 
64 #include "vtkFiltersCoreModule.h" // For export macro
65 #include "vtkSmartPointer.h" // Pointer to locator
66 #include "vtkStaticPointLocator.h" // Locator
68 #include <vector> // API to vtkStaticCleanUnstructuredGrid
69 
70 VTK_ABI_NAMESPACE_BEGIN
71 class vtkCellArray;
72 class vtkCellData;
73 class vtkPointData;
74 class vtkPoints;
75 
77 {
78 public:
80 
85  void PrintSelf(ostream& os, vtkIndent indent) override;
88 
90 
96  vtkSetMacro(ToleranceIsAbsolute, bool);
97  vtkBooleanMacro(ToleranceIsAbsolute, bool);
98  vtkGetMacro(ToleranceIsAbsolute, bool);
100 
102 
106  vtkSetClampMacro(AbsoluteTolerance, double, 0.0, VTK_DOUBLE_MAX);
107  vtkGetMacro(AbsoluteTolerance, double);
109 
111 
116  vtkSetClampMacro(Tolerance, double, 0.0, 1.0);
117  vtkGetMacro(Tolerance, double);
119 
121 
132  vtkSetStringMacro(MergingArray);
133  vtkGetStringMacro(MergingArray);
135 
137 
144  vtkSetMacro(RemoveUnusedPoints, bool);
145  vtkBooleanMacro(RemoveUnusedPoints, bool);
146  vtkGetMacro(RemoveUnusedPoints, bool);
148 
150 
157  vtkSetMacro(ProduceMergeMap, bool);
158  vtkBooleanMacro(ProduceMergeMap, bool);
159  vtkGetMacro(ProduceMergeMap, bool);
161 
163 
171  vtkSetMacro(AveragePointData, bool);
172  vtkBooleanMacro(AveragePointData, bool);
173  vtkGetMacro(AveragePointData, bool);
175 
177 
182  vtkSetMacro(OutputPointsPrecision, int);
183  vtkGetMacro(OutputPointsPrecision, int);
185 
191  vtkGetObjectMacro(Locator, vtkStaticPointLocator);
192 
194  // This filter is difficult to stream. To produce invariant results, the
195  // whole input must be processed at once. This flag allows the user to
196  // select whether strict piece invariance is required. By default it is
197  // on. When off, the filter can stream, but results may change.
198  vtkSetMacro(PieceInvariant, bool);
199  vtkGetMacro(PieceInvariant, bool);
200  vtkBooleanMacro(PieceInvariant, bool);
202 
206  vtkMTimeType GetMTime() override;
207 
208 protected:
210  ~vtkStaticCleanUnstructuredGrid() override = default;
211 
212  // Usual data generation method
215 
217  double Tolerance;
225 
226  // Internal locator for performing point merging
228 
229  // These methods are made static so that vtkStaticCleanPolyData can use them
231  static void MarkPointUses(vtkCellArray* ca, vtkIdType* mergeMap, unsigned char* ptUses);
233  vtkIdType numPts, vtkIdType* pmap, unsigned char* ptUses, std::vector<vtkIdType>& mergeMap);
234  static void CopyPoints(
235  vtkPoints* inPts, vtkPointData* inPD, vtkPoints* outPts, vtkPointData* outPD, vtkIdType* ptMap);
236  static void AveragePoints(vtkPoints* inPts, vtkPointData* inPD, vtkPoints* outPts,
237  vtkPointData* outPD, vtkIdType* ptMap, double tol);
238 
239 private:
241  void operator=(const vtkStaticCleanUnstructuredGrid&) = delete;
242 };
243 
244 VTK_ABI_NAMESPACE_END
245 #endif
object to represent cell connectivity
Definition: vtkCellArray.h:286
represent and manipulate cell attribute data
Definition: vtkCellData.h:141
a simple class to control print indentation
Definition: vtkIndent.h:108
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
represent and manipulate point attribute data
Definition: vtkPointData.h:140
represent and manipulate 3D points
Definition: vtkPoints.h:139
merge duplicate points, and/or remove unused points and/or remove degenerate cells
merge duplicate points, removed unused points, in an vtkUnstructuredGrid
int RequestUpdateExtent(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *) override
This is called by the superclass.
vtkMTimeType GetMTime() override
Get the MTime of this object also considering the locator.
~vtkStaticCleanUnstructuredGrid() override=default
vtkSmartPointer< vtkStaticPointLocator > Locator
void PrintSelf(ostream &os, vtkIndent indent) override
Standard methods for instantiation, obtaining type information, and printing the state of the object.
static void MarkPointUses(vtkCellArray *ca, vtkIdType *mergeMap, unsigned char *ptUses)
static vtkStaticCleanUnstructuredGrid * New()
Standard methods for instantiation, obtaining type information, and printing the state of the object.
static vtkIdType BuildPointMap(vtkIdType numPts, vtkIdType *pmap, unsigned char *ptUses, std::vector< vtkIdType > &mergeMap)
static void AveragePoints(vtkPoints *inPts, vtkPointData *inPD, vtkPoints *outPts, vtkPointData *outPD, vtkIdType *ptMap, double tol)
static void CopyPoints(vtkPoints *inPts, vtkPointData *inPD, vtkPoints *outPts, vtkPointData *outPD, vtkIdType *ptMap)
quickly locate points in 3-space
Superclass for algorithms that produce only unstructured grid as output.
int vtkIdType
Definition: vtkType.h:315
vtkTypeUInt32 vtkMTimeType
Definition: vtkType.h:270
#define VTK_DOUBLE_MAX
Definition: vtkType.h:154