VTK  9.2.6
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 
48 #ifndef vtkMergeCells_h
49 #define vtkMergeCells_h
50 
51 #include "vtkAlgorithm.h" // for vtkAlgorithm::DEFAULT_PRECISION
52 #include "vtkDataSetAttributes.h" // Needed for FieldList
53 #include "vtkFiltersGeneralModule.h" // For export macro
54 #include "vtkObject.h"
55 #include "vtkSmartPointer.h" //fot vtkSmartPointer
56 
57 class vtkCellData;
58 class vtkDataSet;
59 class vtkMergeCellsSTLCloak;
60 class vtkMergePoints;
62 class vtkPointData;
64 
65 class VTKFILTERSGENERAL_EXPORT vtkMergeCells : public vtkObject
66 {
67 public:
68  vtkTypeMacro(vtkMergeCells, vtkObject);
69  void PrintSelf(ostream& os, vtkIndent indent) override;
70 
71  static vtkMergeCells* New();
72 
74 
79  virtual void SetUnstructuredGrid(vtkUnstructuredGrid*);
80  vtkGetObjectMacro(UnstructuredGrid, vtkUnstructuredGrid);
82 
84 
88  vtkSetMacro(TotalNumberOfCells, vtkIdType);
89  vtkGetMacro(TotalNumberOfCells, vtkIdType);
91 
93 
98  vtkSetMacro(TotalNumberOfPoints, vtkIdType);
99  vtkGetMacro(TotalNumberOfPoints, vtkIdType);
101 
103 
109  vtkSetMacro(UseGlobalIds, int);
110  vtkGetMacro(UseGlobalIds, int);
111  vtkBooleanMacro(UseGlobalIds, int);
113 
115 
122  vtkSetClampMacro(PointMergeTolerance, double, 0.0, VTK_DOUBLE_MAX);
123  vtkGetMacro(PointMergeTolerance, double);
125 
127 
131  vtkSetMacro(UseGlobalCellIds, int);
132  vtkGetMacro(UseGlobalCellIds, int);
133  vtkBooleanMacro(UseGlobalCellIds, int);
135 
137 
142  vtkSetMacro(MergeDuplicatePoints, bool);
143  vtkGetMacro(MergeDuplicatePoints, bool);
144  vtkBooleanMacro(MergeDuplicatePoints, bool);
146 
150  void InvalidateCachedLocator();
151 
153 
158  vtkSetMacro(TotalNumberOfDataSets, int);
159  vtkGetMacro(TotalNumberOfDataSets, int);
161 
168  int MergeDataSet(vtkDataSet* set);
169 
171 
176  vtkSetMacro(OutputPointsPrecision, int);
177  vtkGetMacro(OutputPointsPrecision, int);
179 
185  void Finish();
186 
187 protected:
188  vtkMergeCells();
189  ~vtkMergeCells() override;
190 
191  void FreeLists();
192  void StartUGrid(vtkDataSet* set);
193  vtkIdType* MapPointsToIdsUsingGlobalIds(vtkDataSet* set);
194  vtkIdType* MapPointsToIdsUsingLocator(vtkDataSet* set);
195  vtkIdType AddNewCellsUnstructuredGrid(vtkDataSet* set, vtkIdType* idMap);
196  vtkIdType AddNewCellsDataSet(vtkDataSet* set, vtkIdType* idMap);
197 
199 
202 
205 
206  int UseGlobalIds; // point, or node, IDs
207  int UseGlobalCellIds; // cell IDs
208 
211 
212  int OutputPointsPrecision = vtkAlgorithm::DEFAULT_PRECISION;
213 
216 
217  vtkMergeCellsSTLCloak* GlobalIdMap;
218  vtkMergeCellsSTLCloak* GlobalCellIdMap;
219 
222 
224 
225  int NextGrid;
226 
228 
229 private:
230  vtkMergeCells(const vtkMergeCells&) = delete;
231  void operator=(const vtkMergeCells&) = delete;
232 };
233 #endif
abstract base class for most VTK objects
Definition: vtkObject.h:62
represent and manipulate point attribute data
Definition: vtkPointData.h:41
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:165
vtkIdType NumberOfPoints
abstract class to specify dataset behavior
Definition: vtkDataSet.h:62
helps manage arrays from multiple vtkDataSetAttributes.
vtkIdType TotalNumberOfPoints
bool MergeDuplicatePoints
vtkMergeCellsSTLCloak * GlobalCellIdMap
represent and manipulate cell attribute data
Definition: vtkCellData.h:41
Abstract class in support of both point location and point insertion.
int vtkIdType
Definition: vtkType.h:332
vtkIdType NumberOfCells
vtkDataSetAttributes::FieldList * CellList
vtkMergeCellsSTLCloak * GlobalIdMap
merges any number of vtkDataSets back into a single vtkUnstructuredGrid
Definition: vtkMergeCells.h:65
a simple class to control print indentation
Definition: vtkIndent.h:39
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