VTK  9.0.20200705
vtkPolyDataInternals.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkPolyDataInternals.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 
57 #ifndef vtkPolyDataInternals_h
58 #define vtkPolyDataInternals_h
59 
60 #ifndef __VTK_WRAP__ // Don't wrap this class.
61 
62 #include "vtkCommonDataModelModule.h" // For export macro
63 
64 #include "vtkCellType.h"
65 #include "vtkObject.h"
66 #include "vtkType.h"
67 
68 #include <cstdlib> // for std::size_t
69 #include <vector> // for CellMap implementation.
70 
72 {
73 
74 static constexpr vtkTypeUInt64 CELLID_MASK = 0x0fffffffffffffffull;
75 static constexpr vtkTypeUInt64 SHIFTED_TYPE_INDEX_MASK = 0xf000000000000000ull;
76 static constexpr vtkTypeUInt64 TARGET_MASK = 0x3ull << 62;
77 static constexpr vtkTypeUInt64 TYPE_VARIANT_MASK = 0x3ull << 60;
78 
79 // Enumeration of internal cell array targets.
80 enum class Target : vtkTypeUInt64
81 {
82  Verts = 0x0ull << 62,
83  Lines = 0x1ull << 62,
84  Polys = 0x2ull << 62,
85  Strips = 0x3ull << 62,
86 };
87 
88 // Enumeration of type variants.
89 enum class TypeVariant : vtkTypeUInt64
90 {
91  Dead = 0x0ull << 60,
92  Var1 = 0x1ull << 60,
93  Var2 = 0x2ull << 60,
94  Var3 = 0x3ull << 60,
95 };
96 
97 // Lookup table to convert a type index (TaggedCellId.GetTypeIndex()) into
98 // a VTK cell type.
99 // The type index is the highest four bits of the encoded value, eg. the
100 // target and type variant information.
101 static constexpr unsigned char TypeTable[16] = {
102  VTK_EMPTY_CELL, // 0000b | Verts | Dead
103  VTK_VERTEX, // 0001b | Verts | Var1
104  VTK_POLY_VERTEX, // 0010b | Verts | Var2
105  VTK_EMPTY_CELL, // 0011b | Verts | Var3
106  VTK_EMPTY_CELL, // 0100b | Lines | Dead
107  VTK_LINE, // 0101b | Lines | Var1
108  VTK_POLY_LINE, // 0110b | Lines | Var2
109  VTK_EMPTY_CELL, // 0111b | Lines | Var3
110  VTK_EMPTY_CELL, // 1000b | Polys | Dead
111  VTK_TRIANGLE, // 1001b | Polys | Var1
112  VTK_QUAD, // 1010b | Polys | Var2
113  VTK_POLYGON, // 1011b | Polys | Var3
114  VTK_EMPTY_CELL, // 1100b | Strips | Dead
115  VTK_TRIANGLE_STRIP, // 1101b | Strips | Var1
116  VTK_EMPTY_CELL, // 1110b | Strips | Var2
117  VTK_EMPTY_CELL, // 1111b | Strips | Var3
118 };
119 
120 // Convenience method to concatenate a target and type variant into the low
121 // four bits of a single byte. Used to build the TargetVarTable.
122 static constexpr unsigned char GenTargetVar(Target target, TypeVariant var) noexcept
123 {
124  return static_cast<unsigned char>(
125  (static_cast<vtkTypeUInt64>(target) | static_cast<vtkTypeUInt64>(var)) >> 60);
126 }
127 
128 // Lookup table that maps a VTK cell type (eg. VTK_TRIANGLE) into a target +
129 // type variant byte.
130 static constexpr unsigned char TargetVarTable[10] = {
131  GenTargetVar(Target::Verts, TypeVariant::Dead), // 0 | VTK_EMPTY_CELL
132  GenTargetVar(Target::Verts, TypeVariant::Var1), // 1 | VTK_VERTEX
133  GenTargetVar(Target::Verts, TypeVariant::Var2), // 2 | VTK_POLY_VERTEX
134  GenTargetVar(Target::Lines, TypeVariant::Var1), // 3 | VTK_LINE
135  GenTargetVar(Target::Lines, TypeVariant::Var2), // 4 | VTK_POLY_LINE
136  GenTargetVar(Target::Polys, TypeVariant::Var1), // 5 | VTK_TRIANGLE
137  GenTargetVar(Target::Strips, TypeVariant::Var1), // 6 | VTK_TRIANGLE_STRIP
138  GenTargetVar(Target::Polys, TypeVariant::Var3), // 7 | VTK_POLYGON
139  GenTargetVar(Target::Polys, TypeVariant::Var2), // 8 | VTK_PIXEL (treat as quad)
140  GenTargetVar(Target::Polys, TypeVariant::Var2), // 9 | VTK_QUAD
141 };
142 
143 // Thin wrapper around a vtkTypeUInt64 that encodes a target cell array,
144 // cell type, deleted status, and 60-bit cell id.
145 struct VTKCOMMONDATAMODEL_EXPORT TaggedCellId
146 {
147  // Encode a cell id and a VTK cell type (eg. VTK_TRIANGLE) into a
148  // vtkTypeUInt64.
149  static vtkTypeUInt64 Encode(vtkIdType cellId, VTKCellType type) noexcept
150  {
151  const size_t typeIndex = static_cast<size_t>(type);
152  return ((static_cast<vtkTypeUInt64>(cellId) & CELLID_MASK) |
153  (static_cast<vtkTypeUInt64>(TargetVarTable[typeIndex]) << 60));
154  }
155 
156  TaggedCellId() noexcept = default;
157 
158  // Create a TaggedCellId from a cellId and cell type (e.g. VTK_TRIANGLE).
159  TaggedCellId(vtkIdType cellId, VTKCellType cellType) noexcept
160  : Value(Encode(cellId, cellType))
161  {
162  }
163 
164  TaggedCellId(const TaggedCellId&) noexcept = default;
165  TaggedCellId(TaggedCellId&&) noexcept = default;
166  TaggedCellId& operator=(const TaggedCellId&) noexcept = default;
167  TaggedCellId& operator=(TaggedCellId&&) noexcept = default;
168 
169  // Get an enum value describing the internal vtkCellArray target used to
170  // store this cell.
171  Target GetTarget() const noexcept { return static_cast<Target>(this->Value & TARGET_MASK); }
172 
173  // Get the VTK cell type value (eg. VTK_TRIANGLE) as a single byte.
174  unsigned char GetCellType() const noexcept { return TypeTable[this->GetTypeIndex()]; }
175 
176  // Get the cell id used by the target vtkCellArray to store this cell.
177  vtkIdType GetCellId() const noexcept { return static_cast<vtkIdType>(this->Value & CELLID_MASK); }
178 
179  // Update the cell id. Most useful with the CellMap::InsertNextCell(type)
180  // signature.
181  void SetCellId(vtkIdType cellId) noexcept
182  {
183  this->Value &= SHIFTED_TYPE_INDEX_MASK;
184  this->Value |= (static_cast<vtkTypeUInt64>(cellId) & CELLID_MASK);
185  }
186 
187  // Mark this cell as deleted.
188  void MarkDeleted() noexcept { this->Value &= ~TYPE_VARIANT_MASK; }
189 
190  // Returns true if the cell has been deleted.
191  bool IsDeleted() const noexcept { return (this->Value & TYPE_VARIANT_MASK) == 0; }
192 
193 private:
194  vtkTypeUInt64 Value;
195 
196  // These shouldn't be needed outside of this struct. You're probably looking
197  // for GetCellType() instead.
198  TypeVariant GetTypeVariant() const noexcept
199  {
200  return static_cast<TypeVariant>(this->Value & TYPE_VARIANT_MASK);
201  }
202  std::size_t GetTypeIndex() const noexcept { return static_cast<std::size_t>(this->Value >> 60); }
203 };
204 
205 // Thin wrapper around a std::vector<TaggedCellId> to allow shallow copying, etc
206 class VTKCOMMONDATAMODEL_EXPORT CellMap : public vtkObject
207 {
208 public:
209  static CellMap* New();
210  vtkTypeMacro(CellMap, vtkObject);
211 
212  static bool ValidateCellType(VTKCellType cellType) noexcept
213  {
214  // 1-9 excluding 8 (VTK_PIXEL):
215  return cellType > 0 && cellType <= 10 && cellType != VTK_PIXEL;
216  }
217 
218  static bool ValidateCellId(vtkIdType cellId) noexcept
219  {
220  // is positive, won't truncate:
221  return (
222  (static_cast<vtkTypeUInt64>(cellId) & CELLID_MASK) == static_cast<vtkTypeUInt64>(cellId));
223  }
224 
225  void DeepCopy(CellMap* other)
226  {
227  if (other)
228  {
229  this->Map = other->Map;
230  }
231  else
232  {
233  this->Map.clear();
234  }
235  }
236 
237  void SetCapacity(vtkIdType numCells) { this->Map.reserve(static_cast<std::size_t>(numCells)); }
238 
239  TaggedCellId& GetTag(vtkIdType cellId) { return this->Map[static_cast<std::size_t>(cellId)]; }
240 
241  const TaggedCellId& GetTag(vtkIdType cellId) const
242  {
243  return this->Map[static_cast<std::size_t>(cellId)];
244  }
245 
246  // Caller must ValidateCellType first.
247  void InsertNextCell(vtkIdType cellId, VTKCellType cellType)
248  {
249  this->Map.emplace_back(cellId, cellType);
250  }
251 
252  // Caller must ValidateCellType and ValidateCellId first.
253  // useful for reusing the target lookup from cellType and then calling
254  // TaggedCellId::SetCellId later.
256  {
257  this->Map.emplace_back(vtkIdType(0), cellType);
258  return this->Map.back();
259  }
260 
261  vtkIdType GetNumberOfCells() const { return static_cast<vtkIdType>(this->Map.size()); }
262 
263  void Reset() { this->Map.clear(); }
264 
265  void Squeeze()
266  {
267  this->Map.shrink_to_fit(); // yaaaaay c++11
268  }
269 
270  // In rounded-up kibibytes, as is VTK's custom:
271  unsigned long GetActualMemorySize() const
272  {
273  return static_cast<unsigned long>(sizeof(TaggedCellId) * this->Map.capacity() + 1023) / 1024;
274  }
275 
276 protected:
277  CellMap();
278  ~CellMap() override;
279 
280  std::vector<TaggedCellId> Map;
281 
282 private:
283  CellMap(const CellMap&) = delete;
284  CellMap& operator=(const CellMap&) = delete;
285 };
286 
287 } // end namespace vtkPolyData_detail
288 
289 #endif // __VTK_WRAP__
290 #endif // vtkPolyDataInternals.h
291 
292 // VTK-HeaderTest-Exclude: vtkPolyDataInternals.h
vtkPolyData_detail::GenTargetVar
static constexpr unsigned char GenTargetVar(Target target, TypeVariant var) noexcept
Definition: vtkPolyDataInternals.h:122
VTK_TRIANGLE_STRIP
Definition: vtkCellType.h:52
vtkPolyData_detail::TypeVariant::Dead
vtkPolyData_detail::TaggedCellId::Encode
static vtkTypeUInt64 Encode(vtkIdType cellId, VTKCellType type) noexcept
Definition: vtkPolyDataInternals.h:149
vtkX3D::type
Definition: vtkX3D.h:522
vtkIdType
int vtkIdType
Definition: vtkType.h:330
VTK_LINE
Definition: vtkCellType.h:49
VTK_PIXEL
Definition: vtkCellType.h:54
vtkPolyData_detail::Target::Verts
vtkPolyData_detail::Target::Strips
vtkPolyData_detail::CellMap::DeepCopy
void DeepCopy(CellMap *other)
Definition: vtkPolyDataInternals.h:225
vtkPolyData_detail::CellMap::ValidateCellId
static bool ValidateCellId(vtkIdType cellId) noexcept
Definition: vtkPolyDataInternals.h:218
vtkObject
abstract base class for most VTK objects
Definition: vtkObject.h:62
vtkPolyData_detail::CellMap::ValidateCellType
static bool ValidateCellType(VTKCellType cellType) noexcept
Definition: vtkPolyDataInternals.h:212
vtkPolyData_detail::CellMap::GetTag
TaggedCellId & GetTag(vtkIdType cellId)
Definition: vtkPolyDataInternals.h:239
vtkPolyData_detail::TYPE_VARIANT_MASK
static constexpr vtkTypeUInt64 TYPE_VARIANT_MASK
Definition: vtkPolyDataInternals.h:77
vtkPolyData_detail::Target
Target
Definition: vtkPolyDataInternals.h:80
vtkPolyData_detail::TaggedCellId::GetCellId
vtkIdType GetCellId() const noexcept
Definition: vtkPolyDataInternals.h:177
vtkPolyData_detail::Target::Polys
vtkPolyData_detail::CellMap::Reset
void Reset()
Definition: vtkPolyDataInternals.h:263
vtkPolyData_detail::CELLID_MASK
static constexpr vtkTypeUInt64 CELLID_MASK
Definition: vtkPolyDataInternals.h:74
vtkPolyData_detail::TaggedCellId::SetCellId
void SetCellId(vtkIdType cellId) noexcept
Definition: vtkPolyDataInternals.h:181
VTK_POLY_LINE
Definition: vtkCellType.h:50
vtkPolyData_detail
Definition: vtkPolyDataInternals.h:71
VTK_QUAD
Definition: vtkCellType.h:55
vtkType.h
vtkPolyData_detail::Target::Lines
target
boost::graph_traits< vtkGraph * >::vertex_descriptor target(boost::graph_traits< vtkGraph * >::edge_descriptor e, vtkGraph *)
Definition: vtkBoostGraphAdapter.h:965
vtkPolyData_detail::CellMap::InsertNextCell
TaggedCellId & InsertNextCell(VTKCellType cellType)
Definition: vtkPolyDataInternals.h:255
VTK_EMPTY_CELL
Definition: vtkCellType.h:46
vtkPolyData_detail::TaggedCellId::IsDeleted
bool IsDeleted() const noexcept
Definition: vtkPolyDataInternals.h:191
vtkPolyData_detail::CellMap::GetActualMemorySize
unsigned long GetActualMemorySize() const
Definition: vtkPolyDataInternals.h:271
vtkPolyData_detail::TaggedCellId::GetCellType
unsigned char GetCellType() const noexcept
Definition: vtkPolyDataInternals.h:174
vtkCellType.h
VTK_POLY_VERTEX
Definition: vtkCellType.h:48
vtkPolyData_detail::TypeVariant::Var3
vtkPolyData_detail::CellMap::SetCapacity
void SetCapacity(vtkIdType numCells)
Definition: vtkPolyDataInternals.h:237
vtkPolyData_detail::SHIFTED_TYPE_INDEX_MASK
static constexpr vtkTypeUInt64 SHIFTED_TYPE_INDEX_MASK
Definition: vtkPolyDataInternals.h:75
vtkPolyData_detail::CellMap::Map
std::vector< TaggedCellId > Map
Definition: vtkPolyDataInternals.h:280
vtkObject.h
vtkPolyData_detail::CellMap::GetTag
const TaggedCellId & GetTag(vtkIdType cellId) const
Definition: vtkPolyDataInternals.h:241
vtkPolyData_detail::TARGET_MASK
static constexpr vtkTypeUInt64 TARGET_MASK
Definition: vtkPolyDataInternals.h:76
vtkPolyData_detail::TargetVarTable
static constexpr unsigned char TargetVarTable[10]
Definition: vtkPolyDataInternals.h:130
vtkPolyData_detail::TypeVariant::Var2
VTK_POLYGON
Definition: vtkCellType.h:53
vtkPolyData_detail::CellMap
Definition: vtkPolyDataInternals.h:206
VTKCellType
VTKCellType
Definition: vtkCellType.h:43
vtkPolyData_detail::TaggedCellId
Definition: vtkPolyDataInternals.h:145
vtkPolyData_detail::TaggedCellId::MarkDeleted
void MarkDeleted() noexcept
Definition: vtkPolyDataInternals.h:188
vtkPolyData_detail::CellMap::GetNumberOfCells
vtkIdType GetNumberOfCells() const
Definition: vtkPolyDataInternals.h:261
vtkPolyData_detail::CellMap::Squeeze
void Squeeze()
Definition: vtkPolyDataInternals.h:265
VTK_VERTEX
Definition: vtkCellType.h:47
vtkPolyData_detail::TypeVariant
TypeVariant
Definition: vtkPolyDataInternals.h:89
VTK_TRIANGLE
Definition: vtkCellType.h:51
vtkPolyData_detail::TypeTable
static constexpr unsigned char TypeTable[16]
Definition: vtkPolyDataInternals.h:101
vtkPolyData_detail::TypeVariant::Var1
vtkPolyData_detail::CellMap::InsertNextCell
void InsertNextCell(vtkIdType cellId, VTKCellType cellType)
Definition: vtkPolyDataInternals.h:247