VTK  9.0.1
vtkMergeCells.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkMergeCells.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 =========================================================================*/
15 /*----------------------------------------------------------------------------
16  Copyright (c) Sandia Corporation
17  See Copyright.txt or http://www.paraview.org/HTML/Copyright.html for details.
18 ----------------------------------------------------------------------------*/
19 
45 #ifndef vtkMergeCells_h
46 #define vtkMergeCells_h
47 
48 #include "vtkDataSetAttributes.h" // Needed for FieldList
49 #include "vtkFiltersGeneralModule.h" // For export macro
50 #include "vtkObject.h"
51 #include "vtkSmartPointer.h" //fot vtkSmartPointer
52 
53 class vtkCellData;
54 class vtkDataSet;
55 class vtkMergeCellsSTLCloak;
56 class vtkMergePoints;
57 class vtkPointData;
59 
60 class VTKFILTERSGENERAL_EXPORT vtkMergeCells : public vtkObject
61 {
62 public:
63  vtkTypeMacro(vtkMergeCells, vtkObject);
64  void PrintSelf(ostream& os, vtkIndent indent) override;
65 
66  static vtkMergeCells* New();
67 
69 
75  vtkGetObjectMacro(UnstructuredGrid, vtkUnstructuredGrid);
77 
79 
83  vtkSetMacro(TotalNumberOfCells, vtkIdType);
84  vtkGetMacro(TotalNumberOfCells, vtkIdType);
86 
88 
93  vtkSetMacro(TotalNumberOfPoints, vtkIdType);
94  vtkGetMacro(TotalNumberOfPoints, vtkIdType);
96 
98 
104  vtkSetMacro(UseGlobalIds, int);
105  vtkGetMacro(UseGlobalIds, int);
106  vtkBooleanMacro(UseGlobalIds, int);
108 
110 
117  vtkSetClampMacro(PointMergeTolerance, float, 0.0, VTK_FLOAT_MAX);
118  vtkGetMacro(PointMergeTolerance, float);
120 
122 
126  vtkSetMacro(UseGlobalCellIds, int);
127  vtkGetMacro(UseGlobalCellIds, int);
128  vtkBooleanMacro(UseGlobalCellIds, int);
130 
132 
137  vtkSetMacro(MergeDuplicatePoints, bool);
138  vtkGetMacro(MergeDuplicatePoints, bool);
139  vtkBooleanMacro(MergeDuplicatePoints, bool);
141 
146 
148 
153  vtkSetMacro(TotalNumberOfDataSets, int);
154  vtkGetMacro(TotalNumberOfDataSets, int);
156 
164 
170  void Finish();
171 
172 protected:
174  ~vtkMergeCells() override;
175 
176  void FreeLists();
177  void StartUGrid(vtkDataSet* set);
182 
184 
187 
190 
191  int UseGlobalIds; // point, or node, IDs
192  int UseGlobalCellIds; // cell IDs
193 
196 
199 
200  vtkMergeCellsSTLCloak* GlobalIdMap;
201  vtkMergeCellsSTLCloak* GlobalCellIdMap;
202 
205 
207 
208  int NextGrid;
209 
211 
212 private:
213  vtkMergeCells(const vtkMergeCells&) = delete;
214  void operator=(const vtkMergeCells&) = delete;
215 };
216 #endif
represent and manipulate cell attribute data
Definition: vtkCellData.h:33
helps manage arrays from multiple vtkDataSetAttributes.
abstract class to specify dataset behavior
Definition: vtkDataSet.h:57
a simple class to control print indentation
Definition: vtkIndent.h:34
merges any number of vtkDataSets back into a single vtkUnstructuredGrid
Definition: vtkMergeCells.h:61
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)
int TotalNumberOfDataSets
void InvalidateCachedLocator()
Clear the Locator and set it to nullptr.
vtkUnstructuredGrid * UnstructuredGrid
vtkIdType NumberOfPoints
void Finish()
Call Finish() after merging last DataSet to free unneeded memory and to make sure the ugrid's GetNumb...
vtkIdType * MapPointsToIdsUsingGlobalIds(vtkDataSet *set)
float PointMergeTolerance
void StartUGrid(vtkDataSet *set)
static vtkMergeCells * New()
vtkMergeCellsSTLCloak * GlobalCellIdMap
vtkIdType * MapPointsToIdsUsingLocator(vtkDataSet *set)
vtkDataSetAttributes::FieldList * CellList
void FreeLists()
vtkIdType AddNewCellsDataSet(vtkDataSet *set, vtkIdType *idMap)
~vtkMergeCells() override
bool MergeDuplicatePoints
vtkMergeCellsSTLCloak * GlobalIdMap
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkSmartPointer< vtkMergePoints > Locator
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:54
represent and manipulate point attribute data
Definition: vtkPointData.h:32
dataset represents arbitrary combinations of all possible cell types
int vtkIdType
Definition: vtkType.h:338
#define VTK_FLOAT_MAX
Definition: vtkType.h:163