pyfaust.tools module

The module for pyfaust misc tools.

pyfaust.tools.UpdateCholesky(R, P, Pt, index, m)
Args:
R: (scipy.sparse.csr_matrix or a np.ndarray)

A triangular matrix

P: (callable)

A function to implement the matrix product.

Pt: (callable)

A function to implement the matrix transpose product.

index: (list[int])

a list of all non-zero indices of length R.shape[1]+1

m: (int)

dimension of space from which P maps.

Returns:

R is a triangular matrix such that R’*R = Pt(index)*P(index)

pyfaust.tools.UpdateCholeskyFull(R, P, Pt, index, m)

Same as UpdateCholesky() but R is a np.ndarray.

pyfaust.tools.UpdateCholeskySparse(R, P, Pt, index, m)

Same as UpdateCholesky() but R is a scipy.sparse.csr_matrix.

pyfaust.tools.greed_omp_chol(y, D, maxiter=None, tol=0, relerr=True, verbose=False)

Runs the greedy OMP algorithm optimized by Cholesky decomposition.

Args:
y: (np.ndarray)

the vector to approximate by D@x.

D: (np.ndarray or a Faust)

the dictionary as a numpy array or a Faust.

maxiter: (int or NoneType)

the maximum number of iterations of the algorithm. By default (None) it’s y’s dimension: max(y.shape).

tol: (float)

the tolerance error to reach for the algorithm to stop. By default, it’s zero for not stopping on error criterion.

relerr: (bool)

the type of error stopping criterion. Default to True to use relative error, otherwise (False) the absolute error is used.

verbose: (bool)

to enable the verbosity (True).

Returns:

x: the solution of y = D@x (according to the error).

Example:
>>> from pyfaust.tools import omp
>>> from scipy.sparse import random
>>> from pyfaust import rand, seed
>>> import numpy as np
>>> np.random.seed(42) # just for reproducibility
>>> seed(42)
>>> D = rand(1024, 1024)
>>> x0 = random(1024, 1, .01)
>>> y = D @ x0
>>> maxiter = 17
>>> x = omp(y, D, maxiter, tol=10**-16)
Stopping. Exact signal representation found!
>>> # omp() runs at most maxiter iterations until the error tolerance is
>>> # reached
pyfaust.tools.omp(y, D, maxiter=None, tol=0, relerr=True, verbose=False)

Runs the greedy OMP algorithm optimized by Cholesky decomposition.

Args:
y: (np.ndarray)

the vector to approximate by D@x.

D: (np.ndarray or a Faust)

the dictionary as a numpy array or a Faust.

maxiter: (int or NoneType)

the maximum number of iterations of the algorithm. By default (None) it’s y’s dimension: max(y.shape).

tol: (float)

the tolerance error to reach for the algorithm to stop. By default, it’s zero for not stopping on error criterion.

relerr: (bool)

the type of error stopping criterion. Default to True to use relative error, otherwise (False) the absolute error is used.

verbose: (bool)

to enable the verbosity (True).

Returns:

x: the solution of y = D@x (according to the error).

Example:
>>> from pyfaust.tools import omp
>>> from scipy.sparse import random
>>> from pyfaust import rand, seed
>>> import numpy as np
>>> np.random.seed(42) # just for reproducibility
>>> seed(42)
>>> D = rand(1024, 1024)
>>> x0 = random(1024, 1, .01)
>>> y = D @ x0
>>> maxiter = 17
>>> x = omp(y, D, maxiter, tol=10**-16)
Stopping. Exact signal representation found!
>>> # omp() runs at most maxiter iterations until the error tolerance is
>>> # reached