VTK  9.1.0
vtkStaticCellLinksTemplate.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkStaticCellLinksTemplate.h
5 
6  Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7  All rights reserved.
8  See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9 
10  This software is distributed WITHOUT ANY WARRANTY; without even
11  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12  PURPOSE. See the above copyright notice for more information.
13 
14 =========================================================================*/
45 #ifndef vtkStaticCellLinksTemplate_h
46 #define vtkStaticCellLinksTemplate_h
47 
48 class vtkDataSet;
49 class vtkPolyData;
52 class vtkCellArray;
53 
54 #include "vtkAbstractCellLinks.h"
55 
56 template <typename TIds>
58 {
59 public:
61 
67 
71  void Initialize();
72 
78 
83 
88 
93 
97  void SerialBuildLinks(const vtkIdType numPts, const vtkIdType numCells, vtkCellArray* cellArray);
99  const vtkIdType numPts, const vtkIdType numCells, vtkCellArray* cellArray);
100 
102 
105  TIds GetNumberOfCells(vtkIdType ptId) { return (this->Offsets[ptId + 1] - this->Offsets[ptId]); }
106  vtkIdType GetNcells(vtkIdType ptId) { return (this->Offsets[ptId + 1] - this->Offsets[ptId]); }
108 
113  bool MatchesCell(vtkIdType npts, const vtkIdType* pts);
114 
118  TIds* GetCells(vtkIdType ptId) { return (this->Links + this->Offsets[ptId]); }
119 
124  void GetCells(vtkIdType npts, const vtkIdType* pts, vtkIdList* cells);
125 
127 
130  unsigned long GetActualMemorySize();
132  void SelectCells(vtkIdType minMaxDegree[2], unsigned char* cellSelection);
134 
136 
142 
143 protected:
144  // The various templated data members
145  TIds LinksSize;
146  TIds NumPts;
147  TIds NumCells;
148 
149  // These point to the core data structures
150  TIds* Links; // contiguous runs of cell ids
151  TIds* Offsets; // offsets for each point into the links array
152 
153  // Support for execution
154  int Type;
156 
157 private:
159  void operator=(const vtkStaticCellLinksTemplate&) = delete;
160 };
161 
162 #include "vtkStaticCellLinksTemplate.txx"
163 
164 #endif
165 // VTK-HeaderTest-Exclude: vtkStaticCellLinksTemplate.h
vtkStaticCellLinksTemplate::GetActualMemorySize
unsigned long GetActualMemorySize()
Support vtkAbstractCellLinks API.
vtkStaticCellLinksTemplate::NumCells
TIds NumCells
Definition: vtkStaticCellLinksTemplate.h:147
vtkStaticCellLinksTemplate::~vtkStaticCellLinksTemplate
~vtkStaticCellLinksTemplate()
Instantiate and destructor methods.
vtkStaticCellLinksTemplate::GetNumberOfCells
TIds GetNumberOfCells(vtkIdType ptId)
Get the number of cells using the point specified by ptId.
Definition: vtkStaticCellLinksTemplate.h:105
vtkStaticCellLinksTemplate::GetSequentialProcessing
vtkTypeBool GetSequentialProcessing()
Control whether to thread or serial process.
Definition: vtkStaticCellLinksTemplate.h:140
vtkStaticCellLinksTemplate::Type
int Type
Definition: vtkStaticCellLinksTemplate.h:154
vtkIdType
int vtkIdType
Definition: vtkType.h:332
vtkStaticCellLinksTemplate::BuildLinks
void BuildLinks(vtkUnstructuredGrid *ugrid)
Build the link list array for vtkUnstructuredGrid.
vtkStaticCellLinksTemplate::vtkStaticCellLinksTemplate
vtkStaticCellLinksTemplate()
Instantiate and destructor methods.
vtkStaticCellLinksTemplate::Offsets
TIds * Offsets
Definition: vtkStaticCellLinksTemplate.h:151
vtkStaticCellLinksTemplate::GetCells
TIds * GetCells(vtkIdType ptId)
Return a list of cell ids using the point specified by ptId.
Definition: vtkStaticCellLinksTemplate.h:118
vtkStaticCellLinksTemplate::ThreadedBuildLinks
void ThreadedBuildLinks(const vtkIdType numPts, const vtkIdType numCells, vtkCellArray *cellArray)
vtkStaticCellLinksTemplate::LinksSize
TIds LinksSize
Definition: vtkStaticCellLinksTemplate.h:145
vtkStaticCellLinksTemplate::BuildLinks
void BuildLinks(vtkPolyData *pd)
Build the link list array for vtkPolyData.
vtkStaticCellLinksTemplate::SelectCells
void SelectCells(vtkIdType minMaxDegree[2], unsigned char *cellSelection)
Support vtkAbstractCellLinks API.
vtkStaticCellLinksTemplate
object represents upward pointers from points to list of cells using each point (template implementat...
Definition: vtkStaticCellLinksTemplate.h:58
vtkStaticCellLinksTemplate::GetCells
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.
vtkStaticCellLinksTemplate::GetNcells
vtkIdType GetNcells(vtkIdType ptId)
Get the number of cells using the point specified by ptId.
Definition: vtkStaticCellLinksTemplate.h:106
vtkCellArray
object to represent cell connectivity
Definition: vtkCellArray.h:290
vtkIdList
list of point or cell ids
Definition: vtkIdList.h:140
vtkStaticCellLinksTemplate::SetSequentialProcessing
void SetSequentialProcessing(vtkTypeBool seq)
Control whether to thread or serial process.
Definition: vtkStaticCellLinksTemplate.h:139
vtkStaticCellLinksTemplate::BuildLinks
void BuildLinks(vtkExplicitStructuredGrid *esgrid)
Build the link list array for vtkExplicitStructuredGrid.
vtkStaticCellLinksTemplate::BuildLinks
void BuildLinks(vtkDataSet *ds)
Build the link list array for a general dataset.
vtkDataSet
abstract class to specify dataset behavior
Definition: vtkDataSet.h:166
vtkStaticCellLinksTemplate::NumPts
TIds NumPts
Definition: vtkStaticCellLinksTemplate.h:146
vtkStaticCellLinksTemplate::DeepCopy
void DeepCopy(vtkAbstractCellLinks *src)
Support vtkAbstractCellLinks API.
vtkPolyData
concrete dataset represents vertices, lines, polygons, and triangle strips
Definition: vtkPolyData.h:195
vtkExplicitStructuredGrid
structured grid with explicit topology and geometry
Definition: vtkExplicitStructuredGrid.h:97
vtkStaticCellLinksTemplate::SerialBuildLinks
void SerialBuildLinks(const vtkIdType numPts, const vtkIdType numCells, vtkCellArray *cellArray)
Specialized methods for building links from cell array.
vtkUnstructuredGrid
dataset represents arbitrary combinations of all possible cell types
Definition: vtkUnstructuredGrid.h:204
vtkStaticCellLinksTemplate::Initialize
void Initialize()
Make sure any previously created links are cleaned up.
vtkStaticCellLinksTemplate::SequentialProcessing
vtkTypeBool SequentialProcessing
Definition: vtkStaticCellLinksTemplate.h:155
vtkStaticCellLinksTemplate::Links
TIds * Links
Definition: vtkStaticCellLinksTemplate.h:150
vtkTypeBool
int vtkTypeBool
Definition: vtkABI.h:69
vtkStaticCellLinksTemplate::MatchesCell
bool MatchesCell(vtkIdType npts, const vtkIdType *pts)
Indicate whether the point ids provided defines at least one cell, or a portion of a cell.