-  3.39.22
Namespaces | Classes | Functions | Variables
pyfaust Namespace Reference

The FAuST Python Wrapper More...

Namespaces

 demo
 The pyfaust demo module.
 
 fact
 The pyfaust factorization module.
 
 factparams
 The module for the parameterization of FAuST's algorithms (Palm4MSA and Hierarchical Factorization). See also: pyfaust.fact.hierarchical, pyfaust.fact.palm4msa.
 
 poly
 The pyfaust module for polynomial basis as Faust objects.
 
 proj
 This module provides matrix projectors.
 
 tools
 The pyfaust tools module.
 

Classes

class  Faust
 FAuST Python wrapper main class for using multi-layer sparse transforms. More...
 
class  FaustMulMode
 Enumeration class of all matrix chain multiplication methods available to multiply a Faust to a matrix or to compute Faust.toarray(). More...
 

Functions

def version ()
 Returns the FAuST package version. More...
 
def faust_fact (*args, **kwargs)
 This function is a shorthand for pyfaust.fact.hierarchical. More...
 
def license ()
 Prints the FAuST license. More...
 
def norm (F, ord='fro', **kwargs)
 Returns Faust.norm(F, ord) or numpy.linalg.norm(F, ord) depending of F type. More...
 
def dot (A, B, **kwargs)
 Returns Faust.dot(A,B) if A or B is a Faust object, returns numpy.dot(A,B) ortherwise. More...
 
def pinv (F)
 A package function alias for the member function Faust.pinv(). More...
 
def concatenate (_tuple, *args, axis=0, **kwargs)
 A package function alias for the member function Faust.concatenate. More...
 
def hstack (_tuple)
 Concatenates horizontally Faust-s and/or numpy.ndarray objects using Faust.concatenate(). More...
 
def vstack (_tuple)
 Concatenates vertically Faust-s and/or numpy.ndarray arrays using Faust.concatenate(). More...
 
def isFaust (obj)
 Package alias function of Faust.isFaust. More...
 
def wht (n, normed=True, dev="cpu", dtype='float64')
 Constructs a Faust implementing the Walsh-Hadamard Transform (WHT) of order n. More...
 
def bitrev_perm (n)
 Bitreversal permutation. More...
 
def dft (n, normed=True, dev='cpu', diag_opt=False)
 Constructs a Faust F implementing the Discrete Fourier Transform (DFT) of order n. More...
 
def dct (n, normed=True, dev='cpu', dtype='float64')
 Constructs a Faust implementing the Direct Cosine Transform (Type II) Faust of order n. More...
 
def dst (n, normed=True, dev='cpu', dtype='float64')
 Constructs a Faust implementing the Direct Sine Transform (Type II) Faust of order n. More...
 
def circ (c, dev='cpu', diag_opt=False)
 Returns a circulant Faust C defined by the vector c (which is the first column of C.toarray()). More...
 
def anticirc (c, dev='cpu', diag_opt=False)
 Returns an anti-circulant Faust A defined by the vector c (which is the last column of A.toarray()). More...
 
def toeplitz (c, r=None, dev='cpu', diag_opt=False)
 Constructs a toeplitz Faust whose first column is c and first row r. More...
 
def eye (m, n=None, dtype='float64', dev="cpu")
 Faust identity. More...
 
def rand_bsr (num_rows, num_cols, bnrows, bncols, num_factors=None, density=.1, dev='cpu', dtype='float64')
 Generates a random Faust composed only of BSR matrices. More...
 
def rand (num_rows, num_cols, num_factors=None, dim_sizes=None, density=None, fac_type='sparse', per_row=True, dev='cpu', dtype='float64', field=None, seed=0)
 Generates a random Faust. More...
 
def rand_butterfly (n, dtype='float64', dev='cpu', diag_opt=False)
 Constructs a Faust corresponding to the product of log2(n) square factors of size n with butterfly supports and random nonzero coefficients. More...
 
def opt_butterfly_faust (F)
 Optimizes any Faust composed of butterfly factors. More...
 
def enable_gpu_mod (libpaths=None, backend='cuda', silent=False, fatal=False)
 This function loads explicitly the gpu_mod library in memory. More...
 
def is_gpu_mod_enabled ()
 Returns True if the gpu_mod plug-in has been loaded correctly, False otherwise. More...
 
def is_gpu_mod_working ()
 This function returns True if gpu_mod is working properly False otherwise. More...
 
def seed (s)
 (Re)Initializes the pyfaust pseudo-random generator. More...
 
def faust_logo ()
 Generates the FAµST logo and returns it as a Faust. More...
 

Variables

bool WARNING_ON_C_CONTIGUOUS_MATMUL = True
 

Detailed Description

The FAuST Python Wrapper

Function Documentation

◆ anticirc()

def pyfaust.anticirc (   c,
  dev = 'cpu',
  diag_opt = False 
)

Returns an anti-circulant Faust A defined by the vector c (which is the last column of A.toarray()).

Parameters
c(np.ndarray) the vector to define the anti-circulant Faust.
dev(str) the device on which the Faust is created, 'cpu' (default) or 'gpu'.
diag_opt(bool) cf. pyfaust.circ.

Examples

>>> from pyfaust import anticirc
>>> import numpy as np
>>> c = [1, 2, 3, 4, 5, 6, 7, 8]
>>> A = anticirc(c)
>>> A

Faust size 8x8, density 1.5, nnz_sum 96, 6 factor(s):

  • FACTOR 0 (complex) SPARSE, size 8x8, density 0.25, nnz 16
  • FACTOR 1 (complex) SPARSE, size 8x8, density 0.25, nnz 16
  • FACTOR 2 (complex) SPARSE, size 8x8, density 0.25, nnz 16
  • FACTOR 3 (complex) SPARSE, size 8x8, density 0.25, nnz 16
  • FACTOR 4 (complex) SPARSE, size 8x8, density 0.25, nnz 16
  • FACTOR 5 (complex) SPARSE, size 8x8, density 0.25, nnz 16
>>> np.allclose(A.toarray()[:, -1], c)
True
>>> np.real(A.toarray())
array([[2., 3., 4., 5., 6., 7., 8., 1.],
[3., 4., 5., 6., 7., 8., 1., 2.],
[4., 5., 6., 7., 8., 1., 2., 3.],
[5., 6., 7., 8., 1., 2., 3., 4.],
[6., 7., 8., 1., 2., 3., 4., 5.],
[7., 8., 1., 2., 3., 4., 5., 6.],
[8., 1., 2., 3., 4., 5., 6., 7.],
[1., 2., 3., 4., 5., 6., 7., 8.]])
>>> # Look at the density of a larger anticirculant Faust
>>> # it indicates a speedup of the Faust-matrix/vector product
>>> anticirc(np.random.rand(1024)).density()
0.0390625
See also
pyfaust.circ, pyfaust.toeplitz

◆ bitrev_perm()

def pyfaust.bitrev_perm (   n)

Bitreversal permutation.

Parameters
n(int) the size of the permutation, it must be a power of two. P dimensions will be n x n
Returns
P a scipy csr_matrix defining the bit-reversal permutation.
See also
pyfaust.dft

◆ circ()

def pyfaust.circ (   c,
  dev = 'cpu',
  diag_opt = False 
)

Returns a circulant Faust C defined by the vector c (which is the first column of C.toarray()).

Parameters
c(np.ndarray) the vector to define the circulant Faust.
dev(str) the device on which the Faust is created, 'cpu' (default) or 'gpu'.
diag_opt(bool) if True then the returned Faust is optimized using pyfaust.opt_butterfly_faust (because the DFT is used to implement circ).

Examples

>>> from pyfaust import circ
>>> import numpy as np
>>> c = [1, 2, 3, 4, 5, 6, 7, 8]
>>> C = circ(c)
>>> C

Faust size 8x8, density 1.5, nnz_sum 96, 6 factor(s):

  • FACTOR 0 (complex) SPARSE, size 8x8, density 0.25, nnz 16
  • FACTOR 1 (complex) SPARSE, size 8x8, density 0.25, nnz 16
  • FACTOR 2 (complex) SPARSE, size 8x8, density 0.25, nnz 16
  • FACTOR 3 (complex) SPARSE, size 8x8, density 0.25, nnz 16
  • FACTOR 4 (complex) SPARSE, size 8x8, density 0.25, nnz 16
  • FACTOR 5 (complex) SPARSE, size 8x8, density 0.25, nnz 16
>>> np.allclose(C.toarray()[:, 0], c)
True
>>> np.real(C.toarray())
array([[1., 8., 7., 6., 5., 4., 3., 2.],
[2., 1., 8., 7., 6., 5., 4., 3.],
[3., 2., 1., 8., 7., 6., 5., 4.],
[4., 3., 2., 1., 8., 7., 6., 5.],
[5., 4., 3., 2., 1., 8., 7., 6.],
[6., 5., 4., 3., 2., 1., 8., 7.],
[7., 6., 5., 4., 3., 2., 1., 8.],
[8., 7., 6., 5., 4., 3., 2., 1.]])
>>> # Look at the density of a larger circulant Faust
>>> # it indicates a speedup of the Faust-matrix/vector product
>>> circ(np.random.rand(1024)).density()
0.0390625
See also
scipy.linalg.circulant, pyfaust.anticirc, pyfaust.toeplitz

◆ concatenate()

def pyfaust.concatenate (   _tuple,
args,
  axis = 0,
**  kwargs 
)

A package function alias for the member function Faust.concatenate.

Examples

>>> from pyfaust import *
>>> seed(42) # just for reproducibility
>>> F1 = rand(5,50)
>>> F2 = rand(5,50)
>>> concatenate((F1, F2), axis=0)

Faust size 10x50, density 2.17, nnz_sum 1085, 6 factor(s):

  • FACTOR 0 (double) SPARSE, size 10x67, density 0.0746269, nnz 50
  • FACTOR 1 (double) SPARSE, size 67x33, density 0.151515, nnz 335
  • FACTOR 2 (double) SPARSE, size 33x31, density 0.16129, nnz 165
  • FACTOR 3 (double) SPARSE, size 31x56, density 0.0892857, nnz 155
  • FACTOR 4 (double) SPARSE, size 56x100, density 0.05, nnz 280
  • FACTOR 5 (double) SPARSE, size 100x50, density 0.02, nnz 100
See also
numpy.concatenate, Faust.concatenate,

◆ dct()

def pyfaust.dct (   n,
  normed = True,
  dev = 'cpu',
  dtype = 'float64' 
)

Constructs a Faust implementing the Direct Cosine Transform (Type II) Faust of order n.

The analytical formula of DCT II used here is: \(2 \sum_{i=0}^{n-1} x_i cos \left( {\pi k (2i + 1)} \over {2n} \right)\)

Parameters
n(int) the order of the DCT (must be a power of two).
normed(bool) default to True to normalize the DFT Faust as if you called Faust.normalize() and False otherwise.
dev(str) the device on which the Faust is created.
dtype(str) 'float64' (default) or 'float32'.
Returns
The DCT Faust.

Examples

>>> from pyfaust import dct
>>> from scipy.fft import dct as scipy_dct
>>> import numpy as np
>>> DCT8 = dct(8, normed=False)
>>> x = np.arange(8).astype('double')
>>> np.real(DCT8@x)
array([ 56. , -25.76929209, 0. , -2.6938192 ,
0. , -0.80361161, 0. , -0.20280929])
>>> scipy_dct(x)
array([ 56. , -25.76929209, 0. , -2.6938192 ,
0. , -0.80361161, 0. , -0.20280929])
>>> np.allclose(DCT8@x, scipy_dct(x))
True
>>> # check the density with a larger DCT Faust of size 1024
>>> dct(1024).density()
0.076171875
>>> # it is smaller than 1
See also
pyfaust.dft, pyfaust.dst, scipy.fft.dct, pyfaust.fact.butterfly, pyfaust.rand_butterfly, pyfaust.Faust.density

◆ dft()

def pyfaust.dft (   n,
  normed = True,
  dev = 'cpu',
  diag_opt = False 
)

Constructs a Faust F implementing the Discrete Fourier Transform (DFT) of order n.

The factorization corresponds to the butterfly structure of the Cooley-Tukey FFT algorithm. The resulting Faust is complex and has (log2(n)+1) sparse factors. The log2(n) first has 2 nonzeros per row and per column. The last factor is a bit-reversal permutation matrix.

Parameters
n(int) order of the Discrete Fourier Transform (must be a power of two).
normed(bool) default to True to normalize the DFT Faust as if you called Faust.normalize() and False otherwise.
dev(str) device to create the Faust on ('cpu' or 'gpu').
diag_opt(bool) if True then the returned Faust is optimized using pyfaust.opt_butterfly_faust.
Returns
The Faust implementing the DFT of dimension n.

Examples

>>> from pyfaust import dft
>>> dft(1024)

Faust size 1024x1024, density 0.0205078, nnz_sum 21504, 11 factor(s):

  • FACTOR 0 (complex) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 1 (complex) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 2 (complex) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 3 (complex) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 4 (complex) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 5 (complex) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 6 (complex) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 7 (complex) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 8 (complex) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 9 (complex) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 10 (complex) SPARSE, size 1024x1024, density 0.000976562, nnz 1024
>>> dft(1024, normed=True) # is equiv. to next call

Faust size 1024x1024, density 0.0205078, nnz_sum 21504, 11 factor(s):

  • FACTOR 0 (complex) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 1 (complex) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 2 (complex) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 3 (complex) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 4 (complex) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 5 (complex) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 6 (complex) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 7 (complex) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 8 (complex) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 9 (complex) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 10 (complex) SPARSE, size 1024x1024, density 0.000976562, nnz 1024
>>> dft(1024, normed=False).normalize()

Faust size 1024x1024, density 0.0205078, nnz_sum 21504, 11 factor(s):

  • FACTOR 0 (complex) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 1 (complex) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 2 (complex) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 3 (complex) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 4 (complex) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 5 (complex) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 6 (complex) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 7 (complex) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 8 (complex) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 9 (complex) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 10 (complex) SPARSE, size 1024x1024, density 0.000976562, nnz 1024
See also
pyfaust.bitrev, pyfaust.wht, pyfaust.dct, pyfaust.dst, scipy.fft.fft, pyfaust.fact.butterfly, pyfaust.rand_butterfly.

◆ dot()

def pyfaust.dot (   A,
  B,
**  kwargs 
)

Returns Faust.dot(A,B) if A or B is a Faust object, returns numpy.dot(A,B) ortherwise.

See also
Faust.norm,

◆ dst()

def pyfaust.dst (   n,
  normed = True,
  dev = 'cpu',
  dtype = 'float64' 
)

Constructs a Faust implementing the Direct Sine Transform (Type II) Faust of order n.

The analytical formula of DST II used here is: \(2 \sum_{i=0}^{n-1} x_i sin \left( {\pi (k+1) (2i + 1)} \over {2n} \right)\)

Parameters
n(int) the order of the DST (must be a power of two).
normed(bool) default to True to normalize the Hadamard Faust as if you called Faust.normalize() and False otherwise.
dev(str) the device on which the Faust is created.
dtype(str) 'float64' (default) or 'float32'.
Returns
The DST Faust.

Examples

>>> from pyfaust import dst
>>> from scipy.fft import dst as scipy_dst
>>> import numpy as np
>>> DST8 = dst(8, normed=False)
>>> x = np.ones(8)
>>> np.real(DST8@x)
array([1.02516618e+01, 4.93038066e-32, 3.59990489e+00, 3.94430453e-31,
2.40537955e+00, 9.86076132e-32, 2.03918232e+00, 0.00000000e+00])
>>> scipy_dst(x)
array([10.25166179, 0. , 3.59990489, 0. , 2.40537955,
0. , 2.03918232, 0. ])
>>> np.allclose(DST8@x, scipy_dst(x))
True
>>> # check the density with a larger DST Faust of size 1024
>>> dst(1024).density()
0.201171875
>>> # it is smaller than 1
See also
pyfaust.dft, pyfaust.dct, scipy.fft.dst, pyfaust.fact.butterfly, pyfaust.rand_butterfly, pyfaust.Faust.density

◆ enable_gpu_mod()

def pyfaust.enable_gpu_mod (   libpaths = None,
  backend = 'cuda',
  silent = False,
  fatal = False 
)

This function loads explicitly the gpu_mod library in memory.

Normally it's not required to load the library manually, but it could be useful to set a non-default path or to diagnose a loading issue.

Parameters
libpaths(list[str]) the absolute or relative paths where to search the dynamic library (gm) to load. By default, it's none to auto-find the library (if possible).
backend(str) the GPU backend to use, only 'cuda' is available for now.
silent(bool) if True nothing or almost will be displayed on loading, otherwise all messages are visible.

◆ eye()

def pyfaust.eye (   m,
  n = None,
  dtype = 'float64',
  dev = "cpu" 
)

Faust identity.

Parameters
m(int) number of rows,
n(int) number of columns, set to m by default.
dtype(str) the dtype of the identity ('float32', the default 'float64'/'double', or 'complex'/'complex128').

Examples

>>> from pyfaust import eye
>>> eye(5)

Faust size 5x5, density 0.2, nnz_sum 5, 1 factor(s):

  • FACTOR 0 (double) SPARSE, size 5x5, density 0.2, nnz 5 identity matrix flag
>>> eye(5, 4)

Faust size 5x4, density 0.2, nnz_sum 4, 1 factor(s):

  • FACTOR 0 (double) SPARSE, size 5x4, density 0.2, nnz 4 identity matrix flag
>>> eye(5, dtype='complex')

Faust size 5x5, density 0.2, nnz_sum 5, 1 factor(s):

  • FACTOR 0 (complex) SPARSE, size 5x5, density 0.2, nnz 5 identity matrix flag

◆ faust_fact()

def pyfaust.faust_fact ( args,
**  kwargs 
)

This function is a shorthand for pyfaust.fact.hierarchical.

See also
pyfaust.fact.hierarchical

◆ faust_logo()

def pyfaust.faust_logo ( )

Generates the FAµST logo and returns it as a Faust.

Examples

>>> import pyfaust as pf
>>> pf.seed(42)
>>> logo = pf.faust_logo()
>>> logo

Faust size 50x50, density 1.3412, nnz_sum 3353, 5 factor(s):

  • FACTOR 0 (double) DENSE, size 50x50, density 0.2372, nnz 593
  • FACTOR 1 (double) DENSE, size 50x50, density 0.3228, nnz 807
  • FACTOR 2 (double) DENSE, size 50x50, density 0.2432, nnz 608
  • FACTOR 3 (double) DENSE, size 50x50, density 0.3564, nnz 891
  • FACTOR 4 (double) DENSE, size 50x50, density 0.1816, nnz 454
>>> # logo.imshow()
>>> # import matplotlib.pyplot as plt
>>> # plot.show()

◆ hstack()

def pyfaust.hstack (   _tuple)

Concatenates horizontally Faust-s and/or numpy.ndarray objects using Faust.concatenate().

See also numpy.hstack

Note
it could be wiser to encapsulate a Faust in a lazylinop.LazyLinearOp for a lazy concatenation.

Examples

>>> from pyfaust import *
>>> seed(42) # just for reproducibility
>>> F1 = rand(5,50)
>>> F2 = rand(5,50)
>>> hstack((F1, F2))

Faust size 5x100, density 1.99, nnz_sum 995, 6 factor(s):

  • FACTOR 0 (double) SPARSE, size 5x10, density 0.2, nnz 10
  • FACTOR 1 (double) SPARSE, size 10x67, density 0.0746269, nnz 50
  • FACTOR 2 (double) SPARSE, size 67x33, density 0.151515, nnz 335
  • FACTOR 3 (double) SPARSE, size 33x31, density 0.16129, nnz 165
  • FACTOR 4 (double) SPARSE, size 31x56, density 0.0892857, nnz 155
  • FACTOR 5 (double) SPARSE, size 56x100, density 0.05, nnz 280

◆ is_gpu_mod_enabled()

def pyfaust.is_gpu_mod_enabled ( )

Returns True if the gpu_mod plug-in has been loaded correctly, False otherwise.

See also
pyfaust.is_gpu_mod_working

◆ is_gpu_mod_working()

def pyfaust.is_gpu_mod_working ( )

This function returns True if gpu_mod is working properly False otherwise.

is_gpu_mod_working comes as a complement of pyfaust.is_gpu_mod_enabled The latter ensures that gpu_mod shared library/plugin is properly loaded in memory but doesn't ensure that the GPU is available (for example, the NVIDIA driver might not be installed). The former ensures both that the gpu_mod is loaded and the GPU (device 0) is properly available for computing.

◆ isFaust()

def pyfaust.isFaust (   obj)

Package alias function of Faust.isFaust.

See also
Faust.__init__, Faust.isFaust.,

◆ license()

def pyfaust.license ( )

Prints the FAuST license.

◆ norm()

def pyfaust.norm (   F,
  ord = 'fro',
**  kwargs 
)

Returns Faust.norm(F, ord) or numpy.linalg.norm(F, ord) depending of F type.

See also
Faust.norm

◆ opt_butterfly_faust()

def pyfaust.opt_butterfly_faust (   F)

Optimizes any Faust composed of butterfly factors.

The returned Faust will be more efficient if multiplied by a vector or a matrix. This optimization is based on the diagonals of each butterfly factor. Multiplying a butterfly factor B by a vector x (y = B@x) is equivalent to forming two diagonal matrices D1 and D2 from B and compute y' = D1@x + D2 @ x[I] where I is set in the proper order to obtain y' = y.

Parameters
FThe Faust to optimize. If the factors of F are not set according to a butterfly structure, the result is not defined.
Returns
The optimized Faust.
See also
pyfaust.fact.butterfly, pyfaust.dft, pyfaust.rand_butterfly

◆ pinv()

def pyfaust.pinv (   F)

A package function alias for the member function Faust.pinv().

◆ rand()

def pyfaust.rand (   num_rows,
  num_cols,
  num_factors = None,
  dim_sizes = None,
  density = None,
  fac_type = 'sparse',
  per_row = True,
  dev = 'cpu',
  dtype = 'float64',
  field = None,
  seed = 0 
)

Generates a random Faust.

Parameters
num_rowsthe Faust number of rows.
num_colsthe Faust number of columns.
num_factorsif it's an integer it is the number of random factors to set in the Faust. If num_factors is a list or tuple of 2 integers then the number of factors is set randomly between num_factors[0] and num_factors[1] (inclusively). Defaultly, num_factors is None, it means a 5 factors long Faust is generated.
dim_sizesif it's an integer all Faust factors are square matrices (except maybe the first and last ones, depending on num_rows and num_cols). The size of the intermediary square factors is size_dims**2. If it's a list or tuple of 2 integers then the number of rows and columns are both a random number between size_dims[0] and size_dims[1] (inclusively). Note that the first factor number of rows and the last factor number of columns are always fixed (they are the dimension sizes of the Faust: num_rows, num_cols arguments). if dim_sizes is None then dim_sizes is defaultly [num_rows, num_cols].
densitythe approximate density of factors. The default value is such that each factor gets 5 non-zero elements per row, if per_row is True, or per column otherwise. It should be a floating point number greater than 0 and lower or equal to 1. A density of zero is equivalent to the default case.
fac_typethe storage representation of factors. It must be 'sparse', 'dense' or 'mixed'. The latter designates a mix of dense and sparse matrices in the generated Faust (the choice is made according to a uniform distribution).
dtypethe dtype of the Faust ('float32', 'float64' or 'complex').
per_rowif True the factor matrix is constructed per row applying the density to each row. If False the construction is made with the density applied on each column.
devthe device on which to create the Faust ('cpu' or 'gpu').
field(DEPRECATED, use dtype) a str to set the Faust field: 'real' or 'complex'.
seedset PRNG initialization, useful for reproducibility of this function calls, otherwise seed argument shouldn't be used (the PRNG is automatically initialized).
Returns
the random Faust.

Examples

>>> from pyfaust import rand, seed
>>> seed(42)
>>> F = rand(2, 10, density=.5, dtype='complex')
>>> G = rand(10, 20, num_factors=[2, 5], density=.5, fac_type="dense")
>>> F

Faust size 2x10, density 2.6, nnz_sum 52, 5 factor(s):

  • FACTOR 0 (complex) SPARSE, size 2x8, density 0.5, nnz 8
  • FACTOR 1 (complex) SPARSE, size 8x4, density 0.5, nnz 16
  • FACTOR 2 (complex) SPARSE, size 4x5, density 0.4, nnz 8
  • FACTOR 3 (complex) SPARSE, size 5x3, density 0.333333, nnz 5
  • FACTOR 4 (complex) SPARSE, size 3x10, density 0.5, nnz 15
>>> G

Faust size 10x20, density 1.65, nnz_sum 330, 4 factor(s):

  • FACTOR 0 (double) DENSE, size 10x15, density 0.466667, nnz 70
  • FACTOR 1 (double) DENSE, size 15x12, density 0.5, nnz 90
  • FACTOR 2 (double) DENSE, size 12x11, density 0.454545, nnz 60
  • FACTOR 3 (double) DENSE, size 11x20, density 0.5, nnz 110
See also
Faust.__init__, pyfaust.rand_bsr

◆ rand_bsr()

def pyfaust.rand_bsr (   num_rows,
  num_cols,
  bnrows,
  bncols,
  num_factors = None,
  density = .1,
  dev = 'cpu',
  dtype = 'float64' 
)

Generates a random Faust composed only of BSR matrices.

Parameters
num_rows(int) the Faust number of rows.
num_cols(int) the Faust number of columns.
bnrows(int) the nonzero block number of rows (must divide num_rows).
bncols(int) the nonzero block number of columns (must divide num_cols).
num_factors(int or tuple(int,int) or NoneType) If it's an integer it will be the number of random factors to set in the Faust. If num_factors is a tuple of 2 integers then the number of factors will be set randomly between num_factors[0] and num_factors[1] (inclusively). If num_factors is None then 5 factors are generated.
density(float) the Faust factor density (it determines the number of nonzero blocks). It must be between 0 and 1.
dev(str) the device on which the Faust is created, 'cpu' (default) or 'gpu'.
dtype(str) the numpy dtype of the Faust.

Examples

>>> from pyfaust import rand_bsr
>>> rand_bsr(100, 100, 20, 10, num_factors=6)

Faust size 100x100, density 0.6, nnz_sum 6000, 6 factor(s):

  • FACTOR 0 (double) BSR, size 100x100 (blocksize = 20x10), density 0.1, nnz 1000 (nnz blocks: 5)
  • FACTOR 1 (double) BSR, size 100x100 (blocksize = 20x10), density 0.1, nnz 1000 (nnz blocks: 5)
  • FACTOR 2 (double) BSR, size 100x100 (blocksize = 20x10), density 0.1, nnz 1000 (nnz blocks: 5)
  • FACTOR 3 (double) BSR, size 100x100 (blocksize = 20x10), density 0.1, nnz 1000 (nnz blocks: 5)
  • FACTOR 4 (double) BSR, size 100x100 (blocksize = 20x10), density 0.1, nnz 1000 (nnz blocks: 5)
  • FACTOR 5 (double) BSR, size 100x100 (blocksize = 20x10), density 0.1, nnz 1000 (nnz blocks: 5)
See also
bsr_matrix, Faust.__init__, pyfaust.rand

◆ rand_butterfly()

def pyfaust.rand_butterfly (   n,
  dtype = 'float64',
  dev = 'cpu',
  diag_opt = False 
)

Constructs a Faust corresponding to the product of log2(n) square factors of size n with butterfly supports and random nonzero coefficients.

The random coefficients are drawn i.i.d. according to a standard Gaussian (real or complex circular according to the type).

Parameters
norder of the butterfly (must be a power of two).
diag_optif True then the returned Faust is optimized using
pyfaust.opt_butterfly_faust.
dtype'float32', 'float64' or 'complex', the dtype of the random Faust.
dev'cpu' or 'gpu', the device where the Faust is created.
Returns
F, a random butterfly support Faust.
See also
pyfaust.fact.butterfly, pyfaust.dft, pyfaust.opt_butterfly_faust

◆ seed()

def pyfaust.seed (   s)

(Re)Initializes the pyfaust pseudo-random generator.

It is useful to reproduce some code based for example on pyfaust.rand or pyfaust.rand_bsr.

Examples

>>> from pyfaust import rand, seed
>>> seed(42) # just for reproducibility
>>> F = rand(1024, 1024, dim_sizes=[1, 1024], num_factors=5, fac_type='mixed')
>>> # F will always be the same
>>> F

Faust size 1024x1024, density 0.017313, nnz_sum 18154, 5 factor(s):

  • FACTOR 0 (double) SPARSE, size 1024x754, density 0.0066313, nnz 5120
  • FACTOR 1 (double) SPARSE, size 754x386, density 0.0129534, nnz 3770
  • FACTOR 2 (double) DENSE, size 386x1000, density 0.004, nnz 1544
  • FACTOR 3 (double) SPARSE, size 1000x544, density 0.00919118, nnz 5000
  • FACTOR 4 (double) SPARSE, size 544x1024, density 0.00488281, nnz 2720

◆ toeplitz()

def pyfaust.toeplitz (   c,
  r = None,
  dev = 'cpu',
  diag_opt = False 
)

Constructs a toeplitz Faust whose first column is c and first row r.

Parameters
c(np.ndarray) the first column of the toeplitz Faust.
r(np.ndarray) the first row of the toeplitz Faust. If none then r = np.conjugate(c). r[0] is ignored, the first row is always [c[0], r[1:]].
dev(str) the device on which the Faust is created, 'cpu' (default) or 'gpu'.
diag_opt(bool) cf. pyfaust.circ.
Returns
The toeplitz Faust.

Examples

>>> from pyfaust import toeplitz
>>> import numpy as np
>>> c = [1, 2, 3, 4, 5, 6, 7, 8, 9, 10]
>>> T = toeplitz(c)
>>> T

Faust size 10x10, density 5.52, nnz_sum 552, 10 factor(s):

  • FACTOR 0 (complex) SPARSE, size 10x32, density 0.0625, nnz 20
  • FACTOR 1 (complex) SPARSE, size 32x32, density 0.0625, nnz 64
  • FACTOR 2 (complex) SPARSE, size 32x32, density 0.0625, nnz 64
  • FACTOR 3 (complex) SPARSE, size 32x32, density 0.0625, nnz 64
  • FACTOR 4 (complex) SPARSE, size 32x32, density 0.0625, nnz 64
  • FACTOR 5 (complex) SPARSE, size 32x32, density 0.0625, nnz 64
  • FACTOR 6 (complex) SPARSE, size 32x32, density 0.0625, nnz 64
  • FACTOR 7 (complex) SPARSE, size 32x32, density 0.0625, nnz 64
  • FACTOR 8 (complex) SPARSE, size 32x32, density 0.0625, nnz 64
  • FACTOR 9 (complex) SPARSE, size 32x10, density 0.0625, nnz 20
>>> np.allclose(T.toarray()[:, 0], c)
True
>>> np.allclose(T.toarray()[0, :], c)
True
>>> np.real(T.toarray())
array([[ 1., 2., 3., 4., 5., 6., 7., 8., 9., 10.],
[ 2., 1., 2., 3., 4., 5., 6., 7., 8., 9.],
[ 3., 2., 1., 2., 3., 4., 5., 6., 7., 8.],
[ 4., 3., 2., 1., 2., 3., 4., 5., 6., 7.],
[ 5., 4., 3., 2., 1., 2., 3., 4., 5., 6.],
[ 6., 5., 4., 3., 2., 1., 2., 3., 4., 5.],
[ 7., 6., 5., 4., 3., 2., 1., 2., 3., 4.],
[ 8., 7., 6., 5., 4., 3., 2., 1., 2., 3.],
[ 9., 8., 7., 6., 5., 4., 3., 2., 1., 2.],
[10., 9., 8., 7., 6., 5., 4., 3., 2., 1.]])
>>> r = [11, 12, 13, 14, 15, 16, 17, 18, 19, 20]
>>> T2 = toeplitz(c, r)
>>> T2

Faust size 10x10, density 5.52, nnz_sum 552, 10 factor(s):

  • FACTOR 0 (complex) SPARSE, size 10x32, density 0.0625, nnz 20
  • FACTOR 1 (complex) SPARSE, size 32x32, density 0.0625, nnz 64
  • FACTOR 2 (complex) SPARSE, size 32x32, density 0.0625, nnz 64
  • FACTOR 3 (complex) SPARSE, size 32x32, density 0.0625, nnz 64
  • FACTOR 4 (complex) SPARSE, size 32x32, density 0.0625, nnz 64
  • FACTOR 5 (complex) SPARSE, size 32x32, density 0.0625, nnz 64
  • FACTOR 6 (complex) SPARSE, size 32x32, density 0.0625, nnz 64
  • FACTOR 7 (complex) SPARSE, size 32x32, density 0.0625, nnz 64
  • FACTOR 8 (complex) SPARSE, size 32x32, density 0.0625, nnz 64
  • FACTOR 9 (complex) SPARSE, size 32x10, density 0.0625, nnz 20
>>> np.allclose(T2.toarray()[0, :], np.hstack((c[0:1], r[1:])))
True
>>> np.allclose(T2.toarray()[:, 0], c)
True
>>> np.real(T2.toarray())
array([[ 1., 12., 13., 14., 15., 16., 17., 18., 19., 20.],
[ 2., 1., 12., 13., 14., 15., 16., 17., 18., 19.],
[ 3., 2., 1., 12., 13., 14., 15., 16., 17., 18.],
[ 4., 3., 2., 1., 12., 13., 14., 15., 16., 17.],
[ 5., 4., 3., 2., 1., 12., 13., 14., 15., 16.],
[ 6., 5., 4., 3., 2., 1., 12., 13., 14., 15.],
[ 7., 6., 5., 4., 3., 2., 1., 12., 13., 14.],
[ 8., 7., 6., 5., 4., 3., 2., 1., 12., 13.],
[ 9., 8., 7., 6., 5., 4., 3., 2., 1., 12.],
[10., 9., 8., 7., 6., 5., 4., 3., 2., 1.]])
>>> # Look at the density of a larger Toeplitz Faust
>>> # it indicates a speedup of the Faust-matrix/vector product
>>> toeplitz(np.random.rand(1024), np.random.rand(1024)).density()
0.08203125
See also
scipy.linalg.toeplitz, pyfaust.circ, pyfaust.anticirc

◆ version()

def pyfaust.version ( )

Returns the FAuST package version.

◆ vstack()

def pyfaust.vstack (   _tuple)

Concatenates vertically Faust-s and/or numpy.ndarray arrays using Faust.concatenate().

See also numpy.vstack

Note
it could be wiser to encapsulate a Faust in a lazylinop.LazyLinearOp for a lazy concatenation.

Examples

>>> from pyfaust import *
>>> seed(42) # just for reproducibility
>>> F1 = rand(5,50)
>>> F2 = rand(5,50)
>>> vstack((F1, F2))

Faust size 10x50, density 2.17, nnz_sum 1085, 6 factor(s):

  • FACTOR 0 (double) SPARSE, size 10x67, density 0.0746269, nnz 50
  • FACTOR 1 (double) SPARSE, size 67x33, density 0.151515, nnz 335
  • FACTOR 2 (double) SPARSE, size 33x31, density 0.16129, nnz 165
  • FACTOR 3 (double) SPARSE, size 31x56, density 0.0892857, nnz 155
  • FACTOR 4 (double) SPARSE, size 56x100, density 0.05, nnz 280
  • FACTOR 5 (double) SPARSE, size 100x50, density 0.02, nnz 100

◆ wht()

def pyfaust.wht (   n,
  normed = True,
  dev = "cpu",
  dtype = 'float64' 
)

Constructs a Faust implementing the Walsh-Hadamard Transform (WHT) of order n.

The resulting Faust has log2(n) sparse factors of order n, each one having 2 nonzeros per row and per column.

Parameters
n(int) order of the WHT (must be a power of two).
normed(bool) default to True to normalize the Hadamard Faust as if you called Faust.normalize() and False otherwise.
dev(str) device on which to create the Faust ('cpu' or 'gpu').
dtype(str) the Faust dtype, it must be 'float32', 'float64' or 'complex'.
Returns
The Faust implementing the Hadamard transform of dimension n.

Examples

>>> from pyfaust import wht
>>> wht(1024)

Faust size 1024x1024, density 0.0195312, nnz_sum 20480, 10 factor(s):

  • FACTOR 0 (double) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 1 (double) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 2 (double) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 3 (double) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 4 (double) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 5 (double) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 6 (double) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 7 (double) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 8 (double) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 9 (double) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
    >>> wht(1024, normed=True) # is equiv. to next call
    Faust size 1024x1024, density 0.0195312, nnz_sum 20480, 10 factor(s):
  • FACTOR 0 (double) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 1 (double) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 2 (double) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 3 (double) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 4 (double) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 5 (double) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 6 (double) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 7 (double) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 8 (double) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 9 (double) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
>>> wht(1024, normed=False).normalize() # which is less optimized though

Faust size 1024x1024, density 0.0195312, nnz_sum 20480, 10 factor(s):

  • FACTOR 0 (double) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 1 (double) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 2 (double) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 3 (double) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 4 (double) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 5 (double) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 6 (double) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 7 (double) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 8 (double) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 9 (double) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
See also
scipy.linalg.hadamard, pyfaust.dft, pyfaust.fact.butterfly,

Variable Documentation

◆ WARNING_ON_C_CONTIGUOUS_MATMUL

bool pyfaust.WARNING_ON_C_CONTIGUOUS_MATMUL = True
pyfaust.dct
def dct(n, normed=True, dev='cpu', dtype='float64')
Constructs a Faust implementing the Direct Cosine Transform (Type II) Faust of order n.
Definition: __init__.py:4317
pyfaust.seed
def seed(s)
(Re)Initializes the pyfaust pseudo-random generator.
Definition: __init__.py:5277
pyfaust.eye
def eye(m, n=None, dtype='float64', dev="cpu")
Faust identity.
Definition: __init__.py:4743
pyfaust.dst
def dst(n, normed=True, dev='cpu', dtype='float64')
Constructs a Faust implementing the Direct Sine Transform (Type II) Faust of order n.
Definition: __init__.py:4388
pyfaust.dft
def dft(n, normed=True, dev='cpu', diag_opt=False)
Constructs a Faust F implementing the Discrete Fourier Transform (DFT) of order n.
Definition: __init__.py:4257
pyfaust.rand_bsr
def rand_bsr(num_rows, num_cols, bnrows, bncols, num_factors=None, density=.1, dev='cpu', dtype='float64')
Generates a random Faust composed only of BSR matrices.
Definition: __init__.py:4817
pyfaust.rand
def rand(num_rows, num_cols, num_factors=None, dim_sizes=None, density=None, fac_type='sparse', per_row=True, dev='cpu', dtype='float64', field=None, seed=0)
Generates a random Faust.
Definition: __init__.py:4966
pyfaust.anticirc
def anticirc(c, dev='cpu', diag_opt=False)
Returns an anti-circulant Faust A defined by the vector c (which is the last column of A....
Definition: __init__.py:4584
pyfaust.toeplitz
def toeplitz(c, r=None, dev='cpu', diag_opt=False)
Constructs a toeplitz Faust whose first column is c and first row r.
Definition: __init__.py:4694
pyfaust.concatenate
def concatenate(_tuple, *args, axis=0, **kwargs)
A package function alias for the member function Faust.concatenate.
Definition: __init__.py:3950
pyfaust.circ
def circ(c, dev='cpu', diag_opt=False)
Returns a circulant Faust C defined by the vector c (which is the first column of C....
Definition: __init__.py:4498
pyfaust.wht
def wht(n, normed=True, dev="cpu", dtype='float64')
Constructs a Faust implementing the Walsh-Hadamard Transform (WHT) of order n.
Definition: __init__.py:4110
pyfaust.vstack
def vstack(_tuple)
Concatenates vertically Faust-s and/or numpy.ndarray arrays using Faust.concatenate().
Definition: __init__.py:4029
pyfaust.hstack
def hstack(_tuple)
Concatenates horizontally Faust-s and/or numpy.ndarray objects using Faust.concatenate().
Definition: __init__.py:4000