VTK  9.4.20250206
vtkStaticCellLinksTemplate.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
36#ifndef vtkStaticCellLinksTemplate_h
37#define vtkStaticCellLinksTemplate_h
38
39#include "vtkABINamespace.h"
40#include "vtkDeprecation.h" // For VTK_DEPRECATED_IN_9_5_0
41
42#include <memory> // For shared_ptr
43#include <vector> // For vector
44
45VTK_ABI_NAMESPACE_BEGIN
46class vtkDataSet;
47class vtkPolyData;
50class vtkCellArray;
51VTK_ABI_NAMESPACE_END
52
54
55VTK_ABI_NAMESPACE_BEGIN
56template <typename TIds>
58{
59public:
61
67
71 void Initialize();
72
78
83
88
93
95
98 VTK_DEPRECATED_IN_9_5_0("Use BuildLinksFromMultipleArrays instead.")
100 vtkIdType numPts, vtkIdType numCells, std::vector<vtkCellArray*> cellArrays)
101 {
102 this->BuildLinksFromMultipleArrays(numPts, numCells, cellArrays);
103 }
104 VTK_DEPRECATED_IN_9_5_0("Use BuildLinks instead.")
105 void SerialBuildLinks(vtkIdType numPts, vtkIdType numCells, vtkCellArray* cellArray)
106 {
107 this->BuildLinksFromMultipleArrays(numPts, numCells, { cellArray });
108 }
109 VTK_DEPRECATED_IN_9_5_0("Use BuildLinksFromMultipleArrays instead.")
111 vtkIdType numPts, vtkIdType numCells, std::vector<vtkCellArray*> cellArrays)
112 {
113 this->BuildLinksFromMultipleArrays(numPts, numCells, cellArrays);
114 }
115 VTK_DEPRECATED_IN_9_5_0("Use BuildLinks instead.")
116 void ThreadedBuildLinks(vtkIdType numPts, vtkIdType numCells, vtkCellArray* cellArray)
117 {
118 this->BuildLinksFromMultipleArrays(numPts, numCells, { cellArray });
119 }
121 vtkIdType numPts, vtkIdType numCells, std::vector<vtkCellArray*> cellArrays);
122 void BuildLinks(vtkIdType numPts, vtkIdType numCells, vtkCellArray* cellArray)
123 {
124 this->BuildLinksFromMultipleArrays(numPts, numCells, { cellArray });
125 }
127
129
132 TIds GetNumberOfCells(vtkIdType ptId) { return (this->Offsets[ptId + 1] - this->Offsets[ptId]); }
133 vtkIdType GetNcells(vtkIdType ptId) { return (this->Offsets[ptId + 1] - this->Offsets[ptId]); }
135
140 template <typename TGivenIds>
141 bool MatchesCell(TGivenIds npts, const TGivenIds* pts);
142
146 TIds* GetCells(vtkIdType ptId) { return (this->Links + this->Offsets[ptId]); }
147
152 void GetCells(vtkIdType npts, const vtkIdType* pts, vtkIdList* cells);
153
158 TIds GetLinksSize() { return this->LinksSize; }
159
164 TIds GetOffset(vtkIdType ptId) { return this->Offsets[ptId]; }
165
167
170 unsigned long GetActualMemorySize();
171 VTK_DEPRECATED_IN_9_5_0("Use DeepCopy(vtkStaticCellLinksTemplate instead.")
175 void SelectCells(vtkIdType minMaxDegree[2], unsigned char* cellSelection);
177
179
182 VTK_DEPRECATED_IN_9_5_0("No longer used.")
184 VTK_DEPRECATED_IN_9_5_0("No longer used.")
187
188protected:
189 // The various templated data members
191 TIds NumPts;
193
194 // These point to the core data structures
195
196 std::shared_ptr<TIds> LinkSharedPtr; // contiguous runs of cell ids
197 TIds* Links; // Pointer to the links array
198 std::shared_ptr<TIds> OffsetsSharedPtr; // offsets for each point into the links array
199 TIds* Offsets; // Pointer to the offsets array
200
201 // Support for execution
202 int Type;
203 // VTK_DEPRECATED_IN_9_5_0("No longer used.")
205
206private:
208 void operator=(const vtkStaticCellLinksTemplate&) = delete;
209};
210
211VTK_ABI_NAMESPACE_END
212#include "vtkStaticCellLinksTemplate.txx"
213
214#endif
215// VTK-HeaderTest-Exclude: vtkStaticCellLinksTemplate.h
object to represent cell connectivity
abstract class to specify dataset behavior
Definition vtkDataSet.h:165
structured grid with explicit topology and geometry
list of point or cell ids
Definition vtkIdList.h:133
concrete dataset represents vertices, lines, polygons, and triangle strips
object represents upward pointers from points to list of cells using each point (template implementat...
TIds GetLinksSize()
Return the total number of links represented after the links have been built.
void ShallowCopy(vtkStaticCellLinksTemplate *src)
Support vtkAbstractCellLinks API.
vtkStaticCellLinksTemplate()
Instantiate and destructor methods.
void BuildLinks(vtkUnstructuredGrid *ugrid)
Build the link list array for vtkUnstructuredGrid.
vtkIdType GetNcells(vtkIdType ptId)
Get the number of cells using the point specified by ptId.
void SerialBuildLinksFromMultipleArrays(vtkIdType numPts, vtkIdType numCells, std::vector< vtkCellArray * > cellArrays)
Specialized methods for building links from cell array(S).
void DeepCopy(vtkStaticCellLinksTemplate *src)
Support vtkAbstractCellLinks API.
void Initialize()
Make sure any previously created links are cleaned up.
void BuildLinks(vtkExplicitStructuredGrid *esgrid)
Build the link list array for vtkExplicitStructuredGrid.
TIds GetNumberOfCells(vtkIdType ptId)
Get the number of cells using the point specified by ptId.
void ThreadedBuildLinksFromMultipleArrays(vtkIdType numPts, vtkIdType numCells, std::vector< vtkCellArray * > cellArrays)
Specialized methods for building links from cell array(S).
TIds * GetCells(vtkIdType ptId)
Return a list of cell ids using the point specified by ptId.
void BuildLinksFromMultipleArrays(vtkIdType numPts, vtkIdType numCells, std::vector< vtkCellArray * > cellArrays)
Specialized methods for building links from cell array(S).
void BuildLinks(vtkDataSet *ds)
Build the link list array for a general dataset.
void BuildLinks(vtkIdType numPts, vtkIdType numCells, vtkCellArray *cellArray)
Specialized methods for building links from cell array(S).
std::shared_ptr< TIds > OffsetsSharedPtr
~vtkStaticCellLinksTemplate()
Instantiate and destructor methods.
void BuildLinks(vtkPolyData *pd)
Build the link list array for vtkPolyData.
void SetSequentialProcessing(vtkTypeBool seq)
Control whether to thread or serial process.
void DeepCopy(vtkAbstractCellLinks *)
Support vtkAbstractCellLinks API.
unsigned long GetActualMemorySize()
Support vtkAbstractCellLinks API.
bool MatchesCell(TGivenIds npts, const TGivenIds *pts)
Indicate whether the point ids provided defines at least one cell, or a portion of a cell.
void SerialBuildLinks(vtkIdType numPts, vtkIdType numCells, vtkCellArray *cellArray)
Specialized methods for building links from cell array(S).
void GetCells(vtkIdType npts, const vtkIdType *pts, vtkIdList *cells)
Given point ids that define a cell, find the cells that contains all of these point ids.
vtkTypeBool GetSequentialProcessing()
Control whether to thread or serial process.
void SelectCells(vtkIdType minMaxDegree[2], unsigned char *cellSelection)
Support vtkAbstractCellLinks API.
void ThreadedBuildLinks(vtkIdType numPts, vtkIdType numCells, vtkCellArray *cellArray)
Specialized methods for building links from cell array(S).
TIds GetOffset(vtkIdType ptId)
Obtain the offsets into the internal links array.
dataset represents arbitrary combinations of all possible cell types
int vtkTypeBool
Definition vtkABI.h:64
#define VTK_DEPRECATED_IN_9_5_0(reason)
int vtkIdType
Definition vtkType.h:315