VTK  9.6.20260129
vtkCellArray.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
235
236#ifndef vtkCellArray_h
237#define vtkCellArray_h
238
239#include "vtkAbstractCellArray.h"
240#include "vtkCommonDataModelModule.h" // For export macro
241#include "vtkWrappingHints.h" // For VTK_MARSHALMANUAL
242
243#include "vtkAOSDataArrayTemplate.h" // Needed for inline methods
244#include "vtkAffineArray.h" // Needed for inline methods
245#include "vtkCell.h" // Needed for inline methods
246#include "vtkDataArrayRange.h" // Needed for inline methods
247#include "vtkDeprecation.h" // For VTK_DEPRECATED_IN_9_6_0
248#include "vtkFeatures.h" // for VTK_USE_MEMKIND
249#include "vtkSmartPointer.h" // For vtkSmartPointer
250#include "vtkTypeInt32Array.h" // Needed for inline methods
251#include "vtkTypeInt64Array.h" // Needed for inline methods
252#include "vtkTypeList.h" // Needed for ArrayList definition
253
254#include <cassert> // Needed for assert
255#include <initializer_list> // Needed for API
256#include <type_traits> // Needed for std::is_same
257#include <utility> // Needed for std::forward
258
279#define VTK_CELL_ARRAY_V2
280
281VTK_ABI_NAMESPACE_BEGIN
283class vtkIdTypeArray;
284
285class VTKCOMMONDATAMODEL_EXPORT VTK_MARSHALMANUAL vtkCellArray : public vtkAbstractCellArray
286{
287public:
288 using ArrayType32 VTK_DEPRECATED_IN_9_6_0("Use AOSArray32 instead.") = vtkTypeInt32Array;
289 using ArrayType64 VTK_DEPRECATED_IN_9_6_0("Use AOSArray64 instead.") = vtkTypeInt64Array;
294
296
300 static vtkCellArray* New();
302 void PrintSelf(ostream& os, vtkIndent indent) override;
303 void PrintDebug(ostream& os);
305
307
316 "Use vtkArrayDispatch::OffsetsArrays/ConnectivityArrays instead.") =
317 vtkTypeList::Create<vtkTypeInt32Array, vtkTypeInt64Array>;
319
321
333
341 */
342 VTK_DEPRECATED_IN_9_6_0("Use AllocateEstimate or AllocateExact instead.")
343 vtkTypeBool Allocate(vtkIdType sz, vtkIdType vtkNotUsed(ext) = 1000)
344 {
345 return this->AllocateExact(sz, sz) ? 1 : 0;
346 }
347
355 * @sa Squeeze AllocateExact AllocateCopy
356 */
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
380 * @sa Squeeze AllocateEstimate AllocateExact
381 */
382 bool AllocateCopy(vtkCellArray* other)
383 {
384 return this->AllocateExact(other->GetNumberOfCells(), other->GetNumberOfConnectivityIds());
385 }
386
396 bool ResizeExact(vtkIdType numCells, vtkIdType connectivitySize);
397
401 void Initialize() override;
402
406 void Reset();
407
413 void Squeeze();
414
425 bool IsValid();
426
430 vtkIdType GetNumberOfCells() const override { return this->Offsets->GetNumberOfValues() - 1; }
431
436 vtkIdType GetNumberOfOffsets() const override { return this->Offsets->GetNumberOfValues(); }
437
439 * Get the offset (into the connectivity) for a specified cell id.
440 */
441 vtkIdType GetOffset(vtkIdType cellId) override
442 {
443 return static_cast<vtkIdType>(this->Offsets->GetComponent(cellId, 0));
444 }
445
447 * Set the offset (into the connectivity) for a specified cell id.
448 */
449 void SetOffset(vtkIdType cellId, vtkIdType offset)
450 {
451 this->Offsets->SetComponent(cellId, 0, static_cast<double>(offset));
452 }
453
455 * Get the size of the connectivity array that stores the point ids.
456 */
458 {
459 return this->Connectivity->GetNumberOfValues();
460 }
461
468
470
483 void SetData(AOSArray32*, AOSArray32* connectivity);
484 void SetData(AOSArray64*, AOSArray64* connectivity);
485 void SetData(AffineArray32*, AOSArray32* connectivity);
486 void SetData(AffineArray64*, AOSArray64* connectivity);
488
503 bool SetData(vtkDataArray* offsets, vtkDataArray* connectivity);
504
530 Generic,
531 };
532
536 StorageTypes GetStorageType() const noexcept { return this->StorageType; }
537
541 bool IsStorage64Bit() const { return this->StorageType == StorageTypes::Int64; }
542
546 bool IsStorage32Bit() const { return this->StorageType == StorageTypes::Int32; }
547
552
557
559 * @return True if the internal storage is using FixedSize arrays.
560 */
561 bool IsStorageFixedSize() const
562 {
563 return this->IsStorageFixedSize32Bit() || this->IsStorageFixedSize64Bit();
564 }
565
569 bool IsStorageGeneric() const { return this->StorageType == StorageTypes::Generic; }
570
575 * a pointer to vtkIdType can be used.
576 */
577 bool IsStorageShareable() const override
578 {
579 switch (this->StorageType)
580 {
583 return std::is_same_v<vtkTypeInt32, vtkIdType>;
586 return std::is_same_v<vtkTypeInt64, vtkIdType>;
588 default:
589 return false;
590 }
591 }
592
605 void UseFixedSize64BitStorage(vtkIdType cellSize);
608
624 {
625 switch (type)
626 {
628 return this->CanConvertTo32BitStorage();
630 return this->CanConvertTo64BitStorage();
636 default:
637 return true;
638 }
639 }
641
666 {
667 switch (type)
668 {
670 return this->ConvertTo32BitStorage();
672 return this->ConvertTo64BitStorage();
674 return this->ConvertToFixedSize32BitStorage();
676 return this->ConvertToFixedSize64BitStorage();
678 default:
679 return true;
680 }
681 }
683
689 vtkDataArray* GetOffsetsArray() const { return this->Offsets; }
690 VTK_DEPRECATED_IN_9_6_0("Use GetOffsetsAOSArray32() instead.")
691 vtkTypeInt32Array* GetOffsetsArray32() const
692 {
693 return vtkTypeInt32Array::FastDownCast(this->Offsets);
694 }
696 VTK_DEPRECATED_IN_9_6_0("Use GetOffsetsAOSArray64() instead.")
697 vtkTypeInt64Array* GetOffsetsArray64() const
698 {
699 return vtkTypeInt64Array::FastDownCast(this->Offsets);
701 AOSArray64* GetOffsetsAOSArray64() const { return AOSArray64::FastDownCast(this->Offsets); }
702 AffineArray32* GetOffsetsAffineArray32() const
703 {
705 }
706 AffineArray64* GetOffsetsAffineArray64() const
707 {
708 return AffineArray64::FastDownCast(this->Offsets);
709 }
711
719 VTK_DEPRECATED_IN_9_6_0("Use GetConnectivityAOSArray32() instead.")
720 vtkTypeInt32Array* GetConnectivityArray32() const
721 {
722 return vtkTypeInt32Array::FastDownCast(this->Connectivity);
723 }
724 AOSArray32* GetConnectivityAOSArray32() const
725 {
726 return AOSArray32::FastDownCast(this->Connectivity);
728 VTK_DEPRECATED_IN_9_6_0("Use GetConnectivityAOSArray64() instead.")
729 vtkTypeInt64Array* GetConnectivityArray64() const
730 {
731 return vtkTypeInt64Array::FastDownCast(this->Connectivity);
732 }
733 AOSArray64* GetConnectivityAOSArray64() const
734 {
735 return AOSArray64::FastDownCast(this->Connectivity);
736 }
738
747 vtkIdType IsHomogeneous() const override;
748
758 void InitTraversal();
759
774 int GetNextCell(vtkIdType& npts, vtkIdType const*& pts) VTK_SIZEHINT(pts, npts);
775
786 int GetNextCell(vtkIdList* pts);
787
798 inline void GetCellAtId(vtkIdType cellId, vtkIdType& cellSize, vtkIdType const*& cellPoints,
799 vtkIdList* ptIds) VTK_SIZEHINT(cellPoints, cellSize)
800 VTK_EXPECTS(0 <= cellId && cellId < GetNumberOfCells()) override;
801
807 inline void GetCellAtId(vtkIdType cellId, vtkIdList* pts)
808 VTK_EXPECTS(0 <= cellId && cellId < GetNumberOfCells()) override;
809
817 inline void GetCellAtId(vtkIdType cellId, vtkIdType& cellSize, vtkIdType* cellPoints)
818 VTK_SIZEHINT(cellPoints, cellSize)
819 VTK_EXPECTS(0 <= cellId && cellId < GetNumberOfCells()) override;
820
824 vtkIdType GetCellPointAtId(vtkIdType cellId, vtkIdType cellPointIndex) const
825 VTK_EXPECTS(0 <= cellId && cellId < GetNumberOfCells())
826 VTK_EXPECTS(0 <= cellPointIndex && cellPointIndex < this->GetCellSize(cellId));
827
831 vtkIdType GetCellSize(vtkIdType cellId) const override;
832
836 inline vtkIdType InsertNextCell(vtkCell* cell);
837
842 inline vtkIdType InsertNextCell(vtkIdType npts, const vtkIdType* pts) VTK_SIZEHINT(pts, npts);
843
849
855 * `InsertNextCell(int)` overload instead.
856 */
857 vtkIdType InsertNextCell(const std::initializer_list<vtkIdType>& cell)
858 {
859 return this->InsertNextCell(static_cast<vtkIdType>(cell.size()), cell.begin());
860 }
861
868 inline vtkIdType InsertNextCell(int npts);
869
874 void InsertCellPoint(vtkIdType id);
875
880 void UpdateCellCount(int npts);
881
890 void SetTraversalCellId(vtkIdType cellId);
892
896 void ReverseCellAtId(vtkIdType cellId) VTK_EXPECTS(0 <= cellId && cellId < GetNumberOfCells());
897
906 void ReplaceCellAtId(vtkIdType cellId, vtkIdList* list);
907 void ReplaceCellAtId(vtkIdType cellId, vtkIdType cellSize, const vtkIdType* cellPoints)
908 VTK_EXPECTS(0 <= cellId && cellId < GetNumberOfCells()) VTK_SIZEHINT(cellPoints, cellSize);
910
918 void ReplaceCellPointAtId(vtkIdType cellId, vtkIdType cellPointIndex, vtkIdType newPointId);
919
925 * results.
926 */
927 void ReplaceCellAtId(vtkIdType cellId, const std::initializer_list<vtkIdType>& cell)
928 {
929 this->ReplaceCellAtId(cellId, static_cast<vtkIdType>(cell.size()), cell.begin());
930 }
931
936 int GetMaxCellSize() override;
937
941 void DeepCopy(vtkAbstractCellArray* ca) override;
942
946 void ShallowCopy(vtkAbstractCellArray* ca) override;
947
951 void Append(vtkCellArray* src, vtkIdType pointOffset = 0);
952
964
977 void ImportLegacyFormat(const vtkIdType* data, vtkIdType len) VTK_SIZEHINT(data, len);
979
991 void AppendLegacyFormat(vtkIdTypeArray* data, vtkIdType ptOffset = 0);
992 void AppendLegacyFormat(const vtkIdType* data, vtkIdType len, vtkIdType ptOffset = 0)
993 VTK_SIZEHINT(data, len);
995
1004 unsigned long GetActualMemorySize() const;
1005
1006 // The following code is used to support
1007
1008 // The wrappers get understandably confused by some of the template code below
1009#ifndef __VTK_WRAP__
1010 /*
1011 * Utilities class that every dispatch functor used with Dispatch() must inherit
1012 * or optionally use its static methods.
1013 */
1014 struct DispatchUtilities
1016 template <class ArrayT>
1019 template <class OffsetsT>
1020 vtkIdType GetNumberOfCells(OffsetsT* offsets)
1021 {
1022 return offsets->GetNumberOfValues() - 1;
1023 }
1025 template <class ArrayT>
1026 static decltype(vtk::DataArrayValueRange<1, vtkIdType>(std::declval<ArrayT>())) GetRange(
1027 ArrayT* array)
1028 {
1030 }
1032 template <class OffsetsT>
1033 static vtkIdType GetBeginOffset(OffsetsT* offsets, vtkIdType cellId)
1034 {
1035 return static_cast<vtkIdType>(GetRange(offsets)[cellId]);
1036 }
1038 template <class OffsetsT>
1039 static vtkIdType GetEndOffset(OffsetsT* offsets, vtkIdType cellId)
1040 {
1041 return static_cast<vtkIdType>(GetRange(offsets)[cellId + 1]);
1042 }
1044 template <class OffsetsT>
1045 static vtkIdType GetCellSize(OffsetsT* offsets, vtkIdType cellId)
1046 {
1047 auto offsetsRange = GetRange(offsets);
1048 return static_cast<vtkIdType>(offsetsRange[cellId + 1] - offsetsRange[cellId]);
1049 }
1050
1051 template <class OffsetsT, class ConnectivityT>
1052 static decltype(vtk::DataArrayValueRange<1, vtkIdType>(std::declval<ConnectivityT>()))
1053 GetCellRange(OffsetsT* offsets, ConnectivityT* conn, vtkIdType cellId)
1054 {
1055 auto offsetsRange = GetRange(offsets);
1057 conn, offsetsRange[cellId], offsetsRange[cellId + 1]);
1058 }
1059 };
1060
1062
1125 template <typename Functor, typename... Args>
1126 void Dispatch(Functor&& functor, Args&&... args)
1127 {
1128 switch (this->StorageType)
1129 {
1131 functor(static_cast<AOSArray32*>(this->Offsets.Get()),
1132 static_cast<AOSArray32*>(this->Connectivity.Get()), std::forward<Args>(args)...);
1133 break;
1135 functor(static_cast<AOSArray64*>(this->Offsets.Get()),
1136 static_cast<AOSArray64*>(this->Connectivity.Get()), std::forward<Args>(args)...);
1137 break;
1139 functor(static_cast<AffineArray32*>(this->Offsets.Get()),
1140 static_cast<AOSArray32*>(this->Connectivity.Get()), std::forward<Args>(args)...);
1141 break;
1143 functor(static_cast<AffineArray64*>(this->Offsets.Get()),
1144 static_cast<AOSArray64*>(this->Connectivity.Get()), std::forward<Args>(args)...);
1145 break;
1147 default:
1148 functor(this->Offsets.Get(), this->Connectivity.Get(), std::forward<Args>(args)...);
1149 break;
1150 }
1152 template <typename Functor, typename... Args>
1153 void Dispatch(Functor&& functor, Args&&... args) const
1154 {
1155 switch (this->StorageType)
1156 {
1157 case StorageTypes::Int32:
1158 functor(static_cast<AOSArray32*>(this->Offsets.Get()),
1159 static_cast<AOSArray32*>(this->Connectivity.Get()), std::forward<Args>(args)...);
1160 break;
1161 case StorageTypes::Int64:
1162 functor(static_cast<AOSArray64*>(this->Offsets.Get()),
1163 static_cast<AOSArray64*>(this->Connectivity.Get()), std::forward<Args>(args)...);
1164 break;
1165 case StorageTypes::FixedSizeInt32:
1166 functor(static_cast<AffineArray32*>(this->Offsets.Get()),
1167 static_cast<AOSArray32*>(this->Connectivity.Get()), std::forward<Args>(args)...);
1168 break;
1169 case StorageTypes::FixedSizeInt64:
1170 functor(static_cast<AffineArray64*>(this->Offsets.Get()),
1171 static_cast<AOSArray64*>(this->Connectivity.Get()), std::forward<Args>(args)...);
1172 break;
1173 case StorageTypes::Generic:
1174 default:
1175 functor(this->Offsets.Get(), this->Connectivity.Get(), std::forward<Args>(args)...);
1176 break;
1177 }
1178 }
1180
1181 // Holds connectivity and offset arrays of the given ArrayType.
1182 // VTK_DEPRECATED_IN_9_6_0("Use DispatchUtilities")
1183 template <typename ArrayT>
1186 using ArrayType = ArrayT;
1187 using ValueType = typename ArrayType::ValueType;
1188 using CellRangeType = decltype(vtk::DataArrayValueRange<1>(std::declval<ArrayType>()));
1189
1190 // We can't just use is_same here, since binary compatible representations
1191 // (e.g. int and long) are distinct types. Instead, ensure that ValueType
1192 // is a signed integer the same size as vtkIdType.
1193 // If this value is true, ValueType pointers may be safely converted to
1194 // vtkIdType pointers via reinterpret cast.
1195 static constexpr bool ValueTypeIsSameAsIdType = std::is_integral<ValueType>::value &&
1196 std::is_signed<ValueType>::value && (sizeof(ValueType) == sizeof(vtkIdType));
1198 ArrayType* GetOffsets() { return this->Offsets; }
1199 const ArrayType* GetOffsets() const { return this->Offsets; }
1201 ArrayType* GetConnectivity() { return this->Connectivity; }
1202 const ArrayType* GetConnectivity() const { return this->Connectivity; }
1203
1204 vtkIdType GetNumberOfCells() const;
1205
1206 vtkIdType GetBeginOffset(vtkIdType cellId) const;
1207
1208 vtkIdType GetEndOffset(vtkIdType cellId) const;
1209
1210 vtkIdType GetCellSize(vtkIdType cellId) const;
1211
1213
1214 friend class vtkCellArray;
1216 protected:
1217 VisitState()
1218 {
1221 this->Offsets->InsertNextValue(0);
1223 {
1224 this->IsInMemkind = true;
1227 ~VisitState() = default;
1228 void* operator new(size_t nSize)
1229 {
1230 void* r;
1231#ifdef VTK_USE_MEMKIND
1233#else
1234 r = malloc(nSize);
1235#endif
1236 return r;
1237 }
1238 void operator delete(void* p)
1239 {
1240#ifdef VTK_USE_MEMKIND
1241 VisitState* a = static_cast<VisitState*>(p);
1242 if (a->IsInMemkind)
1243 {
1245 }
1246 else
1247 {
1248 free(p);
1249 }
1250#else
1251 free(p);
1252#endif
1257
1258 private:
1259 VisitState(const VisitState&) = delete;
1260 VisitState& operator=(const VisitState&) = delete;
1261 bool IsInMemkind = false;
1262 };
1263
1264private: // Helpers that allow Visit to return a value:
1265 template <typename Functor, typename... Args>
1266 using GetReturnType = decltype(std::declval<Functor>()(
1267 std::declval<VisitState<AOSArray32>&>(), std::declval<Args>()...));
1268
1269 template <typename Functor, typename... Args>
1270 struct ReturnsVoid : std::is_same<GetReturnType<Functor, Args...>, void>
1271 {
1272 };
1273
1274public:
1344 template <typename Functor, typename... Args,
1345 typename = typename std::enable_if<ReturnsVoid<Functor, Args...>::value>::type>
1346 VTK_DEPRECATED_IN_9_6_0("Use Dispatch instead")
1347 void Visit(Functor&& functor, Args&&... args)
1348 {
1349 switch (this->StorageType)
1350 {
1352 {
1356 functor(state, std::forward<Args>(args)...);
1357 break;
1358 }
1360 {
1364 functor(state, std::forward<Args>(args)...);
1365 break;
1366 }
1370 default:
1371 {
1372 vtkWarningMacro("Use Dispatch");
1373 break;
1374 }
1375 }
1376 }
1377
1378 template <typename Functor, typename... Args,
1379 typename = typename std::enable_if<ReturnsVoid<Functor, Args...>::value>::type>
1380 VTK_DEPRECATED_IN_9_6_0("Use Dispatch instead")
1381 void Visit(Functor&& functor, Args&&... args) const
1382 {
1383 switch (this->StorageType)
1384 {
1386 {
1390 functor(state, std::forward<Args>(args)...);
1391 break;
1392 }
1394 {
1398 functor(state, std::forward<Args>(args)...);
1399 break;
1400 }
1404 default:
1405 {
1406 vtkWarningMacro("Use Dispatch");
1407 break;
1408 }
1409 }
1410 }
1411
1412 template <typename Functor, typename... Args,
1413 typename = typename std::enable_if<!ReturnsVoid<Functor, Args...>::value>::type>
1414 VTK_DEPRECATED_IN_9_6_0("Use Dispatch instead")
1415 GetReturnType<Functor, Args...> Visit(Functor&& functor, Args&&... args)
1416 {
1417 switch (this->StorageType)
1418 {
1420 {
1424 return functor(state, std::forward<Args>(args)...);
1425 }
1427 {
1431 return functor(state, std::forward<Args>(args)...);
1432 }
1436 default:
1437 {
1438 vtkWarningMacro("Use Dispatch");
1439 return GetReturnType<Functor, Args...>();
1440 }
1441 }
1442 }
1443 template <typename Functor, typename... Args,
1444 typename = typename std::enable_if<!ReturnsVoid<Functor, Args...>::value>::type>
1445 VTK_DEPRECATED_IN_9_6_0("Use Dispatch instead")
1446 GetReturnType<Functor, Args...> Visit(Functor&& functor, Args&&... args) const
1447 {
1448 switch (this->StorageType)
1449 {
1451 {
1455 return functor(state, std::forward<Args>(args)...);
1456 }
1458 {
1462 return functor(state, std::forward<Args>(args)...);
1463 }
1467 default:
1468 {
1469 vtkWarningMacro("Use Dispatch");
1470 return GetReturnType<Functor, Args...>();
1471 }
1472 }
1473 }
1474
1475#endif // __VTK_WRAP__
1476
1478
1486 static void SetDefaultStorageIs64Bit(bool val) { vtkCellArray::DefaultStorageIs64Bit = val; }
1488
1489 //=================== Begin Legacy Methods ===================================
1490 // These should be deprecated at some point as they are confusing or very slow
1491
1498 VTK_DEPRECATED_IN_9_6_0("This call has no effect.")
1499 virtual void SetNumberOfCells(vtkIdType);
1500
1512 VTK_DEPRECATED_IN_9_6_0("Use AllocateEstimate directly instead.")
1513 vtkIdType EstimateSize(vtkIdType numCells, int maxPtsPerCell);
1514
1523 VTK_DEPRECATED_IN_9_6_0("Method incompatible with current internal storage.")
1525
1532 VTK_DEPRECATED_IN_9_6_0("Method incompatible with current internal storage.")
1534
1544 VTK_DEPRECATED_IN_9_6_0("Use GetCellAtId.")
1545 void GetCell(vtkIdType loc, vtkIdType& npts, const vtkIdType*& pts)
1546 VTK_EXPECTS(0 <= loc && loc < GetNumberOfConnectivityEntries()) VTK_SIZEHINT(pts, npts);
1547
1554 VTK_DEPRECATED_IN_9_6_0("Use GetCellAtId.")
1555 void GetCell(vtkIdType loc, vtkIdList* pts)
1556 VTK_EXPECTS(0 <= loc && loc < GetNumberOfConnectivityEntries());
1557
1564 VTK_DEPRECATED_IN_9_6_0("Use GetNumberOfCells.")
1565 vtkIdType GetInsertLocation(int npts);
1566
1574 VTK_DEPRECATED_IN_9_6_0("Use GetTraversalCellId.")
1576 VTK_DEPRECATED_IN_9_6_0("Use GetTraversalCellId.")
1578 VTK_DEPRECATED_IN_9_6_0("Use SetTraversalCellId.")
1581
1589 VTK_DEPRECATED_IN_9_6_0("Use ReverseCellAtId.")
1590 void ReverseCell(vtkIdType loc) VTK_EXPECTS(0 <= loc && loc < GetNumberOfConnectivityEntries());
1591
1603 VTK_DEPRECATED_IN_9_6_0("Use ReplaceCellAtId.")
1604 void ReplaceCell(vtkIdType loc, int npts, const vtkIdType pts[])
1605 VTK_EXPECTS(0 <= loc && loc < GetNumberOfConnectivityEntries()) VTK_SIZEHINT(pts, npts);
1606
1621 VTK_DEPRECATED_IN_9_6_0("Use ImportLegacyFormat or SetData instead.")
1622 void SetCells(vtkIdType ncells, vtkIdTypeArray* cells);
1623
1635 "Use ExportLegacyFormat, or GetOffsetsArray/GetConnectivityArray instead.")
1637
1638 //=================== End Legacy Methods =====================================
1639
1640 friend class vtkCellArrayIterator;
1642protected:
1643 vtkCellArray();
1644 ~vtkCellArray() override;
1650
1652
1653 static bool DefaultStorageIs64Bit;
1654
1655private:
1656 vtkCellArray(const vtkCellArray&) = delete;
1657 void operator=(const vtkCellArray&) = delete;
1658};
1659
1660template <typename ArrayT>
1662{
1663 return this->Offsets->GetNumberOfValues() - 1;
1664}
1665
1666template <typename ArrayT>
1668{
1669 return static_cast<vtkIdType>(this->Offsets->GetValue(cellId));
1670}
1671
1672template <typename ArrayT>
1674{
1675 return static_cast<vtkIdType>(this->Offsets->GetValue(cellId + 1));
1676}
1677
1678template <typename ArrayT>
1680{
1681 return this->GetEndOffset(cellId) - this->GetBeginOffset(cellId);
1682}
1683
1684template <typename ArrayT>
1691VTK_ABI_NAMESPACE_END
1692
1694{
1695VTK_ABI_NAMESPACE_BEGIN
1696
1698{
1699 // Insert full cell
1700 template <class OffsetsT, class ConnectivityT>
1701 void operator()(OffsetsT* offsets, ConnectivityT* conn, const vtkIdType npts,
1702 const vtkIdType pts[], vtkIdType& cellId)
1703 {
1704 using ValueType = GetAPIType<OffsetsT>;
1705 using OffsetsAccessorType = vtkDataArrayAccessor<OffsetsT>;
1706 using ConnectivityAccessorType = vtkDataArrayAccessor<ConnectivityT>;
1707 OffsetsAccessorType offsetsAccesor(offsets);
1708 ConnectivityAccessorType connAccesor(conn);
1709
1710 cellId = offsets->GetNumberOfValues() - 1;
1711
1712 offsetsAccesor.InsertNext(static_cast<ValueType>(conn->GetNumberOfValues() + npts));
1713
1714 for (vtkIdType i = 0; i < npts; ++i)
1715 {
1716 connAccesor.InsertNext(static_cast<ValueType>(pts[i]));
1717 }
1718 }
1719
1720 // Just update offset table (for incremental API)
1721 template <class OffsetsT, class ConnectivityT>
1722 void operator()(OffsetsT* offsets, ConnectivityT* conn, const vtkIdType npts, vtkIdType& cellId)
1723 {
1724 using ValueType = GetAPIType<OffsetsT>;
1725 using AccessorType = vtkDataArrayAccessor<OffsetsT>;
1726 AccessorType offsetsAccesor(offsets);
1727
1728 cellId = offsets->GetNumberOfValues() - 1;
1729
1730 offsetsAccesor.InsertNext(static_cast<ValueType>(conn->GetNumberOfValues() + npts));
1731 }
1732};
1733
1734// for incremental API:
1736{
1737 template <class OffsetsT, class ConnectivityT>
1738 void operator()(OffsetsT* offsets, ConnectivityT* vtkNotUsed(conn), const vtkIdType npts)
1739 {
1740 using ValueType = GetAPIType<OffsetsT>;
1741
1742 auto offsetsRange = GetRange(offsets);
1743 const ValueType cellBegin = offsetsRange[offsets->GetMaxId() - 1];
1744 offsetsRange[offsets->GetMaxId()] = static_cast<ValueType>(cellBegin + npts);
1745 }
1746};
1747
1749{
1750 template <class OffsetsT, class ConnectivityT>
1752 OffsetsT* offsets, ConnectivityT* vtkNotUsed(conn), vtkIdType cellId, vtkIdType& cellSize)
1753 {
1754 cellSize = GetCellSize(offsets, cellId);
1755 }
1756};
1757
1759{
1760 template <class OffsetsT, class ConnectivityT>
1761 void operator()(OffsetsT* offsets, ConnectivityT* conn, const vtkIdType cellId, vtkIdList* ids)
1762 {
1763 auto offsetsRange = GetRange(offsets);
1764 const auto& beginOffset = offsetsRange[cellId];
1765 const auto& endOffset = offsetsRange[cellId + 1];
1766 const vtkIdType cellSize = static_cast<vtkIdType>(endOffset - beginOffset);
1767 const auto cellConnectivity = GetRange(conn).begin() + beginOffset;
1768
1769 // ValueType differs from vtkIdType, so we have to copy into a temporary buffer:
1770 ids->SetNumberOfIds(cellSize);
1771 vtkIdType* idPtr = ids->GetPointer(0);
1772 for (vtkIdType i = 0; i < cellSize; ++i)
1773 {
1774 idPtr[i] = static_cast<vtkIdType>(cellConnectivity[i]);
1775 }
1776 }
1777
1778 template <class OffsetsT, class ConnectivityT>
1779 void operator()(OffsetsT* offsets, ConnectivityT* conn, const vtkIdType cellId,
1780 vtkIdType& cellSize, vtkIdType* cellPoints)
1781 {
1782 auto offsetsRange = GetRange(offsets);
1783 const auto& beginOffset = offsetsRange[cellId];
1784 const auto& endOffset = offsetsRange[cellId + 1];
1785 cellSize = static_cast<vtkIdType>(endOffset - beginOffset);
1786 const auto cellConnectivity = GetRange(conn).begin() + beginOffset;
1787
1788 // ValueType differs from vtkIdType, so we have to copy into a temporary buffer:
1789 for (vtkIdType i = 0; i < cellSize; ++i)
1790 {
1791 cellPoints[i] = static_cast<vtkIdType>(cellConnectivity[i]);
1792 }
1793 }
1794
1795 // SFINAE helper to check if a Functors's connectivity array's memory can be used as a vtkIdType*.
1796 template <typename ConnectivityT>
1798 {
1799 static constexpr bool value =
1800 std::is_base_of_v<vtkAOSDataArrayTemplate<vtkIdType>, ConnectivityT>;
1801 };
1802
1803 template <class OffsetsT, class ConnectivityT>
1804 typename std::enable_if_t<CanShareConnPtr<ConnectivityT>::value, void> operator()(
1805 OffsetsT* offsets, ConnectivityT* conn, const vtkIdType cellId, vtkIdType& cellSize,
1806 vtkIdType const*& cellPoints, vtkIdList* vtkNotUsed(temp))
1807 {
1808 auto offsetsRange = GetRange(offsets);
1809 const auto& beginOffset = offsetsRange[cellId];
1810 const auto& endOffset = offsetsRange[cellId + 1];
1811 cellSize = static_cast<vtkIdType>(endOffset - beginOffset);
1812 // This is safe, see CanShareConnPtr helper above.
1813 cellPoints = conn->GetPointer(beginOffset);
1814 }
1815
1816 template <class OffsetsT, class ConnectivityT>
1817 typename std::enable_if_t<!CanShareConnPtr<ConnectivityT>::value, void> operator()(
1818 OffsetsT* offsets, ConnectivityT* conn, const vtkIdType cellId, vtkIdType& cellSize,
1819 vtkIdType const*& cellPoints, vtkIdList* temp)
1820 {
1821 auto offsetsRange = GetRange(offsets);
1822 const auto& beginOffset = offsetsRange[cellId];
1823 const auto& endOffset = offsetsRange[cellId + 1];
1824 cellSize = static_cast<vtkIdType>(endOffset - beginOffset);
1825 const auto cellConnectivity = GetRange(conn).begin() + beginOffset;
1826
1827 temp->SetNumberOfIds(cellSize);
1828 vtkIdType* tempPtr = temp->GetPointer(0);
1829 for (vtkIdType i = 0; i < cellSize; ++i)
1830 {
1831 tempPtr[i] = static_cast<vtkIdType>(cellConnectivity[i]);
1832 }
1833
1834 cellPoints = tempPtr;
1835 }
1836};
1837
1839{
1840 template <class OffsetsT, class ConnectivityT>
1841 void operator()(OffsetsT* offsets, ConnectivityT* conn, vtkIdType cellId,
1842 vtkIdType cellPointIndex, vtkIdType& pointId)
1843 {
1844 pointId =
1845 static_cast<vtkIdType>(GetRange(conn)[GetBeginOffset(offsets, cellId) + cellPointIndex]);
1846 }
1847};
1848
1850{
1851 template <class OffsetsT, class ConnectivityT>
1852 void operator()(OffsetsT* offsets, ConnectivityT* conn)
1853 {
1854 using ValueType = GetAPIType<OffsetsT>;
1855 using AccessorType = vtkDataArrayAccessor<OffsetsT>;
1856 offsets->Reset();
1857 conn->Reset();
1858 AccessorType accessor(offsets);
1859 ValueType firstOffset = 0;
1860 accessor.InsertNext(firstOffset);
1861 }
1862};
1863
1865{
1866 template <class OffsetsT, class ConnectivityT>
1867 void operator()(OffsetsT* vtkNotUsed(offsets), ConnectivityT* conn, vtkIdType id)
1868 {
1869 using ValueType = GetAPIType<ConnectivityT>;
1870 using AccessorType = vtkDataArrayAccessor<ConnectivityT>;
1871 AccessorType accessor(conn);
1872 accessor.InsertNext(static_cast<ValueType>(id));
1873 }
1874};
1875
1876VTK_ABI_NAMESPACE_END
1877} // end namespace vtkCellArray_detail
1878
1879VTK_ABI_NAMESPACE_BEGIN
1880//----------------------------------------------------------------------------
1882{
1883 this->TraversalCellId = 0;
1884}
1885
1886//----------------------------------------------------------------------------
1887inline int vtkCellArray::GetNextCell(vtkIdType& npts, vtkIdType const*& pts) VTK_SIZEHINT(pts, npts)
1888{
1889 if (this->TraversalCellId < this->GetNumberOfCells())
1890 {
1891 this->GetCellAtId(this->TraversalCellId, npts, pts);
1892 ++this->TraversalCellId;
1893 return 1;
1894 }
1895
1896 npts = 0;
1897 pts = nullptr;
1898 return 0;
1899}
1900
1901//----------------------------------------------------------------------------
1903{
1905 {
1906 this->GetCellAtId(this->TraversalCellId, pts);
1907 ++this->TraversalCellId;
1908 return 1;
1909 }
1910
1911 pts->Reset();
1912 return 0;
1913}
1914//----------------------------------------------------------------------------
1916{
1917 vtkIdType cellSize;
1918 this->Dispatch(vtkCellArray_detail::GetCellSizeImpl{}, cellId, cellSize);
1919 return cellSize;
1920}
1921
1922//----------------------------------------------------------------------------
1923void vtkCellArray::GetCellAtId(vtkIdType cellId, vtkIdType& cellSize, vtkIdType const*& cellPoints,
1924 vtkIdList* ptIds) VTK_SIZEHINT(cellPoints, cellSize)
1925{
1926 this->Dispatch(vtkCellArray_detail::GetCellAtIdImpl{}, cellId, cellSize, cellPoints, ptIds);
1927}
1928
1929//----------------------------------------------------------------------------
1931{
1932 this->Dispatch(vtkCellArray_detail::GetCellAtIdImpl{}, cellId, pts);
1933}
1934
1935//----------------------------------------------------------------------------
1936void vtkCellArray::GetCellAtId(vtkIdType cellId, vtkIdType& cellSize, vtkIdType* cellPoints)
1937{
1938 this->Dispatch(vtkCellArray_detail::GetCellAtIdImpl{}, cellId, cellSize, cellPoints);
1939}
1940
1941//----------------------------------------------------------------------------
1943{
1944 vtkIdType pointId;
1945 this->Dispatch(vtkCellArray_detail::CellPointAtIdImpl{}, cellId, cellPointIndex, pointId);
1946 return pointId;
1947}
1948
1949//----------------------------------------------------------------------------
1951{
1952 vtkIdType cellId;
1953 this->Dispatch(vtkCellArray_detail::InsertNextCellImpl{}, npts, pts, cellId);
1954 return cellId;
1955}
1956
1957//----------------------------------------------------------------------------
1959{
1960 vtkIdType cellId;
1962 return cellId;
1963}
1964
1965//----------------------------------------------------------------------------
1970
1971//----------------------------------------------------------------------------
1973{
1975}
1976
1977//----------------------------------------------------------------------------
1979{
1980 vtkIdType cellId;
1981 this->Dispatch(
1983 return cellId;
1984}
1985
1986//----------------------------------------------------------------------------
1988{
1989 vtkIdList* pts = cell->GetPointIds();
1990 vtkIdType cellId;
1991 this->Dispatch(
1993 return cellId;
1994}
1995
1996//----------------------------------------------------------------------------
1998{
2000}
2001
2002VTK_ABI_NAMESPACE_END
2003#endif // vtkCellArray.h
Array-Of-Structs implementation of vtkGenericDataArray.
virtual void Initialize()=0
Free any memory and reset to an empty state.
virtual void DeepCopy(vtkAbstractCellArray *ca)=0
Perform a deep copy (no reference counting) of the given cell array.
virtual bool IsStorageShareable() const =0
virtual int GetMaxCellSize()=0
Returns the size of the largest cell.
virtual vtkIdType IsHomogeneous() const =0
Check if all cells have the same number of vertices.
virtual vtkIdType GetNumberOfCells() const =0
Get the number of cells in the array.
virtual vtkIdType GetOffset(vtkIdType cellId)=0
Get the offset (into the connectivity) for a specified cell id.
void GetCellAtId(vtkIdType cellId, vtkIdType &cellSize, vtkIdType const *&cellPoints)
Return the point ids for the cell at cellId.
virtual vtkIdType GetNumberOfConnectivityIds() const =0
Get the size of the connectivity array that stores the point ids.
virtual void ShallowCopy(vtkAbstractCellArray *ca)=0
Shallow copy ca into this cell array.
virtual vtkIdType GetNumberOfOffsets() const =0
Get the number of elements in the offsets array.
Encapsulate traversal logic for vtkCellArray.
vtkIdType GetCellSize(vtkIdType cellId) const override
Return the size of the cell at cellId.
vtkDataArray * GetOffsetsArray() const
Return the array used to store cell offsets.
int GetNextCell(vtkIdType &npts, vtkIdType const *&pts)
vtkIdTypeArray * GetData()
Return the underlying data as a data array.
void UseFixedSizeDefaultStorage(vtkIdType cellSize)
Initialize internal data structures to use 32- or 64-bit storage.
vtkIdType GetNumberOfConnectivityEntries()
Return the size of the array that would be returned from ExportLegacyFormat().
void UseDefaultStorage()
Initialize internal data structures to use 32- or 64-bit storage.
bool AllocateCopy(vtkCellArray *other)
Pre-allocate memory in internal data structures to match the used size of the input vtkCellArray.
void AppendLegacyFormat(vtkIdTypeArray *data, vtkIdType ptOffset=0)
Append an array of data with the legacy vtkCellArray layout, e.g.:
bool IsValid()
Check that internal storage is consistent and in a valid state.
friend class vtkCellArrayIterator
bool CanConvertToFixedSize64BitStorage() const
Check if the existing data can safely be converted to use 32- or 64- bit storage.
void ReplaceCell(vtkIdType loc, int npts, const vtkIdType pts[])
Replace the point ids of the cell at the legacy location with a different list of point ids.
vtkIdType GetTraversalCellId()
Get/Set the current cellId for traversal.
void Visit(Functor &&functor, Args &&... args)
vtkTypeInt32Array ArrayType32
vtkIdType GetNumberOfCells() const override
Get the number of cells in the array.
void Reset()
Reuse list.
vtkAOSDataArrayTemplate< vtkTypeInt32 > AOSArray32
bool CanConvertToFixedSizeDefaultStorage() const
Check if the existing data can safely be converted to use 32- or 64- bit storage.
vtkIdType GetNumberOfConnectivityIds() const override
Get the size of the connectivity array that stores the point ids.
vtkIdType GetTraversalLocation()
Get/Set the current traversal legacy location.
vtkAOSDataArrayTemplate< vtkTypeInt64 > AOSArray64
bool CanConvertToStorageType(StorageTypes type) const
Check if the existing data can safely be converted to use 32- or 64- bit storage.
void SetData(AffineArray32 *, AOSArray32 *connectivity)
Set the internal data arrays to the supplied offsets and connectivity arrays.
bool ConvertToFixedSize64BitStorage()
Convert internal data structures to use 32- or 64-bit storage.
vtkAffineArray< vtkTypeInt64 > AffineArray64
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.
void ReverseCell(vtkIdType loc)
Special method inverts ordering of cell at the specified legacy location.
vtkTypeInt64Array ArrayType64
vtkIdType TraversalCellId
void UseFixedSize64BitStorage(vtkIdType cellSize)
Initialize internal data structures to use 32- or 64-bit storage.
void InitTraversal()
bool IsStorageGeneric() const
void Use64BitStorage()
Initialize internal data structures to use 32- or 64-bit storage.
vtkSmartPointer< vtkDataArray > Offsets
vtkDataArray * GetConnectivityArray() const
Return the array used to store the point ids that define the cells' connectivity.
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.
static bool DefaultStorageIs64Bit
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 GetCellAtId(vtkIdType cellId, vtkIdType &cellSize, vtkIdType const *&cellPoints, vtkIdList *ptIds) override
Return the point ids for the cell at cellId.
void SetData(AffineArray64 *, AOSArray64 *connectivity)
Set the internal data arrays to the supplied offsets and connectivity arrays.
void InsertCellPoint(vtkIdType id)
Used in conjunction with InsertNextCell(npts) to add another point to the list of cells.
unsigned long GetActualMemorySize() const
Return the memory in kibibytes (1024 bytes) consumed by this cell array.
vtkCellArrayIterator * NewIterator()
NewIterator returns a new instance of vtkCellArrayIterator that is initialized to point at the first ...
void SetTraversalLocation(vtkIdType loc)
Get/Set the current traversal legacy location.
void ReplaceCellAtId(vtkIdType cellId, vtkIdList *list)
Replaces the point ids for the specified cell with the supplied list.
void UseFixedSize32BitStorage(vtkIdType cellSize)
Initialize internal data structures to use 32- or 64-bit storage.
void ReverseCellAtId(vtkIdType cellId)
Reverses the order of the point ids for the specified cell.
bool IsStorage32Bit() const
void Squeeze()
Reclaim any extra memory while preserving data.
vtkSmartPointer< vtkDataArray > Connectivity
bool IsStorageFixedSize32Bit() const
bool ConvertTo32BitStorage()
Convert internal data structures to use 32- or 64-bit storage.
vtkTypeBool Allocate(vtkIdType sz, vtkIdType ext=1000)
Allocate memory.
bool ConvertToFixedSizeDefaultStorage()
Convert internal data structures to use 32- or 64-bit storage.
vtkTypeList::Create< vtkTypeInt32Array, vtkTypeInt64Array > StorageArrayList
List of possible array types used for storage.
void SetData(AOSArray32 *, AOSArray32 *connectivity)
Set the internal data arrays to the supplied offsets and connectivity arrays.
bool ConvertToStorageType(StorageTypes type)
Convert internal data structures to use 32- or 64-bit storage.
void SetOffset(vtkIdType cellId, vtkIdType offset)
Set the offset (into the connectivity) for a specified cell id.
bool ConvertToSmallestStorage()
Convert internal data structures to use 32- or 64-bit storage.
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.
bool ConvertToFixedSize32BitStorage()
Convert internal data structures to use 32- or 64-bit storage.
void ImportLegacyFormat(vtkIdTypeArray *data)
Import an array of data with the legacy vtkCellArray layout, e.g.:
vtkAffineArray< vtkTypeInt32 > AffineArray32
void Use32BitStorage()
Initialize internal data structures to use 32- or 64-bit storage.
vtkIdType GetCellPointAtId(vtkIdType cellId, vtkIdType cellPointIndex) const
Return the point id at cellPointIndex for the cell at cellId.
void SetTraversalCellId(vtkIdType cellId)
Get/Set the current cellId for traversal.
bool IsStorageFixedSize64Bit() const
static bool GetDefaultStorageIs64Bit()
Control the default internal storage size.
StorageTypes StorageType
vtkIdType InsertNextCell(vtkCell *cell)
Insert a cell object.
virtual void SetNumberOfCells(vtkIdType)
Set the number of cells in the array.
vtkNew< vtkIdTypeArray > LegacyData
void Dispatch(Functor &&functor, Args &&... args)
AOSArray32 * GetOffsetsAOSArray32() const
Return the array used to store cell offsets.
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.
vtkIdType GetInsertLocation(int npts)
Computes the current legacy insertion location within the internal array.
void UpdateCellCount(int npts)
Used in conjunction with InsertNextCell(int npts) and InsertCellPoint() to update the number of point...
bool IsStorageFixedSize() const
vtkTypeList::Unique< vtkTypeList::Create< vtkAOSDataArrayTemplate< int >, vtkAOSDataArrayTemplate< long >, vtkAOSDataArrayTemplate< long long > > >::Result InputArrayList
List of possible ArrayTypes that are compatible with internal storage.
void SetCells(vtkIdType ncells, vtkIdTypeArray *cells)
Define multiple cells by providing a connectivity list.
static vtkCellArray * New()
Standard methods for instantiation, type information, and printing.
StorageTypes GetStorageType() const noexcept
void ReplaceCellPointAtId(vtkIdType cellId, vtkIdType cellPointIndex, vtkIdType newPointId)
Replaces the pointId at cellPointIndex of a cell with newPointId.
vtkIdType GetSize()
Get the size of the allocated connectivity array.
bool CanConvertToFixedSize32BitStorage() const
Check if the existing data can safely be converted to use 32- or 64- bit storage.
void Append(vtkCellArray *src, vtkIdType pointOffset=0)
Append cells from src into this.
bool IsStorage64Bit() const
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:129
vtkIdList * GetPointIds()
Return the list of point ids defining the cell.
Definition vtkCell.h:222
static vtkDataArray * FastDownCast(vtkAbstractArray *source)
Perform a fast, safe cast from a vtkAbstractArray to a vtkDataArray.
list of point or cell ids
Definition vtkIdList.h:133
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:159
void Reset()
Reset to an empty state but retain previously allocated memory.
Definition vtkIdList.h:244
vtkIdType * GetPointer(vtkIdType i)
Get a pointer to a particular data index.
Definition vtkIdList.h:225
dynamic, self-adjusting array of vtkIdType
static vtkImplicitArray< vtkAffineImplicitBackend< T >, ArrayType > * FastDownCast(vtkAbstractArray *source)
a simple class to control print indentation
Definition vtkIndent.h:108
Allocate and hold a VTK object.
Definition vtkNew.h:167
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...
Hold a reference to a vtkObjectBase instance.
static vtkSmartPointer< T > New()
Create an instance of a VTK object.
typename detail::GetAPITypeImpl< ArrayType, ForceValueTypeForVtkDataArray >::APIType GetAPIType
VTK_ITER_INLINE auto DataArrayValueRange(const ArrayTypePtr &array, ValueIdType start=-1, ValueIdType end=-1) -> typename detail::SelectValueRange< ArrayTypePtr, TupleSize, ForceValueTypeForVtkDataArray >::type
Generate an stl and for-range compatible range of flat AOS iterators from a vtkDataArray.
vtkIdType GetNumberOfCells(OffsetsT *offsets)
static vtkIdType GetCellSize(OffsetsT *offsets, vtkIdType cellId)
static vtkIdType GetEndOffset(OffsetsT *offsets, vtkIdType cellId)
static decltype(vtk::DataArrayValueRange< 1, vtkIdType >(std::declval< ConnectivityT >())) GetCellRange(OffsetsT *offsets, ConnectivityT *conn, vtkIdType cellId)
static vtkIdType GetBeginOffset(OffsetsT *offsets, vtkIdType cellId)
static decltype(vtk::DataArrayValueRange< 1, vtkIdType >(std::declval< ArrayT >())) GetRange(ArrayT *array)
vtk::GetAPIType< ArrayT, vtkIdType > GetAPIType
vtkIdType GetNumberOfCells() const
vtkIdType GetEndOffset(vtkIdType cellId) const
vtkSmartPointer< ArrayType > Offsets
vtkSmartPointer< ArrayType > Connectivity
static constexpr bool ValueTypeIsSameAsIdType
decltype(vtk::DataArrayValueRange< 1 >(std::declval< ArrayType >())) CellRangeType
CellRangeType GetCellRange(vtkIdType cellId)
vtkIdType GetBeginOffset(vtkIdType cellId) const
vtkIdType GetCellSize(vtkIdType cellId) const
typename ArrayType::ValueType ValueType
ArrayType * GetConnectivity()
void operator()(OffsetsT *offsets, ConnectivityT *conn, vtkIdType cellId, vtkIdType cellPointIndex, vtkIdType &pointId)
void operator()(OffsetsT *offsets, ConnectivityT *conn, const vtkIdType cellId, vtkIdList *ids)
void operator()(OffsetsT *offsets, ConnectivityT *conn, const vtkIdType cellId, vtkIdType &cellSize, vtkIdType *cellPoints)
std::enable_if_t< CanShareConnPtr< ConnectivityT >::value, void > operator()(OffsetsT *offsets, ConnectivityT *conn, const vtkIdType cellId, vtkIdType &cellSize, vtkIdType const *&cellPoints, vtkIdList *temp)
std::enable_if_t<!CanShareConnPtr< ConnectivityT >::value, void > operator()(OffsetsT *offsets, ConnectivityT *conn, const vtkIdType cellId, vtkIdType &cellSize, vtkIdType const *&cellPoints, vtkIdList *temp)
void operator()(OffsetsT *offsets, ConnectivityT *conn, vtkIdType cellId, vtkIdType &cellSize)
void operator()(OffsetsT *offsets, ConnectivityT *conn, vtkIdType id)
void operator()(OffsetsT *offsets, ConnectivityT *conn, const vtkIdType npts, vtkIdType &cellId)
void operator()(OffsetsT *offsets, ConnectivityT *conn, const vtkIdType npts, const vtkIdType pts[], vtkIdType &cellId)
void operator()(OffsetsT *offsets, ConnectivityT *conn)
void operator()(OffsetsT *offsets, ConnectivityT *conn, const vtkIdType npts)
Efficient templated access to vtkDataArray.
Remove all duplicate types from TypeList TList, storing the new list in Result.
int vtkTypeBool
Definition vtkABI.h:64
vtkImplicitArray< vtkAffineImplicitBackend< T >, vtkArrayTypes::VTK_AFFINE_ARRAY > vtkAffineArray
A utility alias for wrapping affine functions in implicit arrays.
#define vtkDataArray
STL-compatible iterable ranges that provide access to vtkDataArray elements.
#define VTK_DEPRECATED_IN_9_6_0(reason)
int vtkIdType
Definition vtkType.h:354
#define VTK_SIZEHINT(...)
#define VTK_EXPECTS(x)
#define VTK_MARSHALMANUAL
#define VTK_NEWINSTANCE