WTFIT
vtkZeroSkeleton.h
Go to the documentation of this file.
1 
11 #ifndef _VTK_ZERO_SKELETON_H
12 #define _VTK_ZERO_SKELETON_H
13 
14 // c++ includes
15 #include <algorithm>
16 
17 // base code includes
18 #include <Debug.h>
19 #include <ZeroSkeleton.h>
20 #include <vtkOneSkeleton.h>
21 
22 // VTK includes
23 #include <vtkCellArray.h>
24 #include <vtkDataSet.h>
25 #include <vtkGenericCell.h>
26 #include <vtkIdList.h>
27 #include <vtkPolyData.h>
28 #include <vtkSmartPointer.h>
29 #include <vtkUnstructuredGrid.h>
30 
31 // TODO:
32 // make a VTK filter if needed
33 
34 class vtkZeroSkeleton : public Wrapper{
35 
36  public:
37 
39 
57  int buildVertexEdges(vtkDataSet *input,
58  vector<vector<int> > &vertexEdges,
59  vector<pair<int, int> > *edgeList = NULL,
60  const bool &isTriangulation = false) const;
61 
83  int buildVertexLinks(vtkDataSet *input,
84  vector<vector<long long int > > &vertexLinks,
85  vector<vector<int> > *vertexStars = NULL,
86  const bool &isTriangulation = false) const;
87 
98  int buildVertexNeighbors(vtkDataSet *input,
99  vector<vector<int> > &vertexNeighbors,
100  const bool &isTriangulation = false) const;
101 
112  int buildVertexStars(vtkDataSet *input,
113  vector<vector<int> > &vertexStars,
114  const bool &isTriangulation = false) const;
115 
136  int buildTriangulationVertexLinks(vtkPolyData *input,
137  vector<vector<long long int > > &vertexLinks,
138  vector<vector<int> > *vertexStars = NULL) const{
139 
140  ZeroSkeleton zeroSkeleton;
141  zeroSkeleton.setWrapper(this);
142  return zeroSkeleton.buildVertexLinks(input->GetNumberOfPoints(),
143  input->GetNumberOfCells(), input->GetPolys()->GetPointer(),
144  vertexLinks, vertexStars);
145  }
146 
167  int buildTriangulationVertexLinks(vtkUnstructuredGrid *input,
168  vector<vector<long long int > > &vertexLinks,
169  vector<vector<int> > *vertexStars = NULL) const{
170 
171  ZeroSkeleton zeroSkeleton;
172  zeroSkeleton.setWrapper(this);
173  return zeroSkeleton.buildVertexLinks(input->GetNumberOfPoints(),
174  input->GetNumberOfCells(), input->GetCells()->GetPointer(),
175  vertexLinks, vertexStars);
176  }
177 
193  int buildTriangulationVertexNeighbors(vtkPolyData *input,
194  vector<vector<int> > &vertexNeighbors,
195  vector<pair<int, int> > *edgeList = NULL) const{
196 
197  ZeroSkeleton skeleton;
198  skeleton.setWrapper(this);
199  return skeleton.buildVertexNeighbors(input->GetNumberOfPoints(),
200  input->GetNumberOfCells(), input->GetPolys()->GetPointer(),
201  vertexNeighbors, edgeList);
202  }
203 
219  int buildTriangulationVertexNeighbors(vtkUnstructuredGrid *input,
220  vector<vector<int> > &vertexNeighbors,
221  vector<pair<int, int> > *edgeList = NULL) const{
222 
223  ZeroSkeleton skeleton;
224  skeleton.setWrapper(this);
225  return skeleton.buildVertexNeighbors(input->GetNumberOfPoints(),
226  input->GetNumberOfCells(), input->GetCells()->GetPointer(),
227  vertexNeighbors, edgeList);
228  }
229 
239  int buildTriangulationVertexStars(vtkPolyData *input,
240  vector<vector<int> > &vertexStars) const{
241 
242  ZeroSkeleton zeroSkeleton;
243  return zeroSkeleton.buildVertexStars(input->GetNumberOfPoints(),
244  input->GetNumberOfCells(),
245  input->GetPolys()->GetPointer(),
246  vertexStars);
247  }
248 
258  int buildTriangulationVertexStars(vtkUnstructuredGrid *input,
259  vector<vector<int> > &vertexStars) const{
260 
261  ZeroSkeleton zeroSkeleton;
262  return zeroSkeleton.buildVertexStars(input->GetNumberOfPoints(),
263  input->GetNumberOfCells(),
264  input->GetCells()->GetPointer(),
265  vertexStars);
266  }
267 
268 
269 
270  protected:
271 
272 
273  private:
274 
275  // empty wrapping to VTK for now
276  bool needsToAbort(){ return false;};
277 
278  int updateProgress(const float &progress) {return 0;};
279 };
280 
281 #endif // _VTK_ZERO_SKELETON_H
int buildVertexLinks(vtkDataSet *input, vector< vector< long long int > > &vertexLinks, vector< vector< int > > *vertexStars=NULL, const bool &isTriangulation=false) const
Definition: vtkZeroSkeleton.cpp:43
int buildVertexNeighbors(const int &vertexNumber, const int &cellNumber, const long long int *cellArray, vector< vector< int > > &vertexNeighbors, vector< pair< int, int > > *edgeList=NULL) const
Definition: ZeroSkeleton.cpp:377
ZeroSkeleton processing package.
Definition: ZeroSkeleton.h:25
int buildVertexNeighbors(vtkDataSet *input, vector< vector< int > > &vertexNeighbors, const bool &isTriangulation=false) const
Definition: vtkZeroSkeleton.cpp:213
int buildTriangulationVertexNeighbors(vtkUnstructuredGrid *input, vector< vector< int > > &vertexNeighbors, vector< pair< int, int > > *edgeList=NULL) const
Definition: vtkZeroSkeleton.h:219
int buildTriangulationVertexNeighbors(vtkPolyData *input, vector< vector< int > > &vertexNeighbors, vector< pair< int, int > > *edgeList=NULL) const
Definition: vtkZeroSkeleton.h:193
int buildTriangulationVertexLinks(vtkUnstructuredGrid *input, vector< vector< long long int > > &vertexLinks, vector< vector< int > > *vertexStars=NULL) const
Definition: vtkZeroSkeleton.h:167
Wrapper class to wrap wtfit code.
Definition: Wrapper.h:15
int setWrapper(const Wrapper *wrapper)
Definition: Debug.cpp:76
int buildVertexEdges(vtkDataSet *input, vector< vector< int > > &vertexEdges, vector< pair< int, int > > *edgeList=NULL, const bool &isTriangulation=false) const
Definition: vtkZeroSkeleton.cpp:8
int buildVertexStars(vtkDataSet *input, vector< vector< int > > &vertexStars, const bool &isTriangulation=false) const
Definition: vtkZeroSkeleton.cpp:341
int buildTriangulationVertexLinks(vtkPolyData *input, vector< vector< long long int > > &vertexLinks, vector< vector< int > > *vertexStars=NULL) const
Definition: vtkZeroSkeleton.h:136
int buildTriangulationVertexStars(vtkUnstructuredGrid *input, vector< vector< int > > &vertexStars) const
Definition: vtkZeroSkeleton.h:258
int buildTriangulationVertexStars(vtkPolyData *input, vector< vector< int > > &vertexStars) const
Definition: vtkZeroSkeleton.h:239
vtkZeroSkeleton()
Definition: vtkZeroSkeleton.cpp:3
int buildVertexLinks(const int &vertexNumber, const int &cellNumber, const long long int *cellArray, vector< vector< long long int > > &vertexLinks, vector< vector< int > > *vertexStars=NULL) const
Definition: ZeroSkeleton.cpp:187
VTK helper that processes the 0-skeleton (vertices) of a data-set.
Definition: vtkZeroSkeleton.h:34
int buildVertexStars(const int &vertexNumber, const int &cellNumber, const long long int *cellArray, vector< vector< int > > &vertexStars) const
Definition: ZeroSkeleton.cpp:436