Filters Objects

Main Filters Class

class pygsp.filters.Filter(G, filters=None, **kwargs)[source]

Bases: object

Parent class for all Filters or Filterbanks, contains the shared methods for those classes.

analysis(s, method=None, cheb_order=30, lanczos_order=30, **kwargs)[source]

Operator to analyse a filterbank

Parameters:

s : ndarray

graph signals to analyse

method : string

wether using an exact method, cheby approx or lanczos

cheb_order : int

Order for chebyshev

Returns:

c : ndarray

Transform coefficients

Examples

>>> import numpy as np
>>> from pygsp import graphs, filters
>>> G = graphs.Logo()
>>> MH = filters.MexicanHat(G)
>>> x = np.arange(G.N**2).reshape(G.N, G.N)
>>> co = MH.analysis(x)

[HVG11]

approx(m, N, **kwargs)[source]
can_dual()[source]

Creates a dual graph form a given graph

evaluate(f, *args, **kwargs)
filterbank_bounds(N=999, bounds=None)[source]

Compute approximate frame bounds for a filterbank.

Parameters:

bounds : interval to compute the bound.

Given in an ndarray: np.array([xmin, xnmax]). By default, bounds is None and filtering is bounded by the eigenvalues of G.

N : Number of point for the line search

Default is 999

Returns:

lower : Filterbank lower bound

upper : Filterbank upper bound

filterbank_matrix()[source]

Create the matrix of the filterbank frame.

This function creates the matrix associated to the filterbank g. The size of the matrix is MN x N, where M is the number of filters.

Returns:F : Frame
inverse(c, **kwargs)[source]
plot(**kwargs)[source]

Plot the filter.

See plotting doc.

synthesis(c, order=30, method=None, **kwargs)[source]

Synthesis operator of a filterbank

Parameters:

G : Graph structure.

c : Transform coefficients

method : Select the method ot be used for the computation.

  • ‘exact’ : Exact method using the graph Fourier matrix
  • ‘cheby’ : Chebyshev polynomial approximation
  • ‘lanczos’ : Lanczos approximation

Default : if the Fourier matrix is present: ‘exact’ otherwise ‘cheby’

order : Degree of the Chebyshev approximation

Default is 30

Returns:

signal : sythesis signal

tighten()[source]
wlog_scales(lmin, lmax, Nscales, t1=1, t2=2)[source]

Compute logarithm scales for wavelets

Parameters:

lmin : int

Minimum non-zero eigenvalue

lmax : int

Maximum eigenvalue

Nscales : int

Number of scales

Returns:

s : ndarray

Scale

Abspline

class pygsp.filters.Abspline(G, Nf=6, lpfactor=20, t=None, **kwargs)[source]

Bases: pygsp.filters.filter.Filter

Abspline Filterbank

Inherits its methods from Filters

Parameters:

G : Graph

Nf : int

Number of filters from 0 to lmax (default = 6)

lpfactor : int

Low-pass factor lmin=lmax/lpfactor will be used to determine scales, the scaling function will be created to fill the lowpass gap. (default = 20)

t : ndarray

Vector of scale to be used (Initialized by default at the value of the log scale)

Returns:

out : Abspline

Examples

>>> from pygsp import graphs, filters
>>> G = graphs.Logo()
>>> F = filters.Abspline(G)

Expwin

class pygsp.filters.Expwin(G, bmax=0.2, a=1.0, **kwargs)[source]

Bases: pygsp.filters.filter.Filter

Expwin Filterbank

Inherits its methods from Filters

Parameters:

G : Graph

bmax : float

Maximum relative band (default = 0.2)

a : int

Slope parameter (default = 1)

Returns:

out : Expwin

Examples

>>> from pygsp import graphs, filters
>>> G = graphs.Logo()
>>> F = filters.Expwin(G)

HalfCosine

class pygsp.filters.HalfCosine(G, Nf=6, **kwargs)[source]

Bases: pygsp.filters.filter.Filter

HalfCosine Filterbank

Inherits its methods from Filters

Parameters:

G : Graph

Nf : int

Number of filters from 0 to lmax (default = 6)

Returns

——-

out : HalfCosine

Examples

>>> from pygsp import graphs, filters
>>> G = graphs.Logo()
>>> F = filters.HalfCosine(G)

Itersine

class pygsp.filters.Itersine(G, Nf=6, overlap=2.0, **kwargs)[source]

Bases: pygsp.filters.filter.Filter

Create a itersine filterbanks

This function create a itersine half overlap filterbank of Nf filters Going from 0 to lambda_max

Parameters:

G : Graph

Nf : int

Number of filters from 0 to lmax. (default = 6)

overlap : int

(default = 2)

Returns:

out : Itersine

Examples

>>> from pygsp import graphs, filters
>>> G = graphs.Logo()
>>> F = filters.Itersine(G)

MexicanHat

class pygsp.filters.MexicanHat(G, Nf=6, lpfactor=20, t=None, normalize=False, **kwargs)[source]

Bases: pygsp.filters.filter.Filter

Mexican hat Filterbank

Inherits its methods from Filters

Parameters:

G : Graph

Nf : int

Number of filters from 0 to lmax (default = 6)

lpfactor : int

Low-pass factor lmin=lmax/lpfactor will be used to determine scales, the scaling function will be created to fill the lowpass gap. (default = 20)

t : ndarray

Vector of scale to be used (Initialized by default at the value of the log scale)

normalize : bool

Wether to normalize the wavelet by the factor/sqrt(t). (default = False)

Returns:

out : MexicanHat

Examples

>>> from pygsp import graphs, filters
>>> G = graphs.Logo()
>>> F = filters.MexicanHat(G)

Meyer

class pygsp.filters.Meyer(G, Nf=6, **kwargs)[source]

Bases: pygsp.filters.filter.Filter

Meyer Filterbank

Inherits its methods from Filters

Parameters:

G : Graph

Nf : int

Number of filters from 0 to lmax (default = 6)

Returns:

out : Meyer

Examples

>>> from pygsp import graphs, filters
>>> G = graphs.Logo()
>>> F = filters.Meyer(G)

SimpleTf

class pygsp.filters.SimpleTf(G, Nf=6, t=None, **kwargs)[source]

Bases: pygsp.filters.filter.Filter

SimpleTf Filterbank

Inherits its methods from Filters

Parameters:

G : Graph

Nf : int

Number of filters from 0 to lmax (default = 6)

t : ndarray

Vector of scale to be used (Initialized by default at the value of the log scale)

Returns:

out : SimpleTf

Examples

>>> from pygsp import graphs, filters
>>> G = graphs.Logo()
>>> F = filters.SimpleTf(G)

WarpedTranslates

class pygsp.filters.WarpedTranslates(G, Nf=6, **kwargs)[source]

Bases: pygsp.filters.filter.Filter

Creates a vertex frequency filterbank

Parameters:

G : Graph

Nf : int

Number of filters

Returns:

out : WarpedTranslates

Examples

Not Implemented for now # >>> from pygsp import graphs, filters # >>> G = graphs.Logo() # >>> F = filters.WarpedTranslates(G)

See [SWHV13]

Papadakis

class pygsp.filters.Papadakis(G, a=0.75, **kwargs)[source]

Bases: pygsp.filters.filter.Filter

Papadakis Filterbank

Inherits its methods from Filters

This function create a parseval filterbank of \(2\). The low-pass filter is defined by a function \(f_l(x)\)

\[\begin{split}f_{l}=\begin{cases} 1 & \mbox{if }x\leq a\\ \sqrt{1-\frac{\sin\left(\frac{3\pi}{2a}x\right)}{2}} & \mbox{if }a<x\leq \frac{5a}{3} \\ 0 & \mbox{if }x>\frac{5a}{3} \end{cases}\end{split}\]

The high pass filter is adaptated to obtain a tight frame.

Parameters:

G : Graph

a : float

See equation above for this parameter The spectrum is scaled between 0 and 2 (default = 3/4)

Returns:

out : Papadakis

Examples

>>> from pygsp import graphs, filters
>>> G = graphs.Logo()
>>> F = filters.Papadakis(G)

Regular

class pygsp.filters.Regular(G, d=3, **kwargs)[source]

Bases: pygsp.filters.filter.Filter

Regular Filterbank

Inherits its methods from Filters

This function creates a parseval filterbank \(2\) filters. The low-pass filter is defined by a function \(f_l(x)\) between \(0\) and \(2\). For \(d = 0\).

\[f_{l}= \sin\left( \frac{\pi}{4} x \right)\]

For \(d = 1\)

\[f_{l}= \sin\left( \frac{\pi}{4} \left( 1+ \sin\left(\frac{\pi}{2}(x-1)\right) \right) \right)\]

For \(d = 2\)

\[f_{l}= \sin\left( \frac{\pi}{4} \left( 1+ \sin\left(\frac{\pi}{2} \sin\left(\frac{\pi}{2}(x-1)\right)\right) \right) \right)\]

And so for other degrees \(d\)

The high pass filter is adaptated to obtain a tight frame.

Parameters:

G : Graph

d : float

See equations above for this parameter Degree (default = 3)

Returns:

out : Regular

Examples

>>> from pygsp import graphs, filters
>>> G = graphs.Logo()
>>> F = filters.Regular(G)

Simoncelli

class pygsp.filters.Simoncelli(G, a=0.6666666666666666, **kwargs)[source]

Bases: pygsp.filters.filter.Filter

Simoncelli Filterbank

Inherits its methods from Filters

This function create a parseval filterbank of \(2\). The low-pass filter is defined by a function \(f_l(x)\).

\[\begin{split}f_{l}=\begin{cases} 1 & \mbox{if }x\leq a\\ \cos\left(\frac{\pi}{2}\frac{\log\left(\frac{x}{2}\right)}{\log(2)}\right) & \mbox{if }a<x\leq2a\\ 0 & \mbox{if }x>2a \end{cases}\end{split}\]

The high pass filter is is adaptated to obtain a tight frame.

Parameters:

G : Graph

a : float

See equation above for this parameter The spectrum is scaled between 0 and 2 (default = 2/3)

Returns:

out : Simoncelli

Examples

>>> from pygsp import graphs, filters
>>> G = graphs.Logo()
>>> F = filters.Simoncelli(G)

Held

class pygsp.filters.Held(G, a=0.6666666666666666, **kwargs)[source]

Bases: pygsp.filters.filter.Filter

Held Filterbank

Inherits its methods from Filters

This function create a parseval filterbank of \(2\) filters. The low-pass filter is defined by a function \(f_l(x)\)

\[\begin{split}f_{l}=\begin{cases} 1 & \mbox{if }x\leq a\\ \sin\left(2\pi\mu\left(\frac{x}{8a}\right)\right) & \mbox{if }a<x\leq2a\\ 0 & \mbox{if }x>2a \end{cases}\end{split}\]

with

\[\mu(x) = -1+24x-144*x^2+256*x^3\]

The high pass filter is adaptated to obtain a tight frame.

Parameters:

G : Graph

a : float

See equation above for this parameter The spectrum is scaled between 0 and 2 (default = 2/3)

Returns:

out : Held

Examples

>>> from pygsp import graphs, filters
>>> G = graphs.Logo()
>>> F = filters.Held(G)

Heat

class pygsp.filters.Heat(G, tau=10, normalize=False, **kwargs)[source]

Bases: pygsp.filters.filter.Filter

Heat Filterbank

Inherits its methods from Filters

Parameters:

G : Graph

tau : int or list of ints

Scaling parameter. (default = 10)

normalize : bool

Normalize the kernel (works only if the eigenvalues are present in the graph). (default = 0)

Returns:

out : Heat

Examples

>>> from pygsp import graphs, filters
>>> G = graphs.Logo()
>>> F = filters.Heat(G)