VTK helper that processes the 0-skeleton (vertices) of a data-set.
More...
#include <vtkZeroSkeleton.h>
|
| vtkZeroSkeleton () |
|
int | buildVertexEdges (vtkDataSet *input, vector< vector< int > > &vertexEdges, vector< pair< int, int > > *edgeList=NULL, const bool &isTriangulation=false) const |
|
int | buildVertexLinks (vtkDataSet *input, vector< vector< long long int > > &vertexLinks, vector< vector< int > > *vertexStars=NULL, const bool &isTriangulation=false) const |
|
int | buildVertexNeighbors (vtkDataSet *input, vector< vector< int > > &vertexNeighbors, const bool &isTriangulation=false) const |
|
int | buildVertexStars (vtkDataSet *input, vector< vector< int > > &vertexStars, const bool &isTriangulation=false) const |
|
int | buildTriangulationVertexLinks (vtkPolyData *input, vector< vector< long long int > > &vertexLinks, vector< vector< int > > *vertexStars=NULL) const |
|
int | buildTriangulationVertexLinks (vtkUnstructuredGrid *input, vector< vector< long long int > > &vertexLinks, vector< vector< int > > *vertexStars=NULL) const |
|
int | buildTriangulationVertexNeighbors (vtkPolyData *input, vector< vector< int > > &vertexNeighbors, vector< pair< int, int > > *edgeList=NULL) const |
|
int | buildTriangulationVertexNeighbors (vtkUnstructuredGrid *input, vector< vector< int > > &vertexNeighbors, vector< pair< int, int > > *edgeList=NULL) const |
|
int | buildTriangulationVertexStars (vtkPolyData *input, vector< vector< int > > &vertexStars) const |
|
int | buildTriangulationVertexStars (vtkUnstructuredGrid *input, vector< vector< int > > &vertexStars) const |
|
| Wrapper () |
|
| ~Wrapper () |
|
| Debug () |
|
virtual | ~Debug () |
|
virtual const int | dMsg (ostream &stream, string msg, const int &debugLevel=infoMsg) const |
|
const int | err (const string msg, const int &debugLevel=infoMsg) const |
|
const int | msg (const char *msg, const int &debugLevel=infoMsg) const |
|
virtual const int | setDebugLevel (const int &debugLevel) |
|
int | setThreadNumber (const int threadNumber) |
|
int | setWrapper (const Wrapper *wrapper) |
|
VTK helper that processes the 0-skeleton (vertices) of a data-set.
- Author
- Julien Tierny julie.nosp@m.n.ti.nosp@m.erny@.nosp@m.lip6.nosp@m..fr
- Date
- November 2014.
- See also
- vtkTriangulation
-
wtfit::Triangulation
-
wtfit::ZeroSkeleton
vtkZeroSkeleton::vtkZeroSkeleton |
( |
| ) |
|
int vtkZeroSkeleton::buildTriangulationVertexLinks |
( |
vtkPolyData * |
input, |
|
|
vector< vector< long long int > > & |
vertexLinks, |
|
|
vector< vector< int > > * |
vertexStars = NULL |
|
) |
| const |
|
inline |
Compute the link of each vertex of a triangulation represented by a vtkPolyData object (unspecified behavior if the input mesh is not a valid triangulation).
- Parameters
-
input | Input data-set. |
vertexLinks | Output vertex links. The size of this vector will be equal to the number of vertices in the mesh. Each entry will be a vector listing the simplices of the link of the entry's vertex. In particular, this vector contains, for each simplex, the number of vertices in the simplex (triangles: 3, edges: 2) followed by the corresponding vertex identifiers. |
vertexStars | Optional list of vertex stars (list of 3-dimensional cells connected to each vertex). If NULL, the function will compute this list anyway and free the related memory upon return. If not NULL but pointing to an empty vector, the function will fill this empty vector (useful if this list needs to be used later on by the calling program). If not NULL but pointing to a non-empty vector, this function will use this vector as internal vertex star list. If this vector is not empty but incorrect, the behavior is unspecified. |
- Returns
- Returns 0 upon success, negative values otherwise.
int vtkZeroSkeleton::buildTriangulationVertexLinks |
( |
vtkUnstructuredGrid * |
input, |
|
|
vector< vector< long long int > > & |
vertexLinks, |
|
|
vector< vector< int > > * |
vertexStars = NULL |
|
) |
| const |
|
inline |
Compute the link of each vertex of a triangulation represented by a vtkUnstructuredGrid object (unspecified behavior if the input mesh is not a valid triangulation).
- Parameters
-
input | Input data-set. |
vertexLinks | Output vertex links. The size of this vector will be equal to the number of vertices in the mesh. Each entry will be a vector listing the simplices of the link of the entry's vertex. In particular, this vector contains, for each simplex, the number of vertices in the simplex (triangles: 3, edges: 2) followed by the corresponding vertex identifiers. |
vertexStars | Optional list of vertex stars (list of 3-dimensional cells connected to each vertex). If NULL, the function will compute this list anyway and free the related memory upon return. If not NULL but pointing to an empty vector, the function will fill this empty vector (useful if this list needs to be used later on by the calling program). If not NULL but pointing to a non-empty vector, this function will use this vector as internal vertex star list. If this vector is not empty but incorrect, the behavior is unspecified. |
- Returns
- Returns 0 upon success, negative values otherwise.
int vtkZeroSkeleton::buildTriangulationVertexNeighbors |
( |
vtkPolyData * |
input, |
|
|
vector< vector< int > > & |
vertexNeighbors, |
|
|
vector< pair< int, int > > * |
edgeList = NULL |
|
) |
| const |
|
inline |
Compute the list of neighbors of each vertex of a vtkPolyData. Unspecified behavior if the input mesh is not a valid triangulation).
- Parameters
-
input | Input data-set. |
vertexNeighbors | Output neighbor list. The size of this vector will be equal to the number of vertices in the mesh. Each entry will be a vector listing the vertex identifiers of the entry's vertex' neighbors. |
edgeList | Optional list of edges. If NULL, the function will compute this list anyway and free the related memory upon return. If not NULL but pointing to an empty vector, the function will fill this empty vector (useful if this list needs to be used later on by the calling program). If not NULL but pointing to a non-empty vector, this function will use this vector as internal edge list. If this vector is not empty but incorrect, the behavior is unspecified |
- Returns
- Returns 0 upon success, negative values otherwise.
int vtkZeroSkeleton::buildTriangulationVertexNeighbors |
( |
vtkUnstructuredGrid * |
input, |
|
|
vector< vector< int > > & |
vertexNeighbors, |
|
|
vector< pair< int, int > > * |
edgeList = NULL |
|
) |
| const |
|
inline |
Compute the list of neighbors of each vertex of a vtkUnstructuredGrid. Unspecified behavior if the input mesh is not a valid triangulation).
- Parameters
-
input | Input data-set. |
vertexNeighbors | Output neighbor list. The size of this vector will be equal to the number of vertices in the mesh. Each entry will be a vector listing the vertex identifiers of the entry's vertex' neighbors. |
edgeList | Optional list of edges. If NULL, the function will compute this list anyway and free the related memory upon return. If not NULL but pointing to an empty vector, the function will fill this empty vector (useful if this list needs to be used later on by the calling program). If not NULL but pointing to a non-empty vector, this function will use this vector as internal edge list. If this vector is not empty but incorrect, the behavior is unspecified |
- Returns
- Returns 0 upon success, negative values otherwise.
int vtkZeroSkeleton::buildTriangulationVertexStars |
( |
vtkPolyData * |
input, |
|
|
vector< vector< int > > & |
vertexStars |
|
) |
| const |
|
inline |
Compute the star of each vertex of a triangulation represented by a vtkPolyData object. Unspecified behavior if the input mesh is not a valid triangulation.
- Parameters
-
input | Input data-set. |
vertexStars | Output vertex stars. The size of this vector will be equal to the number of vertices in the mesh. Each entry will be a vector listing the identifiers of the maximum-dimensional cells (3D: tetrahedra, 2D: triangles, etc.) connected to the entry's vertex. |
- Returns
- Returns 0 upon success, negative values otherwise.
int vtkZeroSkeleton::buildTriangulationVertexStars |
( |
vtkUnstructuredGrid * |
input, |
|
|
vector< vector< int > > & |
vertexStars |
|
) |
| const |
|
inline |
Compute the star of each vertex of a triangulation represented by a vtkUnstructuredGrid object. Unspecified behavior if the input mesh is not a valid triangulation.
- Parameters
-
input | Input data-set. |
vertexStars | Output vertex stars. The size of this vector will be equal to the number of vertices in the mesh. Each entry will be a vector listing the identifiers of the maximum-dimensional cells (3D: tetrahedra, 2D: triangles, etc.) connected to the entry's vertex. |
- Returns
- Returns 0 upon success, negative values otherwise.
int vtkZeroSkeleton::buildVertexEdges |
( |
vtkDataSet * |
input, |
|
|
vector< vector< int > > & |
vertexEdges, |
|
|
vector< pair< int, int > > * |
edgeList = NULL , |
|
|
const bool & |
isTriangulation = false |
|
) |
| const |
Compute the list of edges connected to each vertex of a vtkDataSet.
- Parameters
-
input | Input data-set. |
vertexEdges | Output vertex links. The size of this vector will be equal to the number of vertices in the mesh. Each entry will be a vector listing the identifiers of the edges connected to the entry's vertex. |
edgeList | Optional list of edges. If NULL, the function will compute this list anyway and free the related memory upon return. If not NULL but pointing to an empty vector, the function will fill this empty vector (useful if this list needs to be used later on by the calling program). If not NULL but pointing to a non-empty vector, this function will use this vector as internal edge list. If this vector is not empty but incorrect, the behavior is unspecified. |
isTriangulation | Optional flag that speeds up computation if the input mesh is indeed a valid triangulation (unspecified behavior otherwise). |
- Returns
- Returns 0 upon success, negative values otherwise.
int vtkZeroSkeleton::buildVertexLinks |
( |
vtkDataSet * |
input, |
|
|
vector< vector< long long int > > & |
vertexLinks, |
|
|
vector< vector< int > > * |
vertexStars = NULL , |
|
|
const bool & |
isTriangulation = false |
|
) |
| const |
Compute the link of each vertex of a vtkDataSet.
- Parameters
-
input | Input data-set. |
vertexLinks | Output vertex links. The size of this vector will be equal to the number of vertices in the mesh. Each entry will be a vector listing the simplices of the link of the entry's vertex. In particular, this vector contains, for each simplex, the number of vertices in the simplex (triangles: 3, edges: 2) followed by the corresponding vertex identifiers. |
vertexStars | Optional list of vertex stars (list of 3-dimensional cells connected to each vertex). If NULL, the function will compute this list anyway and free the related memory upon return. If not NULL but pointing to an empty vector, the function will fill this empty vector (useful if this list needs to be used later on by the calling program). If not NULL but pointing to a non-empty vector, this function will use this vector as internal vertex star list. If this vector is not empty but incorrect, the behavior is unspecified. |
isTriangulation | Optional flag that speeds up computation if the input mesh is indeed a valid triangulation (unspecified behavior otherwise). |
- Returns
- Returns 0 upon success, negative values otherwise.
int vtkZeroSkeleton::buildVertexNeighbors |
( |
vtkDataSet * |
input, |
|
|
vector< vector< int > > & |
vertexNeighbors, |
|
|
const bool & |
isTriangulation = false |
|
) |
| const |
Compute the list of neighbors of each vertex of a vtkDataSet.
- Parameters
-
input | Input data-set. |
vertexNeighbors | Output neighbor list. The size of this vector will be equal to the number of vertices in the mesh. Each entry will be a vector listing the vertex identifiers of the entry's vertex' neighbors. |
isTriangulation | Optional flag that speeds up computation if the input mesh is indeed a valid triangulation (unspecified behavior otherwise). |
- Returns
- Returns 0 upon success, negative values otherwise.
int vtkZeroSkeleton::buildVertexStars |
( |
vtkDataSet * |
input, |
|
|
vector< vector< int > > & |
vertexStars, |
|
|
const bool & |
isTriangulation = false |
|
) |
| const |
Compute the star of each vertex of a vtkDataSet.
- Parameters
-
input | Input data-set. |
vertexStars | Output vertex stars. The size of this vector will be equal to the number of vertices in the mesh. Each entry will be a vector listing the identifiers of the maximum-dimensional cells (3D: tetrahedra, 2D: triangles, etc.) connected to the entry's vertex. |
isTriangulation | Optional flag that speeds up computation if the input mesh is indeed a valid triangulation (unspecified behavior otherwise). |
- Returns
- Returns 0 upon success, negative values otherwise.
The documentation for this class was generated from the following files: