VTK  9.3.20240423
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
42VTK_ABI_NAMESPACE_BEGIN
43class vtkCellData;
44class vtkDataSet;
45class vtkMergeCellsSTLCloak;
46class vtkMergePoints;
48class vtkPointData;
50
51class VTKFILTERSGENERAL_EXPORT vtkMergeCells : public vtkObject
52{
53public:
54 vtkTypeMacro(vtkMergeCells, vtkObject);
55 void PrintSelf(ostream& os, vtkIndent indent) override;
56
57 static vtkMergeCells* New();
58
60
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
137
139
144 vtkSetMacro(TotalNumberOfDataSets, int);
145 vtkGetMacro(TotalNumberOfDataSets, int);
147
155
157
162 vtkSetMacro(OutputPointsPrecision, int);
163 vtkGetMacro(OutputPointsPrecision, int);
165
171 void Finish();
172
173protected:
175 ~vtkMergeCells() override;
176
177 void FreeLists();
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
212
214
215private:
216 vtkMergeCells(const vtkMergeCells&) = delete;
217 void operator=(const vtkMergeCells&) = delete;
218};
219VTK_ABI_NAMESPACE_END
220#endif
represent and manipulate cell attribute data
helps manage arrays from multiple vtkDataSetAttributes.
abstract class to specify dataset behavior
Definition vtkDataSet.h:166
Abstract class in support of both point location and point insertion.
a simple class to control print indentation
Definition vtkIndent.h:108
merges any number of vtkDataSets back into a single vtkUnstructuredGrid
vtkIdType TotalNumberOfCells
vtkIdType TotalNumberOfPoints
vtkDataSetAttributes::FieldList * PointList
vtkIdType NumberOfCells
int MergeDataSet(vtkDataSet *set)
Provide a DataSet to be merged in to the final UnstructuredGrid.
vtkIdType AddNewCellsUnstructuredGrid(vtkDataSet *set, vtkIdType *idMap)
void InvalidateCachedLocator()
Clear the Locator and set it to nullptr.
vtkUnstructuredGrid * UnstructuredGrid
vtkSmartPointer< vtkIncrementalPointLocator > Locator
vtkIdType NumberOfPoints
void Finish()
Call Finish() after merging last DataSet to free unneeded memory and to make sure the ugrid's GetNumb...
static vtkMergeCells * New()
void StartUGrid(vtkDataSet *set)
vtkMergeCellsSTLCloak * GlobalCellIdMap
vtkDataSetAttributes::FieldList * CellList
void FreeLists()
vtkIdType AddNewCellsDataSet(vtkDataSet *set, vtkIdType *idMap)
vtkIdType * MapPointsToIdsUsingLocator(vtkDataSet *set)
~vtkMergeCells() override
bool MergeDuplicatePoints
vtkMergeCellsSTLCloak * GlobalIdMap
double PointMergeTolerance
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkIdType * MapPointsToIdsUsingGlobalIds(vtkDataSet *set)
virtual void SetUnstructuredGrid(vtkUnstructuredGrid *)
Set the vtkUnstructuredGrid object that will become the union of the DataSets specified in MergeDataS...
merge exactly coincident points
abstract base class for most VTK objects
Definition vtkObject.h:162
represent and manipulate point attribute data
Hold a reference to a vtkObjectBase instance.
dataset represents arbitrary combinations of all possible cell types
int vtkIdType
Definition vtkType.h:315
#define VTK_DOUBLE_MAX
Definition vtkType.h:154