OneSkeleton processing package.
More...
#include <OneSkeleton.h>
|
| OneSkeleton () |
|
| ~OneSkeleton () |
|
int | buildEdgeList (const int &vertexNumber, const int &cellNumber, const long long int *cellArray, vector< pair< int, int > > &edgeList) const |
|
int | buildEdgeLists (const vector< vector< long long int > > &cellArrays, vector< vector< pair< int, int > > > &edgeLists) const |
|
int | buildEdgeStars (const int &vertexNumber, const int &cellNumber, const long long int *cellArray, vector< vector< int > > &starList, vector< pair< int, int > > *edgeList=NULL, vector< vector< int > > *vertexStars=NULL) const |
|
int | buildEdgeSubList (const int &cellNumber, const long long int *cellArray, vector< pair< int, int > > &edgeList) const |
|
| 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) |
|
OneSkeleton processing package.
- Author
- Julien Tierny julie.nosp@m.n.ti.nosp@m.erny@.nosp@m.lip6.nosp@m..fr
- Date
- June 2015.
OneSkeleton is a processing package that handles the 1-skeleton (edges) of a triangulation.
- See also
- Triangulation
-
vtkTriangulation
-
vtkOneSkeleton
OneSkeleton::OneSkeleton |
( |
| ) |
|
OneSkeleton::~OneSkeleton |
( |
| ) |
|
int OneSkeleton::buildEdgeList |
( |
const int & |
vertexNumber, |
|
|
const int & |
cellNumber, |
|
|
const long long int * |
cellArray, |
|
|
vector< pair< int, int > > & |
edgeList |
|
) |
| const |
Compute the list of edges of a valid triangulation.
- Parameters
-
vertexNumber | Number of vertices in the triangulation. |
cellNumber | Number of maximum-dimensional cells in the triangulation (number of tetrahedra in 3D, triangles in 2D, etc.) |
cellArray | Pointer to a contiguous array of cells. Each entry starts by the number of vertices in the cell, followed by the vertex identifiers of the cell. |
edgeList | Output edge list (each entry is an ordered pair of vertex identifiers). |
- Returns
- Returns 0 upon success, negative values otherwise.
int OneSkeleton::buildEdgeLists |
( |
const vector< vector< long long int > > & |
cellArrays, |
|
|
vector< vector< pair< int, int > > > & |
edgeLists |
|
) |
| const |
Compute the list of edges of multiple triangulations.
- Parameters
-
cellArrays | Vector of cells. For each triangulation, each entry starts by the number of vertices in the cell, followed by the vertex identifiers of the cell. |
edgeList | Output edge list (each entry is an ordered pair of vertex identifiers). |
- Returns
- Returns 0 upon success, negative values otherwise.
int OneSkeleton::buildEdgeStars |
( |
const int & |
vertexNumber, |
|
|
const int & |
cellNumber, |
|
|
const long long int * |
cellArray, |
|
|
vector< vector< int > > & |
starList, |
|
|
vector< pair< int, int > > * |
edgeList = NULL , |
|
|
vector< vector< int > > * |
vertexStars = NULL |
|
) |
| const |
Compute the 3-star of all the edges of a triangulation (for each edge, list of the 3-dimensional cells connected to it).
- Parameters
-
vertexNumber | Number of vertices in the triangulation. |
cellNumber | Number of maximum-dimensional cells in the triangulation (number of tetrahedra in 3D, triangles in 2D, etc.) |
cellArray | Pointer to a contiguous array of cells. Each entry starts by the number of vertices in the cell, followed by the vertex identifiers of the cell. |
starList | Output list of 3-stars. The size of this vector will be equal to the number of edges in the mesh. Each entry stores a vector that lists the identifiers of all 3-dimensional cells connected to the entry's edge. |
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. |
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 OneSkeleton::buildEdgeSubList |
( |
const int & |
cellNumber, |
|
|
const long long int * |
cellArray, |
|
|
vector< pair< int, int > > & |
edgeList |
|
) |
| const |
Compute the list of edges of a sub-portion of a valid triangulation.
- Parameters
-
cellNumber | Number of maximum-dimensional cells in the considered subset of the triangulation (number of tetrahedra in 3D, triangles in 2D, etc.) |
cellArray | Pointer to a contiguous array of cells. Each entry starts by the number of vertices in the cell, followed by the vertex identifiers of the cell. |
edgeList | Output edge list (each entry is an ordered pair of vertex identifiers). |
- Returns
- Returns 0 upon success, negative values otherwise.
The documentation for this class was generated from the following files: