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