VTK  9.6.20260219
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 "vtkDataArrayAccessor.h" // Needed for inline methods
247#include "vtkDataArrayRange.h" // Needed for inline methods
248#include "vtkDeprecation.h" // For VTK_DEPRECATED_IN_9_6_0
249#include "vtkFeatures.h" // for VTK_USE_MEMKIND
250#include "vtkSmartPointer.h" // For vtkSmartPointer
251#include "vtkTypeInt32Array.h" // Needed for inline methods
252#include "vtkTypeInt64Array.h" // Needed for inline methods
253#include "vtkTypeList.h" // Needed for ArrayList definition
254
255#include <cassert> // Needed for assert
256#include <initializer_list> // Needed for API
257#include <type_traits> // Needed for std::is_same
258#include <utility> // Needed for std::forward
259
280#define VTK_CELL_ARRAY_V2
281
282VTK_ABI_NAMESPACE_BEGIN
284class vtkIdTypeArray;
285
286class VTKCOMMONDATAMODEL_EXPORT VTK_MARSHALMANUAL vtkCellArray : public vtkAbstractCellArray
287{
288public:
289 using ArrayType32 VTK_DEPRECATED_IN_9_6_0("Use AOSArray32 instead.") = vtkTypeInt32Array;
290 using ArrayType64 VTK_DEPRECATED_IN_9_6_0("Use AOSArray64 instead.") = vtkTypeInt64Array;
295
297
301 static vtkCellArray* New();
303 void PrintSelf(ostream& os, vtkIndent indent) override;
304 void PrintDebug(ostream& os);
306
308
317 "Use vtkArrayDispatch::OffsetsArrays/ConnectivityArrays instead.") =
318 vtkTypeList::Create<vtkTypeInt32Array, vtkTypeInt64Array>;
320
322
334
342 */
343 VTK_DEPRECATED_IN_9_6_0("Use AllocateEstimate or AllocateExact instead.")
344 vtkTypeBool Allocate(vtkIdType sz, vtkIdType vtkNotUsed(ext) = 1000)
345 {
346 return this->AllocateExact(sz, sz) ? 1 : 0;
347 }
348
356 * @sa Squeeze AllocateExact AllocateCopy
357 */
358 bool AllocateEstimate(vtkIdType numCells, vtkIdType maxCellSize)
359 {
360 return this->AllocateExact(numCells, numCells * maxCellSize);
361 }
362
372 bool AllocateExact(vtkIdType numCells, vtkIdType connectivitySize);
373
381 * @sa Squeeze AllocateEstimate AllocateExact
382 */
383 bool AllocateCopy(vtkCellArray* other)
384 {
385 return this->AllocateExact(other->GetNumberOfCells(), other->GetNumberOfConnectivityIds());
386 }
387
397 bool ResizeExact(vtkIdType numCells, vtkIdType connectivitySize);
398
402 void Initialize() override;
403
407 void Reset();
408
414 void Squeeze();
415
426 bool IsValid();
427
431 vtkIdType GetNumberOfCells() const override { return this->Offsets->GetNumberOfValues() - 1; }
432
437 vtkIdType GetNumberOfOffsets() const override { return this->Offsets->GetNumberOfValues(); }
438
440 * Get the offset (into the connectivity) for a specified cell id.
441 */
442 vtkIdType GetOffset(vtkIdType cellId) override
443 {
444 return static_cast<vtkIdType>(this->Offsets->GetComponent(cellId, 0));
445 }
446
448 * Set the offset (into the connectivity) for a specified cell id.
449 */
450 void SetOffset(vtkIdType cellId, vtkIdType offset)
451 {
452 this->Offsets->SetComponent(cellId, 0, static_cast<double>(offset));
453 }
454
456 * Get the size of the connectivity array that stores the point ids.
457 */
459 {
460 return this->Connectivity->GetNumberOfValues();
461 }
462
469
471
484 void SetData(AOSArray32*, AOSArray32* connectivity);
485 void SetData(AOSArray64*, AOSArray64* connectivity);
486 void SetData(AffineArray32*, AOSArray32* connectivity);
487 void SetData(AffineArray64*, AOSArray64* connectivity);
489
504 bool SetData(vtkDataArray* offsets, vtkDataArray* connectivity);
505
531 Generic,
532 };
533
537 StorageTypes GetStorageType() const noexcept { return this->StorageType; }
538
542 bool IsStorage64Bit() const { return this->StorageType == StorageTypes::Int64; }
543
547 bool IsStorage32Bit() const { return this->StorageType == StorageTypes::Int32; }
548
553
558
560 * @return True if the internal storage is using FixedSize arrays.
561 */
562 bool IsStorageFixedSize() const
563 {
564 return this->IsStorageFixedSize32Bit() || this->IsStorageFixedSize64Bit();
565 }
566
570 bool IsStorageGeneric() const { return this->StorageType == StorageTypes::Generic; }
571
576 * a pointer to vtkIdType can be used.
577 */
578 bool IsStorageShareable() const override
579 {
580 switch (this->StorageType)
581 {
584 return std::is_same_v<vtkTypeInt32, vtkIdType>;
587 return std::is_same_v<vtkTypeInt64, vtkIdType>;
589 default:
590 return false;
591 }
592 }
593
606 void UseFixedSize64BitStorage(vtkIdType cellSize);
609
625 {
626 switch (type)
627 {
629 return this->CanConvertTo32BitStorage();
631 return this->CanConvertTo64BitStorage();
637 default:
638 return true;
639 }
640 }
642
667 {
668 switch (type)
669 {
671 return this->ConvertTo32BitStorage();
673 return this->ConvertTo64BitStorage();
675 return this->ConvertToFixedSize32BitStorage();
677 return this->ConvertToFixedSize64BitStorage();
679 default:
680 return true;
681 }
682 }
684
690 vtkDataArray* GetOffsetsArray() const { return this->Offsets; }
691 VTK_DEPRECATED_IN_9_6_0("Use GetOffsetsAOSArray32() instead.")
692 vtkTypeInt32Array* GetOffsetsArray32() const
693 {
694 return vtkTypeInt32Array::FastDownCast(this->Offsets);
695 }
697 VTK_DEPRECATED_IN_9_6_0("Use GetOffsetsAOSArray64() instead.")
698 vtkTypeInt64Array* GetOffsetsArray64() const
699 {
700 return vtkTypeInt64Array::FastDownCast(this->Offsets);
702 AOSArray64* GetOffsetsAOSArray64() const { return AOSArray64::FastDownCast(this->Offsets); }
703 AffineArray32* GetOffsetsAffineArray32() const
704 {
706 }
707 AffineArray64* GetOffsetsAffineArray64() const
708 {
709 return AffineArray64::FastDownCast(this->Offsets);
710 }
712
720 VTK_DEPRECATED_IN_9_6_0("Use GetConnectivityAOSArray32() instead.")
721 vtkTypeInt32Array* GetConnectivityArray32() const
722 {
723 return vtkTypeInt32Array::FastDownCast(this->Connectivity);
724 }
725 AOSArray32* GetConnectivityAOSArray32() const
726 {
727 return AOSArray32::FastDownCast(this->Connectivity);
729 VTK_DEPRECATED_IN_9_6_0("Use GetConnectivityAOSArray64() instead.")
730 vtkTypeInt64Array* GetConnectivityArray64() const
731 {
732 return vtkTypeInt64Array::FastDownCast(this->Connectivity);
733 }
734 AOSArray64* GetConnectivityAOSArray64() const
735 {
736 return AOSArray64::FastDownCast(this->Connectivity);
737 }
739
748 vtkIdType IsHomogeneous() const override;
749
759 void InitTraversal();
760
775 int GetNextCell(vtkIdType& npts, vtkIdType const*& pts) VTK_SIZEHINT(pts, npts);
776
787 int GetNextCell(vtkIdList* pts);
788
799 inline void GetCellAtId(vtkIdType cellId, vtkIdType& cellSize, vtkIdType const*& cellPoints,
800 vtkIdList* ptIds) VTK_SIZEHINT(cellPoints, cellSize)
801 VTK_EXPECTS(0 <= cellId && cellId < GetNumberOfCells()) override;
802
808 inline void GetCellAtId(vtkIdType cellId, vtkIdList* pts)
809 VTK_EXPECTS(0 <= cellId && cellId < GetNumberOfCells()) override;
810
818 inline void GetCellAtId(vtkIdType cellId, vtkIdType& cellSize, vtkIdType* cellPoints)
819 VTK_SIZEHINT(cellPoints, cellSize)
820 VTK_EXPECTS(0 <= cellId && cellId < GetNumberOfCells()) override;
821
825 vtkIdType GetCellPointAtId(vtkIdType cellId, vtkIdType cellPointIndex) const
826 VTK_EXPECTS(0 <= cellId && cellId < GetNumberOfCells())
827 VTK_EXPECTS(0 <= cellPointIndex && cellPointIndex < this->GetCellSize(cellId));
828
832 vtkIdType GetCellSize(vtkIdType cellId) const override;
833
837 inline vtkIdType InsertNextCell(vtkCell* cell);
838
843 inline vtkIdType InsertNextCell(vtkIdType npts, const vtkIdType* pts) VTK_SIZEHINT(pts, npts);
844
850
856 * `InsertNextCell(int)` overload instead.
857 */
858 vtkIdType InsertNextCell(const std::initializer_list<vtkIdType>& cell)
859 {
860 return this->InsertNextCell(static_cast<vtkIdType>(cell.size()), cell.begin());
861 }
862
869 inline vtkIdType InsertNextCell(int npts);
870
875 void InsertCellPoint(vtkIdType id);
876
881 void UpdateCellCount(int npts);
882
891 void SetTraversalCellId(vtkIdType cellId);
893
897 void ReverseCellAtId(vtkIdType cellId) VTK_EXPECTS(0 <= cellId && cellId < GetNumberOfCells());
898
907 void ReplaceCellAtId(vtkIdType cellId, vtkIdList* list);
908 void ReplaceCellAtId(vtkIdType cellId, vtkIdType cellSize, const vtkIdType* cellPoints)
909 VTK_EXPECTS(0 <= cellId && cellId < GetNumberOfCells()) VTK_SIZEHINT(cellPoints, cellSize);
911
919 void ReplaceCellPointAtId(vtkIdType cellId, vtkIdType cellPointIndex, vtkIdType newPointId);
920
926 * results.
927 */
928 void ReplaceCellAtId(vtkIdType cellId, const std::initializer_list<vtkIdType>& cell)
929 {
930 this->ReplaceCellAtId(cellId, static_cast<vtkIdType>(cell.size()), cell.begin());
931 }
932
937 int GetMaxCellSize() override;
938
942 void DeepCopy(vtkAbstractCellArray* ca) override;
943
947 void ShallowCopy(vtkAbstractCellArray* ca) override;
948
952 void Append(vtkCellArray* src, vtkIdType pointOffset = 0);
953
965
978 void ImportLegacyFormat(const vtkIdType* data, vtkIdType len) VTK_SIZEHINT(data, len);
980
992 void AppendLegacyFormat(vtkIdTypeArray* data, vtkIdType ptOffset = 0);
993 void AppendLegacyFormat(const vtkIdType* data, vtkIdType len, vtkIdType ptOffset = 0)
994 VTK_SIZEHINT(data, len);
996
1005 unsigned long GetActualMemorySize() const;
1006
1007 // The following code is used to support
1008
1009 // The wrappers get understandably confused by some of the template code below
1010#ifndef __VTK_WRAP__
1011 /*
1012 * Utilities class that every dispatch functor used with Dispatch() must inherit
1013 * or optionally use its static methods.
1014 */
1015 struct DispatchUtilities
1017 template <class ArrayT>
1020 template <class OffsetsT>
1021 vtkIdType GetNumberOfCells(OffsetsT* offsets)
1022 {
1023 return offsets->GetNumberOfValues() - 1;
1024 }
1026 template <class ArrayT>
1027 static decltype(vtk::DataArrayValueRange<1, vtkIdType>(std::declval<ArrayT>())) GetRange(
1028 ArrayT* array)
1029 {
1031 }
1033 template <class OffsetsT>
1034 static vtkIdType GetBeginOffset(OffsetsT* offsets, vtkIdType cellId)
1035 {
1036 return static_cast<vtkIdType>(GetRange(offsets)[cellId]);
1037 }
1039 template <class OffsetsT>
1040 static vtkIdType GetEndOffset(OffsetsT* offsets, vtkIdType cellId)
1041 {
1042 return static_cast<vtkIdType>(GetRange(offsets)[cellId + 1]);
1043 }
1045 template <class OffsetsT>
1046 static vtkIdType GetCellSize(OffsetsT* offsets, vtkIdType cellId)
1047 {
1048 auto offsetsRange = GetRange(offsets);
1049 return static_cast<vtkIdType>(offsetsRange[cellId + 1] - offsetsRange[cellId]);
1050 }
1051
1052 template <class OffsetsT, class ConnectivityT>
1053 static decltype(vtk::DataArrayValueRange<1, vtkIdType>(std::declval<ConnectivityT>()))
1054 GetCellRange(OffsetsT* offsets, ConnectivityT* conn, vtkIdType cellId)
1055 {
1056 auto offsetsRange = GetRange(offsets);
1058 conn, offsetsRange[cellId], offsetsRange[cellId + 1]);
1059 }
1060 };
1061
1063
1126 template <typename Functor, typename... Args>
1127 void Dispatch(Functor&& functor, Args&&... args)
1128 {
1129 switch (this->StorageType)
1130 {
1132 functor(static_cast<AOSArray32*>(this->Offsets.Get()),
1133 static_cast<AOSArray32*>(this->Connectivity.Get()), std::forward<Args>(args)...);
1134 break;
1136 functor(static_cast<AOSArray64*>(this->Offsets.Get()),
1137 static_cast<AOSArray64*>(this->Connectivity.Get()), std::forward<Args>(args)...);
1138 break;
1140 functor(static_cast<AffineArray32*>(this->Offsets.Get()),
1141 static_cast<AOSArray32*>(this->Connectivity.Get()), std::forward<Args>(args)...);
1142 break;
1144 functor(static_cast<AffineArray64*>(this->Offsets.Get()),
1145 static_cast<AOSArray64*>(this->Connectivity.Get()), std::forward<Args>(args)...);
1146 break;
1148 default:
1149 functor(this->Offsets.Get(), this->Connectivity.Get(), std::forward<Args>(args)...);
1150 break;
1151 }
1153 template <typename Functor, typename... Args>
1154 void Dispatch(Functor&& functor, Args&&... args) const
1155 {
1156 switch (this->StorageType)
1157 {
1158 case StorageTypes::Int32:
1159 functor(static_cast<AOSArray32*>(this->Offsets.Get()),
1160 static_cast<AOSArray32*>(this->Connectivity.Get()), std::forward<Args>(args)...);
1161 break;
1162 case StorageTypes::Int64:
1163 functor(static_cast<AOSArray64*>(this->Offsets.Get()),
1164 static_cast<AOSArray64*>(this->Connectivity.Get()), std::forward<Args>(args)...);
1165 break;
1166 case StorageTypes::FixedSizeInt32:
1167 functor(static_cast<AffineArray32*>(this->Offsets.Get()),
1168 static_cast<AOSArray32*>(this->Connectivity.Get()), std::forward<Args>(args)...);
1169 break;
1170 case StorageTypes::FixedSizeInt64:
1171 functor(static_cast<AffineArray64*>(this->Offsets.Get()),
1172 static_cast<AOSArray64*>(this->Connectivity.Get()), std::forward<Args>(args)...);
1173 break;
1174 case StorageTypes::Generic:
1175 default:
1176 functor(this->Offsets.Get(), this->Connectivity.Get(), std::forward<Args>(args)...);
1177 break;
1178 }
1179 }
1181
1182 // Holds connectivity and offset arrays of the given ArrayType.
1183 // VTK_DEPRECATED_IN_9_6_0("Use DispatchUtilities")
1184 template <typename ArrayT>
1187 using ArrayType = ArrayT;
1188 using ValueType = typename ArrayType::ValueType;
1189 using CellRangeType = decltype(vtk::DataArrayValueRange<1>(std::declval<ArrayType>()));
1190
1191 // We can't just use is_same here, since binary compatible representations
1192 // (e.g. int and long) are distinct types. Instead, ensure that ValueType
1193 // is a signed integer the same size as vtkIdType.
1194 // If this value is true, ValueType pointers may be safely converted to
1195 // vtkIdType pointers via reinterpret cast.
1196 static constexpr bool ValueTypeIsSameAsIdType = std::is_integral<ValueType>::value &&
1197 std::is_signed<ValueType>::value && (sizeof(ValueType) == sizeof(vtkIdType));
1199 ArrayType* GetOffsets() { return this->Offsets; }
1200 const ArrayType* GetOffsets() const { return this->Offsets; }
1202 ArrayType* GetConnectivity() { return this->Connectivity; }
1203 const ArrayType* GetConnectivity() const { return this->Connectivity; }
1204
1205 vtkIdType GetNumberOfCells() const;
1206
1207 vtkIdType GetBeginOffset(vtkIdType cellId) const;
1208
1209 vtkIdType GetEndOffset(vtkIdType cellId) const;
1210
1211 vtkIdType GetCellSize(vtkIdType cellId) const;
1212
1214
1215 friend class vtkCellArray;
1217 protected:
1218 VisitState()
1219 {
1222 this->Offsets->InsertNextValue(0);
1224 {
1225 this->IsInMemkind = true;
1228 ~VisitState() = default;
1229 void* operator new(size_t nSize)
1230 {
1231 void* r;
1232#ifdef VTK_USE_MEMKIND
1234#else
1235 r = malloc(nSize);
1236#endif
1237 return r;
1238 }
1239 void operator delete(void* p)
1240 {
1241#ifdef VTK_USE_MEMKIND
1242 VisitState* a = static_cast<VisitState*>(p);
1243 if (a->IsInMemkind)
1244 {
1246 }
1247 else
1248 {
1249 free(p);
1250 }
1251#else
1252 free(p);
1253#endif
1258
1259 private:
1260 VisitState(const VisitState&) = delete;
1261 VisitState& operator=(const VisitState&) = delete;
1262 bool IsInMemkind = false;
1263 };
1264
1265private: // Helpers that allow Visit to return a value:
1266 template <typename Functor, typename... Args>
1267 using GetReturnType = decltype(std::declval<Functor>()(
1268 std::declval<VisitState<AOSArray32>&>(), std::declval<Args>()...));
1269
1270 template <typename Functor, typename... Args>
1271 struct ReturnsVoid : std::is_same<GetReturnType<Functor, Args...>, void>
1272 {
1273 };
1274
1275public:
1345 template <typename Functor, typename... Args,
1346 typename = typename std::enable_if<ReturnsVoid<Functor, Args...>::value>::type>
1347 VTK_DEPRECATED_IN_9_6_0("Use Dispatch instead")
1348 void Visit(Functor&& functor, Args&&... args)
1349 {
1350 switch (this->StorageType)
1351 {
1353 {
1357 functor(state, std::forward<Args>(args)...);
1358 break;
1359 }
1361 {
1365 functor(state, std::forward<Args>(args)...);
1366 break;
1367 }
1371 default:
1372 {
1373 vtkWarningMacro("Use Dispatch");
1374 break;
1375 }
1376 }
1377 }
1378
1379 template <typename Functor, typename... Args,
1380 typename = typename std::enable_if<ReturnsVoid<Functor, Args...>::value>::type>
1381 VTK_DEPRECATED_IN_9_6_0("Use Dispatch instead")
1382 void Visit(Functor&& functor, Args&&... args) const
1383 {
1384 switch (this->StorageType)
1385 {
1387 {
1391 functor(state, std::forward<Args>(args)...);
1392 break;
1393 }
1395 {
1399 functor(state, std::forward<Args>(args)...);
1400 break;
1401 }
1405 default:
1406 {
1407 vtkWarningMacro("Use Dispatch");
1408 break;
1409 }
1410 }
1411 }
1412
1413 template <typename Functor, typename... Args,
1414 typename = typename std::enable_if<!ReturnsVoid<Functor, Args...>::value>::type>
1415 VTK_DEPRECATED_IN_9_6_0("Use Dispatch instead")
1416 GetReturnType<Functor, Args...> Visit(Functor&& functor, Args&&... args)
1417 {
1418 switch (this->StorageType)
1419 {
1421 {
1425 return functor(state, std::forward<Args>(args)...);
1426 }
1428 {
1432 return functor(state, std::forward<Args>(args)...);
1433 }
1437 default:
1438 {
1439 vtkWarningMacro("Use Dispatch");
1440 return GetReturnType<Functor, Args...>();
1441 }
1442 }
1443 }
1444 template <typename Functor, typename... Args,
1445 typename = typename std::enable_if<!ReturnsVoid<Functor, Args...>::value>::type>
1446 VTK_DEPRECATED_IN_9_6_0("Use Dispatch instead")
1447 GetReturnType<Functor, Args...> Visit(Functor&& functor, Args&&... args) const
1448 {
1449 switch (this->StorageType)
1450 {
1452 {
1456 return functor(state, std::forward<Args>(args)...);
1457 }
1459 {
1463 return functor(state, std::forward<Args>(args)...);
1464 }
1468 default:
1469 {
1470 vtkWarningMacro("Use Dispatch");
1471 return GetReturnType<Functor, Args...>();
1472 }
1473 }
1474 }
1475
1476#endif // __VTK_WRAP__
1477
1479
1487 static void SetDefaultStorageIs64Bit(bool val) { vtkCellArray::DefaultStorageIs64Bit = val; }
1489
1490 //=================== Begin Legacy Methods ===================================
1491 // These should be deprecated at some point as they are confusing or very slow
1492
1499 VTK_DEPRECATED_IN_9_6_0("This call has no effect.")
1500 virtual void SetNumberOfCells(vtkIdType);
1501
1513 VTK_DEPRECATED_IN_9_6_0("Use AllocateEstimate directly instead.")
1514 vtkIdType EstimateSize(vtkIdType numCells, int maxPtsPerCell);
1515
1524 VTK_DEPRECATED_IN_9_6_0("Method incompatible with current internal storage.")
1526
1533 VTK_DEPRECATED_IN_9_6_0("Method incompatible with current internal storage.")
1535
1545 VTK_DEPRECATED_IN_9_6_0("Use GetCellAtId.")
1546 void GetCell(vtkIdType loc, vtkIdType& npts, const vtkIdType*& pts)
1547 VTK_EXPECTS(0 <= loc && loc < GetNumberOfConnectivityEntries()) VTK_SIZEHINT(pts, npts);
1548
1555 VTK_DEPRECATED_IN_9_6_0("Use GetCellAtId.")
1556 void GetCell(vtkIdType loc, vtkIdList* pts)
1557 VTK_EXPECTS(0 <= loc && loc < GetNumberOfConnectivityEntries());
1558
1565 VTK_DEPRECATED_IN_9_6_0("Use GetNumberOfCells.")
1566 vtkIdType GetInsertLocation(int npts);
1567
1575 VTK_DEPRECATED_IN_9_6_0("Use GetTraversalCellId.")
1577 VTK_DEPRECATED_IN_9_6_0("Use GetTraversalCellId.")
1579 VTK_DEPRECATED_IN_9_6_0("Use SetTraversalCellId.")
1582
1590 VTK_DEPRECATED_IN_9_6_0("Use ReverseCellAtId.")
1591 void ReverseCell(vtkIdType loc) VTK_EXPECTS(0 <= loc && loc < GetNumberOfConnectivityEntries());
1592
1604 VTK_DEPRECATED_IN_9_6_0("Use ReplaceCellAtId.")
1605 void ReplaceCell(vtkIdType loc, int npts, const vtkIdType pts[])
1606 VTK_EXPECTS(0 <= loc && loc < GetNumberOfConnectivityEntries()) VTK_SIZEHINT(pts, npts);
1607
1622 VTK_DEPRECATED_IN_9_6_0("Use ImportLegacyFormat or SetData instead.")
1623 void SetCells(vtkIdType ncells, vtkIdTypeArray* cells);
1624
1636 "Use ExportLegacyFormat, or GetOffsetsArray/GetConnectivityArray instead.")
1638
1639 //=================== End Legacy Methods =====================================
1640
1641 friend class vtkCellArrayIterator;
1643protected:
1644 vtkCellArray();
1645 ~vtkCellArray() override;
1651
1653
1654 static bool DefaultStorageIs64Bit;
1655
1656private:
1657 vtkCellArray(const vtkCellArray&) = delete;
1658 void operator=(const vtkCellArray&) = delete;
1659};
1660
1661template <typename ArrayT>
1663{
1664 return this->Offsets->GetNumberOfValues() - 1;
1665}
1666
1667template <typename ArrayT>
1669{
1670 return static_cast<vtkIdType>(this->Offsets->GetValue(cellId));
1671}
1672
1673template <typename ArrayT>
1675{
1676 return static_cast<vtkIdType>(this->Offsets->GetValue(cellId + 1));
1677}
1678
1679template <typename ArrayT>
1681{
1682 return this->GetEndOffset(cellId) - this->GetBeginOffset(cellId);
1683}
1684
1685template <typename ArrayT>
1692VTK_ABI_NAMESPACE_END
1693
1695{
1696VTK_ABI_NAMESPACE_BEGIN
1697
1699{
1700 // Insert full cell
1701 template <class OffsetsT, class ConnectivityT>
1702 void operator()(OffsetsT* offsets, ConnectivityT* conn, const vtkIdType npts,
1703 const vtkIdType pts[], vtkIdType& cellId)
1704 {
1705 using ValueType = GetAPIType<OffsetsT>;
1706 using OffsetsAccessorType = vtkDataArrayAccessor<OffsetsT>;
1707 using ConnectivityAccessorType = vtkDataArrayAccessor<ConnectivityT>;
1708 OffsetsAccessorType offsetsAccesor(offsets);
1709 ConnectivityAccessorType connAccesor(conn);
1710
1711 cellId = offsets->GetNumberOfValues() - 1;
1712
1713 offsetsAccesor.InsertNext(static_cast<ValueType>(conn->GetNumberOfValues() + npts));
1714
1715 for (vtkIdType i = 0; i < npts; ++i)
1716 {
1717 connAccesor.InsertNext(static_cast<ValueType>(pts[i]));
1718 }
1719 }
1720
1721 // Just update offset table (for incremental API)
1722 template <class OffsetsT, class ConnectivityT>
1723 void operator()(OffsetsT* offsets, ConnectivityT* conn, const vtkIdType npts, vtkIdType& cellId)
1724 {
1725 using ValueType = GetAPIType<OffsetsT>;
1726 using AccessorType = vtkDataArrayAccessor<OffsetsT>;
1727 AccessorType offsetsAccesor(offsets);
1728
1729 cellId = offsets->GetNumberOfValues() - 1;
1730
1731 offsetsAccesor.InsertNext(static_cast<ValueType>(conn->GetNumberOfValues() + npts));
1732 }
1733};
1734
1735// for incremental API:
1737{
1738 template <class OffsetsT, class ConnectivityT>
1739 void operator()(OffsetsT* offsets, ConnectivityT* vtkNotUsed(conn), const vtkIdType npts)
1740 {
1741 using ValueType = GetAPIType<OffsetsT>;
1742
1743 auto offsetsRange = GetRange(offsets);
1744 const ValueType cellBegin = offsetsRange[offsets->GetMaxId() - 1];
1745 offsetsRange[offsets->GetMaxId()] = static_cast<ValueType>(cellBegin + npts);
1746 }
1747};
1748
1750{
1751 template <class OffsetsT, class ConnectivityT>
1753 OffsetsT* offsets, ConnectivityT* vtkNotUsed(conn), vtkIdType cellId, vtkIdType& cellSize)
1754 {
1755 cellSize = GetCellSize(offsets, cellId);
1756 }
1757};
1758
1760{
1761 template <class OffsetsT, class ConnectivityT>
1762 void operator()(OffsetsT* offsets, ConnectivityT* conn, const vtkIdType cellId, vtkIdList* ids)
1763 {
1764 auto offsetsRange = GetRange(offsets);
1765 const auto& beginOffset = offsetsRange[cellId];
1766 const auto& endOffset = offsetsRange[cellId + 1];
1767 const vtkIdType cellSize = static_cast<vtkIdType>(endOffset - beginOffset);
1768 const auto cellConnectivity = GetRange(conn).begin() + beginOffset;
1769
1770 // ValueType differs from vtkIdType, so we have to copy into a temporary buffer:
1771 ids->SetNumberOfIds(cellSize);
1772 vtkIdType* idPtr = ids->GetPointer(0);
1773 for (vtkIdType i = 0; i < cellSize; ++i)
1774 {
1775 idPtr[i] = static_cast<vtkIdType>(cellConnectivity[i]);
1776 }
1777 }
1778
1779 template <class OffsetsT, class ConnectivityT>
1780 void operator()(OffsetsT* offsets, ConnectivityT* conn, const vtkIdType cellId,
1781 vtkIdType& cellSize, vtkIdType* cellPoints)
1782 {
1783 auto offsetsRange = GetRange(offsets);
1784 const auto& beginOffset = offsetsRange[cellId];
1785 const auto& endOffset = offsetsRange[cellId + 1];
1786 cellSize = static_cast<vtkIdType>(endOffset - beginOffset);
1787 const auto cellConnectivity = GetRange(conn).begin() + beginOffset;
1788
1789 // ValueType differs from vtkIdType, so we have to copy into a temporary buffer:
1790 for (vtkIdType i = 0; i < cellSize; ++i)
1791 {
1792 cellPoints[i] = static_cast<vtkIdType>(cellConnectivity[i]);
1793 }
1794 }
1795
1796 // SFINAE helper to check if a Functors's connectivity array's memory can be used as a vtkIdType*.
1797 template <typename ConnectivityT>
1799 {
1800 static constexpr bool value =
1801 std::is_base_of_v<vtkAOSDataArrayTemplate<vtkIdType>, ConnectivityT>;
1802 };
1803
1804 template <class OffsetsT, class ConnectivityT>
1805 typename std::enable_if_t<CanShareConnPtr<ConnectivityT>::value, void> operator()(
1806 OffsetsT* offsets, ConnectivityT* conn, const vtkIdType cellId, vtkIdType& cellSize,
1807 vtkIdType const*& cellPoints, vtkIdList* vtkNotUsed(temp))
1808 {
1809 auto offsetsRange = GetRange(offsets);
1810 const auto& beginOffset = offsetsRange[cellId];
1811 const auto& endOffset = offsetsRange[cellId + 1];
1812 cellSize = static_cast<vtkIdType>(endOffset - beginOffset);
1813 // This is safe, see CanShareConnPtr helper above.
1814 cellPoints = conn->GetPointer(beginOffset);
1815 }
1816
1817 template <class OffsetsT, class ConnectivityT>
1818 typename std::enable_if_t<!CanShareConnPtr<ConnectivityT>::value, void> operator()(
1819 OffsetsT* offsets, ConnectivityT* conn, const vtkIdType cellId, vtkIdType& cellSize,
1820 vtkIdType const*& cellPoints, vtkIdList* temp)
1821 {
1822 auto offsetsRange = GetRange(offsets);
1823 const auto& beginOffset = offsetsRange[cellId];
1824 const auto& endOffset = offsetsRange[cellId + 1];
1825 cellSize = static_cast<vtkIdType>(endOffset - beginOffset);
1826 const auto cellConnectivity = GetRange(conn).begin() + beginOffset;
1827
1828 temp->SetNumberOfIds(cellSize);
1829 vtkIdType* tempPtr = temp->GetPointer(0);
1830 for (vtkIdType i = 0; i < cellSize; ++i)
1831 {
1832 tempPtr[i] = static_cast<vtkIdType>(cellConnectivity[i]);
1833 }
1834
1835 cellPoints = tempPtr;
1836 }
1837};
1838
1840{
1841 template <class OffsetsT, class ConnectivityT>
1842 void operator()(OffsetsT* offsets, ConnectivityT* conn, vtkIdType cellId,
1843 vtkIdType cellPointIndex, vtkIdType& pointId)
1844 {
1845 pointId =
1846 static_cast<vtkIdType>(GetRange(conn)[GetBeginOffset(offsets, cellId) + cellPointIndex]);
1847 }
1848};
1849
1851{
1852 template <class OffsetsT, class ConnectivityT>
1853 void operator()(OffsetsT* offsets, ConnectivityT* conn)
1854 {
1855 using ValueType = GetAPIType<OffsetsT>;
1856 using AccessorType = vtkDataArrayAccessor<OffsetsT>;
1857 offsets->Reset();
1858 conn->Reset();
1859 AccessorType accessor(offsets);
1860 ValueType firstOffset = 0;
1861 accessor.InsertNext(firstOffset);
1862 }
1863};
1864
1866{
1867 template <class OffsetsT, class ConnectivityT>
1868 void operator()(OffsetsT* vtkNotUsed(offsets), ConnectivityT* conn, vtkIdType id)
1869 {
1870 using ValueType = GetAPIType<ConnectivityT>;
1871 using AccessorType = vtkDataArrayAccessor<ConnectivityT>;
1872 AccessorType accessor(conn);
1873 accessor.InsertNext(static_cast<ValueType>(id));
1874 }
1875};
1876
1877VTK_ABI_NAMESPACE_END
1878} // end namespace vtkCellArray_detail
1879
1880VTK_ABI_NAMESPACE_BEGIN
1881//----------------------------------------------------------------------------
1883{
1884 this->TraversalCellId = 0;
1885}
1886
1887//----------------------------------------------------------------------------
1888inline int vtkCellArray::GetNextCell(vtkIdType& npts, vtkIdType const*& pts) VTK_SIZEHINT(pts, npts)
1889{
1890 if (this->TraversalCellId < this->GetNumberOfCells())
1891 {
1892 this->GetCellAtId(this->TraversalCellId, npts, pts);
1893 ++this->TraversalCellId;
1894 return 1;
1895 }
1896
1897 npts = 0;
1898 pts = nullptr;
1899 return 0;
1900}
1901
1902//----------------------------------------------------------------------------
1904{
1906 {
1907 this->GetCellAtId(this->TraversalCellId, pts);
1908 ++this->TraversalCellId;
1909 return 1;
1910 }
1911
1912 pts->Reset();
1913 return 0;
1914}
1915//----------------------------------------------------------------------------
1917{
1918 vtkIdType cellSize;
1919 this->Dispatch(vtkCellArray_detail::GetCellSizeImpl{}, cellId, cellSize);
1920 return cellSize;
1921}
1922
1923//----------------------------------------------------------------------------
1924void vtkCellArray::GetCellAtId(vtkIdType cellId, vtkIdType& cellSize, vtkIdType const*& cellPoints,
1925 vtkIdList* ptIds) VTK_SIZEHINT(cellPoints, cellSize)
1926{
1927 this->Dispatch(vtkCellArray_detail::GetCellAtIdImpl{}, cellId, cellSize, cellPoints, ptIds);
1928}
1929
1930//----------------------------------------------------------------------------
1932{
1933 this->Dispatch(vtkCellArray_detail::GetCellAtIdImpl{}, cellId, pts);
1934}
1935
1936//----------------------------------------------------------------------------
1937void vtkCellArray::GetCellAtId(vtkIdType cellId, vtkIdType& cellSize, vtkIdType* cellPoints)
1938{
1939 this->Dispatch(vtkCellArray_detail::GetCellAtIdImpl{}, cellId, cellSize, cellPoints);
1940}
1941
1942//----------------------------------------------------------------------------
1944{
1945 vtkIdType pointId;
1946 this->Dispatch(vtkCellArray_detail::CellPointAtIdImpl{}, cellId, cellPointIndex, pointId);
1947 return pointId;
1948}
1949
1950//----------------------------------------------------------------------------
1952{
1953 vtkIdType cellId;
1954 this->Dispatch(vtkCellArray_detail::InsertNextCellImpl{}, npts, pts, cellId);
1955 return cellId;
1956}
1957
1958//----------------------------------------------------------------------------
1960{
1961 vtkIdType cellId;
1963 return cellId;
1964}
1965
1966//----------------------------------------------------------------------------
1971
1972//----------------------------------------------------------------------------
1974{
1976}
1977
1978//----------------------------------------------------------------------------
1980{
1981 vtkIdType cellId;
1982 this->Dispatch(
1984 return cellId;
1985}
1986
1987//----------------------------------------------------------------------------
1989{
1990 vtkIdList* pts = cell->GetPointIds();
1991 vtkIdType cellId;
1992 this->Dispatch(
1994 return cellId;
1995}
1996
1997//----------------------------------------------------------------------------
1999{
2001}
2002
2003VTK_ABI_NAMESPACE_END
2004#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:363
#define VTK_SIZEHINT(...)
#define VTK_EXPECTS(x)
#define VTK_MARSHALMANUAL
#define VTK_NEWINSTANCE