VTK  9.3.0
vtkCellTypes.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-License-Identifier: BSD-3-Clause
39#ifndef vtkCellTypes_h
40#define vtkCellTypes_h
41
42#include "vtkCommonDataModelModule.h" // For export macro
43#include "vtkObject.h"
44
45#include "vtkCellType.h" // Needed for inline methods
46#include "vtkDeprecation.h" // For VTK_DEPRECATED_IN_9_2_0
47#include "vtkIdTypeArray.h" // Needed for inline methods
48#include "vtkSmartPointer.h" // Needed for internals
49#include "vtkUnsignedCharArray.h" // Needed for inline methods
50
51VTK_ABI_NAMESPACE_BEGIN
52class vtkIntArray;
53
54class VTKCOMMONDATAMODEL_EXPORT vtkCellTypes : public vtkObject
55{
56public:
57 static vtkCellTypes* New();
58 vtkTypeMacro(vtkCellTypes, vtkObject);
59 void PrintSelf(ostream& os, vtkIndent indent) override;
60
64 int Allocate(vtkIdType sz = 512, vtkIdType ext = 1000);
65
69 void InsertCell(vtkIdType id, unsigned char type, vtkIdType loc);
70
74 vtkIdType InsertNextCell(unsigned char type, vtkIdType loc);
75
81 VTK_DEPRECATED_IN_9_2_0("Please use version without cellLocations.")
82 void SetCellTypes(
83 vtkIdType ncells, vtkUnsignedCharArray* cellTypes, vtkIdTypeArray* cellLocations);
84
88 void SetCellTypes(vtkIdType ncells, vtkUnsignedCharArray* cellTypes);
89
90 VTK_DEPRECATED_IN_9_2_0("Please use version without cellLocations.")
91 void SetCellTypes(vtkIdType ncells, vtkUnsignedCharArray* cellTypes, vtkIntArray* cellLocations);
92
99 VTK_DEPRECATED_IN_9_2_0("The Location API will disappear.")
100 vtkIdType GetCellLocation(vtkIdType cellId) { return this->LocationArray->GetValue(cellId); }
101
105 void DeleteCell(vtkIdType cellId) { this->TypeArray->SetValue(cellId, VTK_EMPTY_CELL); }
106
110 vtkIdType GetNumberOfTypes() { return (this->MaxId + 1); }
111
115 int IsType(unsigned char type);
116
120 vtkIdType InsertNextType(unsigned char type) { return this->InsertNextCell(type, -1); }
121
125 unsigned char GetCellType(vtkIdType cellId) { return this->TypeArray->GetValue(cellId); }
126
130 void Squeeze();
131
135 void Reset();
136
145 unsigned long GetActualMemorySize();
146
152
157 static const char* GetClassNameFromTypeId(int typeId);
158
163 static int GetTypeIdFromClassName(const char* classname);
164
171 static int IsLinear(unsigned char type);
172
176 static int GetDimension(unsigned char type);
177
179
182 vtkUnsignedCharArray* GetCellTypesArray() { return this->TypeArray; }
183 vtkIdTypeArray* GetCellLocationsArray() { return this->LocationArray; }
185
186protected:
188 ~vtkCellTypes() override = default;
189
191
192 // DEPRECATION_IN_9_2_0 Note for whoever is in deprecation duties:
193 // The attribute LocationArray needs to be deleted, and any code in this class that would fail
194 // compiling because of its removal deleted as well.
195 vtkSmartPointer<vtkIdTypeArray> LocationArray; // pointer to array of offsets
196
197 vtkIdType MaxId; // maximum index inserted thus far
198
199private:
200 vtkCellTypes(const vtkCellTypes&) = delete;
201 void operator=(const vtkCellTypes&) = delete;
202};
203
204//----------------------------------------------------------------------------
205inline int vtkCellTypes::IsType(unsigned char type)
206{
207 vtkIdType numTypes = this->GetNumberOfTypes();
208
209 for (vtkIdType i = 0; i < numTypes; i++)
210 {
211 if (type == this->GetCellType(i))
212 {
213 return 1;
214 }
215 }
216 return 0;
217}
218
219//-----------------------------------------------------------------------------
220inline int vtkCellTypes::IsLinear(unsigned char type)
221{
222 return ((type <= 20) || (type == VTK_CONVEX_POINT_SET) || (type == VTK_POLYHEDRON));
223}
224
225VTK_ABI_NAMESPACE_END
226#endif
object provides direct access to cells in vtkCellArray and type information
vtkIdType InsertNextCell(unsigned char type, vtkIdType loc)
Add a cell to the object in the next available slot.
int Allocate(vtkIdType sz=512, vtkIdType ext=1000)
Allocate memory for this array.
void InsertCell(vtkIdType id, unsigned char type, vtkIdType loc)
Add a cell at specified id.
vtkIdType MaxId
void Squeeze()
Reclaim any extra memory.
void Reset()
Initialize object without releasing memory.
vtkIdType InsertNextType(unsigned char type)
Add the type specified to the end of the list.
~vtkCellTypes() override=default
vtkSmartPointer< vtkIdTypeArray > LocationArray
static vtkCellTypes * New()
static int IsLinear(unsigned char type)
This convenience method is a fast check to determine if a cell type represents a linear or nonlinear ...
unsigned long GetActualMemorySize()
Return the memory in kibibytes (1024 bytes) consumed by this cell type array.
vtkIdType GetNumberOfTypes()
Return the number of types in the list.
vtkSmartPointer< vtkUnsignedCharArray > TypeArray
void DeleteCell(vtkIdType cellId)
Delete cell by setting to nullptr cell type.
static const char * GetClassNameFromTypeId(int typeId)
Given an int (as defined in vtkCellType.h) identifier for a class return it's classname.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
static int GetDimension(unsigned char type)
Get the dimension of a cell.
int IsType(unsigned char type)
Return 1 if type specified is contained in list; 0 otherwise.
unsigned char GetCellType(vtkIdType cellId)
Return the type of cell.
void DeepCopy(vtkCellTypes *src)
Standard DeepCopy method.
static int GetTypeIdFromClassName(const char *classname)
Given a data object classname, return it's int identified (as defined in vtkCellType....
vtkUnsignedCharArray * GetCellTypesArray()
Methods for obtaining the arrays representing types and locations.
vtkIdTypeArray * GetCellLocationsArray()
Methods for obtaining the arrays representing types and locations.
dynamic, self-adjusting array of vtkIdType
a simple class to control print indentation
Definition vtkIndent.h:38
dynamic, self-adjusting array of int
Definition vtkIntArray.h:44
abstract base class for most VTK objects
Definition vtkObject.h:58
Hold a reference to a vtkObjectBase instance.
dynamic, self-adjusting array of unsigned char
@ VTK_EMPTY_CELL
Definition vtkCellType.h:37
@ VTK_POLYHEDRON
Definition vtkCellType.h:80
@ VTK_CONVEX_POINT_SET
Definition vtkCellType.h:77
#define VTK_DEPRECATED_IN_9_2_0(reason)
int vtkIdType
Definition vtkType.h:315