VTK  9.6.1
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
126
127#ifndef vtkCellArray_h
128#define vtkCellArray_h
129
130#include "vtkAbstractCellArray.h"
131#include "vtkCommonDataModelModule.h" // For export macro
132#include "vtkWrappingHints.h" // For VTK_MARSHALMANUAL
133
134#include "vtkAOSDataArrayTemplate.h" // Needed for inline methods
135#include "vtkAffineArray.h" // Needed for inline methods
136#include "vtkCell.h" // Needed for inline methods
137#include "vtkDataArrayRange.h" // Needed for inline methods
138#include "vtkDeprecation.h" // For VTK_DEPRECATED_IN_9_6_0
139#include "vtkFeatures.h" // for VTK_USE_MEMKIND
140#include "vtkSmartPointer.h" // For vtkSmartPointer
141#include "vtkTypeInt32Array.h" // Needed for inline methods
142#include "vtkTypeInt64Array.h" // Needed for inline methods
143#include "vtkTypeList.h" // Needed for ArrayList definition
144
145#include <cassert> // Needed for assert
146#include <initializer_list> // Needed for API
147#include <type_traits> // Needed for std::is_same
148#include <utility> // Needed for std::forward
149
170#define VTK_CELL_ARRAY_V2
171
172VTK_ABI_NAMESPACE_BEGIN
174class vtkIdTypeArray;
175
176class VTKCOMMONDATAMODEL_EXPORT VTK_MARSHALMANUAL vtkCellArray : public vtkAbstractCellArray
177{
178public:
179 using ArrayType32 VTK_DEPRECATED_IN_9_6_0("Use AOSArray32 instead.") = vtkTypeInt32Array;
180 using ArrayType64 VTK_DEPRECATED_IN_9_6_0("Use AOSArray64 instead.") = vtkTypeInt64Array;
185
187
191 static vtkCellArray* New();
193 void PrintSelf(ostream& os, vtkIndent indent) override;
194 void PrintDebug(ostream& os);
196
198
207 "Use vtkArrayDispatch::OffsetsArrays/ConnectivityArrays instead.") =
208 vtkTypeList::Create<vtkTypeInt32Array, vtkTypeInt64Array>;
210
212
224
232 */
233 VTK_DEPRECATED_IN_9_6_0("Use AllocateEstimate or AllocateExact instead.")
234 vtkTypeBool Allocate(vtkIdType sz, vtkIdType vtkNotUsed(ext) = 1000)
235 {
236 return this->AllocateExact(sz, sz) ? 1 : 0;
237 }
238
246 * @sa Squeeze AllocateExact AllocateCopy
247 */
248 bool AllocateEstimate(vtkIdType numCells, vtkIdType maxCellSize)
249 {
250 return this->AllocateExact(numCells, numCells * maxCellSize);
251 }
252
262 bool AllocateExact(vtkIdType numCells, vtkIdType connectivitySize);
263
271 * @sa Squeeze AllocateEstimate AllocateExact
272 */
273 bool AllocateCopy(vtkCellArray* other)
274 {
275 return this->AllocateExact(other->GetNumberOfCells(), other->GetNumberOfConnectivityIds());
276 }
277
287 bool ResizeExact(vtkIdType numCells, vtkIdType connectivitySize);
288
292 void Initialize() override;
293
297 void Reset();
298
304 void Squeeze();
305
316 bool IsValid();
317
321 vtkIdType GetNumberOfCells() const override { return this->Offsets->GetNumberOfValues() - 1; }
322
327 vtkIdType GetNumberOfOffsets() const override { return this->Offsets->GetNumberOfValues(); }
328
330 * Get the offset (into the connectivity) for a specified cell id.
331 */
332 vtkIdType GetOffset(vtkIdType cellId) override
333 {
334 return static_cast<vtkIdType>(this->Offsets->GetComponent(cellId, 0));
335 }
336
338 * Set the offset (into the connectivity) for a specified cell id.
339 */
340 void SetOffset(vtkIdType cellId, vtkIdType offset)
341 {
342 this->Offsets->SetComponent(cellId, 0, static_cast<double>(offset));
343 }
344
346 * Get the size of the connectivity array that stores the point ids.
347 */
349 {
350 return this->Connectivity->GetNumberOfValues();
351 }
352
359
361
374 void SetData(AOSArray32*, AOSArray32* connectivity);
375 void SetData(AOSArray64*, AOSArray64* connectivity);
376 void SetData(AffineArray32*, AOSArray32* connectivity);
377 void SetData(AffineArray64*, AOSArray64* connectivity);
379
394 bool SetData(vtkDataArray* offsets, vtkDataArray* connectivity);
395
421 Generic,
422 };
423
427 StorageTypes GetStorageType() const noexcept { return this->StorageType; }
428
432 bool IsStorage64Bit() const { return this->StorageType == StorageTypes::Int64; }
433
437 bool IsStorage32Bit() const { return this->StorageType == StorageTypes::Int32; }
438
443
448
450 * @return True if the internal storage is using FixedSize arrays.
451 */
452 bool IsStorageFixedSize() const
453 {
454 return this->IsStorageFixedSize32Bit() || this->IsStorageFixedSize64Bit();
455 }
456
460 bool IsStorageGeneric() const { return this->StorageType == StorageTypes::Generic; }
461
466 * a pointer to vtkIdType can be used.
467 */
468 bool IsStorageShareable() const override
469 {
470 switch (this->StorageType)
471 {
474 return std::is_same_v<vtkTypeInt32, vtkIdType>;
477 return std::is_same_v<vtkTypeInt64, vtkIdType>;
479 default:
480 return false;
481 }
482 }
483
496 void UseFixedSize64BitStorage(vtkIdType cellSize);
499
515 {
516 switch (type)
517 {
519 return this->CanConvertTo32BitStorage();
521 return this->CanConvertTo64BitStorage();
527 default:
528 return true;
529 }
530 }
532
557 {
558 switch (type)
559 {
561 return this->ConvertTo32BitStorage();
563 return this->ConvertTo64BitStorage();
565 return this->ConvertToFixedSize32BitStorage();
567 return this->ConvertToFixedSize64BitStorage();
569 default:
570 return true;
571 }
572 }
574
580 vtkDataArray* GetOffsetsArray() const { return this->Offsets; }
581 VTK_DEPRECATED_IN_9_6_0("Use GetOffsetsAOSArray32() instead.")
582 vtkTypeInt32Array* GetOffsetsArray32() const
583 {
584 return vtkTypeInt32Array::FastDownCast(this->Offsets);
585 }
587 VTK_DEPRECATED_IN_9_6_0("Use GetOffsetsAOSArray64() instead.")
588 vtkTypeInt64Array* GetOffsetsArray64() const
589 {
590 return vtkTypeInt64Array::FastDownCast(this->Offsets);
592 AOSArray64* GetOffsetsAOSArray64() const { return AOSArray64::FastDownCast(this->Offsets); }
593 AffineArray32* GetOffsetsAffineArray32() const
594 {
596 }
597 AffineArray64* GetOffsetsAffineArray64() const
598 {
599 return AffineArray64::FastDownCast(this->Offsets);
600 }
602
610 VTK_DEPRECATED_IN_9_6_0("Use GetConnectivityAOSArray32() instead.")
611 vtkTypeInt32Array* GetConnectivityArray32() const
612 {
613 return vtkTypeInt32Array::FastDownCast(this->Connectivity);
614 }
615 AOSArray32* GetConnectivityAOSArray32() const
616 {
617 return AOSArray32::FastDownCast(this->Connectivity);
619 VTK_DEPRECATED_IN_9_6_0("Use GetConnectivityAOSArray64() instead.")
620 vtkTypeInt64Array* GetConnectivityArray64() const
621 {
622 return vtkTypeInt64Array::FastDownCast(this->Connectivity);
623 }
624 AOSArray64* GetConnectivityAOSArray64() const
625 {
626 return AOSArray64::FastDownCast(this->Connectivity);
627 }
629
638 vtkIdType IsHomogeneous() const override;
639
649 void InitTraversal();
650
665 int GetNextCell(vtkIdType& npts, vtkIdType const*& pts) VTK_SIZEHINT(pts, npts);
666
677 int GetNextCell(vtkIdList* pts);
678
689 void GetCellAtId(vtkIdType cellId, vtkIdType& cellSize, vtkIdType const*& cellPoints,
690 vtkIdList* ptIds) VTK_SIZEHINT(cellPoints, cellSize)
691 VTK_EXPECTS(0 <= cellId && cellId < GetNumberOfCells()) override;
692
698 void GetCellAtId(vtkIdType cellId, vtkIdList* pts)
699 VTK_EXPECTS(0 <= cellId && cellId < GetNumberOfCells()) override;
700
708 void GetCellAtId(vtkIdType cellId, vtkIdType& cellSize, vtkIdType* cellPoints) VTK_SIZEHINT(
709 cellPoints, cellSize) VTK_EXPECTS(0 <= cellId && cellId < GetNumberOfCells()) override;
710
714 vtkIdType GetCellPointAtId(vtkIdType cellId, vtkIdType cellPointIndex) const
715 VTK_EXPECTS(0 <= cellId && cellId < GetNumberOfCells())
716 VTK_EXPECTS(0 <= cellPointIndex && cellPointIndex < this->GetCellSize(cellId));
717
721 vtkIdType GetCellSize(vtkIdType cellId) const override;
722
727
732 vtkIdType InsertNextCell(vtkIdType npts, const vtkIdType* pts) VTK_SIZEHINT(pts, npts);
733
739
745 * `InsertNextCell(int)` overload instead.
746 */
747 vtkIdType InsertNextCell(const std::initializer_list<vtkIdType>& cell)
748 {
749 return this->InsertNextCell(static_cast<vtkIdType>(cell.size()), cell.begin());
750 }
751
758 vtkIdType InsertNextCell(int npts);
759
764 void InsertCellPoint(vtkIdType id);
765
770 void UpdateCellCount(int npts);
771
780 void SetTraversalCellId(vtkIdType cellId);
782
786 void ReverseCellAtId(vtkIdType cellId) VTK_EXPECTS(0 <= cellId && cellId < GetNumberOfCells());
787
796 void ReplaceCellAtId(vtkIdType cellId, vtkIdList* list);
797 void ReplaceCellAtId(vtkIdType cellId, vtkIdType cellSize, const vtkIdType* cellPoints)
798 VTK_EXPECTS(0 <= cellId && cellId < GetNumberOfCells()) VTK_SIZEHINT(cellPoints, cellSize);
800
808 void ReplaceCellPointAtId(vtkIdType cellId, vtkIdType cellPointIndex, vtkIdType newPointId);
809
815 * results.
816 */
817 void ReplaceCellAtId(vtkIdType cellId, const std::initializer_list<vtkIdType>& cell)
818 {
819 this->ReplaceCellAtId(cellId, static_cast<vtkIdType>(cell.size()), cell.begin());
820 }
821
826 int GetMaxCellSize() override;
827
831 void DeepCopy(vtkAbstractCellArray* ca) override;
832
836 void ShallowCopy(vtkAbstractCellArray* ca) override;
837
841 void Append(vtkCellArray* src, vtkIdType pointOffset = 0);
842
854
867 void ImportLegacyFormat(const vtkIdType* data, vtkIdType len) VTK_SIZEHINT(data, len);
869
881 void AppendLegacyFormat(vtkIdTypeArray* data, vtkIdType ptOffset = 0);
882 void AppendLegacyFormat(const vtkIdType* data, vtkIdType len, vtkIdType ptOffset = 0)
883 VTK_SIZEHINT(data, len);
885
894 unsigned long GetActualMemorySize() const;
895
896 // The following code is used to support
897
898 // The wrappers get understandably confused by some of the template code below
899#ifndef __VTK_WRAP__
900 /*
901 * Utilities class that every dispatch functor used with Dispatch() must inherit
902 * or optionally use its static methods.
903 */
904 struct DispatchUtilities
906 template <class ArrayT>
909 template <class OffsetsT>
910 vtkIdType GetNumberOfCells(OffsetsT* offsets)
911 {
912 return offsets->GetNumberOfValues() - 1;
913 }
915 template <class ArrayT>
916 static decltype(vtk::DataArrayValueRange<1, vtkIdType>(std::declval<ArrayT>())) GetRange(
917 ArrayT* array)
918 {
920 }
922 template <class OffsetsT>
923 static vtkIdType GetBeginOffset(OffsetsT* offsets, vtkIdType cellId)
924 {
925 return static_cast<vtkIdType>(GetRange(offsets)[cellId]);
926 }
928 template <class OffsetsT>
929 static vtkIdType GetEndOffset(OffsetsT* offsets, vtkIdType cellId)
930 {
931 return static_cast<vtkIdType>(GetRange(offsets)[cellId + 1]);
932 }
934 template <class OffsetsT>
935 static vtkIdType GetCellSize(OffsetsT* offsets, vtkIdType cellId)
936 {
937 auto offsetsRange = GetRange(offsets);
938 return static_cast<vtkIdType>(offsetsRange[cellId + 1] - offsetsRange[cellId]);
939 }
940
941 template <class OffsetsT, class ConnectivityT>
942 static decltype(vtk::DataArrayValueRange<1, vtkIdType>(std::declval<ConnectivityT>()))
943 GetCellRange(OffsetsT* offsets, ConnectivityT* conn, vtkIdType cellId)
944 {
945 auto offsetsRange = GetRange(offsets);
947 conn, offsetsRange[cellId], offsetsRange[cellId + 1]);
948 }
949 };
950
952
1015 template <typename Functor, typename... Args>
1016 void Dispatch(Functor&& functor, Args&&... args)
1017 {
1018 switch (this->StorageType)
1019 {
1021 functor(static_cast<AOSArray32*>(this->Offsets.Get()),
1022 static_cast<AOSArray32*>(this->Connectivity.Get()), std::forward<Args>(args)...);
1023 break;
1025 functor(static_cast<AOSArray64*>(this->Offsets.Get()),
1026 static_cast<AOSArray64*>(this->Connectivity.Get()), std::forward<Args>(args)...);
1027 break;
1029 functor(static_cast<AffineArray32*>(this->Offsets.Get()),
1030 static_cast<AOSArray32*>(this->Connectivity.Get()), std::forward<Args>(args)...);
1031 break;
1033 functor(static_cast<AffineArray64*>(this->Offsets.Get()),
1034 static_cast<AOSArray64*>(this->Connectivity.Get()), std::forward<Args>(args)...);
1035 break;
1037 default:
1038 functor(this->Offsets.Get(), this->Connectivity.Get(), std::forward<Args>(args)...);
1039 break;
1040 }
1042 template <typename Functor, typename... Args>
1043 void Dispatch(Functor&& functor, Args&&... args) const
1044 {
1045 switch (this->StorageType)
1046 {
1047 case StorageTypes::Int32:
1048 functor(static_cast<AOSArray32*>(this->Offsets.Get()),
1049 static_cast<AOSArray32*>(this->Connectivity.Get()), std::forward<Args>(args)...);
1050 break;
1051 case StorageTypes::Int64:
1052 functor(static_cast<AOSArray64*>(this->Offsets.Get()),
1053 static_cast<AOSArray64*>(this->Connectivity.Get()), std::forward<Args>(args)...);
1054 break;
1055 case StorageTypes::FixedSizeInt32:
1056 functor(static_cast<AffineArray32*>(this->Offsets.Get()),
1057 static_cast<AOSArray32*>(this->Connectivity.Get()), std::forward<Args>(args)...);
1058 break;
1059 case StorageTypes::FixedSizeInt64:
1060 functor(static_cast<AffineArray64*>(this->Offsets.Get()),
1061 static_cast<AOSArray64*>(this->Connectivity.Get()), std::forward<Args>(args)...);
1062 break;
1063 case StorageTypes::Generic:
1064 default:
1065 functor(this->Offsets.Get(), this->Connectivity.Get(), std::forward<Args>(args)...);
1066 break;
1067 }
1068 }
1070
1071 // Holds connectivity and offset arrays of the given ArrayType.
1072 // VTK_DEPRECATED_IN_9_6_0("Use DispatchUtilities")
1073 template <typename ArrayT>
1076 using ArrayType = ArrayT;
1077 using ValueType = typename ArrayType::ValueType;
1078 using CellRangeType = decltype(vtk::DataArrayValueRange<1>(std::declval<ArrayType>()));
1079
1080 // We can't just use is_same here, since binary compatible representations
1081 // (e.g. int and long) are distinct types. Instead, ensure that ValueType
1082 // is a signed integer the same size as vtkIdType.
1083 // If this value is true, ValueType pointers may be safely converted to
1084 // vtkIdType pointers via reinterpret cast.
1085 static constexpr bool ValueTypeIsSameAsIdType = std::is_integral<ValueType>::value &&
1086 std::is_signed<ValueType>::value && (sizeof(ValueType) == sizeof(vtkIdType));
1088 ArrayType* GetOffsets() { return this->Offsets; }
1089 const ArrayType* GetOffsets() const { return this->Offsets; }
1091 ArrayType* GetConnectivity() { return this->Connectivity; }
1092 const ArrayType* GetConnectivity() const { return this->Connectivity; }
1093
1094 vtkIdType GetNumberOfCells() const;
1095
1096 vtkIdType GetBeginOffset(vtkIdType cellId) const;
1097
1098 vtkIdType GetEndOffset(vtkIdType cellId) const;
1099
1100 vtkIdType GetCellSize(vtkIdType cellId) const;
1101
1103
1104 friend class vtkCellArray;
1106 protected:
1107 VisitState()
1108 {
1111 this->Offsets->InsertNextValue(0);
1113 {
1114 this->IsInMemkind = true;
1117 ~VisitState() = default;
1118 void* operator new(size_t nSize)
1119 {
1120 void* r;
1121#ifdef VTK_USE_MEMKIND
1123#else
1124 r = malloc(nSize);
1125#endif
1126 return r;
1127 }
1128 void operator delete(void* p)
1129 {
1130#ifdef VTK_USE_MEMKIND
1131 VisitState* a = static_cast<VisitState*>(p);
1132 if (a->IsInMemkind)
1133 {
1135 }
1136 else
1137 {
1138 free(p);
1139 }
1140#else
1141 free(p);
1142#endif
1147
1148 private:
1149 VisitState(const VisitState&) = delete;
1150 VisitState& operator=(const VisitState&) = delete;
1151 bool IsInMemkind = false;
1152 };
1153
1154private: // Helpers that allow Visit to return a value:
1155 template <typename Functor, typename... Args>
1156 using GetReturnType = decltype(std::declval<Functor>()(
1157 std::declval<VisitState<AOSArray32>&>(), std::declval<Args>()...));
1158
1159 template <typename Functor, typename... Args>
1160 struct ReturnsVoid : std::is_same<GetReturnType<Functor, Args...>, void>
1161 {
1162 };
1163
1164public:
1234 template <typename Functor, typename... Args,
1235 typename = typename std::enable_if<ReturnsVoid<Functor, Args...>::value>::type>
1236 VTK_DEPRECATED_IN_9_6_0("Use Dispatch instead")
1237 void Visit(Functor&& functor, Args&&... args)
1238 {
1239 switch (this->StorageType)
1240 {
1242 {
1246 functor(state, std::forward<Args>(args)...);
1247 break;
1248 }
1250 {
1254 functor(state, std::forward<Args>(args)...);
1255 break;
1256 }
1260 default:
1261 {
1262 vtkWarningMacro("Use Dispatch");
1263 break;
1264 }
1265 }
1266 }
1267
1268 template <typename Functor, typename... Args,
1269 typename = typename std::enable_if<ReturnsVoid<Functor, Args...>::value>::type>
1270 VTK_DEPRECATED_IN_9_6_0("Use Dispatch instead")
1271 void Visit(Functor&& functor, Args&&... args) const
1272 {
1273 switch (this->StorageType)
1274 {
1276 {
1280 functor(state, std::forward<Args>(args)...);
1281 break;
1282 }
1284 {
1288 functor(state, std::forward<Args>(args)...);
1289 break;
1290 }
1294 default:
1295 {
1296 vtkWarningMacro("Use Dispatch");
1297 break;
1298 }
1299 }
1300 }
1301
1302 template <typename Functor, typename... Args,
1303 typename = typename std::enable_if<!ReturnsVoid<Functor, Args...>::value>::type>
1304 VTK_DEPRECATED_IN_9_6_0("Use Dispatch instead")
1305 GetReturnType<Functor, Args...> Visit(Functor&& functor, Args&&... args)
1306 {
1307 switch (this->StorageType)
1308 {
1310 {
1314 return functor(state, std::forward<Args>(args)...);
1315 }
1317 {
1321 return functor(state, std::forward<Args>(args)...);
1322 }
1326 default:
1327 {
1328 vtkWarningMacro("Use Dispatch");
1329 return GetReturnType<Functor, Args...>();
1330 }
1331 }
1332 }
1333 template <typename Functor, typename... Args,
1334 typename = typename std::enable_if<!ReturnsVoid<Functor, Args...>::value>::type>
1335 VTK_DEPRECATED_IN_9_6_0("Use Dispatch instead")
1336 GetReturnType<Functor, Args...> Visit(Functor&& functor, Args&&... args) const
1337 {
1338 switch (this->StorageType)
1339 {
1341 {
1345 return functor(state, std::forward<Args>(args)...);
1346 }
1348 {
1352 return functor(state, std::forward<Args>(args)...);
1353 }
1357 default:
1358 {
1359 vtkWarningMacro("Use Dispatch");
1360 return GetReturnType<Functor, Args...>();
1361 }
1362 }
1363 }
1364
1365#endif // __VTK_WRAP__
1366
1368
1376 static void SetDefaultStorageIs64Bit(bool val) { vtkCellArray::DefaultStorageIs64Bit = val; }
1378
1379 //=================== Begin Legacy Methods ===================================
1380 // These should be deprecated at some point as they are confusing or very slow
1381
1388 VTK_DEPRECATED_IN_9_6_0("This call has no effect.")
1389 virtual void SetNumberOfCells(vtkIdType);
1390
1402 VTK_DEPRECATED_IN_9_6_0("Use AllocateEstimate directly instead.")
1403 vtkIdType EstimateSize(vtkIdType numCells, int maxPtsPerCell);
1404
1413 VTK_DEPRECATED_IN_9_6_0("Method incompatible with current internal storage.")
1415
1422 VTK_DEPRECATED_IN_9_6_0("Method incompatible with current internal storage.")
1424
1434 VTK_DEPRECATED_IN_9_6_0("Use GetCellAtId.")
1435 void GetCell(vtkIdType loc, vtkIdType& npts, const vtkIdType*& pts)
1436 VTK_EXPECTS(0 <= loc && loc < GetNumberOfConnectivityEntries()) VTK_SIZEHINT(pts, npts);
1437
1444 VTK_DEPRECATED_IN_9_6_0("Use GetCellAtId.")
1445 void GetCell(vtkIdType loc, vtkIdList* pts)
1446 VTK_EXPECTS(0 <= loc && loc < GetNumberOfConnectivityEntries());
1447
1454 VTK_DEPRECATED_IN_9_6_0("Use GetNumberOfCells.")
1455 vtkIdType GetInsertLocation(int npts);
1456
1464 VTK_DEPRECATED_IN_9_6_0("Use GetTraversalCellId.")
1466 VTK_DEPRECATED_IN_9_6_0("Use GetTraversalCellId.")
1468 VTK_DEPRECATED_IN_9_6_0("Use SetTraversalCellId.")
1471
1479 VTK_DEPRECATED_IN_9_6_0("Use ReverseCellAtId.")
1480 void ReverseCell(vtkIdType loc) VTK_EXPECTS(0 <= loc && loc < GetNumberOfConnectivityEntries());
1481
1493 VTK_DEPRECATED_IN_9_6_0("Use ReplaceCellAtId.")
1494 void ReplaceCell(vtkIdType loc, int npts, const vtkIdType pts[])
1495 VTK_EXPECTS(0 <= loc && loc < GetNumberOfConnectivityEntries()) VTK_SIZEHINT(pts, npts);
1496
1511 VTK_DEPRECATED_IN_9_6_0("Use ImportLegacyFormat or SetData instead.")
1512 void SetCells(vtkIdType ncells, vtkIdTypeArray* cells);
1513
1525 "Use ExportLegacyFormat, or GetOffsetsArray/GetConnectivityArray instead.")
1527
1528 //=================== End Legacy Methods =====================================
1529
1530 friend class vtkCellArrayIterator;
1532protected:
1533 vtkCellArray();
1534 ~vtkCellArray() override;
1540
1542
1543 static bool DefaultStorageIs64Bit;
1544
1545private:
1546 vtkCellArray(const vtkCellArray&) = delete;
1547 void operator=(const vtkCellArray&) = delete;
1548};
1549
1550template <typename ArrayT>
1552{
1553 return this->Offsets->GetNumberOfValues() - 1;
1554}
1555
1556template <typename ArrayT>
1558{
1559 return static_cast<vtkIdType>(this->Offsets->GetValue(cellId));
1560}
1561
1562template <typename ArrayT>
1564{
1565 return static_cast<vtkIdType>(this->Offsets->GetValue(cellId + 1));
1566}
1567
1568template <typename ArrayT>
1570{
1571 return this->GetEndOffset(cellId) - this->GetBeginOffset(cellId);
1572}
1573
1574template <typename ArrayT>
1581VTK_ABI_NAMESPACE_END
1582
1584{
1585VTK_ABI_NAMESPACE_BEGIN
1586
1588{
1589 // Insert full cell
1590 template <class OffsetsT, class ConnectivityT>
1591 void operator()(OffsetsT* offsets, ConnectivityT* conn, const vtkIdType npts,
1592 const vtkIdType pts[], vtkIdType& cellId)
1593 {
1594 using ValueType = GetAPIType<OffsetsT>;
1595 using OffsetsAccessorType = vtkDataArrayAccessor<OffsetsT>;
1596 using ConnectivityAccessorType = vtkDataArrayAccessor<ConnectivityT>;
1597 OffsetsAccessorType offsetsAccesor(offsets);
1598 ConnectivityAccessorType connAccesor(conn);
1599
1600 cellId = offsets->GetNumberOfValues() - 1;
1601
1602 offsetsAccesor.InsertNext(static_cast<ValueType>(conn->GetNumberOfValues() + npts));
1603
1604 for (vtkIdType i = 0; i < npts; ++i)
1605 {
1606 connAccesor.InsertNext(static_cast<ValueType>(pts[i]));
1607 }
1608 }
1609
1610 // Just update offset table (for incremental API)
1611 template <class OffsetsT, class ConnectivityT>
1612 void operator()(OffsetsT* offsets, ConnectivityT* conn, const vtkIdType npts, vtkIdType& cellId)
1613 {
1614 using ValueType = GetAPIType<OffsetsT>;
1615 using AccessorType = vtkDataArrayAccessor<OffsetsT>;
1616 AccessorType offsetsAccesor(offsets);
1617
1618 cellId = offsets->GetNumberOfValues() - 1;
1619
1620 offsetsAccesor.InsertNext(static_cast<ValueType>(conn->GetNumberOfValues() + npts));
1621 }
1622};
1623
1624// for incremental API:
1626{
1627 template <class OffsetsT, class ConnectivityT>
1628 void operator()(OffsetsT* offsets, ConnectivityT* vtkNotUsed(conn), const vtkIdType npts)
1629 {
1630 using ValueType = GetAPIType<OffsetsT>;
1631
1632 auto offsetsRange = GetRange(offsets);
1633 const ValueType cellBegin = offsetsRange[offsets->GetMaxId() - 1];
1634 offsetsRange[offsets->GetMaxId()] = static_cast<ValueType>(cellBegin + npts);
1635 }
1636};
1637
1639{
1640 template <class OffsetsT, class ConnectivityT>
1642 OffsetsT* offsets, ConnectivityT* vtkNotUsed(conn), vtkIdType cellId, vtkIdType& cellSize)
1643 {
1644 cellSize = GetCellSize(offsets, cellId);
1645 }
1646};
1647
1649{
1650 template <class OffsetsT, class ConnectivityT>
1651 void operator()(OffsetsT* offsets, ConnectivityT* conn, const vtkIdType cellId, vtkIdList* ids)
1652 {
1653 auto offsetsRange = GetRange(offsets);
1654 const auto& beginOffset = offsetsRange[cellId];
1655 const auto& endOffset = offsetsRange[cellId + 1];
1656 const vtkIdType cellSize = static_cast<vtkIdType>(endOffset - beginOffset);
1657 const auto cellConnectivity = GetRange(conn).begin() + beginOffset;
1658
1659 // ValueType differs from vtkIdType, so we have to copy into a temporary buffer:
1660 ids->SetNumberOfIds(cellSize);
1661 vtkIdType* idPtr = ids->GetPointer(0);
1662 for (vtkIdType i = 0; i < cellSize; ++i)
1663 {
1664 idPtr[i] = static_cast<vtkIdType>(cellConnectivity[i]);
1665 }
1666 }
1667
1668 template <class OffsetsT, class ConnectivityT>
1669 void operator()(OffsetsT* offsets, ConnectivityT* conn, const vtkIdType cellId,
1670 vtkIdType& cellSize, vtkIdType* cellPoints)
1671 {
1672 auto offsetsRange = GetRange(offsets);
1673 const auto& beginOffset = offsetsRange[cellId];
1674 const auto& endOffset = offsetsRange[cellId + 1];
1675 cellSize = static_cast<vtkIdType>(endOffset - beginOffset);
1676 const auto cellConnectivity = GetRange(conn).begin() + beginOffset;
1677
1678 // ValueType differs from vtkIdType, so we have to copy into a temporary buffer:
1679 for (vtkIdType i = 0; i < cellSize; ++i)
1680 {
1681 cellPoints[i] = static_cast<vtkIdType>(cellConnectivity[i]);
1682 }
1683 }
1684
1685 // SFINAE helper to check if a Functors's connectivity array's memory can be used as a vtkIdType*.
1686 template <typename ConnectivityT>
1688 {
1689 static constexpr bool value =
1690 std::is_base_of_v<vtkAOSDataArrayTemplate<vtkIdType>, ConnectivityT>;
1691 };
1692
1693 template <class OffsetsT, class ConnectivityT>
1694 typename std::enable_if_t<CanShareConnPtr<ConnectivityT>::value, void> operator()(
1695 OffsetsT* offsets, ConnectivityT* conn, const vtkIdType cellId, vtkIdType& cellSize,
1696 vtkIdType const*& cellPoints, vtkIdList* vtkNotUsed(temp))
1697 {
1698 auto offsetsRange = GetRange(offsets);
1699 const auto& beginOffset = offsetsRange[cellId];
1700 const auto& endOffset = offsetsRange[cellId + 1];
1701 cellSize = static_cast<vtkIdType>(endOffset - beginOffset);
1702 // This is safe, see CanShareConnPtr helper above.
1703 cellPoints = conn->GetPointer(beginOffset);
1704 }
1705
1706 template <class OffsetsT, class ConnectivityT>
1707 typename std::enable_if_t<!CanShareConnPtr<ConnectivityT>::value, void> operator()(
1708 OffsetsT* offsets, ConnectivityT* conn, const vtkIdType cellId, vtkIdType& cellSize,
1709 vtkIdType const*& cellPoints, vtkIdList* temp)
1710 {
1711 auto offsetsRange = GetRange(offsets);
1712 const auto& beginOffset = offsetsRange[cellId];
1713 const auto& endOffset = offsetsRange[cellId + 1];
1714 cellSize = static_cast<vtkIdType>(endOffset - beginOffset);
1715 const auto cellConnectivity = GetRange(conn).begin() + beginOffset;
1716
1717 temp->SetNumberOfIds(cellSize);
1718 vtkIdType* tempPtr = temp->GetPointer(0);
1719 for (vtkIdType i = 0; i < cellSize; ++i)
1720 {
1721 tempPtr[i] = static_cast<vtkIdType>(cellConnectivity[i]);
1722 }
1723
1724 cellPoints = tempPtr;
1725 }
1726};
1727
1729{
1730 template <class OffsetsT, class ConnectivityT>
1731 void operator()(OffsetsT* offsets, ConnectivityT* conn, vtkIdType cellId,
1732 vtkIdType cellPointIndex, vtkIdType& pointId)
1733 {
1734 pointId =
1735 static_cast<vtkIdType>(GetRange(conn)[GetBeginOffset(offsets, cellId) + cellPointIndex]);
1736 }
1737};
1738
1740{
1741 template <class OffsetsT, class ConnectivityT>
1742 void operator()(OffsetsT* offsets, ConnectivityT* conn)
1743 {
1744 using ValueType = GetAPIType<OffsetsT>;
1745 using AccessorType = vtkDataArrayAccessor<OffsetsT>;
1746 offsets->Reset();
1747 conn->Reset();
1748 AccessorType accessor(offsets);
1749 ValueType firstOffset = 0;
1750 accessor.InsertNext(firstOffset);
1751 }
1752};
1753
1755{
1756 template <class OffsetsT, class ConnectivityT>
1757 void operator()(OffsetsT* vtkNotUsed(offsets), ConnectivityT* conn, vtkIdType id)
1758 {
1759 using ValueType = GetAPIType<ConnectivityT>;
1760 using AccessorType = vtkDataArrayAccessor<ConnectivityT>;
1761 AccessorType accessor(conn);
1762 accessor.InsertNext(static_cast<ValueType>(id));
1763 }
1764};
1765
1766VTK_ABI_NAMESPACE_END
1767} // end namespace vtkCellArray_detail
1768
1769VTK_ABI_NAMESPACE_BEGIN
1770//----------------------------------------------------------------------------
1772{
1773 this->TraversalCellId = 0;
1774}
1775
1776//----------------------------------------------------------------------------
1777inline int vtkCellArray::GetNextCell(vtkIdType& npts, vtkIdType const*& pts) VTK_SIZEHINT(pts, npts)
1778{
1779 if (this->TraversalCellId < this->GetNumberOfCells())
1780 {
1781 this->GetCellAtId(this->TraversalCellId, npts, pts);
1782 ++this->TraversalCellId;
1783 return 1;
1784 }
1785
1786 npts = 0;
1787 pts = nullptr;
1788 return 0;
1789}
1790
1791//----------------------------------------------------------------------------
1793{
1795 {
1796 this->GetCellAtId(this->TraversalCellId, pts);
1797 ++this->TraversalCellId;
1798 return 1;
1799 }
1800
1801 pts->Reset();
1802 return 0;
1803}
1804//----------------------------------------------------------------------------
1806{
1807 vtkIdType cellSize;
1808 this->Dispatch(vtkCellArray_detail::GetCellSizeImpl{}, cellId, cellSize);
1809 return cellSize;
1810}
1811
1812//----------------------------------------------------------------------------
1813inline void vtkCellArray::GetCellAtId(vtkIdType cellId, vtkIdType& cellSize,
1814 vtkIdType const*& cellPoints, vtkIdList* ptIds) VTK_SIZEHINT(cellPoints, cellSize)
1815{
1816 this->Dispatch(vtkCellArray_detail::GetCellAtIdImpl{}, cellId, cellSize, cellPoints, ptIds);
1817}
1818
1819//----------------------------------------------------------------------------
1821{
1822 this->Dispatch(vtkCellArray_detail::GetCellAtIdImpl{}, cellId, pts);
1823}
1824
1825//----------------------------------------------------------------------------
1826inline void vtkCellArray::GetCellAtId(vtkIdType cellId, vtkIdType& cellSize, vtkIdType* cellPoints)
1827{
1828 this->Dispatch(vtkCellArray_detail::GetCellAtIdImpl{}, cellId, cellSize, cellPoints);
1829}
1830
1831//----------------------------------------------------------------------------
1833{
1834 vtkIdType pointId;
1835 this->Dispatch(vtkCellArray_detail::CellPointAtIdImpl{}, cellId, cellPointIndex, pointId);
1836 return pointId;
1837}
1838
1839//----------------------------------------------------------------------------
1841 VTK_SIZEHINT(pts, npts)
1842{
1843 vtkIdType cellId;
1844 this->Dispatch(vtkCellArray_detail::InsertNextCellImpl{}, npts, pts, cellId);
1845 return cellId;
1846}
1847
1848//----------------------------------------------------------------------------
1850{
1851 vtkIdType cellId;
1853 return cellId;
1854}
1855
1856//----------------------------------------------------------------------------
1861
1862//----------------------------------------------------------------------------
1864{
1866}
1867
1868//----------------------------------------------------------------------------
1870{
1871 vtkIdType cellId;
1872 this->Dispatch(
1874 return cellId;
1875}
1876
1877//----------------------------------------------------------------------------
1879{
1880 vtkIdList* pts = cell->GetPointIds();
1881 vtkIdType cellId;
1882 this->Dispatch(
1884 return cellId;
1885}
1886
1887//----------------------------------------------------------------------------
1889{
1891}
1892
1893VTK_ABI_NAMESPACE_END
1894#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.
vtkTypeList::Unique< vtkTypeList::Create< vtkAOSDataArrayTemplate< int >, vtkAOSDataArrayTemplate< long >, vtkAOSDataArrayTemplate< long long > > >::Result InputArrayList
List of possible ArrayTypes that are compatible with internal storage.
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.
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.
vtkTypeList::Create< vtkTypeInt32Array, vtkTypeInt64Array > StorageArrayList
List of possible array types used for 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
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:50
vtkIdList * GetPointIds()
Return the list of point ids defining the cell.
Definition vtkCell.h:143
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:24
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:50
void Reset()
Reset to an empty state but retain previously allocated memory.
Definition vtkIdList.h:135
vtkIdType * GetPointer(vtkIdType i)
Get a pointer to a particular data index.
Definition vtkIdList.h:116
dynamic, self-adjusting array of vtkIdType
static vtkImplicitArray< vtkAffineImplicitBackend< T >, ArrayType > * FastDownCast(vtkAbstractArray *source)
a simple class to control print indentation
Definition vtkIndent.h:29
Allocate and hold a VTK object.
Definition vtkNew.h:58
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.
void Reset()
Initializes the field list to empty.
#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