VTK  9.5.20251217
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
259VTK_ABI_NAMESPACE_BEGIN
261class vtkIdTypeArray;
262
263class VTKCOMMONDATAMODEL_EXPORT VTK_MARSHALMANUAL vtkCellArray : public vtkAbstractCellArray
264{
265public:
266 using ArrayType32 VTK_DEPRECATED_IN_9_6_0("Use AOSArray32 instead.") = vtkTypeInt32Array;
267 using ArrayType64 VTK_DEPRECATED_IN_9_6_0("Use AOSArray64 instead.") = vtkTypeInt64Array;
272
274
278 static vtkCellArray* New();
280 void PrintSelf(ostream& os, vtkIndent indent) override;
281 void PrintDebug(ostream& os);
283
285
294 "Use vtkArrayDispatch::OffsetsArrays/ConnectivityArrays instead.") =
295 vtkTypeList::Create<vtkTypeInt32Array, vtkTypeInt64Array>;
297
299
311
319 */
320 VTK_DEPRECATED_IN_9_6_0("Use AllocateEstimate or AllocateExact instead.")
321 vtkTypeBool Allocate(vtkIdType sz, vtkIdType vtkNotUsed(ext) = 1000)
322 {
323 return this->AllocateExact(sz, sz) ? 1 : 0;
324 }
325
333 * @sa Squeeze AllocateExact AllocateCopy
334 */
335 bool AllocateEstimate(vtkIdType numCells, vtkIdType maxCellSize)
336 {
337 return this->AllocateExact(numCells, numCells * maxCellSize);
338 }
339
349 bool AllocateExact(vtkIdType numCells, vtkIdType connectivitySize);
350
358 * @sa Squeeze AllocateEstimate AllocateExact
359 */
360 bool AllocateCopy(vtkCellArray* other)
361 {
362 return this->AllocateExact(other->GetNumberOfCells(), other->GetNumberOfConnectivityIds());
363 }
364
374 bool ResizeExact(vtkIdType numCells, vtkIdType connectivitySize);
375
379 void Initialize() override;
380
384 void Reset();
385
391 void Squeeze();
392
403 bool IsValid();
404
408 vtkIdType GetNumberOfCells() const override { return this->Offsets->GetNumberOfValues() - 1; }
409
414 vtkIdType GetNumberOfOffsets() const override { return this->Offsets->GetNumberOfValues(); }
415
417 * Get the offset (into the connectivity) for a specified cell id.
418 */
419 vtkIdType GetOffset(vtkIdType cellId) override
420 {
421 return static_cast<vtkIdType>(this->Offsets->GetComponent(cellId, 0));
422 }
423
425 * Set the offset (into the connectivity) for a specified cell id.
426 */
427 void SetOffset(vtkIdType cellId, vtkIdType offset)
428 {
429 this->Offsets->SetComponent(cellId, 0, static_cast<double>(offset));
430 }
431
433 * Get the size of the connectivity array that stores the point ids.
434 */
436 {
437 return this->Connectivity->GetNumberOfValues();
438 }
439
446
448
461 void SetData(AOSArray32*, AOSArray32* connectivity);
462 void SetData(AOSArray64*, AOSArray64* connectivity);
463 void SetData(AffineArray32*, AOSArray32* connectivity);
464 void SetData(AffineArray64*, AOSArray64* connectivity);
466
481 bool SetData(vtkDataArray* offsets, vtkDataArray* connectivity);
482
508 Generic,
509 };
510
514 StorageTypes GetStorageType() const noexcept { return this->StorageType; }
515
519 bool IsStorage64Bit() const { return this->StorageType == StorageTypes::Int64; }
520
524 bool IsStorage32Bit() const { return this->StorageType == StorageTypes::Int32; }
525
530
535
537 * @return True if the internal storage is using FixedSize arrays.
538 */
539 bool IsStorageFixedSize() const
540 {
541 return this->IsStorageFixedSize32Bit() || this->IsStorageFixedSize64Bit();
542 }
543
547 bool IsStorageGeneric() const { return this->StorageType == StorageTypes::Generic; }
548
553 * a pointer to vtkIdType can be used.
554 */
555 bool IsStorageShareable() const override
556 {
557 switch (this->StorageType)
558 {
561 return std::is_same_v<vtkTypeInt32, vtkIdType>;
564 return std::is_same_v<vtkTypeInt64, vtkIdType>;
566 default:
567 return false;
568 }
569 }
570
583 void UseFixedSize64BitStorage(vtkIdType cellSize);
586
602 {
603 switch (type)
604 {
606 return this->CanConvertTo32BitStorage();
608 return this->CanConvertTo64BitStorage();
614 default:
615 return true;
616 }
617 }
619
644 {
645 switch (type)
646 {
648 return this->ConvertTo32BitStorage();
650 return this->ConvertTo64BitStorage();
652 return this->ConvertToFixedSize32BitStorage();
654 return this->ConvertToFixedSize64BitStorage();
656 default:
657 return true;
658 }
659 }
661
667 vtkDataArray* GetOffsetsArray() const { return this->Offsets; }
668 VTK_DEPRECATED_IN_9_6_0("Use GetOffsetsAOSArray32() instead.")
669 vtkTypeInt32Array* GetOffsetsArray32() const
670 {
671 return vtkTypeInt32Array::FastDownCast(this->Offsets);
672 }
674 VTK_DEPRECATED_IN_9_6_0("Use GetOffsetsAOSArray64() instead.")
675 vtkTypeInt64Array* GetOffsetsArray64() const
676 {
677 return vtkTypeInt64Array::FastDownCast(this->Offsets);
679 AOSArray64* GetOffsetsAOSArray64() const { return AOSArray64::FastDownCast(this->Offsets); }
680 AffineArray32* GetOffsetsAffineArray32() const
681 {
683 }
684 AffineArray64* GetOffsetsAffineArray64() const
685 {
686 return AffineArray64::FastDownCast(this->Offsets);
687 }
689
697 VTK_DEPRECATED_IN_9_6_0("Use GetConnectivityAOSArray32() instead.")
698 vtkTypeInt32Array* GetConnectivityArray32() const
699 {
700 return vtkTypeInt32Array::FastDownCast(this->Connectivity);
701 }
702 AOSArray32* GetConnectivityAOSArray32() const
703 {
704 return AOSArray32::FastDownCast(this->Connectivity);
706 VTK_DEPRECATED_IN_9_6_0("Use GetConnectivityAOSArray64() instead.")
707 vtkTypeInt64Array* GetConnectivityArray64() const
708 {
709 return vtkTypeInt64Array::FastDownCast(this->Connectivity);
710 }
711 AOSArray64* GetConnectivityAOSArray64() const
712 {
713 return AOSArray64::FastDownCast(this->Connectivity);
714 }
716
725 vtkIdType IsHomogeneous() const override;
726
736 void InitTraversal();
737
752 int GetNextCell(vtkIdType& npts, vtkIdType const*& pts) VTK_SIZEHINT(pts, npts);
753
764 int GetNextCell(vtkIdList* pts);
765
776 void GetCellAtId(vtkIdType cellId, vtkIdType& cellSize, vtkIdType const*& cellPoints,
777 vtkIdList* ptIds) VTK_SIZEHINT(cellPoints, cellSize)
778 VTK_EXPECTS(0 <= cellId && cellId < GetNumberOfCells()) override;
779
785 void GetCellAtId(vtkIdType cellId, vtkIdList* pts)
786 VTK_EXPECTS(0 <= cellId && cellId < GetNumberOfCells()) override;
787
795 void GetCellAtId(vtkIdType cellId, vtkIdType& cellSize, vtkIdType* cellPoints) VTK_SIZEHINT(
796 cellPoints, cellSize) VTK_EXPECTS(0 <= cellId && cellId < GetNumberOfCells()) override;
797
801 vtkIdType GetCellPointAtId(vtkIdType cellId, vtkIdType cellPointIndex) const
802 VTK_EXPECTS(0 <= cellId && cellId < GetNumberOfCells())
803 VTK_EXPECTS(0 <= cellPointIndex && cellPointIndex < this->GetCellSize(cellId));
804
808 vtkIdType GetCellSize(vtkIdType cellId) const override;
809
814
819 vtkIdType InsertNextCell(vtkIdType npts, const vtkIdType* pts) VTK_SIZEHINT(pts, npts);
820
826
832 * `InsertNextCell(int)` overload instead.
833 */
834 vtkIdType InsertNextCell(const std::initializer_list<vtkIdType>& cell)
835 {
836 return this->InsertNextCell(static_cast<vtkIdType>(cell.size()), cell.begin());
837 }
838
845 vtkIdType InsertNextCell(int npts);
846
851 void InsertCellPoint(vtkIdType id);
852
857 void UpdateCellCount(int npts);
858
867 void SetTraversalCellId(vtkIdType cellId);
869
873 void ReverseCellAtId(vtkIdType cellId) VTK_EXPECTS(0 <= cellId && cellId < GetNumberOfCells());
874
883 void ReplaceCellAtId(vtkIdType cellId, vtkIdList* list);
884 void ReplaceCellAtId(vtkIdType cellId, vtkIdType cellSize, const vtkIdType* cellPoints)
885 VTK_EXPECTS(0 <= cellId && cellId < GetNumberOfCells()) VTK_SIZEHINT(cellPoints, cellSize);
887
895 void ReplaceCellPointAtId(vtkIdType cellId, vtkIdType cellPointIndex, vtkIdType newPointId);
896
902 * results.
903 */
904 void ReplaceCellAtId(vtkIdType cellId, const std::initializer_list<vtkIdType>& cell)
905 {
906 this->ReplaceCellAtId(cellId, static_cast<vtkIdType>(cell.size()), cell.begin());
907 }
908
913 int GetMaxCellSize() override;
914
918 void DeepCopy(vtkAbstractCellArray* ca) override;
919
923 void ShallowCopy(vtkAbstractCellArray* ca) override;
924
928 void Append(vtkCellArray* src, vtkIdType pointOffset = 0);
929
941
954 void ImportLegacyFormat(const vtkIdType* data, vtkIdType len) VTK_SIZEHINT(data, len);
956
968 void AppendLegacyFormat(vtkIdTypeArray* data, vtkIdType ptOffset = 0);
969 void AppendLegacyFormat(const vtkIdType* data, vtkIdType len, vtkIdType ptOffset = 0)
970 VTK_SIZEHINT(data, len);
972
981 unsigned long GetActualMemorySize() const;
982
983 // The following code is used to support
984
985 // The wrappers get understandably confused by some of the template code below
986#ifndef __VTK_WRAP__
987 /*
988 * Utilities class that every dispatch functor used with Dispatch() must inherit
989 * or optionally use its static methods.
990 */
991 struct DispatchUtilities
993 template <class ArrayT>
996 template <class OffsetsT>
997 vtkIdType GetNumberOfCells(OffsetsT* offsets)
998 {
999 return offsets->GetNumberOfValues() - 1;
1000 }
1002 template <class ArrayT>
1003 static decltype(vtk::DataArrayValueRange<1, vtkIdType>(std::declval<ArrayT>())) GetRange(
1004 ArrayT* array)
1005 {
1007 }
1009 template <class OffsetsT>
1010 static vtkIdType GetBeginOffset(OffsetsT* offsets, vtkIdType cellId)
1011 {
1012 return static_cast<vtkIdType>(GetRange(offsets)[cellId]);
1013 }
1015 template <class OffsetsT>
1016 static vtkIdType GetEndOffset(OffsetsT* offsets, vtkIdType cellId)
1017 {
1018 return static_cast<vtkIdType>(GetRange(offsets)[cellId + 1]);
1019 }
1021 template <class OffsetsT>
1022 static vtkIdType GetCellSize(OffsetsT* offsets, vtkIdType cellId)
1023 {
1024 auto offsetsRange = GetRange(offsets);
1025 return static_cast<vtkIdType>(offsetsRange[cellId + 1] - offsetsRange[cellId]);
1026 }
1027
1028 template <class OffsetsT, class ConnectivityT>
1029 static decltype(vtk::DataArrayValueRange<1, vtkIdType>(std::declval<ConnectivityT>()))
1030 GetCellRange(OffsetsT* offsets, ConnectivityT* conn, vtkIdType cellId)
1031 {
1032 auto offsetsRange = GetRange(offsets);
1034 conn, offsetsRange[cellId], offsetsRange[cellId + 1]);
1035 }
1036 };
1037
1039
1102 template <typename Functor, typename... Args>
1103 void Dispatch(Functor&& functor, Args&&... args)
1104 {
1105 switch (this->StorageType)
1106 {
1108 functor(static_cast<AOSArray32*>(this->Offsets.Get()),
1109 static_cast<AOSArray32*>(this->Connectivity.Get()), std::forward<Args>(args)...);
1110 break;
1112 functor(static_cast<AOSArray64*>(this->Offsets.Get()),
1113 static_cast<AOSArray64*>(this->Connectivity.Get()), std::forward<Args>(args)...);
1114 break;
1116 functor(static_cast<AffineArray32*>(this->Offsets.Get()),
1117 static_cast<AOSArray32*>(this->Connectivity.Get()), std::forward<Args>(args)...);
1118 break;
1120 functor(static_cast<AffineArray64*>(this->Offsets.Get()),
1121 static_cast<AOSArray64*>(this->Connectivity.Get()), std::forward<Args>(args)...);
1122 break;
1124 default:
1125 functor(this->Offsets.Get(), this->Connectivity.Get(), std::forward<Args>(args)...);
1126 break;
1127 }
1129 template <typename Functor, typename... Args>
1130 void Dispatch(Functor&& functor, Args&&... args) const
1131 {
1132 switch (this->StorageType)
1133 {
1134 case StorageTypes::Int32:
1135 functor(static_cast<AOSArray32*>(this->Offsets.Get()),
1136 static_cast<AOSArray32*>(this->Connectivity.Get()), std::forward<Args>(args)...);
1137 break;
1138 case StorageTypes::Int64:
1139 functor(static_cast<AOSArray64*>(this->Offsets.Get()),
1140 static_cast<AOSArray64*>(this->Connectivity.Get()), std::forward<Args>(args)...);
1141 break;
1142 case StorageTypes::FixedSizeInt32:
1143 functor(static_cast<AffineArray32*>(this->Offsets.Get()),
1144 static_cast<AOSArray32*>(this->Connectivity.Get()), std::forward<Args>(args)...);
1145 break;
1146 case StorageTypes::FixedSizeInt64:
1147 functor(static_cast<AffineArray64*>(this->Offsets.Get()),
1148 static_cast<AOSArray64*>(this->Connectivity.Get()), std::forward<Args>(args)...);
1149 break;
1150 case StorageTypes::Generic:
1151 default:
1152 functor(this->Offsets.Get(), this->Connectivity.Get(), std::forward<Args>(args)...);
1153 break;
1154 }
1155 }
1157
1158 // Holds connectivity and offset arrays of the given ArrayType.
1159 // VTK_DEPRECATED_IN_9_6_0("Use DispatchUtilities")
1160 template <typename ArrayT>
1163 using ArrayType = ArrayT;
1164 using ValueType = typename ArrayType::ValueType;
1165 using CellRangeType = decltype(vtk::DataArrayValueRange<1>(std::declval<ArrayType>()));
1166
1167 // We can't just use is_same here, since binary compatible representations
1168 // (e.g. int and long) are distinct types. Instead, ensure that ValueType
1169 // is a signed integer the same size as vtkIdType.
1170 // If this value is true, ValueType pointers may be safely converted to
1171 // vtkIdType pointers via reinterpret cast.
1172 static constexpr bool ValueTypeIsSameAsIdType = std::is_integral<ValueType>::value &&
1173 std::is_signed<ValueType>::value && (sizeof(ValueType) == sizeof(vtkIdType));
1175 ArrayType* GetOffsets() { return this->Offsets; }
1176 const ArrayType* GetOffsets() const { return this->Offsets; }
1178 ArrayType* GetConnectivity() { return this->Connectivity; }
1179 const ArrayType* GetConnectivity() const { return this->Connectivity; }
1180
1181 vtkIdType GetNumberOfCells() const;
1182
1183 vtkIdType GetBeginOffset(vtkIdType cellId) const;
1184
1185 vtkIdType GetEndOffset(vtkIdType cellId) const;
1186
1187 vtkIdType GetCellSize(vtkIdType cellId) const;
1188
1190
1191 friend class vtkCellArray;
1193 protected:
1194 VisitState()
1195 {
1198 this->Offsets->InsertNextValue(0);
1200 {
1201 this->IsInMemkind = true;
1204 ~VisitState() = default;
1205 void* operator new(size_t nSize)
1206 {
1207 void* r;
1208#ifdef VTK_USE_MEMKIND
1210#else
1211 r = malloc(nSize);
1212#endif
1213 return r;
1214 }
1215 void operator delete(void* p)
1216 {
1217#ifdef VTK_USE_MEMKIND
1218 VisitState* a = static_cast<VisitState*>(p);
1219 if (a->IsInMemkind)
1220 {
1222 }
1223 else
1224 {
1225 free(p);
1226 }
1227#else
1228 free(p);
1229#endif
1234
1235 private:
1236 VisitState(const VisitState&) = delete;
1237 VisitState& operator=(const VisitState&) = delete;
1238 bool IsInMemkind = false;
1239 };
1240
1241private: // Helpers that allow Visit to return a value:
1242 template <typename Functor, typename... Args>
1243 using GetReturnType = decltype(std::declval<Functor>()(
1244 std::declval<VisitState<AOSArray32>&>(), std::declval<Args>()...));
1245
1246 template <typename Functor, typename... Args>
1247 struct ReturnsVoid : std::is_same<GetReturnType<Functor, Args...>, void>
1248 {
1249 };
1250
1251public:
1321 template <typename Functor, typename... Args,
1322 typename = typename std::enable_if<ReturnsVoid<Functor, Args...>::value>::type>
1323 VTK_DEPRECATED_IN_9_6_0("Use Dispatch instead")
1324 void Visit(Functor&& functor, Args&&... args)
1325 {
1326 switch (this->StorageType)
1327 {
1329 {
1333 functor(state, std::forward<Args>(args)...);
1334 break;
1335 }
1337 {
1341 functor(state, std::forward<Args>(args)...);
1342 break;
1343 }
1347 default:
1348 {
1349 vtkWarningMacro("Use Dispatch");
1350 break;
1351 }
1352 }
1353 }
1354
1355 template <typename Functor, typename... Args,
1356 typename = typename std::enable_if<ReturnsVoid<Functor, Args...>::value>::type>
1357 VTK_DEPRECATED_IN_9_6_0("Use Dispatch instead")
1358 void Visit(Functor&& functor, Args&&... args) const
1359 {
1360 switch (this->StorageType)
1361 {
1363 {
1367 functor(state, std::forward<Args>(args)...);
1368 break;
1369 }
1371 {
1375 functor(state, std::forward<Args>(args)...);
1376 break;
1377 }
1381 default:
1382 {
1383 vtkWarningMacro("Use Dispatch");
1384 break;
1385 }
1386 }
1387 }
1388
1389 template <typename Functor, typename... Args,
1390 typename = typename std::enable_if<!ReturnsVoid<Functor, Args...>::value>::type>
1391 VTK_DEPRECATED_IN_9_6_0("Use Dispatch instead")
1392 GetReturnType<Functor, Args...> Visit(Functor&& functor, Args&&... args)
1393 {
1394 switch (this->StorageType)
1395 {
1397 {
1401 return functor(state, std::forward<Args>(args)...);
1402 }
1404 {
1408 return functor(state, std::forward<Args>(args)...);
1409 }
1413 default:
1414 {
1415 vtkWarningMacro("Use Dispatch");
1416 return GetReturnType<Functor, Args...>();
1417 }
1418 }
1419 }
1420 template <typename Functor, typename... Args,
1421 typename = typename std::enable_if<!ReturnsVoid<Functor, Args...>::value>::type>
1422 VTK_DEPRECATED_IN_9_6_0("Use Dispatch instead")
1423 GetReturnType<Functor, Args...> Visit(Functor&& functor, Args&&... args) const
1424 {
1425 switch (this->StorageType)
1426 {
1428 {
1432 return functor(state, std::forward<Args>(args)...);
1433 }
1435 {
1439 return functor(state, std::forward<Args>(args)...);
1440 }
1444 default:
1445 {
1446 vtkWarningMacro("Use Dispatch");
1447 return GetReturnType<Functor, Args...>();
1448 }
1449 }
1450 }
1451
1452#endif // __VTK_WRAP__
1453
1455
1463 static void SetDefaultStorageIs64Bit(bool val) { vtkCellArray::DefaultStorageIs64Bit = val; }
1465
1466 //=================== Begin Legacy Methods ===================================
1467 // These should be deprecated at some point as they are confusing or very slow
1468
1475 VTK_DEPRECATED_IN_9_6_0("This call has no effect.")
1476 virtual void SetNumberOfCells(vtkIdType);
1477
1489 VTK_DEPRECATED_IN_9_6_0("Use AllocateEstimate directly instead.")
1490 vtkIdType EstimateSize(vtkIdType numCells, int maxPtsPerCell);
1491
1500 VTK_DEPRECATED_IN_9_6_0("Method incompatible with current internal storage.")
1502
1509 VTK_DEPRECATED_IN_9_6_0("Method incompatible with current internal storage.")
1511
1521 VTK_DEPRECATED_IN_9_6_0("Use GetCellAtId.")
1522 void GetCell(vtkIdType loc, vtkIdType& npts, const vtkIdType*& pts)
1523 VTK_EXPECTS(0 <= loc && loc < GetNumberOfConnectivityEntries()) VTK_SIZEHINT(pts, npts);
1524
1531 VTK_DEPRECATED_IN_9_6_0("Use GetCellAtId.")
1532 void GetCell(vtkIdType loc, vtkIdList* pts)
1533 VTK_EXPECTS(0 <= loc && loc < GetNumberOfConnectivityEntries());
1534
1541 VTK_DEPRECATED_IN_9_6_0("Use GetNumberOfCells.")
1542 vtkIdType GetInsertLocation(int npts);
1543
1551 VTK_DEPRECATED_IN_9_6_0("Use GetTraversalCellId.")
1553 VTK_DEPRECATED_IN_9_6_0("Use GetTraversalCellId.")
1555 VTK_DEPRECATED_IN_9_6_0("Use SetTraversalCellId.")
1558
1566 VTK_DEPRECATED_IN_9_6_0("Use ReverseCellAtId.")
1567 void ReverseCell(vtkIdType loc) VTK_EXPECTS(0 <= loc && loc < GetNumberOfConnectivityEntries());
1568
1580 VTK_DEPRECATED_IN_9_6_0("Use ReplaceCellAtId.")
1581 void ReplaceCell(vtkIdType loc, int npts, const vtkIdType pts[])
1582 VTK_EXPECTS(0 <= loc && loc < GetNumberOfConnectivityEntries()) VTK_SIZEHINT(pts, npts);
1583
1598 VTK_DEPRECATED_IN_9_6_0("Use ImportLegacyFormat or SetData instead.")
1599 void SetCells(vtkIdType ncells, vtkIdTypeArray* cells);
1600
1612 "Use ExportLegacyFormat, or GetOffsetsArray/GetConnectivityArray instead.")
1614
1615 //=================== End Legacy Methods =====================================
1616
1617 friend class vtkCellArrayIterator;
1619protected:
1620 vtkCellArray();
1621 ~vtkCellArray() override;
1627
1629
1630 static bool DefaultStorageIs64Bit;
1631
1632private:
1633 vtkCellArray(const vtkCellArray&) = delete;
1634 void operator=(const vtkCellArray&) = delete;
1635};
1636
1637template <typename ArrayT>
1639{
1640 return this->Offsets->GetNumberOfValues() - 1;
1641}
1642
1643template <typename ArrayT>
1645{
1646 return static_cast<vtkIdType>(this->Offsets->GetValue(cellId));
1647}
1648
1649template <typename ArrayT>
1651{
1652 return static_cast<vtkIdType>(this->Offsets->GetValue(cellId + 1));
1653}
1654
1655template <typename ArrayT>
1657{
1658 return this->GetEndOffset(cellId) - this->GetBeginOffset(cellId);
1659}
1660
1661template <typename ArrayT>
1668VTK_ABI_NAMESPACE_END
1669
1671{
1672VTK_ABI_NAMESPACE_BEGIN
1673
1675{
1676 // Insert full cell
1677 template <class OffsetsT, class ConnectivityT>
1678 void operator()(OffsetsT* offsets, ConnectivityT* conn, const vtkIdType npts,
1679 const vtkIdType pts[], vtkIdType& cellId)
1680 {
1681 using ValueType = GetAPIType<OffsetsT>;
1682 using OffsetsAccessorType = vtkDataArrayAccessor<OffsetsT>;
1683 using ConnectivityAccessorType = vtkDataArrayAccessor<ConnectivityT>;
1684 OffsetsAccessorType offsetsAccesor(offsets);
1685 ConnectivityAccessorType connAccesor(conn);
1686
1687 cellId = offsets->GetNumberOfValues() - 1;
1688
1689 offsetsAccesor.InsertNext(static_cast<ValueType>(conn->GetNumberOfValues() + npts));
1690
1691 for (vtkIdType i = 0; i < npts; ++i)
1692 {
1693 connAccesor.InsertNext(static_cast<ValueType>(pts[i]));
1694 }
1695 }
1696
1697 // Just update offset table (for incremental API)
1698 template <class OffsetsT, class ConnectivityT>
1699 void operator()(OffsetsT* offsets, ConnectivityT* conn, const vtkIdType npts, vtkIdType& cellId)
1700 {
1701 using ValueType = GetAPIType<OffsetsT>;
1702 using AccessorType = vtkDataArrayAccessor<OffsetsT>;
1703 AccessorType offsetsAccesor(offsets);
1704
1705 cellId = offsets->GetNumberOfValues() - 1;
1706
1707 offsetsAccesor.InsertNext(static_cast<ValueType>(conn->GetNumberOfValues() + npts));
1708 }
1709};
1710
1711// for incremental API:
1713{
1714 template <class OffsetsT, class ConnectivityT>
1715 void operator()(OffsetsT* offsets, ConnectivityT* vtkNotUsed(conn), const vtkIdType npts)
1716 {
1717 using ValueType = GetAPIType<OffsetsT>;
1718
1719 auto offsetsRange = GetRange(offsets);
1720 const ValueType cellBegin = offsetsRange[offsets->GetMaxId() - 1];
1721 offsetsRange[offsets->GetMaxId()] = static_cast<ValueType>(cellBegin + npts);
1722 }
1723};
1724
1726{
1727 template <class OffsetsT, class ConnectivityT>
1729 OffsetsT* offsets, ConnectivityT* vtkNotUsed(conn), vtkIdType cellId, vtkIdType& cellSize)
1730 {
1731 cellSize = GetCellSize(offsets, cellId);
1732 }
1733};
1734
1736{
1737 template <class OffsetsT, class ConnectivityT>
1738 void operator()(OffsetsT* offsets, ConnectivityT* conn, const vtkIdType cellId, vtkIdList* ids)
1739 {
1740 auto offsetsRange = GetRange(offsets);
1741 const auto& beginOffset = offsetsRange[cellId];
1742 const auto& endOffset = offsetsRange[cellId + 1];
1743 const vtkIdType cellSize = static_cast<vtkIdType>(endOffset - beginOffset);
1744 const auto cellConnectivity = GetRange(conn).begin() + beginOffset;
1745
1746 // ValueType differs from vtkIdType, so we have to copy into a temporary buffer:
1747 ids->SetNumberOfIds(cellSize);
1748 vtkIdType* idPtr = ids->GetPointer(0);
1749 for (vtkIdType i = 0; i < cellSize; ++i)
1750 {
1751 idPtr[i] = static_cast<vtkIdType>(cellConnectivity[i]);
1752 }
1753 }
1754
1755 template <class OffsetsT, class ConnectivityT>
1756 void operator()(OffsetsT* offsets, ConnectivityT* conn, const vtkIdType cellId,
1757 vtkIdType& cellSize, vtkIdType* cellPoints)
1758 {
1759 auto offsetsRange = GetRange(offsets);
1760 const auto& beginOffset = offsetsRange[cellId];
1761 const auto& endOffset = offsetsRange[cellId + 1];
1762 cellSize = static_cast<vtkIdType>(endOffset - beginOffset);
1763 const auto cellConnectivity = GetRange(conn).begin() + beginOffset;
1764
1765 // ValueType differs from vtkIdType, so we have to copy into a temporary buffer:
1766 for (vtkIdType i = 0; i < cellSize; ++i)
1767 {
1768 cellPoints[i] = static_cast<vtkIdType>(cellConnectivity[i]);
1769 }
1770 }
1771
1772 // SFINAE helper to check if a Functors's connectivity array's memory can be used as a vtkIdType*.
1773 template <typename ConnectivityT>
1775 {
1776 static constexpr bool value =
1777 std::is_base_of_v<vtkAOSDataArrayTemplate<vtkIdType>, ConnectivityT>;
1778 };
1779
1780 template <class OffsetsT, class ConnectivityT>
1781 typename std::enable_if_t<CanShareConnPtr<ConnectivityT>::value, void> operator()(
1782 OffsetsT* offsets, ConnectivityT* conn, const vtkIdType cellId, vtkIdType& cellSize,
1783 vtkIdType const*& cellPoints, vtkIdList* vtkNotUsed(temp))
1784 {
1785 auto offsetsRange = GetRange(offsets);
1786 const auto& beginOffset = offsetsRange[cellId];
1787 const auto& endOffset = offsetsRange[cellId + 1];
1788 cellSize = static_cast<vtkIdType>(endOffset - beginOffset);
1789 // This is safe, see CanShareConnPtr helper above.
1790 cellPoints = conn->GetPointer(beginOffset);
1791 }
1792
1793 template <class OffsetsT, class ConnectivityT>
1794 typename std::enable_if_t<!CanShareConnPtr<ConnectivityT>::value, void> operator()(
1795 OffsetsT* offsets, ConnectivityT* conn, const vtkIdType cellId, vtkIdType& cellSize,
1796 vtkIdType const*& cellPoints, vtkIdList* temp)
1797 {
1798 auto offsetsRange = GetRange(offsets);
1799 const auto& beginOffset = offsetsRange[cellId];
1800 const auto& endOffset = offsetsRange[cellId + 1];
1801 cellSize = static_cast<vtkIdType>(endOffset - beginOffset);
1802 const auto cellConnectivity = GetRange(conn).begin() + beginOffset;
1803
1804 temp->SetNumberOfIds(cellSize);
1805 vtkIdType* tempPtr = temp->GetPointer(0);
1806 for (vtkIdType i = 0; i < cellSize; ++i)
1807 {
1808 tempPtr[i] = static_cast<vtkIdType>(cellConnectivity[i]);
1809 }
1810
1811 cellPoints = tempPtr;
1812 }
1813};
1814
1816{
1817 template <class OffsetsT, class ConnectivityT>
1818 void operator()(OffsetsT* offsets, ConnectivityT* conn, vtkIdType cellId,
1819 vtkIdType cellPointIndex, vtkIdType& pointId)
1820 {
1821 pointId =
1822 static_cast<vtkIdType>(GetRange(conn)[GetBeginOffset(offsets, cellId) + cellPointIndex]);
1823 }
1824};
1825
1827{
1828 template <class OffsetsT, class ConnectivityT>
1829 void operator()(OffsetsT* offsets, ConnectivityT* conn)
1830 {
1831 using ValueType = GetAPIType<OffsetsT>;
1832 using AccessorType = vtkDataArrayAccessor<OffsetsT>;
1833 offsets->Reset();
1834 conn->Reset();
1835 AccessorType accessor(offsets);
1836 ValueType firstOffset = 0;
1837 accessor.InsertNext(firstOffset);
1838 }
1839};
1840
1842{
1843 template <class OffsetsT, class ConnectivityT>
1844 void operator()(OffsetsT* vtkNotUsed(offsets), ConnectivityT* conn, vtkIdType id)
1845 {
1846 using ValueType = GetAPIType<ConnectivityT>;
1847 using AccessorType = vtkDataArrayAccessor<ConnectivityT>;
1848 AccessorType accessor(conn);
1849 accessor.InsertNext(static_cast<ValueType>(id));
1850 }
1851};
1852
1853VTK_ABI_NAMESPACE_END
1854} // end namespace vtkCellArray_detail
1855
1856VTK_ABI_NAMESPACE_BEGIN
1857//----------------------------------------------------------------------------
1859{
1860 this->TraversalCellId = 0;
1861}
1862
1863//----------------------------------------------------------------------------
1864inline int vtkCellArray::GetNextCell(vtkIdType& npts, vtkIdType const*& pts) VTK_SIZEHINT(pts, npts)
1865{
1866 if (this->TraversalCellId < this->GetNumberOfCells())
1867 {
1868 this->GetCellAtId(this->TraversalCellId, npts, pts);
1869 ++this->TraversalCellId;
1870 return 1;
1871 }
1872
1873 npts = 0;
1874 pts = nullptr;
1875 return 0;
1876}
1877
1878//----------------------------------------------------------------------------
1880{
1882 {
1883 this->GetCellAtId(this->TraversalCellId, pts);
1884 ++this->TraversalCellId;
1885 return 1;
1886 }
1887
1888 pts->Reset();
1889 return 0;
1890}
1891//----------------------------------------------------------------------------
1893{
1894 vtkIdType cellSize;
1895 this->Dispatch(vtkCellArray_detail::GetCellSizeImpl{}, cellId, cellSize);
1896 return cellSize;
1897}
1898
1899//----------------------------------------------------------------------------
1900inline void vtkCellArray::GetCellAtId(vtkIdType cellId, vtkIdType& cellSize,
1901 vtkIdType const*& cellPoints, vtkIdList* ptIds) VTK_SIZEHINT(cellPoints, cellSize)
1902{
1903 this->Dispatch(vtkCellArray_detail::GetCellAtIdImpl{}, cellId, cellSize, cellPoints, ptIds);
1904}
1905
1906//----------------------------------------------------------------------------
1908{
1909 this->Dispatch(vtkCellArray_detail::GetCellAtIdImpl{}, cellId, pts);
1910}
1911
1912//----------------------------------------------------------------------------
1913inline void vtkCellArray::GetCellAtId(vtkIdType cellId, vtkIdType& cellSize, vtkIdType* cellPoints)
1914{
1915 this->Dispatch(vtkCellArray_detail::GetCellAtIdImpl{}, cellId, cellSize, cellPoints);
1916}
1917
1918//----------------------------------------------------------------------------
1920{
1921 vtkIdType pointId;
1922 this->Dispatch(vtkCellArray_detail::CellPointAtIdImpl{}, cellId, cellPointIndex, pointId);
1923 return pointId;
1924}
1925
1926//----------------------------------------------------------------------------
1928 VTK_SIZEHINT(pts, npts)
1929{
1930 vtkIdType cellId;
1931 this->Dispatch(vtkCellArray_detail::InsertNextCellImpl{}, npts, pts, cellId);
1932 return cellId;
1933}
1934
1935//----------------------------------------------------------------------------
1937{
1938 vtkIdType cellId;
1940 return cellId;
1941}
1942
1943//----------------------------------------------------------------------------
1948
1949//----------------------------------------------------------------------------
1951{
1953}
1954
1955//----------------------------------------------------------------------------
1957{
1958 vtkIdType cellId;
1959 this->Dispatch(
1961 return cellId;
1962}
1963
1964//----------------------------------------------------------------------------
1966{
1967 vtkIdList* pts = cell->GetPointIds();
1968 vtkIdType cellId;
1969 this->Dispatch(
1971 return cellId;
1972}
1973
1974//----------------------------------------------------------------------------
1976{
1978}
1979
1980VTK_ABI_NAMESPACE_END
1981#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:368
#define VTK_SIZEHINT(...)
#define VTK_EXPECTS(x)
#define VTK_MARSHALMANUAL
#define VTK_NEWINSTANCE