Skip to content

Commit

Permalink
get edge list of directed graphs
Browse files Browse the repository at this point in the history
  • Loading branch information
mdeff committed Nov 30, 2018
1 parent de21533 commit 49fa17a
Show file tree
Hide file tree
Showing 4 changed files with 76 additions and 50 deletions.
8 changes: 4 additions & 4 deletions pygsp/graphs/difference.py
Original file line number Diff line number Diff line change
Expand Up @@ -54,19 +54,19 @@ def compute_differential_operator(self):
"""

v_in, v_out, weights = self.get_edge_list()
sources, targets, weights = self.get_edge_list()

n = self.n_edges
rows = np.concatenate([v_in, v_out])
rows = np.concatenate([sources, targets])
columns = np.concatenate([np.arange(n), np.arange(n)])
values = np.empty(2*n)

if self.lap_type == 'combinatorial':
values[:n] = np.sqrt(weights)
values[n:] = -values[:n]
elif self.lap_type == 'normalized':
values[:n] = np.sqrt(weights / self.dw[v_in])
values[n:] = -np.sqrt(weights / self.dw[v_out])
values[:n] = np.sqrt(weights / self.dw[sources])
values[n:] = -np.sqrt(weights / self.dw[targets])
else:
raise ValueError('Unknown lap_type {}'.format(self.lap_type))

Expand Down
71 changes: 52 additions & 19 deletions pygsp/graphs/graph.py
Original file line number Diff line number Diff line change
Expand Up @@ -646,37 +646,70 @@ def estimate_lmax(self, recompute=False):
def get_edge_list(self):
r"""Return an edge list, an alternative representation of the graph.
The weighted adjacency matrix is the canonical form used in this
package to represent a graph as it is the easiest to work with when
considering spectral methods.
Each edge :math:`e_k = (v_i, v_j) \in \mathcal{E}` from :math:`v_i` to
:math:`v_j` is associated with the weight :math:`W[i, j]`. For each
edge :math:`e_k`, the method returns :math:`(i, j, W[i, j])` as
`(sources[k], targets[k], weights[k])`, with :math:`i \in [0,
|\mathcal{V}|-1], j \in [0, |\mathcal{V}|-1], k \in [0,
|\mathcal{E}|-1]`.
Returns
-------
v_in : vector of int
v_out : vector of int
sources : vector of int
Source node indices.
targets : vector of int
Target node indices.
weights : vector of float
Edge weights.
Notes
-----
The weighted adjacency matrix is the canonical form used in this
package to represent a graph as it is the easiest to work with when
considering spectral methods.
Edge orientation (i.e., which node is the source or the target) is
arbitrary for undirected graphs.
The implementation uses the low triangular part of the adjacency
matrix, hence :math:`j \leq i \ \forall k`.
Examples
--------
>>> G = graphs.Logo()
>>> v_in, v_out, weights = G.get_edge_list()
>>> v_in.shape, v_out.shape, weights.shape
((3131,), (3131,), (3131,))
"""
Edge list of a directed graph.
if self.is_directed():
raise NotImplementedError('Directed graphs not supported yet.')
>>> adjacency = np.array([
... [0, 3, 0],
... [3, 0, 4],
... [0, 0, 0],
... ])
>>> graph = graphs.Graph(adjacency)
>>> sources, targets, weights = graph.get_edge_list()
>>> print(sources, targets, weights)
[0 1 1] [1 0 2] [3 3 4]
else:
v_in, v_out = sparse.tril(self.W).nonzero()
weights = self.W[v_in, v_out]
weights = np.asarray(weights).squeeze()
Edge list of an undirected graph.
>>> adjacency = np.array([
... [0, 3, 0],
... [3, 0, 4],
... [0, 4, 0],
... ])
>>> graph = graphs.Graph(adjacency)
>>> sources, targets, weights = graph.get_edge_list()
>>> print(sources, targets, weights)
[1 2] [0 1] [3 4]
"""

W = self.W if self.is_directed() else sparse.tril(self.W)

# TODO G.ind_edges = sub2ind(size(G.W), G.v_in, G.v_out)
sources, targets = W.nonzero()
weights = self.W[sources, targets]
weights = np.asarray(weights).squeeze()

assert self.Ne == v_in.size == v_out.size == weights.size
return v_in, v_out, weights
assert self.n_edges == sources.size == targets.size == weights.size
return sources, targets, weights

def plot(self, vertex_color=None, vertex_size=None, highlight=[],
edges=None, edge_color=None, edge_width=None,
Expand Down
31 changes: 11 additions & 20 deletions pygsp/plotting.py
Original file line number Diff line number Diff line change
Expand Up @@ -346,6 +346,11 @@ def _plot_graph(G, vertex_color, vertex_size, highlight,
ax : :class:`matplotlib.axes.Axes`
The axes the plot belongs to. Only with the matplotlib backend.
Notes
-----
The orientation of directed edges is not shown. If edges exist in both
directions, they will be drawn on top of each other.
Examples
--------
>>> import matplotlib
Expand Down Expand Up @@ -448,8 +453,6 @@ def is_single_color(color):
if title is None:
title = G.__repr__(limit=4)

G = _handle_directed(G)

if backend == 'pyqtgraph':
if vertex_color is None:
_qtg_plot_graph(G, edges=edges, vertex_size=vertex_size,
Expand All @@ -476,10 +479,10 @@ def _plt_plot_graph(G, vertex_color, vertex_size, highlight,

if edges and (G.coords.ndim != 1): # No edges for 1D plots.

v_in, v_out, _ = G.get_edge_list()
sources, targets, _ = G.get_edge_list()
edges = [
G.coords[v_in],
G.coords[v_out],
G.coords[sources],
G.coords[targets],
]
edges = np.stack(edges, axis=1)

Expand Down Expand Up @@ -660,28 +663,16 @@ def _plot_spectrogram(G, node_idx):

def _get_coords(G, edge_list=False):

v_in, v_out, _ = G.get_edge_list()
sources, targets, _ = G.get_edge_list()

if edge_list:
return np.stack((v_in, v_out), axis=1)
return np.stack((sources, targets), axis=1)

coords = [np.stack((G.coords[v_in, d], G.coords[v_out, d]), axis=0)
coords = [np.stack((G.coords[sources, d], G.coords[targets, d]), axis=0)
for d in range(G.coords.shape[1])]

if G.coords.shape[1] == 2:
return coords

elif G.coords.shape[1] == 3:
return [coord.reshape(-1, order='F') for coord in coords]


def _handle_directed(G):
# FIXME: plot edge direction. For now we just symmetrize the weight matrix.
if not G.is_directed():
return G
else:
from pygsp import graphs
G2 = graphs.Graph(utils.symmetrize(G.W))
G2.coords = G.coords
G2.plotting = G.plotting
return G2
16 changes: 9 additions & 7 deletions pygsp/tests/test_graphs.py
Original file line number Diff line number Diff line change
Expand Up @@ -137,12 +137,14 @@ def test_fourier_transform(self):
np.testing.assert_allclose(s, s_star)

def test_edge_list(self):
G = graphs.StochasticBlockModel(N=100, directed=False)
v_in, v_out, weights = G.get_edge_list()
self.assertEqual(G.W[v_in[42], v_out[42]], weights[42])

G = graphs.StochasticBlockModel(N=100, directed=True)
self.assertRaises(NotImplementedError, G.get_edge_list)
for directed in [False, True]:
G = graphs.ErdosRenyi(100, directed=directed)
sources, targets, weights = G.get_edge_list()
if not directed:
self.assertTrue(np.all(sources >= targets))
edges = np.arange(G.n_edges)
np.testing.assert_equal(G.W[sources[edges], targets[edges]],
weights[edges][np.newaxis, :])

def test_differential_operator(self):
G = graphs.StochasticBlockModel(N=100, directed=False)
Expand All @@ -153,7 +155,7 @@ def test_differential_operator(self):
np.testing.assert_allclose(L.toarray(), G.L.toarray())

G = graphs.StochasticBlockModel(N=100, directed=True)
self.assertRaises(NotImplementedError, G.compute_differential_operator)
# self.assertRaises(NotImplementedError, G.compute_differential_operator)

def test_difference(self):
for lap_type in ['combinatorial', 'normalized']:
Expand Down

0 comments on commit 49fa17a

Please sign in to comment.