VTK  9.3.0
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
41VTK_ABI_NAMESPACE_BEGIN
42class vtkDataSet;
43class vtkPolyData;
46class vtkCellArray;
47VTK_ABI_NAMESPACE_END
48
50
51VTK_ABI_NAMESPACE_BEGIN
52template <typename TIds>
54{
55public:
57
63
67 void Initialize();
68
74
79
84
89
93 void SerialBuildLinks(vtkIdType numPts, vtkIdType numCells, vtkCellArray* cellArray);
94 void ThreadedBuildLinks(vtkIdType numPts, vtkIdType numCells, vtkCellArray* cellArray);
95
97
100 TIds GetNumberOfCells(vtkIdType ptId) { return (this->Offsets[ptId + 1] - this->Offsets[ptId]); }
101 vtkIdType GetNcells(vtkIdType ptId) { return (this->Offsets[ptId + 1] - this->Offsets[ptId]); }
103
108 template <typename TGivenIds>
109 bool MatchesCell(TGivenIds npts, const TGivenIds* pts);
110
114 TIds* GetCells(vtkIdType ptId) { return (this->Links + this->Offsets[ptId]); }
115
120 void GetCells(vtkIdType npts, const vtkIdType* pts, vtkIdList* cells);
121
126 TIds GetLinksSize() { return this->LinksSize; }
127
132 TIds GetOffset(vtkIdType ptId) { return this->Offsets[ptId]; }
133
135
138 unsigned long GetActualMemorySize();
140 void SelectCells(vtkIdType minMaxDegree[2], unsigned char* cellSelection);
142
144
150
151protected:
152 // The various templated data members
154 TIds NumPts;
156
157 // These point to the core data structures
158 TIds* Links; // contiguous runs of cell ids
159 TIds* Offsets; // offsets for each point into the links array
160
161 // Support for execution
162 int Type;
164
165private:
167 void operator=(const vtkStaticCellLinksTemplate&) = delete;
168};
169
170VTK_ABI_NAMESPACE_END
171#include "vtkStaticCellLinksTemplate.txx"
172
173#endif
174// VTK-HeaderTest-Exclude: vtkStaticCellLinksTemplate.h
object to represent cell connectivity
abstract class to specify dataset behavior
Definition vtkDataSet.h:53
structured grid with explicit topology and geometry
list of point or cell ids
Definition vtkIdList.h:23
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition vtkPolyData.h:80
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.
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 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.
TIds * GetCells(vtkIdType ptId)
Return a list of cell ids using the point specified by ptId.
void BuildLinks(vtkDataSet *ds)
Build the link list array for a general dataset.
~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.
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.
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.
void DeepCopy(vtkAbstractCellLinks *src)
Support vtkAbstractCellLinks API.
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)
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
int vtkIdType
Definition vtkType.h:315