VTK  9.2.20230603
vtkCellArray.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkCellArray.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 =========================================================================*/
250 #ifndef vtkCellArray_h
251 #define vtkCellArray_h
252 
253 #include "vtkCommonDataModelModule.h" // For export macro
254 #include "vtkObject.h"
255 
256 #include "vtkAOSDataArrayTemplate.h" // Needed for inline methods
257 #include "vtkCell.h" // Needed for inline methods
258 #include "vtkDataArrayRange.h" // Needed for inline methods
259 #include "vtkFeatures.h" // for VTK_USE_MEMKIND
260 #include "vtkSmartPointer.h" // For vtkSmartPointer
261 #include "vtkTypeInt32Array.h" // Needed for inline methods
262 #include "vtkTypeInt64Array.h" // Needed for inline methods
263 #include "vtkTypeList.h" // Needed for ArrayList definition
264 
265 #include <cassert> // for assert
266 #include <initializer_list> // for API
267 #include <type_traits> // for std::is_same
268 #include <utility> // for std::forward
269 
290 #define VTK_CELL_ARRAY_V2
291 
292 VTK_ABI_NAMESPACE_BEGIN
294 class vtkIdTypeArray;
295 
296 class VTKCOMMONDATAMODEL_EXPORT vtkCellArray : public vtkObject
297 {
298 public:
299  using ArrayType32 = vtkTypeInt32Array;
300  using ArrayType64 = vtkTypeInt64Array;
301 
303 
307  static vtkCellArray* New();
308  vtkTypeMacro(vtkCellArray, vtkObject);
309  void PrintSelf(ostream& os, vtkIndent indent) override;
310  void PrintDebug(ostream& os);
312 
321  using StorageArrayList = vtkTypeList::Create<ArrayType32, ArrayType64>;
322 
334 
343  vtkTypeBool Allocate(vtkIdType sz, vtkIdType vtkNotUsed(ext) = 1000)
344  {
345  return this->AllocateExact(sz, sz) ? 1 : 0;
346  }
347 
357  bool AllocateEstimate(vtkIdType numCells, vtkIdType maxCellSize)
358  {
359  return this->AllocateExact(numCells, numCells * maxCellSize);
360  }
361 
371  bool AllocateExact(vtkIdType numCells, vtkIdType connectivitySize);
372 
383  {
384  return this->AllocateExact(other->GetNumberOfCells(), other->GetNumberOfConnectivityIds());
385  }
386 
396  bool ResizeExact(vtkIdType numCells, vtkIdType connectivitySize);
397 
401  void Initialize();
402 
406  void Reset();
407 
413  void Squeeze();
414 
425  bool IsValid();
426 
431  {
432  if (this->Storage.Is64Bit())
433  {
434  return this->Storage.GetArrays64().Offsets->GetNumberOfValues() - 1;
435  }
436  else
437  {
438  return this->Storage.GetArrays32().Offsets->GetNumberOfValues() - 1;
439  }
440  }
441 
447  {
448  if (this->Storage.Is64Bit())
449  {
450  return this->Storage.GetArrays64().Offsets->GetNumberOfValues();
451  }
452  else
453  {
454  return this->Storage.GetArrays32().Offsets->GetNumberOfValues();
455  }
456  }
457 
462  {
463  if (this->Storage.Is64Bit())
464  {
465  return this->Storage.GetArrays64().Offsets->GetValue(cellId);
466  }
467  else
468  {
469  return this->Storage.GetArrays32().Offsets->GetValue(cellId);
470  }
471  }
472 
480  {
481  if (this->Storage.Is64Bit())
482  {
483  return this->Storage.GetArrays64().Connectivity->GetNumberOfValues();
484  }
485  else
486  {
487  return this->Storage.GetArrays32().Connectivity->GetNumberOfValues();
488  }
489  }
490 
497 
508  void SetData(vtkIdTypeArray* offsets, vtkIdTypeArray* connectivity);
511  void SetData(
513  void SetData(vtkTypeInt32Array* offsets, vtkTypeInt32Array* connectivity);
514  void SetData(vtkTypeInt64Array* offsets, vtkTypeInt64Array* connectivity);
529  bool SetData(vtkDataArray* offsets, vtkDataArray* connectivity);
530 
544  bool SetData(vtkIdType cellSize, vtkDataArray* connectivity);
545 
550  bool IsStorage64Bit() const { return this->Storage.Is64Bit(); }
551 
558  bool IsStorageShareable() const
559  {
560  if (this->Storage.Is64Bit())
561  {
563  }
564  else
565  {
567  }
568  }
569 
624  {
625  if (this->Storage.Is64Bit())
626  {
627  return this->GetOffsetsArray64();
628  }
629  else
630  {
631  return this->GetOffsetsArray32();
632  }
633  }
645  {
646  if (this->Storage.Is64Bit())
647  {
648  return this->GetConnectivityArray64();
649  }
650  else
651  {
652  return this->GetConnectivityArray32();
653  }
654  }
668 
678  void InitTraversal();
679 
694  int GetNextCell(vtkIdType& npts, vtkIdType const*& pts) VTK_SIZEHINT(pts, npts);
695 
706  int GetNextCell(vtkIdList* pts);
707 
718  void GetCellAtId(vtkIdType cellId, vtkIdType& cellSize, vtkIdType const*& cellPoints)
719  VTK_SIZEHINT(cellPoints, cellSize) VTK_EXPECTS(0 <= cellId && cellId < GetNumberOfCells());
720 
730  void GetCellAtId(
731  vtkIdType cellId, vtkIdType& cellSize, vtkIdType const*& cellPoints, vtkIdList* ptIds)
732  VTK_SIZEHINT(cellPoints, cellSize) VTK_EXPECTS(0 <= cellId && cellId < GetNumberOfCells());
733 
739  void GetCellAtId(vtkIdType cellId, vtkIdList* pts)
740  VTK_EXPECTS(0 <= cellId && cellId < GetNumberOfCells());
741 
745  vtkIdType GetCellSize(vtkIdType cellId) const;
746 
750  vtkIdType InsertNextCell(vtkCell* cell);
751 
756  vtkIdType InsertNextCell(vtkIdType npts, const vtkIdType* pts) VTK_SIZEHINT(pts, npts);
757 
762  vtkIdType InsertNextCell(vtkIdList* pts);
763 
771  vtkIdType InsertNextCell(const std::initializer_list<vtkIdType>& cell)
772  {
773  return this->InsertNextCell(static_cast<vtkIdType>(cell.size()), cell.begin());
774  }
775 
782  vtkIdType InsertNextCell(int npts);
783 
788  void InsertCellPoint(vtkIdType id);
789 
794  void UpdateCellCount(int npts);
795 
810  void ReverseCellAtId(vtkIdType cellId) VTK_EXPECTS(0 <= cellId && cellId < GetNumberOfCells());
811 
820  void ReplaceCellAtId(vtkIdType cellId, vtkIdList* list);
821  void ReplaceCellAtId(vtkIdType cellId, vtkIdType cellSize, const vtkIdType* cellPoints)
822  VTK_EXPECTS(0 <= cellId && cellId < GetNumberOfCells()) VTK_SIZEHINT(cellPoints, cellSize);
832  void ReplaceCellPointAtId(vtkIdType cellId, vtkIdType cellPointIndex, vtkIdType newPointId);
833 
841  void ReplaceCellAtId(vtkIdType cellId, const std::initializer_list<vtkIdType>& cell)
842  {
843  return this->ReplaceCellAtId(cellId, static_cast<vtkIdType>(cell.size()), cell.begin());
844  }
845 
851 
856 
861 
865  void Append(vtkCellArray* src, vtkIdType pointOffset = 0);
866 
878 
906  void AppendLegacyFormat(const vtkIdType* data, vtkIdType len, vtkIdType ptOffset = 0)
907  VTK_SIZEHINT(data, len);
918  unsigned long GetActualMemorySize() const;
919 
920  // The following code is used to support
921 
922  // The wrappers get understandably confused by some of the template code below
923 #ifndef __VTK_WRAP__
924 
925  // Holds connectivity and offset arrays of the given ArrayType.
926  template <typename ArrayT>
927  struct VisitState
928  {
929  using ArrayType = ArrayT;
930  using ValueType = typename ArrayType::ValueType;
931  using CellRangeType = decltype(vtk::DataArrayValueRange<1>(std::declval<ArrayType>()));
932 
933  // We can't just use is_same here, since binary compatible representations
934  // (e.g. int and long) are distinct types. Instead, ensure that ValueType
935  // is a signed integer the same size as vtkIdType.
936  // If this value is true, ValueType pointers may be safely converted to
937  // vtkIdType pointers via reinterpret cast.
938  static constexpr bool ValueTypeIsSameAsIdType = std::is_integral<ValueType>::value &&
939  std::is_signed<ValueType>::value && (sizeof(ValueType) == sizeof(vtkIdType));
940 
941  ArrayType* GetOffsets() { return this->Offsets; }
942  const ArrayType* GetOffsets() const { return this->Offsets; }
943 
944  ArrayType* GetConnectivity() { return this->Connectivity; }
945  const ArrayType* GetConnectivity() const { return this->Connectivity; }
946 
948 
950 
952 
954 
956 
957  friend class vtkCellArray;
958 
959  protected:
961  {
962  this->Connectivity = vtkSmartPointer<ArrayType>::New();
963  this->Offsets = vtkSmartPointer<ArrayType>::New();
964  this->Offsets->InsertNextValue(0);
966  {
967  this->IsInMemkind = true;
968  }
969  }
970  ~VisitState() = default;
971  void* operator new(size_t nSize)
972  {
973  void* r;
974 #ifdef VTK_USE_MEMKIND
976 #else
977  r = malloc(nSize);
978 #endif
979  return r;
980  }
981  void operator delete(void* p)
982  {
983 #ifdef VTK_USE_MEMKIND
984  VisitState* a = static_cast<VisitState*>(p);
985  if (a->IsInMemkind)
986  {
988  }
989  else
990  {
991  free(p);
992  }
993 #else
994  free(p);
995 #endif
996  }
997 
1000 
1001  private:
1002  VisitState(const VisitState&) = delete;
1003  VisitState& operator=(const VisitState&) = delete;
1004  bool IsInMemkind = false;
1005  };
1006 
1007 private: // Helpers that allow Visit to return a value:
1008  template <typename Functor, typename... Args>
1009  using GetReturnType = decltype(
1010  std::declval<Functor>()(std::declval<VisitState<ArrayType32>&>(), std::declval<Args>()...));
1011 
1012  template <typename Functor, typename... Args>
1013  struct ReturnsVoid : std::is_same<GetReturnType<Functor, Args...>, void>
1014  {
1015  };
1016 
1017 public:
1087  template <typename Functor, typename... Args,
1088  typename = typename std::enable_if<ReturnsVoid<Functor, Args...>::value>::type>
1089  void Visit(Functor&& functor, Args&&... args)
1090  {
1091  if (this->Storage.Is64Bit())
1092  {
1093  // If you get an error on the next line, a call to Visit(functor, Args...)
1094  // is being called with arguments that do not match the functor's call
1095  // signature. See the Visit documentation for details.
1096  functor(this->Storage.GetArrays64(), std::forward<Args>(args)...);
1097  }
1098  else
1099  {
1100  // If you get an error on the next line, a call to Visit(functor, Args...)
1101  // is being called with arguments that do not match the functor's call
1102  // signature. See the Visit documentation for details.
1103  functor(this->Storage.GetArrays32(), std::forward<Args>(args)...);
1104  }
1105  }
1106 
1107  template <typename Functor, typename... Args,
1108  typename = typename std::enable_if<ReturnsVoid<Functor, Args...>::value>::type>
1109  void Visit(Functor&& functor, Args&&... args) const
1110  {
1111  if (this->Storage.Is64Bit())
1112  {
1113  // If you get an error on the next line, a call to Visit(functor, Args...)
1114  // is being called with arguments that do not match the functor's call
1115  // signature. See the Visit documentation for details.
1116  functor(this->Storage.GetArrays64(), std::forward<Args>(args)...);
1117  }
1118  else
1119  {
1120  // If you get an error on the next line, a call to Visit(functor, Args...)
1121  // is being called with arguments that do not match the functor's call
1122  // signature. See the Visit documentation for details.
1123  functor(this->Storage.GetArrays32(), std::forward<Args>(args)...);
1124  }
1125  }
1126 
1127  template <typename Functor, typename... Args,
1128  typename = typename std::enable_if<!ReturnsVoid<Functor, Args...>::value>::type>
1129  GetReturnType<Functor, Args...> Visit(Functor&& functor, Args&&... args)
1130  {
1131  if (this->Storage.Is64Bit())
1132  {
1133  // If you get an error on the next line, a call to Visit(functor, Args...)
1134  // is being called with arguments that do not match the functor's call
1135  // signature. See the Visit documentation for details.
1136  return functor(this->Storage.GetArrays64(), std::forward<Args>(args)...);
1137  }
1138  else
1139  {
1140  // If you get an error on the next line, a call to Visit(functor, Args...)
1141  // is being called with arguments that do not match the functor's call
1142  // signature. See the Visit documentation for details.
1143  return functor(this->Storage.GetArrays32(), std::forward<Args>(args)...);
1144  }
1145  }
1146  template <typename Functor, typename... Args,
1147  typename = typename std::enable_if<!ReturnsVoid<Functor, Args...>::value>::type>
1148  GetReturnType<Functor, Args...> Visit(Functor&& functor, Args&&... args) const
1149  {
1150  if (this->Storage.Is64Bit())
1151  {
1152  // If you get an error on the next line, a call to Visit(functor, Args...)
1153  // is being called with arguments that do not match the functor's call
1154  // signature. See the Visit documentation for details.
1155  return functor(this->Storage.GetArrays64(), std::forward<Args>(args)...);
1156  }
1157  else
1158  {
1159  // If you get an error on the next line, a call to Visit(functor, Args...)
1160  // is being called with arguments that do not match the functor's call
1161  // signature. See the Visit documentation for details.
1162  return functor(this->Storage.GetArrays32(), std::forward<Args>(args)...);
1163  }
1164  }
1165 
1168 #endif // __VTK_WRAP__
1169 
1170  //=================== Begin Legacy Methods ===================================
1171  // These should be deprecated at some point as they are confusing or very slow
1172 
1180 
1192  vtkIdType EstimateSize(vtkIdType numCells, int maxPtsPerCell);
1193 
1203 
1211 
1221  void GetCell(vtkIdType loc, vtkIdType& npts, const vtkIdType*& pts)
1222  VTK_EXPECTS(0 <= loc && loc < GetNumberOfConnectivityEntries()) VTK_SIZEHINT(pts, npts);
1223 
1230  void GetCell(vtkIdType loc, vtkIdList* pts)
1231  VTK_EXPECTS(0 <= loc && loc < GetNumberOfConnectivityEntries());
1232 
1239  vtkIdType GetInsertLocation(int npts);
1240 
1248  vtkIdType GetTraversalLocation();
1249  vtkIdType GetTraversalLocation(vtkIdType npts);
1250  void SetTraversalLocation(vtkIdType loc);
1260  void ReverseCell(vtkIdType loc) VTK_EXPECTS(0 <= loc && loc < GetNumberOfConnectivityEntries());
1261 
1273  void ReplaceCell(vtkIdType loc, int npts, const vtkIdType pts[])
1274  VTK_EXPECTS(0 <= loc && loc < GetNumberOfConnectivityEntries()) VTK_SIZEHINT(pts, npts);
1275 
1290  void SetCells(vtkIdType ncells, vtkIdTypeArray* cells);
1291 
1303 
1304  //=================== End Legacy Methods =====================================
1305 
1306  friend class vtkCellArrayIterator;
1307 
1308 protected:
1310  ~vtkCellArray() override;
1311 
1312  // Encapsulates storage of the internal arrays as a discriminated union
1313  // between 32-bit and 64-bit storage.
1314  struct Storage
1315  {
1316  // Union type that switches 32 and 64 bit array storage
1317  union ArraySwitch {
1318  ArraySwitch() = default; // handled by Storage
1319  ~ArraySwitch() = default; // handle by Storage
1322  };
1323 
1325  {
1326 #ifdef VTK_USE_MEMKIND
1327  this->Arrays =
1329 #else
1330  this->Arrays = new ArraySwitch;
1331 #endif
1332 
1333  // Default to the compile-time setting:
1334 #ifdef VTK_USE_64BIT_IDS
1335 
1336  this->Arrays->Int64 = new VisitState<ArrayType64>;
1337  this->StorageIs64Bit = true;
1338 
1339 #else // VTK_USE_64BIT_IDS
1340 
1341  this->Arrays->Int32 = new VisitState<ArrayType32>;
1342  this->StorageIs64Bit = false;
1343 
1344 #endif // VTK_USE_64BIT_IDS
1345 #ifdef VTK_USE_MEMKIND
1347  {
1348  this->IsInMemkind = true;
1349  }
1350 #else
1351  (void)this->IsInMemkind; // comp warning workaround
1352 #endif
1353  }
1354 
1356  {
1357  if (this->StorageIs64Bit)
1358  {
1359  this->Arrays->Int64->~VisitState();
1360  delete this->Arrays->Int64;
1361  }
1362  else
1363  {
1364  this->Arrays->Int32->~VisitState();
1365  delete this->Arrays->Int32;
1366  }
1367 #ifdef VTK_USE_MEMKIND
1368  if (this->IsInMemkind)
1369  {
1371  }
1372  else
1373  {
1374  free(this->Arrays);
1375  }
1376 #else
1377  delete this->Arrays;
1378 #endif
1379  }
1380 
1381  // Switch the internal arrays to be 32-bit. Any old data is lost. Returns
1382  // true if the storage changes.
1384  {
1385  if (!this->StorageIs64Bit)
1386  {
1387  return false;
1388  }
1389 
1390  this->Arrays->Int64->~VisitState();
1391  delete this->Arrays->Int64;
1392  this->Arrays->Int32 = new VisitState<ArrayType32>;
1393  this->StorageIs64Bit = false;
1394 
1395  return true;
1396  }
1397 
1398  // Switch the internal arrays to be 64-bit. Any old data is lost. Returns
1399  // true if the storage changes.
1401  {
1402  if (this->StorageIs64Bit)
1403  {
1404  return false;
1405  }
1406 
1407  this->Arrays->Int32->~VisitState();
1408  delete this->Arrays->Int32;
1409  this->Arrays->Int64 = new VisitState<ArrayType64>;
1410  this->StorageIs64Bit = true;
1411 
1412  return true;
1413  }
1414 
1415  // Returns true if the storage is currently configured to be 64 bit.
1416  bool Is64Bit() const { return this->StorageIs64Bit; }
1417 
1418  // Get the VisitState for 32-bit arrays
1420  {
1421  assert(!this->StorageIs64Bit);
1422  return *this->Arrays->Int32;
1423  }
1424 
1426  {
1427  assert(!this->StorageIs64Bit);
1428  return *this->Arrays->Int32;
1429  }
1430 
1431  // Get the VisitState for 64-bit arrays
1433  {
1434  assert(this->StorageIs64Bit);
1435  return *this->Arrays->Int64;
1436  }
1437 
1439  {
1440  assert(this->StorageIs64Bit);
1441  return *this->Arrays->Int64;
1442  }
1443 
1444  private:
1445  // Access restricted to ensure proper union construction/destruction thru
1446  // API.
1447  ArraySwitch* Arrays;
1448  bool StorageIs64Bit;
1449  bool IsInMemkind = false;
1450  };
1451 
1454  vtkIdType TraversalCellId{ 0 };
1455 
1457 
1458 private:
1459  vtkCellArray(const vtkCellArray&) = delete;
1460  void operator=(const vtkCellArray&) = delete;
1461 };
1462 
1463 template <typename ArrayT>
1465 {
1466  return this->Offsets->GetNumberOfValues() - 1;
1467 }
1468 
1469 template <typename ArrayT>
1471 {
1472  return static_cast<vtkIdType>(this->Offsets->GetValue(cellId));
1473 }
1474 
1475 template <typename ArrayT>
1477 {
1478  return static_cast<vtkIdType>(this->Offsets->GetValue(cellId + 1));
1479 }
1480 
1481 template <typename ArrayT>
1483 {
1484  return this->GetEndOffset(cellId) - this->GetBeginOffset(cellId);
1485 }
1486 
1487 template <typename ArrayT>
1490 {
1491  return vtk::DataArrayValueRange<1>(
1492  this->GetConnectivity(), this->GetBeginOffset(cellId), this->GetEndOffset(cellId));
1493 }
1494 VTK_ABI_NAMESPACE_END
1495 
1497 {
1498 VTK_ABI_NAMESPACE_BEGIN
1499 
1501 {
1502  // Insert full cell
1503  template <typename CellStateT>
1504  vtkIdType operator()(CellStateT& state, const vtkIdType npts, const vtkIdType pts[])
1505  {
1506  using ValueType = typename CellStateT::ValueType;
1507  auto* conn = state.GetConnectivity();
1508  auto* offsets = state.GetOffsets();
1509 
1510  const vtkIdType cellId = offsets->GetNumberOfValues() - 1;
1511 
1512  offsets->InsertNextValue(static_cast<ValueType>(conn->GetNumberOfValues() + npts));
1513 
1514  for (vtkIdType i = 0; i < npts; ++i)
1515  {
1516  conn->InsertNextValue(static_cast<ValueType>(pts[i]));
1517  }
1518 
1519  return cellId;
1520  }
1521 
1522  // Just update offset table (for incremental API)
1523  template <typename CellStateT>
1524  vtkIdType operator()(CellStateT& state, const vtkIdType npts)
1525  {
1526  using ValueType = typename CellStateT::ValueType;
1527  auto* conn = state.GetConnectivity();
1528  auto* offsets = state.GetOffsets();
1529 
1530  const vtkIdType cellId = offsets->GetNumberOfValues() - 1;
1531 
1532  offsets->InsertNextValue(static_cast<ValueType>(conn->GetNumberOfValues() + npts));
1533 
1534  return cellId;
1535  }
1536 };
1537 
1538 // for incremental API:
1540 {
1541  template <typename CellStateT>
1542  void operator()(CellStateT& state, const vtkIdType npts)
1543  {
1544  using ValueType = typename CellStateT::ValueType;
1545 
1546  auto* offsets = state.GetOffsets();
1547  const ValueType cellBegin = offsets->GetValue(offsets->GetMaxId() - 1);
1548  offsets->SetValue(offsets->GetMaxId(), static_cast<ValueType>(cellBegin + npts));
1549  }
1550 };
1551 
1553 {
1554  template <typename CellStateT>
1555  vtkIdType operator()(CellStateT& state, vtkIdType cellId)
1556  {
1557  return state.GetCellSize(cellId);
1558  }
1559 };
1560 
1562 {
1563  template <typename CellStateT>
1564  void operator()(CellStateT& state, const vtkIdType cellId, vtkIdList* ids)
1565  {
1566  using ValueType = typename CellStateT::ValueType;
1567 
1568  const vtkIdType beginOffset = state.GetBeginOffset(cellId);
1569  const vtkIdType endOffset = state.GetEndOffset(cellId);
1570  const vtkIdType cellSize = endOffset - beginOffset;
1571  const auto cellConnectivity = state.GetConnectivity()->GetPointer(beginOffset);
1572 
1573  // ValueType differs from vtkIdType, so we have to copy into a temporary buffer:
1574  ids->SetNumberOfIds(cellSize);
1575  vtkIdType* idPtr = ids->GetPointer(0);
1576  for (ValueType i = 0; i < cellSize; ++i)
1577  {
1578  idPtr[i] = static_cast<vtkIdType>(cellConnectivity[i]);
1579  }
1580  }
1581 
1582  // SFINAE helper to check if a VisitState's connectivity array's memory
1583  // can be used as a vtkIdType*.
1584  template <typename CellStateT>
1586  {
1587  private:
1588  using ValueType = typename CellStateT::ValueType;
1589  using ArrayType = typename CellStateT::ArrayType;
1591  static constexpr bool ValueTypeCompat = CellStateT::ValueTypeIsSameAsIdType;
1592  static constexpr bool ArrayTypeCompat = std::is_base_of<AOSArrayType, ArrayType>::value;
1593 
1594  public:
1595  static constexpr bool value = ValueTypeCompat && ArrayTypeCompat;
1596  };
1597 
1598  template <typename CellStateT>
1600  CellStateT& state, const vtkIdType cellId, vtkIdType& cellSize, vtkIdType const*& cellPoints,
1601  vtkIdList* vtkNotUsed(temp))
1602  {
1603  const vtkIdType beginOffset = state.GetBeginOffset(cellId);
1604  const vtkIdType endOffset = state.GetEndOffset(cellId);
1605  cellSize = endOffset - beginOffset;
1606  // This is safe, see CanShareConnPtr helper above.
1607  cellPoints = reinterpret_cast<vtkIdType*>(state.GetConnectivity()->GetPointer(beginOffset));
1608  }
1609 
1610  template <typename CellStateT>
1612  CellStateT& state, const vtkIdType cellId, vtkIdType& cellSize, vtkIdType const*& cellPoints,
1613  vtkIdList* temp)
1614  {
1615  using ValueType = typename CellStateT::ValueType;
1616 
1617  const vtkIdType beginOffset = state.GetBeginOffset(cellId);
1618  const vtkIdType endOffset = state.GetEndOffset(cellId);
1619  cellSize = endOffset - beginOffset;
1620  const ValueType* cellConnectivity = state.GetConnectivity()->GetPointer(beginOffset);
1621 
1622  // ValueType differs from vtkIdType, so we have to copy into a temporary buffer:
1623  temp->SetNumberOfIds(cellSize);
1624  vtkIdType* tempPtr = temp->GetPointer(0);
1625  for (vtkIdType i = 0; i < cellSize; ++i)
1626  {
1627  tempPtr[i] = static_cast<vtkIdType>(cellConnectivity[i]);
1628  }
1629 
1630  cellPoints = temp->GetPointer(0);
1631  }
1632 };
1633 
1635 {
1636  template <typename CellStateT>
1637  void operator()(CellStateT& state)
1638  {
1639  state.GetOffsets()->Reset();
1640  state.GetConnectivity()->Reset();
1641  state.GetOffsets()->InsertNextValue(0);
1642  }
1643 };
1644 
1645 VTK_ABI_NAMESPACE_END
1646 } // end namespace vtkCellArray_detail
1647 
1648 VTK_ABI_NAMESPACE_BEGIN
1649 //----------------------------------------------------------------------------
1651 {
1652  this->TraversalCellId = 0;
1653 }
1654 
1655 //----------------------------------------------------------------------------
1656 inline int vtkCellArray::GetNextCell(vtkIdType& npts, vtkIdType const*& pts) VTK_SIZEHINT(pts, npts)
1657 {
1658  if (this->TraversalCellId < this->GetNumberOfCells())
1659  {
1660  this->GetCellAtId(this->TraversalCellId, npts, pts);
1661  ++this->TraversalCellId;
1662  return 1;
1663  }
1664 
1665  npts = 0;
1666  pts = nullptr;
1667  return 0;
1668 }
1669 
1670 //----------------------------------------------------------------------------
1672 {
1673  if (this->TraversalCellId < this->GetNumberOfCells())
1674  {
1675  this->GetCellAtId(this->TraversalCellId, pts);
1676  ++this->TraversalCellId;
1677  return 1;
1678  }
1679 
1680  pts->Reset();
1681  return 0;
1682 }
1683 //----------------------------------------------------------------------------
1685 {
1686  return this->Visit(vtkCellArray_detail::GetCellSizeImpl{}, cellId);
1687 }
1688 
1689 //----------------------------------------------------------------------------
1690 inline void vtkCellArray::GetCellAtId(vtkIdType cellId, vtkIdType& cellSize,
1691  vtkIdType const*& cellPoints) VTK_SIZEHINT(cellPoints, cellSize)
1692 {
1693  this->Visit(vtkCellArray_detail::GetCellAtIdImpl{}, cellId, cellSize, cellPoints, this->TempCell);
1694 }
1695 
1696 //----------------------------------------------------------------------------
1697 inline void vtkCellArray::GetCellAtId(vtkIdType cellId, vtkIdType& cellSize,
1698  vtkIdType const*& cellPoints, vtkIdList* ptIds) VTK_SIZEHINT(cellPoints, cellSize)
1699 {
1700  this->Visit(vtkCellArray_detail::GetCellAtIdImpl{}, cellId, cellSize, cellPoints, ptIds);
1701 }
1702 
1703 //----------------------------------------------------------------------------
1705 {
1706  this->Visit(vtkCellArray_detail::GetCellAtIdImpl{}, cellId, pts);
1707 }
1708 
1709 //----------------------------------------------------------------------------
1711  VTK_SIZEHINT(pts, npts)
1712 {
1713  return this->Visit(vtkCellArray_detail::InsertNextCellImpl{}, npts, pts);
1714 }
1715 
1716 //----------------------------------------------------------------------------
1718 {
1719  return this->Visit(vtkCellArray_detail::InsertNextCellImpl{}, npts);
1720 }
1721 
1722 //----------------------------------------------------------------------------
1724 {
1725  if (this->Storage.Is64Bit())
1726  {
1727  using ValueType = typename ArrayType64::ValueType;
1728  this->Storage.GetArrays64().Connectivity->InsertNextValue(static_cast<ValueType>(id));
1729  }
1730  else
1731  {
1732  using ValueType = typename ArrayType32::ValueType;
1733  this->Storage.GetArrays32().Connectivity->InsertNextValue(static_cast<ValueType>(id));
1734  }
1735 }
1736 
1737 //----------------------------------------------------------------------------
1738 inline void vtkCellArray::UpdateCellCount(int npts)
1739 {
1741 }
1742 
1743 //----------------------------------------------------------------------------
1745 {
1746  return this->Visit(
1748 }
1749 
1750 //----------------------------------------------------------------------------
1752 {
1753  vtkIdList* pts = cell->GetPointIds();
1754  return this->Visit(
1756 }
1757 
1758 //----------------------------------------------------------------------------
1759 inline void vtkCellArray::Reset()
1760 {
1762 }
1763 
1764 VTK_ABI_NAMESPACE_END
1765 #endif // vtkCellArray.h
Encapsulate traversal logic for vtkCellArray.
object to represent cell connectivity
Definition: vtkCellArray.h:297
void SetData(vtkAOSDataArrayTemplate< int > *offsets, vtkAOSDataArrayTemplate< int > *connectivity)
Set the internal data arrays to the supplied offsets and connectivity arrays.
int GetNextCell(vtkIdType &npts, vtkIdType const *&pts)
vtkIdType GetNumberOfConnectivityEntries()
Return the size of the array that would be returned from ExportLegacyFormat().
void SetData(vtkTypeInt32Array *offsets, vtkTypeInt32Array *connectivity)
Set the internal data arrays to the supplied offsets and connectivity arrays.
vtkTypeBool Allocate(vtkIdType sz, vtkIdType vtkNotUsed(ext)=1000)
Allocate memory.
Definition: vtkCellArray.h:343
void UseDefaultStorage()
Initialize internal data structures to use 32- or 64-bit storage.
Storage Storage
vtkIdType GetOffset(vtkIdType cellId)
Get the offset (into the connectivity) for a specified cell id.
Definition: vtkCellArray.h:461
vtkIdType GetCellSize(vtkIdType cellId) const
Return the size of the cell at cellId.
bool AllocateCopy(vtkCellArray *other)
Pre-allocate memory in internal data structures to match the used size of the input vtkCellArray.
Definition: vtkCellArray.h:382
void AppendLegacyFormat(vtkIdTypeArray *data, vtkIdType ptOffset=0)
Append an array of data with the legacy vtkCellArray layout, e.g.
bool IsStorageShareable() const
Definition: vtkCellArray.h:558
ArrayType32 * GetOffsetsArray32()
Return the array used to store cell offsets.
Definition: vtkCellArray.h:634
bool IsValid()
Check that internal storage is consistent and in a valid state.
void AppendLegacyFormat(const vtkIdType *data, vtkIdType len, vtkIdType ptOffset=0)
Append an array of data with the legacy vtkCellArray layout, e.g.
bool SetData(vtkIdType cellSize, vtkDataArray *connectivity)
Sets the internal arrays to the supported connectivity array with an offsets array automatically gene...
vtkIdType GetTraversalCellId()
Get/Set the current cellId for traversal.
void Visit(Functor &&functor, Args &&... args)
vtkTypeInt32Array ArrayType32
Definition: vtkCellArray.h:299
void Reset()
Reuse list.
ArrayType64 * GetOffsetsArray64()
Return the array used to store cell offsets.
Definition: vtkCellArray.h:635
void SetData(vtkAOSDataArrayTemplate< long long > *offsets, vtkAOSDataArrayTemplate< long long > *connectivity)
Set the internal data arrays to the supplied offsets and connectivity arrays.
GetReturnType< Functor, Args... > Visit(Functor &&functor, Args &&... args)
bool AllocateExact(vtkIdType numCells, vtkIdType connectivitySize)
Pre-allocate memory in internal data structures.
bool ResizeExact(vtkIdType numCells, vtkIdType connectivitySize)
ResizeExact() resizes the internal structures to hold numCells total cell offsets and connectivitySiz...
vtkIdType EstimateSize(vtkIdType numCells, int maxPtsPerCell)
Utility routines help manage memory of cell array.
bool AllocateEstimate(vtkIdType numCells, vtkIdType maxCellSize)
Pre-allocate memory in internal data structures.
Definition: vtkCellArray.h:357
vtkIdType GetNumberOfOffsets() const
Get the number of elements in the offsets array.
Definition: vtkCellArray.h:446
vtkTypeInt64Array ArrayType64
Definition: vtkCellArray.h:300
vtkIdType TraversalCellId
void InitTraversal()
void Use64BitStorage()
Initialize internal data structures to use 32- or 64-bit storage.
bool ConvertToDefaultStorage()
Convert internal data structures to use 32- or 64-bit storage.
bool CanConvertToDefaultStorage() const
Check if the existing data can safely be converted to use 32- or 64- bit storage.
void GetCell(vtkIdType loc, vtkIdType &npts, const vtkIdType *&pts)
Internal method used to retrieve a cell given a legacy offset location.
bool CanConvertTo32BitStorage() const
Check if the existing data can safely be converted to use 32- or 64- bit storage.
void InsertCellPoint(vtkIdType id)
Used in conjunction with InsertNextCell(npts) to add another point to the list of cells.
void GetCellAtId(vtkIdType cellId, vtkIdType &cellSize, vtkIdType const *&cellPoints)
Return the point ids for the cell at cellId.
void ReplaceCellAtId(vtkIdType cellId, vtkIdList *list)
Replaces the point ids for the specified cell with the supplied list.
void ReverseCellAtId(vtkIdType cellId)
Reverses the order of the point ids for the specified cell.
bool SetData(vtkDataArray *offsets, vtkDataArray *connectivity)
Sets the internal arrays to the supplied offsets and connectivity arrays.
void Squeeze()
Reclaim any extra memory while preserving data.
void SetData(vtkIdTypeArray *offsets, vtkIdTypeArray *connectivity)
Set the internal data arrays to the supplied offsets and connectivity arrays.
static vtkCellArray * New()
Standard methods for instantiation, type information, and printing.
ArrayType32 * GetConnectivityArray32()
Return the array used to store the point ids that define the cells' connectivity.
Definition: vtkCellArray.h:655
ArrayType64 * GetConnectivityArray64()
Return the array used to store the point ids that define the cells' connectivity.
Definition: vtkCellArray.h:656
bool ConvertTo32BitStorage()
Convert internal data structures to use 32- or 64-bit storage.
vtkIdType IsHomogeneous()
Check if all cells have the same number of vertices.
void Visit(Functor &&functor, Args &&... args) const
int GetMaxCellSize()
Returns the size of the largest cell.
bool ConvertToSmallestStorage()
Convert internal data structures to use 32- or 64-bit storage.
void ShallowCopy(vtkCellArray *ca)
Shallow copy ca into this cell array.
void ExportLegacyFormat(vtkIdTypeArray *data)
Fill data with the old-style vtkCellArray data layout, e.g.
bool ConvertTo64BitStorage()
Convert internal data structures to use 32- or 64-bit storage.
void SetData(vtkAOSDataArrayTemplate< long > *offsets, vtkAOSDataArrayTemplate< long > *connectivity)
Set the internal data arrays to the supplied offsets and connectivity arrays.
void ImportLegacyFormat(const vtkIdType *data, vtkIdType len)
Import an array of data with the legacy vtkCellArray layout, e.g.
void ImportLegacyFormat(vtkIdTypeArray *data)
Import an array of data with the legacy vtkCellArray layout, e.g.
void Use32BitStorage()
Initialize internal data structures to use 32- or 64-bit storage.
void Initialize()
Free any memory and reset to an empty state.
void DeepCopy(vtkCellArray *ca)
Perform a deep copy (no reference counting) of the given cell array.
typename vtkTypeList::Unique< vtkTypeList::Create< vtkAOSDataArrayTemplate< int >, vtkAOSDataArrayTemplate< long >, vtkAOSDataArrayTemplate< long long > >>::Result InputArrayList
List of possible ArrayTypes that are compatible with internal storage.
Definition: vtkCellArray.h:333
vtkDataArray * GetConnectivityArray()
Return the array used to store the point ids that define the cells' connectivity.
Definition: vtkCellArray.h:644
void SetTraversalCellId(vtkIdType cellId)
Get/Set the current cellId for traversal.
void SetData(vtkTypeInt64Array *offsets, vtkTypeInt64Array *connectivity)
Set the internal data arrays to the supplied offsets and connectivity arrays.
vtkIdType GetNumberOfCells() const
Get the number of cells in the array.
Definition: vtkCellArray.h:430
vtkIdType InsertNextCell(vtkCell *cell)
Insert a cell object.
virtual void SetNumberOfCells(vtkIdType)
Set the number of cells in the array.
vtkIdType GetNumberOfConnectivityIds() const
Get the size of the connectivity array that stores the point ids.
Definition: vtkCellArray.h:479
vtkNew< vtkIdTypeArray > LegacyData
void PrintSelf(ostream &os, vtkIndent indent) override
Standard methods for instantiation, type information, and printing.
void PrintDebug(ostream &os)
Standard methods for instantiation, type information, and printing.
vtkCellArrayIterator * NewIterator()
NewIterator returns a new instance of vtkCellArrayIterator that is initialized to point at the first ...
void UpdateCellCount(int npts)
Used in conjunction with InsertNextCell(int npts) and InsertCellPoint() to update the number of point...
GetReturnType< Functor, Args... > Visit(Functor &&functor, Args &&... args) const
vtkNew< vtkIdList > TempCell
vtkIdType GetSize()
Get the size of the allocated connectivity array.
vtkDataArray * GetOffsetsArray()
Return the array used to store cell offsets.
Definition: vtkCellArray.h:623
vtkTypeList::Create< ArrayType32, ArrayType64 > StorageArrayList
List of possible array types used for storage.
Definition: vtkCellArray.h:321
vtkIdType InsertNextCell(const std::initializer_list< vtkIdType > &cell)
Overload that allows InsertNextCell({0, 1, 2}) syntax.
Definition: vtkCellArray.h:771
void Append(vtkCellArray *src, vtkIdType pointOffset=0)
Append cells from src into this.
bool IsStorage64Bit() const
Definition: vtkCellArray.h:550
void ReplaceCellAtId(vtkIdType cellId, vtkIdType cellSize, const vtkIdType *cellPoints)
Replaces the point ids for the specified cell with the supplied list.
bool CanConvertTo64BitStorage() const
Check if the existing data can safely be converted to use 32- or 64- bit storage.
abstract class to specify cell behavior
Definition: vtkCell.h:151
vtkIdList * GetPointIds()
Return the list of point ids defining the cell.
Definition: vtkCell.h:246
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:166
list of point or cell ids
Definition: vtkIdList.h:144
void SetNumberOfIds(vtkIdType number)
Specify the number of ids for this object to hold.
vtkIdType GetNumberOfIds() const noexcept
Return the number of id's in the list.
Definition: vtkIdList.h:170
void Reset()
Reset to an empty state but retain previously allocated memory.
Definition: vtkIdList.h:255
vtkIdType * GetPointer(vtkIdType i)
Get a pointer to a particular data index.
Definition: vtkIdList.h:236
dynamic, self-adjusting array of vtkIdType
a simple class to control print indentation
Definition: vtkIndent.h:120
static vtkMallocingFunction GetCurrentMallocFunction()
static vtkFreeingFunction GetAlternateFreeFunction()
static bool GetUsingMemkind()
A global state flag that controls whether vtkObjects are constructed in the usual way (the default) o...
abstract base class for most VTK objects
Definition: vtkObject.h:83
static vtkSmartPointer< T > New()
Create an instance of a VTK object.
vtkSmartPointer< vtkDataArray > GetData(const Ioss::GroupingEntity *entity, const std::string &fieldname, Ioss::Transform *transform=nullptr, Cache *cache=nullptr, const std::string &cachekey=std::string())
Returns a VTK array for a given field (fieldname) on the chosen block (or set) entity.
vtkSmartPointer< vtkCellArray > GetConnectivity(Ioss::GroupingEntity *group_entity, int &vtk_topology_type, Cache *cache=nullptr)
Read connectivity information from the group_entity.
@ value
Definition: vtkX3D.h:232
@ type
Definition: vtkX3D.h:528
@ data
Definition: vtkX3D.h:327
VisitState< ArrayType32 > & GetArrays32()
const VisitState< ArrayType64 > & GetArrays64() const
const VisitState< ArrayType32 > & GetArrays32() const
VisitState< ArrayType64 > & GetArrays64()
ArrayType * GetOffsets()
Definition: vtkCellArray.h:941
vtkIdType GetNumberOfCells() const
vtkIdType GetEndOffset(vtkIdType cellId) const
vtkSmartPointer< ArrayType > Offsets
Definition: vtkCellArray.h:999
vtkSmartPointer< ArrayType > Connectivity
Definition: vtkCellArray.h:998
ArrayType * GetConnectivity()
Definition: vtkCellArray.h:944
decltype(vtk::DataArrayValueRange< 1 >(std::declval< ArrayType >())) CellRangeType
Definition: vtkCellArray.h:931
CellRangeType GetCellRange(vtkIdType cellId)
const ArrayType * GetOffsets() const
Definition: vtkCellArray.h:942
vtkIdType GetBeginOffset(vtkIdType cellId) const
vtkIdType GetCellSize(vtkIdType cellId) const
const ArrayType * GetConnectivity() const
Definition: vtkCellArray.h:945
typename ArrayType::ValueType ValueType
Definition: vtkCellArray.h:930
void operator()(CellStateT &state, const vtkIdType cellId, vtkIdList *ids)
std::enable_if<!CanShareConnPtr< CellStateT >::value, void >::type operator()(CellStateT &state, const vtkIdType cellId, vtkIdType &cellSize, vtkIdType const *&cellPoints, vtkIdList *temp)
std::enable_if< CanShareConnPtr< CellStateT >::value, void >::type operator()(CellStateT &state, const vtkIdType cellId, vtkIdType &cellSize, vtkIdType const *&cellPoints, vtkIdList *vtkNotUsed(temp))
vtkIdType operator()(CellStateT &state, vtkIdType cellId)
vtkIdType operator()(CellStateT &state, const vtkIdType npts, const vtkIdType pts[])
vtkIdType operator()(CellStateT &state, const vtkIdType npts)
void operator()(CellStateT &state)
void operator()(CellStateT &state, const vtkIdType npts)
Remove all duplicate types from TypeList TList, storing the new list in Result.
Definition: vtkTypeList.h:128
VisitState< ArrayType32 > * Int32
VisitState< ArrayType64 > * Int64
int vtkTypeBool
Definition: vtkABI.h:71
STL-compatible iterable ranges that provide access to vtkDataArray elements.
int vtkIdType
Definition: vtkType.h:327
#define VTK_SIZEHINT(...)
#define VTK_EXPECTS(x)
#define VTK_NEWINSTANCE