VTK  9.2.6
vtkStaticCleanUnstructuredGrid.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkStaticCleanUnstructuredGrid.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 =========================================================================*/
73 #ifndef vtkStaticCleanUnstructuredGrid_h
74 #define vtkStaticCleanUnstructuredGrid_h
75 
76 #include "vtkFiltersCoreModule.h" // For export macro
77 #include "vtkSmartPointer.h" // Pointer to locator
78 #include "vtkStaticPointLocator.h" // Locator
80 #include <vector> // API to vtkStaticCleanUnstructuredGrid
81 
82 class vtkCellArray;
83 class vtkCellData;
84 class vtkPointData;
85 class vtkPoints;
86 
88 {
89 public:
91 
96  void PrintSelf(ostream& os, vtkIndent indent) override;
99 
101 
107  vtkSetMacro(ToleranceIsAbsolute, bool);
108  vtkBooleanMacro(ToleranceIsAbsolute, bool);
109  vtkGetMacro(ToleranceIsAbsolute, bool);
111 
113 
117  vtkSetClampMacro(AbsoluteTolerance, double, 0.0, VTK_DOUBLE_MAX);
118  vtkGetMacro(AbsoluteTolerance, double);
120 
122 
127  vtkSetClampMacro(Tolerance, double, 0.0, 1.0);
128  vtkGetMacro(Tolerance, double);
130 
132 
143  vtkSetStringMacro(MergingArray);
144  vtkGetStringMacro(MergingArray);
146 
148 
155  vtkSetMacro(RemoveUnusedPoints, bool);
156  vtkBooleanMacro(RemoveUnusedPoints, bool);
157  vtkGetMacro(RemoveUnusedPoints, bool);
159 
161 
168  vtkSetMacro(ProduceMergeMap, bool);
169  vtkBooleanMacro(ProduceMergeMap, bool);
170  vtkGetMacro(ProduceMergeMap, bool);
172 
174 
182  vtkSetMacro(AveragePointData, bool);
183  vtkBooleanMacro(AveragePointData, bool);
184  vtkGetMacro(AveragePointData, bool);
186 
188 
193  vtkSetMacro(OutputPointsPrecision, int);
194  vtkGetMacro(OutputPointsPrecision, int);
196 
202  vtkGetObjectMacro(Locator, vtkStaticPointLocator);
203 
205  // This filter is difficult to stream. To produce invariant results, the
206  // whole input must be processed at once. This flag allows the user to
207  // select whether strict piece invariance is required. By default it is
208  // on. When off, the filter can stream, but results may change.
209  vtkSetMacro(PieceInvariant, bool);
210  vtkGetMacro(PieceInvariant, bool);
211  vtkBooleanMacro(PieceInvariant, bool);
213 
217  vtkMTimeType GetMTime() override;
218 
219 protected:
221  ~vtkStaticCleanUnstructuredGrid() override = default;
222 
223  // Usual data generation method
226 
228  double Tolerance;
236 
237  // Internal locator for performing point merging
239 
240  // These methods are made static so that vtkStaticCleanPolyData can use them
242  static void MarkPointUses(vtkCellArray* ca, vtkIdType* mergeMap, unsigned char* ptUses);
244  vtkIdType numPts, vtkIdType* pmap, unsigned char* ptUses, std::vector<vtkIdType>& mergeMap);
245  static void CopyPoints(
246  vtkPoints* inPts, vtkPointData* inPD, vtkPoints* outPts, vtkPointData* outPD, vtkIdType* ptMap);
247  static void AveragePoints(vtkPoints* inPts, vtkPointData* inPD, vtkPoints* outPts,
248  vtkPointData* outPD, vtkIdType* ptMap, double tol);
249 
250 private:
252  void operator=(const vtkStaticCleanUnstructuredGrid&) = delete;
253 };
254 
255 #endif
object to represent cell connectivity
Definition: vtkCellArray.h:296
represent and manipulate cell attribute data
Definition: vtkCellData.h:151
a simple class to control print indentation
Definition: vtkIndent.h:119
Store zero or more vtkInformation instances.
Store vtkAlgorithm input/output information.
represent and manipulate point attribute data
Definition: vtkPointData.h:151
represent and manipulate 3D points
Definition: vtkPoints.h:149
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:332
vtkTypeUInt32 vtkMTimeType
Definition: vtkType.h:287
#define VTK_DOUBLE_MAX
Definition: vtkType.h:165