PhysicsBasedAnimationToolkit 0.0.10
Cross-platform C++20 library of algorithms and data structures commonly used in computer graphics research on physically-based simulation.
Loading...
Searching...
No Matches
pbat::graph Namespace Reference

Classes

struct  BreadthFirstSearch
 
struct  DepthFirstSearch
 
class  DisjointSetUnionFind
 Disjoint Set Union-Find data structure. More...
 
struct  PartitioningOptions
 Options for graph partitioning. More...
 
struct  WeightedEdgeTraits
 Traits for WeightedEdge. More...
 

Typedefs

template<class TWeight = Scalar, class TIndex = Index>
using WeightedEdge = Eigen::Triplet<TWeight, TIndex>
 Weighted edge (wrapper around Eigen triplet type)
 

Enumerations

enum class  EGreedyColorSelectionStrategy { LeastUsed , FirstAvailable }
 Enumeration of color selection strategies for graph coloring algorithms. More...
 
enum class  EGreedyColorOrderingStrategy { Natural , SmallestDegree , LargestDegree }
 Enumeration of vertex traversal ordering strategies for graph coloring algorithms. More...
 
enum class  EMeshDualGraphOptions : std::int32_t { VertexAdjacent = 0b001 , EdgeAdjacent = 0b010 , FaceAdjacent = 0b100 , All = 0b111 }
 Types of dual graph adjacencies.
 

Functions

template<class TWeightedEdgeIterator, class TWeightedEdge = typename std::iterator_traits<TWeightedEdgeIterator>::value_type, class TScalar = typename WeightedEdgeTraits<TWeightedEdge>::ScalarType, class TIndex = typename WeightedEdgeTraits<TWeightedEdge>::IndexType>
auto AdjacencyMatrixFromEdges (TWeightedEdgeIterator begin, TWeightedEdgeIterator end, TIndex m=TIndex(-1), TIndex n=TIndex(-1)) -> Eigen::SparseMatrix< TScalar, Eigen::RowMajor, TIndex >
 Construct adjacency matrix from edge/triplet list.
 
template<class TDerivedA, class TScalar = typename TDerivedA::Scalar, std::integral TIndex = typename TDerivedA::StorageIndex>
auto AdjacencyMatrixPrefix (Eigen::SparseCompressedBase< TDerivedA > const &A)
 Non-owning wrapper around the offset pointers of a compressed sparse matrix.
 
template<class TDerivedA, class TScalar = typename TDerivedA::Scalar, std::integral TIndex = typename TDerivedA::StorageIndex>
auto AdjacencyMatrixIndices (Eigen::SparseCompressedBase< TDerivedA > const &A)
 Non-owning wrapper around the indices of a compressed sparse matrix.
 
template<class TDerivedA, class TScalar = typename TDerivedA::Scalar, std::integral TIndex = typename TDerivedA::StorageIndex>
auto AdjacencyMatrixWeights (Eigen::SparseCompressedBase< TDerivedA > const &A)
 Non-owning wrapper around the weights of a compressed sparse matrix.
 
template<class TDerivedP, std::integral TIndex = typename TDerivedP::Scalar>
auto MapToAdjacency (Eigen::DenseBase< TDerivedP > const &p, TIndex n=TIndex(-1)) -> std::tuple< Eigen::Vector< TIndex, Eigen::Dynamic >, Eigen::Vector< TIndex, Eigen::Dynamic > >
 Construct adjacency list in compressed sparse format from a map p s.t. p(i) is the index of the vertex adjacent to i.
 
template<class TDerivedA, class TScalar = typename TDerivedA::Scalar, std::integral TIndex = typename TDerivedA::StorageIndex>
auto MatrixToAdjacency (Eigen::SparseCompressedBase< TDerivedA > const &A)
 Obtain the offset and indices arrays of an input adjacency matrix.
 
template<class TDerivedA, class TScalar = typename TDerivedA::Scalar, std::integral TIndex = typename TDerivedA::StorageIndex>
auto MatrixToWeightedAdjacency (Eigen::SparseCompressedBase< TDerivedA > const &A)
 Construct adjacency list in compressed sparse format from an input adjacency matrix.
 
template<class TIndex = Index>
auto ListOfListsToAdjacency (std::vector< std::vector< TIndex > > const &lil)
 Construct adjacency list in compressed sparse format from an input adjacency list in list of lists format.
 
template<class TDerivedPtr, class TDerivedAdj, class TIndex = typename TDerivedPtr::Scalar, class Func>
void ForEachEdge (Eigen::DenseBase< TDerivedPtr > &ptr, Eigen::DenseBase< TDerivedAdj > &adj, Func &&f)
 Edge iteration over the adjacency list in compressed sparse format.
 
template<class TDerivedPtr, class TDerivedAdj, class TIndex = typename TDerivedPtr::Scalar, class Func>
void RemoveEdges (Eigen::DenseBase< TDerivedPtr > &ptr, Eigen::DenseBase< TDerivedAdj > &adj, Func fShouldDeleteEdge)
 In-place removal of edges from the adjacency list in compressed sparse format.
 
template<class TDerivedPtr, class TDerivedAdj, int NC = 128, std::integral TIndex = typename TDerivedPtr::Scalar>
auto GreedyColor (Eigen::DenseBase< TDerivedPtr > const &ptr, Eigen::DenseBase< TDerivedAdj > const &adj, EGreedyColorOrderingStrategy eOrderingStrategy=EGreedyColorOrderingStrategy::LargestDegree, EGreedyColorSelectionStrategy eSelectionStrategy=EGreedyColorSelectionStrategy::LeastUsed) -> Eigen::Vector< TIndex, Eigen::Dynamic >
 Greedy graph coloring algorithm.
 
template<common::CIndex TIndex, class TDerivedP, class TDerivedAdj>
TIndex ConnectedComponents (Eigen::DenseBase< TDerivedP > const &ptr, Eigen::DenseBase< TDerivedAdj > const &adj, Eigen::Ref< Eigen::Vector< TIndex, Eigen::Dynamic > > components, DepthFirstSearch< TIndex > &dfs)
 Compute connected components of a graph using depth-first search.
 
template<common::CIndex TIndex, class TDerivedP, class TDerivedAdj>
TIndex ConnectedComponents (Eigen::DenseBase< TDerivedP > const &ptr, Eigen::DenseBase< TDerivedAdj > const &adj, Eigen::Ref< Eigen::Vector< TIndex, Eigen::Dynamic > > components, BreadthFirstSearch< TIndex > &bfs)
 Compute connected components of a graph using breadth-first search.
 
template<class TDerivedE, class TDerivedW, common::CIndex TIndex = typename TDerivedE::Scalar, common::CArithmetic TScalar = typename TDerivedW::Scalar>
auto MeshAdjacencyMatrix (Eigen::DenseBase< TDerivedE > const &E, Eigen::DenseBase< TDerivedW > const &w, TIndex nNodes=TIndex(-1), bool bVertexToElement=false, bool bHasDuplicates=false) -> Eigen::SparseMatrix< TScalar, Eigen::ColMajor, TIndex >
 Construct adjacency matrix from mesh.
 
template<class TDerivedE, common::CIndex TIndex = typename TDerivedE::Scalar>
auto MeshAdjacencyMatrix (Eigen::DenseBase< TDerivedE > const &E, TIndex nNodes=TIndex(-1), bool bVertexToElement=false, bool bHasDuplicates=false) -> Eigen::SparseMatrix< TIndex, Eigen::ColMajor, TIndex >
 Construct adjacency matrix from mesh.
 
template<class TDerivedE, common::CIndex TIndex = typename TDerivedE::Scalar>
auto MeshPrimalGraph (Eigen::DenseBase< TDerivedE > const &E, TIndex nNodes=TIndex(-1)) -> Eigen::SparseMatrix< TIndex, Eigen::ColMajor, TIndex >
 Construct primal graph of input mesh, i.e. the graph of adjacent vertices.
 
template<class TDerivedE, common::CIndex TIndex = typename TDerivedE::Scalar>
auto MeshDualGraph (Eigen::DenseBase< TDerivedE > const &E, TIndex nNodes=TIndex(-1), EMeshDualGraphOptions opts=EMeshDualGraphOptions::All) -> Eigen::SparseMatrix< TIndex, Eigen::ColMajor, TIndex >
 Construct dual graph of input mesh, i.e. the graph of adjacent elements.
 
IndexVectorX Partition (Eigen::Ref< IndexVectorX const > const &ptr, Eigen::Ref< IndexVectorX const > const &adj, Eigen::Ref< IndexVectorX const > const &wadj, Index nPartitions, PartitioningOptions opts=PartitioningOptions{})
 Partition input graph.
 

Detailed Description

Graph algorithms and data structures

Typedef Documentation

◆ WeightedEdge

template<class TWeight = Scalar, class TIndex = Index>
using pbat::graph::WeightedEdge = Eigen::Triplet<TWeight, TIndex>

Weighted edge (wrapper around Eigen triplet type)

Template Parameters
TWeightScalar type of the graph edge weights
TIndexIndex type of the graph vertices

Enumeration Type Documentation

◆ EGreedyColorOrderingStrategy

Enumeration of vertex traversal ordering strategies for graph coloring algorithms.

Enumerator
Natural 

Natural ordering of the vertices (i.e. [0,n-1])

SmallestDegree 

Always visit the vertex with the smallest degree next.

LargestDegree 

Always visit the vertex with the largest degree next.

◆ EGreedyColorSelectionStrategy

Enumeration of color selection strategies for graph coloring algorithms.

Enumerator
LeastUsed 

Select the least used color from the color palette.

FirstAvailable 

Select the first available color from the color palette.

Function Documentation

◆ AdjacencyMatrixFromEdges()

template<class TWeightedEdgeIterator, class TWeightedEdge = typename std::iterator_traits<TWeightedEdgeIterator>::value_type, class TScalar = typename WeightedEdgeTraits<TWeightedEdge>::ScalarType, class TIndex = typename WeightedEdgeTraits<TWeightedEdge>::IndexType>
auto pbat::graph::AdjacencyMatrixFromEdges ( TWeightedEdgeIterator begin,
TWeightedEdgeIterator end,
TIndex m = TIndex(-1),
TIndex n = TIndex(-1) ) -> Eigen::SparseMatrix<TScalar, Eigen::RowMajor, TIndex>

Construct adjacency matrix from edge/triplet list.

Template Parameters
TWeightedEdgeIteratorIterator type of the edge list
TWeightedEdgeType of the edge
TScalarScalar type of the graph edge weights
TIndexIndex type of the graph vertices
Parameters
beginIterator to the beginning of the edge list
endIterator to the end of the edge list
mNumber of rows of the adjacency matrix. If not provided (i.e. m < 0), it is inferred from the edge list.
nNumber of columns of the adjacency matrix. If not provided (i.e. n < 0), it is inferred from the edge list.
Returns
Adjacency matrix

◆ AdjacencyMatrixIndices()

template<class TDerivedA, class TScalar = typename TDerivedA::Scalar, std::integral TIndex = typename TDerivedA::StorageIndex>
auto pbat::graph::AdjacencyMatrixIndices ( Eigen::SparseCompressedBase< TDerivedA > const & A)

Non-owning wrapper around the indices of a compressed sparse matrix.

Template Parameters
TDerivedAType of input adjacency matrix
TScalarScalar type of the graph edge weights
TIndexIndex type of the graph vertices
Parameters
AInput adjacency matrix
Returns
Non-owning wrapper around the indices of the adjacency matrix

◆ AdjacencyMatrixPrefix()

template<class TDerivedA, class TScalar = typename TDerivedA::Scalar, std::integral TIndex = typename TDerivedA::StorageIndex>
auto pbat::graph::AdjacencyMatrixPrefix ( Eigen::SparseCompressedBase< TDerivedA > const & A)

Non-owning wrapper around the offset pointers of a compressed sparse matrix.

Template Parameters
TDerivedAType of input adjacency matrix
TScalarScalar type of the graph edge weights
TIndexIndex type of the graph vertices
Parameters
AInput adjacency matrix
Returns
Non-owning wrapper around the offset pointers of the adjacency matrix

◆ AdjacencyMatrixWeights()

template<class TDerivedA, class TScalar = typename TDerivedA::Scalar, std::integral TIndex = typename TDerivedA::StorageIndex>
auto pbat::graph::AdjacencyMatrixWeights ( Eigen::SparseCompressedBase< TDerivedA > const & A)

Non-owning wrapper around the weights of a compressed sparse matrix.

Template Parameters
TDerivedAType of input adjacency matrix
TScalarScalar type of the graph edge weights
TIndexIndex type of the graph vertices
Parameters
AInput adjacency matrix
Returns
Non-owning wrapper around the weights of the adjacency matrix

◆ ConnectedComponents() [1/2]

template<common::CIndex TIndex, class TDerivedP, class TDerivedAdj>
TIndex pbat::graph::ConnectedComponents ( Eigen::DenseBase< TDerivedP > const & ptr,
Eigen::DenseBase< TDerivedAdj > const & adj,
Eigen::Ref< Eigen::Vector< TIndex, Eigen::Dynamic > > components,
BreadthFirstSearch< TIndex > & bfs )

Compute connected components of a graph using breadth-first search.

Template Parameters
TIndexIndex type used in the graph
TDerivedPType of the pointer vector (adjacency list start indices for each vertex)
TDerivedAdjType of the adjacency list (vector of vertex indices)
Parameters
ptrPointer to the start of each vertex's adjacency list
adjAdjacency list of the graph
componentsOutput vector to store component labels for each vertex
bfsBreadth-first search object
Returns
Number of connected components found in the graph
Precondition
components[u] < 0 for all vertices u in the graph

◆ ConnectedComponents() [2/2]

template<common::CIndex TIndex, class TDerivedP, class TDerivedAdj>
TIndex pbat::graph::ConnectedComponents ( Eigen::DenseBase< TDerivedP > const & ptr,
Eigen::DenseBase< TDerivedAdj > const & adj,
Eigen::Ref< Eigen::Vector< TIndex, Eigen::Dynamic > > components,
DepthFirstSearch< TIndex > & dfs )

Compute connected components of a graph using depth-first search.

Template Parameters
TIndexIndex type used in the graph
TDerivedPType of the pointer vector (adjacency list start indices for each vertex)
TDerivedAdjType of the adjacency list (vector of vertex indices)
Parameters
ptrPointer to the start of each vertex's adjacency list
adjAdjacency list of the graph
componentsOutput vector to store component labels for each vertex
dfsDepth-first search object
Returns
Number of connected components found in the graph
Precondition
components[u] < 0 for all vertices u in the graph

◆ ForEachEdge()

template<class TDerivedPtr, class TDerivedAdj, class TIndex = typename TDerivedPtr::Scalar, class Func>
void pbat::graph::ForEachEdge ( Eigen::DenseBase< TDerivedPtr > & ptr,
Eigen::DenseBase< TDerivedAdj > & adj,
Func && f )

Edge iteration over the adjacency list in compressed sparse format.

Template Parameters
TDerivedPtrEigen dense expression of the offset pointers of the adjacency list
TDerivedAdjEigen dense expression of the indices of the adjacency list
TIndexIndex type of the graph vertices
FuncCallable type of the edge iteration function
Parameters
ptrOffset pointers of the adjacency list
adjIndices of the adjacency list
fCallable function to be applied to each edge with signature void(TIndex i, TIndex j, TIndex k)

◆ GreedyColor()

template<class TDerivedPtr, class TDerivedAdj, int NC = 128, std::integral TIndex = typename TDerivedPtr::Scalar>
auto pbat::graph::GreedyColor ( Eigen::DenseBase< TDerivedPtr > const & ptr,
Eigen::DenseBase< TDerivedAdj > const & adj,
EGreedyColorOrderingStrategy eOrderingStrategy = EGreedyColorOrderingStrategy::LargestDegree,
EGreedyColorSelectionStrategy eSelectionStrategy = EGreedyColorSelectionStrategy::LeastUsed ) -> Eigen::Vector<TIndex, Eigen::Dynamic>

Greedy graph coloring algorithm.

Template Parameters
TDerivedPtrEigen dense expression for offset pointers of the adjacency list
TDerivedAdjEigen dense expression for indices of the adjacency list
NCColor palette capacity
TIndexIndex type for vertices
Parameters
ptrOffset pointers array of the adjacency list
adjIndices array of the adjacency list
eOrderingStrategyVertex visiting order strategy
eSelectionStrategyColor selection strategy
Returns
|# vertices| array mapping vertices to their associated color

◆ ListOfListsToAdjacency()

template<class TIndex = Index>
auto pbat::graph::ListOfListsToAdjacency ( std::vector< std::vector< TIndex > > const & lil)

Construct adjacency list in compressed sparse format from an input adjacency list in list of lists format.

Template Parameters
TIndexIndex type of the graph vertices
Parameters
lilInput adjacency list in list of lists format
Returns
Tuple of the offset pointers and indices of the adjacency list

◆ MapToAdjacency()

template<class TDerivedP, std::integral TIndex = typename TDerivedP::Scalar>
auto pbat::graph::MapToAdjacency ( Eigen::DenseBase< TDerivedP > const & p,
TIndex n = TIndex(-1) ) -> std::tuple<Eigen::Vector<TIndex, Eigen::Dynamic>, Eigen::Vector<TIndex, Eigen::Dynamic>>

Construct adjacency list in compressed sparse format from a map p s.t. p(i) is the index of the vertex adjacent to i.

Template Parameters
TDerivedPType of input map
TScalarScalar type of the graph edge weights
Parameters
pInput map
nNumber of vertices in the graph. If not provided (i.e. n < 0), it is inferred from the map.
Returns
Tuple of the offset pointers and indices of the adjacency list

◆ MatrixToAdjacency()

template<class TDerivedA, class TScalar = typename TDerivedA::Scalar, std::integral TIndex = typename TDerivedA::StorageIndex>
auto pbat::graph::MatrixToAdjacency ( Eigen::SparseCompressedBase< TDerivedA > const & A)

Obtain the offset and indices arrays of an input adjacency matrix.

Template Parameters
TDerivedAType of input adjacency matrix
TScalarScalar type of the graph edge weights
TIndexIndex type of the graph vertices
Parameters
AInput adjacency matrix
Returns
Tuple of the offset pointers and indices of the adjacency matrix

◆ MatrixToWeightedAdjacency()

template<class TDerivedA, class TScalar = typename TDerivedA::Scalar, std::integral TIndex = typename TDerivedA::StorageIndex>
auto pbat::graph::MatrixToWeightedAdjacency ( Eigen::SparseCompressedBase< TDerivedA > const & A)

Construct adjacency list in compressed sparse format from an input adjacency matrix.

Template Parameters
TDerivedAType of input adjacency matrix
TScalarScalar type of the graph edge weights
TIndexIndex type of the graph vertices
Parameters
AInput adjacency matrix
Returns
Tuple of the offset pointers, indices, and weights of the adjacency matrix

◆ MeshAdjacencyMatrix() [1/2]

template<class TDerivedE, class TDerivedW, common::CIndex TIndex = typename TDerivedE::Scalar, common::CArithmetic TScalar = typename TDerivedW::Scalar>
auto pbat::graph::MeshAdjacencyMatrix ( Eigen::DenseBase< TDerivedE > const & E,
Eigen::DenseBase< TDerivedW > const & w,
TIndex nNodes = TIndex(-1),
bool bVertexToElement = false,
bool bHasDuplicates = false ) -> Eigen::SparseMatrix<TScalar, Eigen::ColMajor, TIndex>

Construct adjacency matrix from mesh.

Template Parameters
TDerivedEEigen dense expression for mesh elements
TDerivedWEigen dense expression for element-vertex weights
TIndexType of indices used in element array
TScalarType of weights
Parameters
E|# nodes per element|x|# elements| array of element indices
w|# nodes per element|x|# elements| array of element-vertex weights
nNodesNumber of nodes in the mesh. If nNodes < 1, the number of nodes is inferred from E.
bVertexToElementIf true, the adjacency matrix maps vertices to elements, rather than elements to vertices
bHasDuplicatesIf true, duplicate entries in the input mesh will be handled
Returns
Adjacency matrix of requested mesh connectivity
Precondition
E.rows() == w.rows() and E.cols() == w.cols()

◆ MeshAdjacencyMatrix() [2/2]

template<class TDerivedE, common::CIndex TIndex = typename TDerivedE::Scalar>
auto pbat::graph::MeshAdjacencyMatrix ( Eigen::DenseBase< TDerivedE > const & E,
TIndex nNodes = TIndex(-1),
bool bVertexToElement = false,
bool bHasDuplicates = false ) -> Eigen::SparseMatrix<TIndex, Eigen::ColMajor, TIndex>

Construct adjacency matrix from mesh.

Template Parameters
TDerivedEEigen dense expression for mesh elements
TIndexType of indices used in element array
Parameters
E|# nodes per element|x|# elements| array of element indices
nNodesNumber of nodes in the mesh. If nNodes < 1, the number of nodes is inferred from E.
bVertexToElementIf true, the adjacency matrix maps vertices to elements, rather than
bHasDuplicatesIf true, duplicate entries in the input mesh will be handled
Returns
Adjacency matrix of requested mesh connectivity

◆ MeshDualGraph()

template<class TDerivedE, common::CIndex TIndex = typename TDerivedE::Scalar>
auto pbat::graph::MeshDualGraph ( Eigen::DenseBase< TDerivedE > const & E,
TIndex nNodes = TIndex(-1),
EMeshDualGraphOptions opts = EMeshDualGraphOptions::All ) -> Eigen::SparseMatrix<TIndex, Eigen::ColMajor, TIndex>

Construct dual graph of input mesh, i.e. the graph of adjacent elements.

Template Parameters
TDerivedEEigen dense expression for mesh elements
TIndexType of indices used in element array
Parameters
E|# nodes per element|x|# elements| array of element indices
nNodesNumber of nodes in the mesh. If nNodes < 1, the number of nodes is inferred from E.
optsAdjacency types to keep in the dual graph
Returns
Dual graph of the input mesh

◆ MeshPrimalGraph()

template<class TDerivedE, common::CIndex TIndex = typename TDerivedE::Scalar>
auto pbat::graph::MeshPrimalGraph ( Eigen::DenseBase< TDerivedE > const & E,
TIndex nNodes = TIndex(-1) ) -> Eigen::SparseMatrix<TIndex, Eigen::ColMajor, TIndex>

Construct primal graph of input mesh, i.e. the graph of adjacent vertices.

Template Parameters
TDerivedEEigen dense expression for mesh elements
TIndexType of indices used in element array
Parameters
E|# nodes per element|x|# elements| array of element indices
nNodesNumber of nodes in the mesh. If nNodes < 1, the number of nodes is inferred from E.
Returns
Primal graph of the input mesh

◆ Partition()

PBAT_API IndexVectorX pbat::graph::Partition ( Eigen::Ref< IndexVectorX const > const & ptr,
Eigen::Ref< IndexVectorX const > const & adj,
Eigen::Ref< IndexVectorX const > const & wadj,
Index nPartitions,
PartitioningOptions opts = PartitioningOptions{} )

Partition input graph.

Internally delegates to METIS [7]

Parameters
ptrOffset pointers of adjacency list
adjIndices of adjacency list
wadjEdge weights of adjacency list
nPartitionsNumber of desired partitions
optsPartitioning options
Returns
|# vertices| array p of partition indices (i.e. p[i] = partition of vertex i)

◆ RemoveEdges()

template<class TDerivedPtr, class TDerivedAdj, class TIndex = typename TDerivedPtr::Scalar, class Func>
void pbat::graph::RemoveEdges ( Eigen::DenseBase< TDerivedPtr > & ptr,
Eigen::DenseBase< TDerivedAdj > & adj,
Func fShouldDeleteEdge )

In-place removal of edges from the adjacency list in compressed sparse format.

Template Parameters
TDerivedPtrEigen dense expression of the offset pointers of the adjacency list
TDerivedAdjEigen dense expression of the indices of the adjacency list
TIndexIndex type of the graph vertices
FuncCallable type of the edge removal function with signature bool(TIndex i, TIndex j)
Parameters
ptrOffset pointers of the adjacency list
adjIndices of the adjacency list
fShouldDeleteEdgeCallable function to determine if an edge should be removed with signature bool(TIndex i, TIndex j)