Reduction

The pygsp.reduction module implements functionalities for the reduction of graphs’ vertex set while keeping the graph structure.

tree_multiresolution(G, Nlevel[, …]) Compute a multiresolution of trees
graph_multiresolution(G, levels[, sparsify, …]) Compute a pyramid of graphs (by Kron reduction).
kron_reduction(G, ind) Compute the Kron reduction.
pyramid_analysis(Gs, f, **kwargs) Compute the graph pyramid transform coefficients.
pyramid_synthesis(Gs, cap, pe[, order]) Synthesize a signal from its pyramid coefficients.
interpolate(G, f_subsampled, keep_inds[, …]) Interpolate a graph signal.
graph_sparsify(M, epsilon[, maxiter]) Sparsify a graph (with Spielman-Srivastava).
pygsp.reduction.graph_multiresolution(G, levels, sparsify=True, sparsify_eps=None, downsampling_method='largest_eigenvector', reduction_method='kron', compute_full_eigen=False, reg_eps=0.005)[source]

Compute a pyramid of graphs (by Kron reduction).

‘graph_multiresolution(G,levels)’ computes a multiresolution of graph by repeatedly downsampling and performing graph reduction. The default downsampling method is the largest eigenvector method based on the polarity of the components of the eigenvector associated with the largest graph Laplacian eigenvalue. The default graph reduction method is Kron reduction followed by a graph sparsification step. param is a structure of optional parameters.

Parameters:

G : Graph structure

The graph to reduce.

levels : int

Number of level of decomposition

lambd : float

Stability parameter. It adds self loop to the graph to give the algorithm some stability (default = 0.025). [UNUSED?!]

sparsify : bool

To perform a spectral sparsification step immediately after the graph reduction (default is True).

sparsify_eps : float

Parameter epsilon used in the spectral sparsification (default is min(10/sqrt(G.N),.3)).

downsampling_method: string

The graph downsampling method (default is ‘largest_eigenvector’).

reduction_method : string

The graph reduction method (default is ‘kron’)

compute_full_eigen : bool

To also compute the graph Laplacian eigenvalues and eigenvectors for every graph in the multiresolution sequence (default is False).

reg_eps : float

The regularized graph Laplacian is \(\bar{L}=L+\epsilon I\). A smaller epsilon may lead to better regularization, but will also require a higher order Chebyshev approximation. (default is 0.005)

Returns:

Gs : list

A list of graph layers.

Examples

>>> from pygsp import reduction
>>> levels = 5
>>> G = graphs.Sensor(N=512)
>>> G.compute_fourier_basis()
>>> Gs = reduction.graph_multiresolution(G, levels, sparsify=False)
>>> for idx in range(levels):
...     Gs[idx].plotting['plot_name'] = 'Reduction level: {}'.format(idx)
...     Gs[idx].plot()
pygsp.reduction.graph_sparsify(M, epsilon, maxiter=10)[source]

Sparsify a graph (with Spielman-Srivastava).

Parameters:

M : Graph or sparse matrix

Graph structure or a Laplacian matrix

epsilon : int

Sparsification parameter

Returns:

Mnew : Graph or sparse matrix

New graph structure or sparse matrix

Notes

Epsilon should be between 1/sqrt(N) and 1

References

See [SS11], [Rud99] and [RV07]. for more informations

Examples

>>> from pygsp import reduction
>>> G = graphs.Sensor(256, Nc=20, distribute=True)
>>> epsilon = 0.4
>>> G2 = reduction.graph_sparsify(G, epsilon)
pygsp.reduction.interpolate(G, f_subsampled, keep_inds, order=100, reg_eps=0.005, **kwargs)[source]

Interpolate a graph signal.

Parameters:

G : Graph

f_subsampled : ndarray

A graph signal on the graph G.

keep_inds : ndarray

List of indices on which the signal is sampled.

order : int

Degree of the Chebyshev approximation (default = 100).

reg_eps : float

The regularized graph Laplacian is $bar{L}=L+epsilon I$. A smaller epsilon may lead to better regularization, but will also require a higher order Chebyshev approximation.

Returns:

f_interpolated : ndarray

Interpolated graph signal on the full vertex set of G.

References

See [Pes09]

pygsp.reduction.kron_reduction(G, ind)[source]

Compute the Kron reduction.

This function perform the Kron reduction of the weight matrix in the graph G, with boundary nodes labeled by ind. This function will create a new graph with a weight matrix Wnew that contain only boundary nodes and is computed as the Schur complement of the original matrix with respect to the selected indices.

Parameters:

G : Graph or sparse matrix

Graph structure or weight matrix

ind : list

indices of the nodes to keep

Returns:

Gnew : Graph or sparse matrix

New graph structure or weight matrix

References

See [DB13]

pygsp.reduction.pyramid_analysis(Gs, f, **kwargs)[source]

Compute the graph pyramid transform coefficients.

Parameters:

Gs : list of graphs

A multiresolution sequence of graph structures.

f : ndarray

Graph signal to analyze.

h_filters : list

A list of filter that will be used for the analysis and sythesis operator. If only one filter is given, it will be used for all levels. Default is h(x) = 1 / (2x+1)

Returns:

ca : ndarray

Coarse approximation at each level

pe : ndarray

Prediction error at each level

h_filters : list

Graph spectral filters applied

References

See [SFV13] and [Pes09].

pygsp.reduction.pyramid_synthesis(Gs, cap, pe, order=30, **kwargs)[source]

Synthesize a signal from its pyramid coefficients.

Parameters:

Gs : Array of Graphs

A multiresolution sequence of graph structures.

cap : ndarray

Coarsest approximation of the original signal.

pe : ndarray

Prediction error at each level.

use_exact : bool

To use exact graph spectral filtering instead of the Chebyshev approximation.

order : int

Degree of the Chebyshev approximation (default=30).

least_squares : bool

To use the least squares synthesis (default=False).

h_filters : ndarray

The filters used in the analysis operator. These are required for least squares synthesis, but not for the direct synthesis method.

use_landweber : bool

To use the Landweber iteration approximation in the least squares synthesis.

reg_eps : float

Interpolation parameter.

landweber_its : int

Number of iterations in the Landweber approximation for least squares synthesis.

landweber_tau : float

Parameter for the Landweber iteration.

Returns:

reconstruction : ndarray

The reconstructed signal.

ca : ndarray

Coarse approximations at each level

pygsp.reduction.tree_multiresolution(G, Nlevel, reduction_method='resistance_distance', compute_full_eigen=False, root=None)[source]

Compute a multiresolution of trees

Parameters:

G : Graph

Graph structure of a tree.

Nlevel : Number of times to downsample and coarsen the tree

root : int

The index of the root of the tree. (default = 1)

reduction_method : str

The graph reduction method (default = ‘resistance_distance’)

compute_full_eigen : bool

To also compute the graph Laplacian eigenvalues for every tree in the sequence

Returns:

Gs : ndarray

Ndarray, with each element containing a graph structure represent a reduced tree.

subsampled_vertex_indices : ndarray

Indices of the vertices of the previous tree that are kept for the subsequent tree.