VTK  9.3.1
vtkMergeCells.h
Go to the documentation of this file.
1 // SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2 // SPDX-FileCopyrightText: Copyright (c) Sandia Corporation
3 // SPDX-License-Identifier: BSD-3-Clause
4 
33 #ifndef vtkMergeCells_h
34 #define vtkMergeCells_h
35 
36 #include "vtkAlgorithm.h" // for vtkAlgorithm::DEFAULT_PRECISION
37 #include "vtkDataSetAttributes.h" // Needed for FieldList
38 #include "vtkFiltersGeneralModule.h" // For export macro
39 #include "vtkObject.h"
40 #include "vtkSmartPointer.h" // for vtkSmartPointer
41 
42 VTK_ABI_NAMESPACE_BEGIN
43 class vtkCellData;
44 class vtkDataSet;
45 class vtkMergeCellsSTLCloak;
46 class vtkMergePoints;
48 class vtkPointData;
50 
51 class VTKFILTERSGENERAL_EXPORT vtkMergeCells : public vtkObject
52 {
53 public:
54  vtkTypeMacro(vtkMergeCells, vtkObject);
55  void PrintSelf(ostream& os, vtkIndent indent) override;
56 
57  static vtkMergeCells* New();
58 
60 
65  virtual void SetUnstructuredGrid(vtkUnstructuredGrid*);
66  vtkGetObjectMacro(UnstructuredGrid, vtkUnstructuredGrid);
68 
70 
74  vtkSetMacro(TotalNumberOfCells, vtkIdType);
75  vtkGetMacro(TotalNumberOfCells, vtkIdType);
77 
79 
84  vtkSetMacro(TotalNumberOfPoints, vtkIdType);
85  vtkGetMacro(TotalNumberOfPoints, vtkIdType);
87 
89 
95  vtkSetMacro(UseGlobalIds, int);
96  vtkGetMacro(UseGlobalIds, int);
97  vtkBooleanMacro(UseGlobalIds, int);
99 
101 
108  vtkSetClampMacro(PointMergeTolerance, double, 0.0, VTK_DOUBLE_MAX);
109  vtkGetMacro(PointMergeTolerance, double);
111 
113 
117  vtkSetMacro(UseGlobalCellIds, int);
118  vtkGetMacro(UseGlobalCellIds, int);
119  vtkBooleanMacro(UseGlobalCellIds, int);
121 
123 
128  vtkSetMacro(MergeDuplicatePoints, bool);
129  vtkGetMacro(MergeDuplicatePoints, bool);
130  vtkBooleanMacro(MergeDuplicatePoints, bool);
132 
136  void InvalidateCachedLocator();
137 
139 
144  vtkSetMacro(TotalNumberOfDataSets, int);
145  vtkGetMacro(TotalNumberOfDataSets, int);
147 
154  int MergeDataSet(vtkDataSet* set);
155 
157 
162  vtkSetMacro(OutputPointsPrecision, int);
163  vtkGetMacro(OutputPointsPrecision, int);
165 
171  void Finish();
172 
173 protected:
174  vtkMergeCells();
175  ~vtkMergeCells() override;
176 
177  void FreeLists();
178  void StartUGrid(vtkDataSet* set);
179  vtkIdType* MapPointsToIdsUsingGlobalIds(vtkDataSet* set);
180  vtkIdType* MapPointsToIdsUsingLocator(vtkDataSet* set);
181  vtkIdType AddNewCellsUnstructuredGrid(vtkDataSet* set, vtkIdType* idMap);
182  vtkIdType AddNewCellsDataSet(vtkDataSet* set, vtkIdType* idMap);
183 
185 
188 
191 
192  int UseGlobalIds; // point, or node, IDs
193  int UseGlobalCellIds; // cell IDs
194 
197 
198  int OutputPointsPrecision = vtkAlgorithm::DEFAULT_PRECISION;
199 
202 
203  vtkMergeCellsSTLCloak* GlobalIdMap;
204  vtkMergeCellsSTLCloak* GlobalCellIdMap;
205 
208 
210 
211  int NextGrid;
212 
214 
215 private:
216  vtkMergeCells(const vtkMergeCells&) = delete;
217  void operator=(const vtkMergeCells&) = delete;
218 };
219 VTK_ABI_NAMESPACE_END
220 #endif
abstract base class for most VTK objects
Definition: vtkObject.h:51
represent and manipulate point attribute data
Definition: vtkPointData.h:29
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
#define VTK_DOUBLE_MAX
Definition: vtkType.h:154
vtkIdType NumberOfPoints
abstract class to specify dataset behavior
Definition: vtkDataSet.h:52
helps manage arrays from multiple vtkDataSetAttributes.
vtkIdType TotalNumberOfPoints
bool MergeDuplicatePoints
vtkMergeCellsSTLCloak * GlobalCellIdMap
represent and manipulate cell attribute data
Definition: vtkCellData.h:30
Abstract class in support of both point location and point insertion.
int vtkIdType
Definition: vtkType.h:315
vtkIdType NumberOfCells
vtkDataSetAttributes::FieldList * CellList
vtkMergeCellsSTLCloak * GlobalIdMap
merges any number of vtkDataSets back into a single vtkUnstructuredGrid
Definition: vtkMergeCells.h:51
a simple class to control print indentation
Definition: vtkIndent.h:28
merge exactly coincident points
dataset represents arbitrary combinations of all possible cell types
int TotalNumberOfDataSets
vtkIdType TotalNumberOfCells
double PointMergeTolerance
vtkSmartPointer< vtkIncrementalPointLocator > Locator
static vtkObject * New()
Create an object with Debug turned off, modified time initialized to zero, and reference counting on...
vtkDataSetAttributes::FieldList * PointList
vtkUnstructuredGrid * UnstructuredGrid