VTK  9.3.20240328
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_4_0
41 
42 #include <memory> // For shared_ptr
43 #include <vector> // For vector
44 
45 VTK_ABI_NAMESPACE_BEGIN
46 class vtkDataSet;
47 class vtkPolyData;
50 class vtkCellArray;
51 VTK_ABI_NAMESPACE_END
52 
53 #include "vtkAbstractCellLinks.h"
54 
55 VTK_ABI_NAMESPACE_BEGIN
56 template <typename TIds>
58 {
59 public:
61 
67 
71  void Initialize();
72 
78 
83 
88 
93 
95 
99  vtkIdType numPts, vtkIdType numCells, std::vector<vtkCellArray*> cellArrays);
100  void SerialBuildLinks(vtkIdType numPts, vtkIdType numCells, vtkCellArray* cellArray)
101  {
102  this->SerialBuildLinksFromMultipleArrays(numPts, numCells, { cellArray });
103  }
105  vtkIdType numPts, vtkIdType numCells, std::vector<vtkCellArray*> cellArrays);
106  void ThreadedBuildLinks(vtkIdType numPts, vtkIdType numCells, vtkCellArray* cellArray)
107  {
108  this->ThreadedBuildLinksFromMultipleArrays(numPts, numCells, { cellArray });
109  }
111 
113 
116  TIds GetNumberOfCells(vtkIdType ptId) { return (this->Offsets[ptId + 1] - this->Offsets[ptId]); }
117  vtkIdType GetNcells(vtkIdType ptId) { return (this->Offsets[ptId + 1] - this->Offsets[ptId]); }
119 
124  template <typename TGivenIds>
125  bool MatchesCell(TGivenIds npts, const TGivenIds* pts);
126 
130  TIds* GetCells(vtkIdType ptId) { return (this->Links + this->Offsets[ptId]); }
131 
136  void GetCells(vtkIdType npts, const vtkIdType* pts, vtkIdList* cells);
137 
142  TIds GetLinksSize() { return this->LinksSize; }
143 
148  TIds GetOffset(vtkIdType ptId) { return this->Offsets[ptId]; }
149 
151 
154  unsigned long GetActualMemorySize();
155  VTK_DEPRECATED_IN_9_4_0("Use DeepCopy(vtkStaticCellLinksTemplate instead.")
159  void SelectCells(vtkIdType minMaxDegree[2], unsigned char* cellSelection);
161 
163 
169 
170 protected:
171  // The various templated data members
172  TIds LinksSize;
173  TIds NumPts;
174  TIds NumCells;
175 
176  // These point to the core data structures
177 
178  std::shared_ptr<TIds> LinkSharedPtr; // contiguous runs of cell ids
179  TIds* Links; // Pointer to the links array
180  std::shared_ptr<TIds> OffsetsSharedPtr; // offsets for each point into the links array
181  TIds* Offsets; // Pointer to the offsets array
182 
183  // Support for execution
184  int Type;
186 
187 private:
189  void operator=(const vtkStaticCellLinksTemplate&) = delete;
190 };
191 
192 VTK_ABI_NAMESPACE_END
193 #include "vtkStaticCellLinksTemplate.txx"
194 
195 #endif
196 // VTK-HeaderTest-Exclude: vtkStaticCellLinksTemplate.h
object to represent cell connectivity
Definition: vtkCellArray.h:285
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:132
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:180
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).
void BuildLinks(vtkDataSet *ds)
Build the link list array for a general dataset.
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).
std::shared_ptr< TIds > LinkSharedPtr
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.
TIds * GetCells(vtkIdType ptId)
Return a list of cell ids using the point specified by ptId.
dataset represents arbitrary combinations of all possible cell types
int vtkTypeBool
Definition: vtkABI.h:64
#define VTK_DEPRECATED_IN_9_4_0(reason)
int vtkIdType
Definition: vtkType.h:315