VTK  9.1.0
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=========================================================================*/
125#ifndef vtkCellTypes_h
126#define vtkCellTypes_h
127
128#include "vtkCommonDataModelModule.h" // For export macro
129#include "vtkObject.h"
130
131#include "vtkCellType.h" // Needed for VTK_EMPTY_CELL
132#include "vtkIdTypeArray.h" // Needed for inline methods
133#include "vtkIntArray.h" // Needed for inline methods
134#include "vtkUnsignedCharArray.h" // Needed for inline methods
135
136class VTKCOMMONDATAMODEL_EXPORT vtkCellTypes : public vtkObject
137{
138public:
139 static vtkCellTypes* New();
140 vtkTypeMacro(vtkCellTypes, vtkObject);
141 void PrintSelf(ostream& os, vtkIndent indent) override;
142
146 int Allocate(vtkIdType sz = 512, vtkIdType ext = 1000);
147
151 void InsertCell(vtkIdType id, unsigned char type, vtkIdType loc);
152
156 vtkIdType InsertNextCell(unsigned char type, vtkIdType loc);
157
162 vtkIdType ncells, vtkUnsignedCharArray* cellTypes, vtkIdTypeArray* cellLocations);
163
168 void SetCellTypes(vtkIdType ncells, vtkUnsignedCharArray* cellTypes, vtkIntArray* cellLocations);
169
173 vtkIdType GetCellLocation(vtkIdType cellId) { return this->LocationArray->GetValue(cellId); }
174
178 void DeleteCell(vtkIdType cellId) { this->TypeArray->SetValue(cellId, VTK_EMPTY_CELL); }
179
183 vtkIdType GetNumberOfTypes() { return (this->MaxId + 1); }
184
188 int IsType(unsigned char type);
189
193 vtkIdType InsertNextType(unsigned char type) { return this->InsertNextCell(type, -1); }
194
198 unsigned char GetCellType(vtkIdType cellId) { return this->TypeArray->GetValue(cellId); }
199
203 void Squeeze();
204
208 void Reset();
209
218 unsigned long GetActualMemorySize();
219
225
230 static const char* GetClassNameFromTypeId(int typeId);
231
236 static int GetTypeIdFromClassName(const char* classname);
237
244 static int IsLinear(unsigned char type);
245
247
250 vtkUnsignedCharArray* GetCellTypesArray() { return this->TypeArray; }
251 vtkIdTypeArray* GetCellLocationsArray() { return this->LocationArray; }
253
254protected:
256 ~vtkCellTypes() override;
257
258 vtkUnsignedCharArray* TypeArray; // pointer to types array
259 vtkIdTypeArray* LocationArray; // pointer to array of offsets
260 vtkIdType Size; // allocated size of data
261 vtkIdType MaxId; // maximum index inserted thus far
262 vtkIdType Extend; // grow array by this point
263
264private:
265 vtkCellTypes(const vtkCellTypes&) = delete;
266 void operator=(const vtkCellTypes&) = delete;
267};
268
269//----------------------------------------------------------------------------
270inline int vtkCellTypes::IsType(unsigned char type)
271{
272 vtkIdType numTypes = this->GetNumberOfTypes();
273
274 for (vtkIdType i = 0; i < numTypes; i++)
275 {
276 if (type == this->GetCellType(i))
277 {
278 return 1;
279 }
280 }
281 return 0;
282}
283
284//-----------------------------------------------------------------------------
285inline int vtkCellTypes::IsLinear(unsigned char type)
286{
287 return ((type <= 20) || (type == VTK_CONVEX_POINT_SET) || (type == VTK_POLYHEDRON));
288}
289
290#endif
object provides direct access to cells in vtkCellArray and type information
Definition: vtkCellTypes.h:137
vtkIdType InsertNextCell(unsigned char type, vtkIdType loc)
Add a cell to the object in the next available slot.
vtkIdTypeArray * LocationArray
Definition: vtkCellTypes.h:259
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
Definition: vtkCellTypes.h:261
vtkIdType GetCellLocation(vtkIdType cellId)
Return the location of the cell in the associated vtkCellArray.
Definition: vtkCellTypes.h:173
void Squeeze()
Reclaim any extra memory.
void SetCellTypes(vtkIdType ncells, vtkUnsignedCharArray *cellTypes, vtkIdTypeArray *cellLocations)
Specify a group of cell types.
void Reset()
Initialize object without releasing memory.
vtkIdType InsertNextType(unsigned char type)
Add the type specified to the end of the list.
Definition: vtkCellTypes.h:193
~vtkCellTypes() override
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 ...
Definition: vtkCellTypes.h:285
void SetCellTypes(vtkIdType ncells, vtkUnsignedCharArray *cellTypes, vtkIntArray *cellLocations)
Specify a group of cell types.
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.
Definition: vtkCellTypes.h:183
void DeleteCell(vtkIdType cellId)
Delete cell by setting to nullptr cell type.
Definition: vtkCellTypes.h:178
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.
vtkUnsignedCharArray * TypeArray
Definition: vtkCellTypes.h:258
vtkIdType Extend
Definition: vtkCellTypes.h:262
int IsType(unsigned char type)
Return 1 if type specified is contained in list; 0 otherwise.
Definition: vtkCellTypes.h:270
unsigned char GetCellType(vtkIdType cellId)
Return the type of cell.
Definition: vtkCellTypes.h:198
void DeepCopy(vtkCellTypes *src)
Standard DeepCopy method.
vtkIdType Size
Definition: vtkCellTypes.h:260
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.
Definition: vtkCellTypes.h:250
vtkIdTypeArray * GetCellLocationsArray()
Methods for obtaining the arrays representing types and locations.
Definition: vtkCellTypes.h:251
dynamic, self-adjusting array of vtkIdType
a simple class to control print indentation
Definition: vtkIndent.h:113
dynamic, self-adjusting array of int
Definition: vtkIntArray.h:149
abstract base class for most VTK objects
Definition: vtkObject.h:82
dynamic, self-adjusting array of unsigned char
@ VTK_EMPTY_CELL
Definition: vtkCellType.h:85
@ VTK_POLYHEDRON
Definition: vtkCellType.h:128
@ VTK_CONVEX_POINT_SET
Definition: vtkCellType.h:125
int vtkIdType
Definition: vtkType.h:332