incidence

Contents

incidence#

graph_tool.spectral.incidence(g, vindex=None, eindex=None, operator=False, csr=True)[source]#

Return the incidence matrix of the graph.

Parameters:
gGraph

Graph to be used.

vindexVertexPropertyMap (optional, default: None)

Vertex property map specifying the row indices. If not provided, the internal vertex index is used.

eindexEdgePropertyMap (optional, default: None)

Edge property map specifying the column indices. If not provided, the internal edge index is used.

operatorbool (optional, default: False)

If True, a scipy.sparse.linalg.LinearOperator subclass is returned, instead of a sparse matrix.

csrbool (optional, default: True)

If True, and operator is False, a scipy.sparse.csr_matrix sparse matrix is returned, otherwise a scipy.sparse.coo_matrix is returned instead.

Returns:
acsr_matrix or IncidenceOperator

The (sparse) incidence matrix.

Notes

For undirected graphs, the incidence matrix is defined as

\[\begin{split}b_{i,j} = \begin{cases} 1 & \text{if vertex } v_i \text{and edge } e_j \text{ are incident}, \\ 0 & \text{otherwise} \end{cases}\end{split}\]

For directed graphs, the definition is

\[\begin{split}b_{i,j} = \begin{cases} 1 & \text{if edge } e_j \text{ enters vertex } v_i, \\ -1 & \text{if edge } e_j \text{ leaves vertex } v_i, \\ 0 & \text{otherwise} \end{cases}\end{split}\]

LinearOperator vs. sparse matrices

For many linear algebra computations it is more efficient to pass operator=True to this function. In this case, it will return a scipy.sparse.linalg.LinearOperator subclass, which implements matrix-vector and matrix-matrix multiplication, and is sufficient for the sparse linear algebra operations available in the scipy module scipy.sparse.linalg. This avoids copying the whole graph as a sparse matrix, and performs the multiplication operations in parallel (if enabled during compilation) — see note below.

Parallel implementation.

If enabled during compilation, this algorithm will run in parallel using OpenMP. See the parallel algorithms section for information about how to control several aspects of parallelization.

(The above is only applicable if operator == True, and when the object returned is used to perform matrix-vector or matrix-matrix multiplications.)

References

Examples

>>> g = gt.collection.data["polblogs"]
>>> B = gt.incidence(g, operator=True)
>>> N = g.num_vertices()
>>> s1 = scipy.sparse.linalg.svds(B, k=N//2, which="LM", return_singular_vectors=False)
>>> s2 = scipy.sparse.linalg.svds(B, k=N-N//2, which="SM", return_singular_vectors=False)
>>> s = np.concatenate((s1, s2))
>>> s.sort()
>>> figure(figsize=(8, 2))
<...>
>>> plot(s, "s")
[...]
>>> xlabel(r"$i$")
Text(...)
>>> ylabel(r"$\lambda_i$")
Text(...)
>>> tight_layout()
>>> savefig("polblogs-indidence-svd.svg")
../_images/polblogs-indidence-svd.svg

Incidence singular values for the political blogs network.#