TTK
Loading...
Searching...
No Matches
Public Member Functions | List of all members
ttk::TwoSkeleton Class Reference

TwoSkeleton processing package. More...

#include <TwoSkeleton.h>

Inheritance diagram for ttk::TwoSkeleton:
ttk::Debug ttk::BaseClass

Public Member Functions

 TwoSkeleton ()
 
int buildCellNeighborsFromEdges (const CellArray &cellArray, FlatJaggedArray &cellNeighbors, const FlatJaggedArray &edgeStars) const
 
int buildCellNeighborsFromVertices (const SimplexId &vertexNumber, const CellArray &cellArray, FlatJaggedArray &cellNeighbors, FlatJaggedArray *vertexStars=nullptr) const
 
int buildEdgeTriangles (const SimplexId &vertexNumber, const CellArray &cellArray, FlatJaggedArray &edgeTriangleList, const std::vector< std::array< SimplexId, 2 > > &edgeList, std::vector< std::array< SimplexId, 3 > > *triangleEdgeList=nullptr) const
 
int buildTriangleList (const SimplexId &vertexNumber, const CellArray &cellArray, std::vector< std::array< SimplexId, 3 > > *triangleList=nullptr, FlatJaggedArray *triangleStars=nullptr, std::vector< std::array< SimplexId, 4 > > *cellTriangleList=nullptr) const
 
int buildTriangleEdgeList (const SimplexId &vertexNumber, const CellArray &cellArray, std::vector< std::array< SimplexId, 3 > > &triangleEdgeList, const std::vector< std::array< SimplexId, 2 > > &edgeList, FlatJaggedArray *vertexEdgeList=nullptr, std::vector< std::array< SimplexId, 3 > > *triangleList=nullptr, FlatJaggedArray *triangleStarList=nullptr, std::vector< std::array< SimplexId, 4 > > *cellTriangleList=nullptr) const
 
int buildTriangleLinks (const std::vector< std::array< SimplexId, 3 > > &triangleList, const FlatJaggedArray &triangleStars, const CellArray &cellArray, FlatJaggedArray &triangleLinks) const
 
int buildVertexTriangles (const SimplexId &vertexNumber, const std::vector< std::array< SimplexId, 3 > > &triangleList, FlatJaggedArray &vertexTriangles) const
 
- Public Member Functions inherited from ttk::Debug
 Debug ()
 
 ~Debug () override
 
virtual int setDebugLevel (const int &debugLevel)
 
int setWrapper (const Wrapper *wrapper) override
 
int printMsg (const std::string &msg, const debug::Priority &priority=debug::Priority::INFO, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cout) const
 
int printMsg (const std::vector< std::string > &msgs, const debug::Priority &priority=debug::Priority::INFO, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cout) const
 
int printErr (const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
 
int printWrn (const std::string &msg, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cerr) const
 
int printMsg (const std::string &msg, const double &progress, const double &time, const int &threads, const double &memory, const debug::LineMode &lineMode=debug::LineMode::NEW, const debug::Priority &priority=debug::Priority::PERFORMANCE, std::ostream &stream=std::cout) const
 
int printMsg (const std::string &msg, const double &progress, const double &time, const debug::LineMode &lineMode=debug::LineMode::NEW, const debug::Priority &priority=debug::Priority::PERFORMANCE, std::ostream &stream=std::cout) const
 
int printMsg (const std::string &msg, const double &progress, const double &time, const int &threads, const debug::LineMode &lineMode=debug::LineMode::NEW, const debug::Priority &priority=debug::Priority::PERFORMANCE, std::ostream &stream=std::cout) const
 
int printMsg (const std::string &msg, const double &progress, const debug::LineMode &lineMode=debug::LineMode::NEW, const debug::Priority &priority=debug::Priority::PERFORMANCE, std::ostream &stream=std::cout) const
 
int printMsg (const std::string &msg, const double &progress, const debug::Priority &priority, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cout) const
 
int printMsg (const std::vector< std::vector< std::string > > &rows, const debug::Priority &priority=debug::Priority::INFO, const bool hasHeader=true, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cout) const
 
int printMsg (const debug::Separator &separator, const debug::LineMode &lineMode=debug::LineMode::NEW, const debug::Priority &priority=debug::Priority::INFO, std::ostream &stream=std::cout) const
 
int printMsg (const debug::Separator &separator, const debug::Priority &priority, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cout) const
 
int printMsg (const std::string &msg, const debug::Separator &separator, const debug::LineMode &lineMode=debug::LineMode::NEW, const debug::Priority &priority=debug::Priority::INFO, std::ostream &stream=std::cout) const
 
void setDebugMsgPrefix (const std::string &prefix)
 
- Public Member Functions inherited from ttk::BaseClass
 BaseClass ()
 
virtual ~BaseClass ()=default
 
int getThreadNumber () const
 
virtual int setThreadNumber (const int threadNumber)
 

Additional Inherited Members

- Protected Member Functions inherited from ttk::Debug
int printMsgInternal (const std::string &msg, const std::string &right, const std::string &filler, const debug::Priority &priority=debug::Priority::INFO, const debug::LineMode &lineMode=debug::LineMode::NEW, std::ostream &stream=std::cout) const
 
int printMsgInternal (const std::string &msg, const debug::Priority &priority, const debug::LineMode &lineMode, std::ostream &stream=std::cout) const
 
int welcomeMsg (std::ostream &stream)
 
- Protected Attributes inherited from ttk::Debug
int debugLevel_
 
std::string debugMsgPrefix_
 
std::string debugMsgNamePrefix_
 
- Protected Attributes inherited from ttk::BaseClass
bool lastObject_
 
int threadNumber_
 
Wrapperwrapper_
 
- Static Protected Attributes inherited from ttk::Debug
static COMMON_EXPORTS debug::LineMode lastLineMode = ttk::debug::LineMode::NEW
 

Detailed Description

TwoSkeleton processing package.

Author
Julien Tierny julie.nosp@m.n.ti.nosp@m.erny@.nosp@m.lip6.nosp@m..fr
Date
August 2015.

TwoSkeleton is a processing package that handles the 2-skeleton (triangles) of a triangulation.

See also
Triangulation
ttkTriangulation

Definition at line 24 of file TwoSkeleton.h.

Constructor & Destructor Documentation

◆ TwoSkeleton()

TwoSkeleton::TwoSkeleton ( )

Definition at line 6 of file TwoSkeleton.cpp.

Member Function Documentation

◆ buildCellNeighborsFromEdges()

int TwoSkeleton::buildCellNeighborsFromEdges ( const CellArray cellArray,
FlatJaggedArray cellNeighbors,
const FlatJaggedArray edgeStars 
) const

Compute the list of cell-neighbors of each cell of a 2D triangulation (unspecified behavior if the input mesh is not a triangulation).

Parameters
cellArrayCell container allowing to retrieve the vertices ids of each cell.
cellNeighborsOutput neighbor list. The size of this std::vector will be equal to the number of cells in the mesh. Each entry will be a std::vector listing the cell identifiers of the entry's cell's neighbors.
edgeStarsArray of edge stars (list of 2-dimensional cells connected to each edge).
Returns
Returns 0 upon success, negative values otherwise.

Definition at line 10 of file TwoSkeleton.cpp.

◆ buildCellNeighborsFromVertices()

int TwoSkeleton::buildCellNeighborsFromVertices ( const SimplexId vertexNumber,
const CellArray cellArray,
FlatJaggedArray cellNeighbors,
FlatJaggedArray vertexStars = nullptr 
) const

Compute the list of cell-neighbors of each cell of a 2D triangulation (unspecified behavior if the input mesh is not a triangulation).

Parameters
vertexNumberNumber of vertices in the triangulation.
cellArrayCell container allowing to retrieve the vertices ids of each cell.
cellNeighborsOutput neighbor list. The size of this std::vector will be equal to the number of cells in the mesh. Each entry will be a std::vector listing the cell identifiers of the entry's cell's neighbors.
vertexStarsOptional list of vertex stars (list of 2-dimensional cells connected to each vertex). If nullptr, the function will compute this list anyway and free the related memory upon return. If not nullptr but pointing to an empty std::vector, the function will fill this empty std::vector (useful if this list needs to be used later on by the calling program). If not nullptr but pointing to a non-empty std::vector, this function will use this std::vector as internal vertex star list. If this std::vector is not empty but incorrect, the behavior is unspecified.
Returns
Returns 0 upon success, negative values otherwise.

Definition at line 65 of file TwoSkeleton.cpp.

◆ buildEdgeTriangles()

int TwoSkeleton::buildEdgeTriangles ( const SimplexId vertexNumber,
const CellArray cellArray,
FlatJaggedArray edgeTriangleList,
const std::vector< std::array< SimplexId, 2 > > &  edgeList,
std::vector< std::array< SimplexId, 3 > > *  triangleEdgeList = nullptr 
) const

Compute the list of triangles connected to each edge for 3D triangulations (unspecified behavior if the input mesh is not a triangulation).

Parameters
vertexNumberNumber of vertices in the triangulation.
cellArrayCell container allowing to retrieve the vertices ids of each cell.
edgeTriangleListOutput edge triangle list. The size of this std::vector will be equal to the number of edges in the triangulation. Each entry will be a std::vector listing the triangle identifiers for each triangle connected to the entry's edge.
edgeListEdge list (list of std::pairs of vertex identifiers). If this std::vector is not empty but incorrect, the behavior is unspecified.
triangleEdgeListOptional output triangle edge list (list of std::vectors of edges identifiers per triangle). If nullptr, the function will compute this list anyway and free the related memory upon return. If not nullptr but pointing to an empty std::vector, the function will fill this empty std::vector (useful if this list needs to be used later on by the calling program). If not nullptr but pointing to a non-empty std::vector, this function will use this std::vector as internal edge list. If this std::vector is not empty but incorrect, the behavior is unspecified.

Definition at line 159 of file TwoSkeleton.cpp.

◆ buildTriangleEdgeList()

int TwoSkeleton::buildTriangleEdgeList ( const SimplexId vertexNumber,
const CellArray cellArray,
std::vector< std::array< SimplexId, 3 > > &  triangleEdgeList,
const std::vector< std::array< SimplexId, 2 > > &  edgeList,
FlatJaggedArray vertexEdgeList = nullptr,
std::vector< std::array< SimplexId, 3 > > *  triangleList = nullptr,
FlatJaggedArray triangleStarList = nullptr,
std::vector< std::array< SimplexId, 4 > > *  cellTriangleList = nullptr 
) const

Compute the list of edges connected to each triangle for 3D triangulations (unspecified behavior if the input mesh is not a 3D triangulation).

Parameters
vertexNumberNumber of vertices in the triangulation.
cellArrayCell container allowing to retrieve the vertices ids of each cell.
triangleEdgeListOutput triangle edge list. The size of this std::vector will be equal to the number of triangles in the triangulation. Each entry will be a std::vector listing the edge identifiers for each edge connected to the entry's triangle.
edgeListEdge list (list of std::pairs of vertex identifiers). If this std::vector is not empty but incorrect, the behavior is unspecified.
vertexEdgeListOptional output vertex edge list (list of edge identifiers for each vertex). If nullptr, the function will compute this list anyway and free the related memory upon return. If not nullptr but pointing to an empty std::vector, the function will fill this empty std::vector (useful if this list needs to be used later on by the calling program). If not nullptr but pointing to a non-empty std::vector, this function will use this std::vector as internal vertex edge list. If this std::vector is not empty but incorrect, the behavior is unspecified.
triangleListOptional output triangle list (list of std::vectors of vertex identifiers). If nullptr, the function will compute this list anyway and free the related memory upon return. If not nullptr but pointing to an empty std::vector, the function will fill this empty std::vector (useful if this list needs to be used later on by the calling program). If not nullptr but pointing to a non-empty std::vector, this function will use this std::vector as internal triangle list. If this std::vector is not empty but incorrect, the behavior is unspecified.
triangleStarListOptional output triangle star list (list of tetrahedron identifiers for each triangle). If nullptr, the function will compute this list anyway and free the related memory upon return. If not nullptr but pointing to an empty std::vector, the function will fill this empty std::vector (useful if this list needs to be used later on by the calling program). If not nullptr but pointing to a non-empty std::vector, this function will use this std::vector as internal triangle star list. If this std::vector is not empty but incorrect, the behavior is unspecified.
cellTriangleListOptional output cell triangle list (list of triangle identifiers for each tetrahedron). If nullptr, the function will compute this list anyway and free the related memory upon return. If not nullptr but pointing to an empty std::vector, the function will fill this empty std::vector (useful if this list needs to be used later on by the calling program). If not nullptr but pointing to a non-empty std::vector, this function will use this std::vector as internal cell triangle list. If this std::vector is not empty but incorrect, the behavior is unspecified.
Returns
Returns 0 upon success, negative values otherwise.

Definition at line 400 of file TwoSkeleton.cpp.

◆ buildTriangleLinks()

int TwoSkeleton::buildTriangleLinks ( const std::vector< std::array< SimplexId, 3 > > &  triangleList,
const FlatJaggedArray triangleStars,
const CellArray cellArray,
FlatJaggedArray triangleLinks 
) const

Compute the links of triangles in a 3D triangulation.

Parameters
triangleListInput triangle list. The number of entries of this list is equal to the number of triangles in the triangulation. Each entry lists the vertex identifiers of the corresponding triangle.
triangleStarsInput triangle star list. The number of entries of this list is equal to the number of triangles in the triangulation. Each entry lists the identifiers of the tetrahedra which are the co-faces of the corresponding triangle.
cellArrayCell container allowing to retrieve the vertices ids of each cell.
triangleLinksOutput triangle link list. The number of entries of this list is equal to the number of triangles in the triangulation. Each entry lists the identifiers of the vertices in the link of the corresponding triangle.
Returns
Returns 0 upon success, negative values otherwise.

Definition at line 476 of file TwoSkeleton.cpp.

◆ buildTriangleList()

int TwoSkeleton::buildTriangleList ( const SimplexId vertexNumber,
const CellArray cellArray,
std::vector< std::array< SimplexId, 3 > > *  triangleList = nullptr,
FlatJaggedArray triangleStars = nullptr,
std::vector< std::array< SimplexId, 4 > > *  cellTriangleList = nullptr 
) const

Compute the list of triangles of a triangulation represented by a vtkUnstructuredGrid object. Unspecified behavior if the input mesh is not a valid triangulation.

Parameters
vertexNumberNumber of vertices in the triangulation.
cellArrayCell container allowing to retrieve the vertices ids of each cell.
triangleListOptional output triangle list (each entry is the ordered std::vector of the vertex identifiers of the entry's triangle).
triangleStarsOptional output for triangle tet-adjacency (for each triangle, list of its adjacent tetrahedra).
cellTriangleListOptional list of triangles per tetrahedron cell.
Returns
Returns 0 upon success, negative values otherwise.

Definition at line 229 of file TwoSkeleton.cpp.

◆ buildVertexTriangles()

int TwoSkeleton::buildVertexTriangles ( const SimplexId vertexNumber,
const std::vector< std::array< SimplexId, 3 > > &  triangleList,
FlatJaggedArray vertexTriangles 
) const

Compute the list of triangles connected to each vertex for 3D triangulations (unspecified behavior if the input mesh is not a triangulation).

Parameters
vertexNumberNumber of vertices in the triangulation.
triangleListInput triangle list (list of std::vectors of vertex identifiers).
vertexTrianglesOutput vertex triangle list (list of std::vectors of triangle identifiers).

Definition at line 538 of file TwoSkeleton.cpp.


The documentation for this class was generated from the following files: