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
30#ifndef vtkCellTypes_h
31#define vtkCellTypes_h
32
33#include "vtkCommonDataModelModule.h" // For export macro
34#include "vtkObject.h"
35
36#include "vtkCellType.h" // Needed for inline methods
37#include "vtkDeprecation.h" // For VTK_DEPRECATED_IN_9_2_0
38#include "vtkIdTypeArray.h" // Needed for inline methods
39#include "vtkSmartPointer.h" // Needed for internals
40#include "vtkUnsignedCharArray.h" // Needed for inline methods
41
42VTK_ABI_NAMESPACE_BEGIN
43class vtkIntArray;
44
45class VTKCOMMONDATAMODEL_EXPORT vtkCellTypes : public vtkObject
46{
47public:
48 static vtkCellTypes* New();
49 vtkTypeMacro(vtkCellTypes, vtkObject);
50 void PrintSelf(ostream& os, vtkIndent indent) override;
51
55 int Allocate(vtkIdType sz = 512, vtkIdType ext = 1000);
56
60 void InsertCell(vtkIdType id, unsigned char type, vtkIdType loc);
61
65 vtkIdType InsertNextCell(unsigned char type, vtkIdType loc);
66
72 VTK_DEPRECATED_IN_9_2_0("Please use version without cellLocations.")
73 void SetCellTypes(
74 vtkIdType ncells, vtkUnsignedCharArray* cellTypes, vtkIdTypeArray* cellLocations);
75
79 void SetCellTypes(vtkIdType ncells, vtkUnsignedCharArray* cellTypes);
80
81 VTK_DEPRECATED_IN_9_2_0("Please use version without cellLocations.")
82 void SetCellTypes(vtkIdType ncells, vtkUnsignedCharArray* cellTypes, vtkIntArray* cellLocations);
83
90 VTK_DEPRECATED_IN_9_2_0("The Location API will disappear.")
91 vtkIdType GetCellLocation(vtkIdType cellId) { return this->LocationArray->GetValue(cellId); }
92
96 void DeleteCell(vtkIdType cellId) { this->TypeArray->SetValue(cellId, VTK_EMPTY_CELL); }
97
101 vtkIdType GetNumberOfTypes() { return (this->MaxId + 1); }
102
106 int IsType(unsigned char type);
107
111 vtkIdType InsertNextType(unsigned char type) { return this->InsertNextCell(type, -1); }
112
116 unsigned char GetCellType(vtkIdType cellId) { return this->TypeArray->GetValue(cellId); }
117
121 void Squeeze();
122
126 void Reset();
127
136 unsigned long GetActualMemorySize();
137
143
148 static const char* GetClassNameFromTypeId(int typeId);
149
154 static int GetTypeIdFromClassName(const char* classname);
155
162 static int IsLinear(unsigned char type);
163
167 static int GetDimension(unsigned char type);
168
170
173 vtkUnsignedCharArray* GetCellTypesArray() { return this->TypeArray; }
174 vtkIdTypeArray* GetCellLocationsArray() { return this->LocationArray; }
176
177protected:
179 ~vtkCellTypes() override = default;
180
182
183 // DEPRECATION_IN_9_2_0 Note for whoever is in deprecation duties:
184 // The attribute LocationArray needs to be deleted, and any code in this class that would fail
185 // compiling because of its removal deleted as well.
186 vtkSmartPointer<vtkIdTypeArray> LocationArray; // pointer to array of offsets
187
188 vtkIdType MaxId; // maximum index inserted thus far
189
190private:
191 vtkCellTypes(const vtkCellTypes&) = delete;
192 void operator=(const vtkCellTypes&) = delete;
193};
194
195//----------------------------------------------------------------------------
196inline int vtkCellTypes::IsType(unsigned char type)
197{
198 vtkIdType numTypes = this->GetNumberOfTypes();
199
200 for (vtkIdType i = 0; i < numTypes; i++)
201 {
202 if (type == this->GetCellType(i))
203 {
204 return 1;
205 }
206 }
207 return 0;
208}
209
210//-----------------------------------------------------------------------------
211inline int vtkCellTypes::IsLinear(unsigned char type)
212{
213 return ((type <= 20) || (type == VTK_CONVEX_POINT_SET) || (type == VTK_POLYHEDRON));
214}
215
216VTK_ABI_NAMESPACE_END
217#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:29
dynamic, self-adjusting array of int
Definition vtkIntArray.h:35
abstract base class for most VTK objects
Definition vtkObject.h:49
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