VTK  9.5.20250910
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 "vtkWrappingHints.h" // For VTK_MARSHALAUTO
69#include <vector> // API to vtkStaticCleanUnstructuredGrid
70
71VTK_ABI_NAMESPACE_BEGIN
72class vtkCellArray;
73class vtkCellData;
74class vtkPointData;
75class vtkPoints;
76
79{
80public:
82
87 void PrintSelf(ostream& os, vtkIndent indent) override;
90
92
98 vtkSetMacro(ToleranceIsAbsolute, bool);
99 vtkBooleanMacro(ToleranceIsAbsolute, bool);
100 vtkGetMacro(ToleranceIsAbsolute, bool);
102
104
108 vtkSetClampMacro(AbsoluteTolerance, double, 0.0, VTK_DOUBLE_MAX);
109 vtkGetMacro(AbsoluteTolerance, double);
111
113
118 vtkSetClampMacro(Tolerance, double, 0.0, 1.0);
119 vtkGetMacro(Tolerance, double);
121
123
134 vtkSetStringMacro(MergingArray);
135 vtkGetStringMacro(MergingArray);
137
139
146 vtkSetMacro(RemoveUnusedPoints, bool);
147 vtkBooleanMacro(RemoveUnusedPoints, bool);
148 vtkGetMacro(RemoveUnusedPoints, bool);
150
152
159 vtkSetMacro(ProduceMergeMap, bool);
160 vtkBooleanMacro(ProduceMergeMap, bool);
161 vtkGetMacro(ProduceMergeMap, bool);
163
165
173 vtkSetMacro(AveragePointData, bool);
174 vtkBooleanMacro(AveragePointData, bool);
175 vtkGetMacro(AveragePointData, bool);
177
179
184 vtkSetMacro(OutputPointsPrecision, int);
185 vtkGetMacro(OutputPointsPrecision, int);
187
193 vtkGetObjectMacro(Locator, vtkStaticPointLocator);
194
196 // This filter is difficult to stream. To produce invariant results, the
197 // whole input must be processed at once. This flag allows the user to
198 // select whether strict piece invariance is required. By default it is
199 // on. When off, the filter can stream, but results may change.
200 vtkSetMacro(PieceInvariant, bool);
201 vtkGetMacro(PieceInvariant, bool);
202 vtkBooleanMacro(PieceInvariant, bool);
204
209
210 // These methods are made static so that other filters like vtkStaticCleanPolyData can use them
211 static void MarkPointUses(vtkCellArray* ca, vtkIdType* mergeMap, unsigned char* ptUses);
213 vtkIdType numPts, vtkIdType* pmap, unsigned char* ptUses, std::vector<vtkIdType>& mergeMap);
214 static void CopyPoints(
215 vtkPoints* inPts, vtkPointData* inPD, vtkPoints* outPts, vtkPointData* outPD, vtkIdType* ptMap);
216 static void AveragePoints(vtkPoints* inPts, vtkPointData* inPD, vtkPoints* outPts,
217 vtkPointData* outPD, vtkIdType* ptMap, double tol);
218
219protected:
222
223 // Usual data generation method
226
228 double Tolerance;
236
237 // Internal locator for performing point merging
239
240private:
242 void operator=(const vtkStaticCleanUnstructuredGrid&) = delete;
243};
244
245VTK_ABI_NAMESPACE_END
246#endif
object to represent cell connectivity
represent and manipulate cell attribute data
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
represent and manipulate 3D points
Definition vtkPoints.h:139
Hold a reference to a vtkObjectBase instance.
merge duplicate points, removed unused points, in an vtkUnstructuredGrid
static vtkStaticCleanUnstructuredGrid * New()
Standard methods for instantiation, obtaining type information, and printing the state of the object.
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 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:332
vtkTypeUInt32 vtkMTimeType
Definition vtkType.h:287
#define VTK_DOUBLE_MAX
Definition vtkType.h:171
#define VTK_MARSHALAUTO