VTK  9.2.6
vtkCellTypes.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkCellTypes.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 =========================================================================*/
42 #ifndef vtkCellTypes_h
43 #define vtkCellTypes_h
44 
45 #include "vtkCommonDataModelModule.h" // For export macro
46 #include "vtkObject.h"
47 
48 #include "vtkCellType.h" // Needed for inline methods
49 #include "vtkDeprecation.h" // Needed for deprecation
50 #include "vtkIdTypeArray.h" // Needed for inline methods
51 #include "vtkSmartPointer.h" // Needed for internals
52 #include "vtkUnsignedCharArray.h" // Needed for inline methods
53 
54 class vtkIntArray;
55 
56 class VTKCOMMONDATAMODEL_EXPORT vtkCellTypes : public vtkObject
57 {
58 public:
59  static vtkCellTypes* New();
60  vtkTypeMacro(vtkCellTypes, vtkObject);
61  void PrintSelf(ostream& os, vtkIndent indent) override;
62 
66  int Allocate(vtkIdType sz = 512, vtkIdType ext = 1000);
67 
71  void InsertCell(vtkIdType id, unsigned char type, vtkIdType loc);
72 
76  vtkIdType InsertNextCell(unsigned char type, vtkIdType loc);
77 
83  VTK_DEPRECATED_IN_9_2_0("Please use version without cellLocations.")
84  void SetCellTypes(
85  vtkIdType ncells, vtkUnsignedCharArray* cellTypes, vtkIdTypeArray* cellLocations);
86 
90  void SetCellTypes(vtkIdType ncells, vtkUnsignedCharArray* cellTypes);
91 
92  VTK_DEPRECATED_IN_9_2_0("Please use version without cellLocations.")
93  void SetCellTypes(vtkIdType ncells, vtkUnsignedCharArray* cellTypes, vtkIntArray* cellLocations);
94 
101  VTK_DEPRECATED_IN_9_2_0("The Location API will disappear.")
102  vtkIdType GetCellLocation(vtkIdType cellId) { return this->LocationArray->GetValue(cellId); }
103 
107  void DeleteCell(vtkIdType cellId) { this->TypeArray->SetValue(cellId, VTK_EMPTY_CELL); }
108 
112  vtkIdType GetNumberOfTypes() { return (this->MaxId + 1); }
113 
117  int IsType(unsigned char type);
118 
122  vtkIdType InsertNextType(unsigned char type) { return this->InsertNextCell(type, -1); }
123 
127  unsigned char GetCellType(vtkIdType cellId) { return this->TypeArray->GetValue(cellId); }
128 
132  void Squeeze();
133 
137  void Reset();
138 
147  unsigned long GetActualMemorySize();
148 
153  void DeepCopy(vtkCellTypes* src);
154 
159  static const char* GetClassNameFromTypeId(int typeId);
160 
165  static int GetTypeIdFromClassName(const char* classname);
166 
173  static int IsLinear(unsigned char type);
174 
178  static int GetDimension(unsigned char type);
179 
181 
184  vtkUnsignedCharArray* GetCellTypesArray() { return this->TypeArray; }
185  vtkIdTypeArray* GetCellLocationsArray() { return this->LocationArray; }
187 
188 protected:
189  vtkCellTypes();
190  ~vtkCellTypes() override = default;
191 
193 
194  // DEPRECATION_IN_9_2_0 Note for whoever is in deprecation duties:
195  // The attribute LocationArray needs to be deleted, and any code in this class that would fail
196  // compiling because of its removal deleted as well.
197  vtkSmartPointer<vtkIdTypeArray> LocationArray; // pointer to array of offsets
198 
199  vtkIdType MaxId; // maximum index inserted thus far
200 
201 private:
202  vtkCellTypes(const vtkCellTypes&) = delete;
203  void operator=(const vtkCellTypes&) = delete;
204 };
205 
206 //----------------------------------------------------------------------------
207 inline int vtkCellTypes::IsType(unsigned char type)
208 {
209  vtkIdType numTypes = this->GetNumberOfTypes();
210 
211  for (vtkIdType i = 0; i < numTypes; i++)
212  {
213  if (type == this->GetCellType(i))
214  {
215  return 1;
216  }
217  }
218  return 0;
219 }
220 
221 //-----------------------------------------------------------------------------
222 inline int vtkCellTypes::IsLinear(unsigned char type)
223 {
224  return ((type <= 20) || (type == VTK_CONVEX_POINT_SET) || (type == VTK_POLYHEDRON));
225 }
226 
227 #endif
unsigned char GetCellType(vtkIdType cellId)
Return the type of cell.
Definition: vtkCellTypes.h:127
abstract base class for most VTK objects
Definition: vtkObject.h:62
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
void DeleteCell(vtkIdType cellId)
Delete cell by setting to nullptr cell type.
Definition: vtkCellTypes.h:107
int IsType(unsigned char type)
Return 1 if type specified is contained in list; 0 otherwise.
Definition: vtkCellTypes.h:207
#define VTK_DEPRECATED_IN_9_2_0(reason)
vtkUnsignedCharArray * GetCellTypesArray()
Methods for obtaining the arrays representing types and locations.
Definition: vtkCellTypes.h:184
dynamic, self-adjusting array of vtkIdType
int vtkIdType
Definition: vtkType.h:332
vtkIdType MaxId
Definition: vtkCellTypes.h:199
vtkIdTypeArray * GetCellLocationsArray()
Methods for obtaining the arrays representing types and locations.
Definition: vtkCellTypes.h:185
dynamic, self-adjusting array of int
Definition: vtkIntArray.h:45
a simple class to control print indentation
Definition: vtkIndent.h:39
vtkIdType InsertNextType(unsigned char type)
Add the type specified to the end of the list.
Definition: vtkCellTypes.h:122
vtkSmartPointer< vtkIdTypeArray > LocationArray
Definition: vtkCellTypes.h:197
static int IsLinear(unsigned char type)
This convenience method is a fast check to determine if a cell type represents a linear or nonlinear ...
Definition: vtkCellTypes.h:222
dynamic, self-adjusting array of unsigned char
static vtkObject * New()
Create an object with Debug turned off, modified time initialized to zero, and reference counting on...
object provides direct access to cells in vtkCellArray and type information
Definition: vtkCellTypes.h:56
vtkIdType GetNumberOfTypes()
Return the number of types in the list.
Definition: vtkCellTypes.h:112
vtkSmartPointer< vtkUnsignedCharArray > TypeArray
Definition: vtkCellTypes.h:192