VTK  9.1.0
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 "vtkAlgorithm.h" // for vtkAlgorithm::DEFAULT_PRECISION
49 #include "vtkDataSetAttributes.h" // Needed for FieldList
50 #include "vtkFiltersGeneralModule.h" // For export macro
51 #include "vtkObject.h"
52 #include "vtkSmartPointer.h" //fot vtkSmartPointer
53 
54 class vtkCellData;
55 class vtkDataSet;
56 class vtkMergeCellsSTLCloak;
57 class vtkMergePoints;
59 class vtkPointData;
61 
62 class VTKFILTERSGENERAL_EXPORT vtkMergeCells : public vtkObject
63 {
64 public:
65  vtkTypeMacro(vtkMergeCells, vtkObject);
66  void PrintSelf(ostream& os, vtkIndent indent) override;
67 
68  static vtkMergeCells* New();
69 
71 
77  vtkGetObjectMacro(UnstructuredGrid, vtkUnstructuredGrid);
79 
81 
85  vtkSetMacro(TotalNumberOfCells, vtkIdType);
86  vtkGetMacro(TotalNumberOfCells, vtkIdType);
88 
90 
95  vtkSetMacro(TotalNumberOfPoints, vtkIdType);
96  vtkGetMacro(TotalNumberOfPoints, vtkIdType);
98 
100 
106  vtkSetMacro(UseGlobalIds, int);
107  vtkGetMacro(UseGlobalIds, int);
108  vtkBooleanMacro(UseGlobalIds, int);
110 
112 
119  vtkSetClampMacro(PointMergeTolerance, double, 0.0, VTK_DOUBLE_MAX);
120  vtkGetMacro(PointMergeTolerance, double);
122 
124 
128  vtkSetMacro(UseGlobalCellIds, int);
129  vtkGetMacro(UseGlobalCellIds, int);
130  vtkBooleanMacro(UseGlobalCellIds, int);
132 
134 
139  vtkSetMacro(MergeDuplicatePoints, bool);
140  vtkGetMacro(MergeDuplicatePoints, bool);
141  vtkBooleanMacro(MergeDuplicatePoints, bool);
143 
148 
150 
155  vtkSetMacro(TotalNumberOfDataSets, int);
156  vtkGetMacro(TotalNumberOfDataSets, int);
158 
166 
168 
173  vtkSetMacro(OutputPointsPrecision, int);
174  vtkGetMacro(OutputPointsPrecision, int);
176 
182  void Finish();
183 
184 protected:
186  ~vtkMergeCells() override;
187 
188  void FreeLists();
189  void StartUGrid(vtkDataSet* set);
194 
196 
199 
202 
203  int UseGlobalIds; // point, or node, IDs
204  int UseGlobalCellIds; // cell IDs
205 
208 
209  int OutputPointsPrecision = vtkAlgorithm::DEFAULT_PRECISION;
210 
213 
214  vtkMergeCellsSTLCloak* GlobalIdMap;
215  vtkMergeCellsSTLCloak* GlobalCellIdMap;
216 
219 
221 
222  int NextGrid;
223 
225 
226 private:
227  vtkMergeCells(const vtkMergeCells&) = delete;
228  void operator=(const vtkMergeCells&) = delete;
229 };
230 #endif
vtkMergeCells::TotalNumberOfCells
vtkIdType TotalNumberOfCells
Definition: vtkMergeCells.h:197
vtkMergeCells::InputIsUGrid
char InputIsUGrid
Definition: vtkMergeCells.h:211
vtkMergeCells::UseGlobalCellIds
int UseGlobalCellIds
Definition: vtkMergeCells.h:204
vtkMergeCells::CellList
vtkDataSetAttributes::FieldList * CellList
Definition: vtkMergeCells.h:218
vtkMergeCells::NextGrid
int NextGrid
Definition: vtkMergeCells.h:222
vtkMergeCells::~vtkMergeCells
~vtkMergeCells() override
vtkPointData
represent and manipulate point attribute data
Definition: vtkPointData.h:142
vtkMergeCells::vtkMergeCells
vtkMergeCells()
vtkIdType
int vtkIdType
Definition: vtkType.h:332
vtkMergeCells::PointMergeTolerance
double PointMergeTolerance
Definition: vtkMergeCells.h:206
vtkMergeCells::GlobalIdMap
vtkMergeCellsSTLCloak * GlobalIdMap
Definition: vtkMergeCells.h:214
vtkMergeCells::TotalNumberOfDataSets
int TotalNumberOfDataSets
Definition: vtkMergeCells.h:195
vtkMergeCells::UseGlobalIds
int UseGlobalIds
Definition: vtkMergeCells.h:203
vtkSmartPointer< vtkIncrementalPointLocator >
vtkObject
abstract base class for most VTK objects
Definition: vtkObject.h:82
vtkMergeCells::InputIsPointSet
char InputIsPointSet
Definition: vtkMergeCells.h:212
vtkMergeCells::MergeDataSet
int MergeDataSet(vtkDataSet *set)
Provide a DataSet to be merged in to the final UnstructuredGrid.
vtkMergeCells::TotalNumberOfPoints
vtkIdType TotalNumberOfPoints
Definition: vtkMergeCells.h:198
vtkAlgorithm::DEFAULT_PRECISION
@ DEFAULT_PRECISION
Definition: vtkAlgorithm.h:153
vtkAlgorithm.h
vtkMergeCells::SetUnstructuredGrid
virtual void SetUnstructuredGrid(vtkUnstructuredGrid *)
Set the vtkUnstructuredGrid object that will become the union of the DataSets specified in MergeDataS...
vtkMergeCells::UnstructuredGrid
vtkUnstructuredGrid * UnstructuredGrid
Definition: vtkMergeCells.h:220
vtkMergeCells::NumberOfCells
vtkIdType NumberOfCells
Definition: vtkMergeCells.h:200
vtkDataSetAttributesFieldList
helps manage arrays from multiple vtkDataSetAttributes.
Definition: vtkDataSetAttributesFieldList.h:70
vtkCellData
represent and manipulate cell attribute data
Definition: vtkCellData.h:142
vtkMergeCells::PrintSelf
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
vtkIndent
a simple class to control print indentation
Definition: vtkIndent.h:113
vtkMergeCells
merges any number of vtkDataSets back into a single vtkUnstructuredGrid
Definition: vtkMergeCells.h:63
vtkSmartPointer.h
vtkIncrementalPointLocator
Abstract class in support of both point location and point insertion.
Definition: vtkIncrementalPointLocator.h:52
vtkMergeCells::NumberOfPoints
vtkIdType NumberOfPoints
Definition: vtkMergeCells.h:201
vtkMergeCells::FreeLists
void FreeLists()
vtkObject.h
vtkDataSet
abstract class to specify dataset behavior
Definition: vtkDataSet.h:166
vtkMergeCells::AddNewCellsDataSet
vtkIdType AddNewCellsDataSet(vtkDataSet *set, vtkIdType *idMap)
vtkMergeCells::MapPointsToIdsUsingLocator
vtkIdType * MapPointsToIdsUsingLocator(vtkDataSet *set)
vtkMergeCells::Locator
vtkSmartPointer< vtkIncrementalPointLocator > Locator
Definition: vtkMergeCells.h:224
vtkMergeCells::New
static vtkMergeCells * New()
vtkMergeCells::Finish
void Finish()
Call Finish() after merging last DataSet to free unneeded memory and to make sure the ugrid's GetNumb...
vtkMergeCells::MergeDuplicatePoints
bool MergeDuplicatePoints
Definition: vtkMergeCells.h:207
vtkUnstructuredGrid
dataset represents arbitrary combinations of all possible cell types
Definition: vtkUnstructuredGrid.h:204
vtkMergeCells::InvalidateCachedLocator
void InvalidateCachedLocator()
Clear the Locator and set it to nullptr.
vtkMergeCells::GlobalCellIdMap
vtkMergeCellsSTLCloak * GlobalCellIdMap
Definition: vtkMergeCells.h:215
vtkMergeCells::MapPointsToIdsUsingGlobalIds
vtkIdType * MapPointsToIdsUsingGlobalIds(vtkDataSet *set)
VTK_DOUBLE_MAX
#define VTK_DOUBLE_MAX
Definition: vtkType.h:165
vtkMergeCells::AddNewCellsUnstructuredGrid
vtkIdType AddNewCellsUnstructuredGrid(vtkDataSet *set, vtkIdType *idMap)
vtkDataSetAttributes.h
vtkMergeCells::StartUGrid
void StartUGrid(vtkDataSet *set)
vtkMergePoints
merge exactly coincident points
Definition: vtkMergePoints.h:63
vtkMergeCells::PointList
vtkDataSetAttributes::FieldList * PointList
Definition: vtkMergeCells.h:217