PhysicsBasedAnimationToolkit 0.0.10
Cross-platform C++20 library of algorithms and data structures commonly used in computer graphics research on physically-based simulation.
|
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. | |
Graph algorithms and data structures
using pbat::graph::WeightedEdge = Eigen::Triplet<TWeight, TIndex> |
|
strong |
|
strong |
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.
TWeightedEdgeIterator | Iterator type of the edge list |
TWeightedEdge | Type of the edge |
TScalar | Scalar type of the graph edge weights |
TIndex | Index type of the graph vertices |
begin | Iterator to the beginning of the edge list |
end | Iterator to the end of the edge list |
m | Number of rows of the adjacency matrix. If not provided (i.e. m < 0), it is inferred from the edge list. |
n | Number of columns of the adjacency matrix. If not provided (i.e. n < 0), it is inferred from the edge list. |
auto pbat::graph::AdjacencyMatrixIndices | ( | Eigen::SparseCompressedBase< TDerivedA > const & | A | ) |
Non-owning wrapper around the indices of a compressed sparse matrix.
TDerivedA | Type of input adjacency matrix |
TScalar | Scalar type of the graph edge weights |
TIndex | Index type of the graph vertices |
A | Input adjacency matrix |
auto pbat::graph::AdjacencyMatrixPrefix | ( | Eigen::SparseCompressedBase< TDerivedA > const & | A | ) |
Non-owning wrapper around the offset pointers of a compressed sparse matrix.
TDerivedA | Type of input adjacency matrix |
TScalar | Scalar type of the graph edge weights |
TIndex | Index type of the graph vertices |
A | Input adjacency matrix |
auto pbat::graph::AdjacencyMatrixWeights | ( | Eigen::SparseCompressedBase< TDerivedA > const & | A | ) |
Non-owning wrapper around the weights of a compressed sparse matrix.
TDerivedA | Type of input adjacency matrix |
TScalar | Scalar type of the graph edge weights |
TIndex | Index type of the graph vertices |
A | Input adjacency matrix |
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.
TIndex | Index type used in the graph |
TDerivedP | Type of the pointer vector (adjacency list start indices for each vertex) |
TDerivedAdj | Type of the adjacency list (vector of vertex indices) |
ptr | Pointer to the start of each vertex's adjacency list |
adj | Adjacency list of the graph |
components | Output vector to store component labels for each vertex |
bfs | Breadth-first search object |
components[u] < 0
for all vertices u
in the graph 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.
TIndex | Index type used in the graph |
TDerivedP | Type of the pointer vector (adjacency list start indices for each vertex) |
TDerivedAdj | Type of the adjacency list (vector of vertex indices) |
ptr | Pointer to the start of each vertex's adjacency list |
adj | Adjacency list of the graph |
components | Output vector to store component labels for each vertex |
dfs | Depth-first search object |
components[u] < 0
for all vertices u
in the graph void pbat::graph::ForEachEdge | ( | Eigen::DenseBase< TDerivedPtr > & | ptr, |
Eigen::DenseBase< TDerivedAdj > & | adj, | ||
Func && | f ) |
Edge iteration over the adjacency list in compressed sparse format.
TDerivedPtr | Eigen dense expression of the offset pointers of the adjacency list |
TDerivedAdj | Eigen dense expression of the indices of the adjacency list |
TIndex | Index type of the graph vertices |
Func | Callable type of the edge iteration function |
ptr | Offset pointers of the adjacency list |
adj | Indices of the adjacency list |
f | Callable function to be applied to each edge with signature void(TIndex i, TIndex j,
TIndex k) |
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.
TDerivedPtr | Eigen dense expression for offset pointers of the adjacency list |
TDerivedAdj | Eigen dense expression for indices of the adjacency list |
NC | Color palette capacity |
TIndex | Index type for vertices |
ptr | Offset pointers array of the adjacency list |
adj | Indices array of the adjacency list |
eOrderingStrategy | Vertex visiting order strategy |
eSelectionStrategy | Color selection strategy |
|# vertices|
array mapping vertices to their associated color 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.
TIndex | Index type of the graph vertices |
lil | Input adjacency list in list of lists format |
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.
TDerivedP | Type of input map |
TScalar | Scalar type of the graph edge weights |
p | Input map |
n | Number of vertices in the graph. If not provided (i.e. n < 0), it is inferred from the map. |
auto pbat::graph::MatrixToAdjacency | ( | Eigen::SparseCompressedBase< TDerivedA > const & | A | ) |
Obtain the offset and indices arrays of an input adjacency matrix.
TDerivedA | Type of input adjacency matrix |
TScalar | Scalar type of the graph edge weights |
TIndex | Index type of the graph vertices |
A | Input adjacency matrix |
auto pbat::graph::MatrixToWeightedAdjacency | ( | Eigen::SparseCompressedBase< TDerivedA > const & | A | ) |
Construct adjacency list in compressed sparse format from an input adjacency matrix.
TDerivedA | Type of input adjacency matrix |
TScalar | Scalar type of the graph edge weights |
TIndex | Index type of the graph vertices |
A | Input adjacency matrix |
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.
TDerivedE | Eigen dense expression for mesh elements |
TDerivedW | Eigen dense expression for element-vertex weights |
TIndex | Type of indices used in element array |
TScalar | Type of weights |
E | |# nodes per element|x|# elements| array of element indices |
w | |# nodes per element|x|# elements| array of element-vertex weights |
nNodes | Number of nodes in the mesh. If nNodes < 1 , the number of nodes is inferred from E. |
bVertexToElement | If true, the adjacency matrix maps vertices to elements, rather than elements to vertices |
bHasDuplicates | If true, duplicate entries in the input mesh will be handled |
E.rows() == w.rows()
and E.cols() == w.cols()
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.
TDerivedE | Eigen dense expression for mesh elements |
TIndex | Type of indices used in element array |
E | |# nodes per element|x|# elements| array of element indices |
nNodes | Number of nodes in the mesh. If nNodes < 1 , the number of nodes is inferred from E. |
bVertexToElement | If true, the adjacency matrix maps vertices to elements, rather than |
bHasDuplicates | If true, duplicate entries in the input mesh will be handled |
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.
TDerivedE | Eigen dense expression for mesh elements |
TIndex | Type of indices used in element array |
E | |# nodes per element|x|# elements| array of element indices |
nNodes | Number of nodes in the mesh. If nNodes < 1 , the number of nodes is inferred from E. |
opts | Adjacency types to keep in the dual graph |
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.
TDerivedE | Eigen dense expression for mesh elements |
TIndex | Type of indices used in element array |
E | |# nodes per element|x|# elements| array of element indices |
nNodes | Number of nodes in the mesh. If nNodes < 1 , the number of nodes is inferred from E. |
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]
ptr | Offset pointers of adjacency list |
adj | Indices of adjacency list |
wadj | Edge weights of adjacency list |
nPartitions | Number of desired partitions |
opts | Partitioning options |
|# vertices|
array p of partition indices (i.e. p[i] = partition of vertex i) 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.
TDerivedPtr | Eigen dense expression of the offset pointers of the adjacency list |
TDerivedAdj | Eigen dense expression of the indices of the adjacency list |
TIndex | Index type of the graph vertices |
Func | Callable type of the edge removal function with signature bool(TIndex i, TIndex j) |
ptr | Offset pointers of the adjacency list |
adj | Indices of the adjacency list |
fShouldDeleteEdge | Callable function to determine if an edge should be removed with signature bool(TIndex i, TIndex j) |