pyfaust core module

The pyfaust core module which provides the Faust class and many functions to generate Fausts.

class pyfaust.Faust(factors=None, filepath=None, **kwargs)

pyfaust main class for using multi-layer sparse transforms.

This class provides a NumPy-like interface for operations with FAuST data structures, which correspond ideally to matrices that can be written exactly as the product of sparse matrices.

The Faust class is designed to allow fast matrix-vector multiplications together with reduced memory storage compared to what would be obtained by manipulating directly the corresponding dense matrix.

A particular example is the matrix associated to the discrete Fourier transform, which can be represented exactly as a Faust, leading to a fast and compact implementation (see pyfaust.dft())

Although sparse matrices are more interesting for optimization it’s not forbidden to define a Faust as a product of dense matrices or a mix of dense and sparse matrices.

The matrices composing the Faust product, also called the factors, are defined on complex or real fields. Hence a Faust can be a complex Faust or a real Faust.

Several Python built-ins have been overloaded to ensure that a Faust is almost handled as a native NumPy array.

The main exception is that contrary to a NumPy array a Faust is immutable. It means that you cannot modify elements of a Faust using the assignment operator = as made with a NumPy array (e.g. M[i, j] = 2). That limitation is the reason why the Python built-in __setitem__() is not implemented in this class. Note however that you can optionally contravene the immuability in certain functions (e.g. with the inplace argument of the functions Faust.swap_rows(), Faust.swap_cols(), Faust.optimize_time()).

Other noticeable limitations:

  • Elementwise multiplication, addition/subtraction are available, but performing elementwise operations between two Fausts is discouraged,

  • One cannot reshape a Faust.

A last but not least caveat is that Faust doesn’t support the NumPy universal functions (ufuncs) except if the contrary is specified in the API doc. for a particular function.

Mainly for convenience and test purposes, a Faust can be converted into the corresponding full matrix using the function Faust.toarray().

Warning

using Faust.toarray() is discouraged except for test purposes, as it loses the main potential interests of the FAuST structure: compressed memory storage and faster matrix-vector multiplication compared to its equivalent full matrix representation.

In this documentation, the expression ‘full matrix’ designates the array Faust.toarray() obtained by the multiplication of the Faust factors.

Note

it could be wiser to encapsulate a Faust in a lazylinop.LazyLinOp for a totally lazy paradigm on all available operations.

List of functions that are memory costly:

For more information about FAuST take a look at https://faust.inria.fr.

property H

Returns the conjugate transpose of F.

This function is intended to be used as a property (see the examples).

Returns:

a Faust object H implementing the conjugate transpose of F.toarray() such that: H.toarray() == F.toarray().H

Examples:
>>> from pyfaust import rand
>>> F = rand(50, 50, field='complex')
>>> H1 = F.H
>>> H2 = F.transpose()
>>> H2 = H2.conj()
>>> (H1.toarray() == H2.toarray()).all()
True
property T

Returns the transpose of F.

This function is intended to be used as a property (see the examples).

Args:

F: the Faust object.

Returns:

a Faust object implementing the transpose of F.toarray(), such that: <code>F.T.toarray() == F.toarray().T</code>

Examples:
>>> from pyfaust import rand
>>> F = rand(10, 23)
>>> tF = F.T
>>> tF.shape
(23, 10)
__add__(*args, **kwargs)

Sums F to one or a sequence of variables (Faust objects, arrays or scalars).

Note

This method overloads the Python function/operator +.

Args:
args: (list[Faust or np.ndarray or scipy.sparse.csr_matrix])

the list of variables to sum all together with F. These can be Faust objects, numpy arrays (and scipy.sparse.csr_matrix) or scalars. Faust and arrays/matrices must be of compatible sizes.

Returns:

the sum as a Faust object.

Example:
>>> import pyfaust as pf
>>> from numpy.linalg import norm
>>> pf.seed(42) # just for reproducibility
>>> F = pf.rand(10, 12)
>>> F
Faust size 10x12, density 2.04167, nnz_sum 245, 5 factor(s):
- FACTOR 0 (double) SPARSE, size 10x12, density 0.333333, nnz 40
- FACTOR 1 (double) SPARSE, size 12x10, density 0.5, nnz 60
- FACTOR 2 (double) SPARSE, size 10x11, density 0.454545, nnz 50
- FACTOR 3 (double) SPARSE, size 11x10, density 0.5, nnz 55
- FACTOR 4 (double) SPARSE, size 10x12, density 0.333333, nnz 40
>>> G = pf.rand(10, 12)
>>> G
Faust size 10x12, density 2.025, nnz_sum 243, 5 factor(s):
- FACTOR 0 (double) SPARSE, size 10x11, density 0.454545, nnz 50
- FACTOR 1 (double) SPARSE, size 11x10, density 0.5, nnz 55
- FACTOR 2 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 3 (double) SPARSE, size 10x12, density 0.333333, nnz 40
- FACTOR 4 (double) SPARSE, size 12x12, density 0.333333, nnz 48
>>> F+G
Faust size 10x12, density 4.43333, nnz_sum 532, 7 factor(s):
- FACTOR 0 (double) SPARSE, size 10x20, density 0.1, nnz 20
- FACTOR 1 (double) SPARSE, size 20x23, density 0.195652, nnz 90
- FACTOR 2 (double) SPARSE, size 23x20, density 0.25, nnz 115
- FACTOR 3 (double) SPARSE, size 20x21, density 0.238095, nnz 100
- FACTOR 4 (double) SPARSE, size 21x22, density 0.205628, nnz 95
- FACTOR 5 (double) SPARSE, size 22x24, density 0.166667, nnz 88
- FACTOR 6 (double) SPARSE, size 24x12, density 0.0833333, nnz 24
>>> norm((F+G).toarray() - F.toarray() - G.toarray())
5.975125014346798e-15
>>> F+2
Faust size 10x12, density 2.84167, nnz_sum 341, 7 factor(s):
- FACTOR 0 (double) SPARSE, size 10x20, density 0.1, nnz 20
- FACTOR 1 (double) SPARSE, size 20x22, density 0.113636, nnz 50
- FACTOR 2 (double) SPARSE, size 22x20, density 0.159091, nnz 70
- FACTOR 3 (double) SPARSE, size 20x21, density 0.142857, nnz 60
- FACTOR 4 (double) SPARSE, size 21x11, density 0.281385, nnz 65
- FACTOR 5 (double) SPARSE, size 11x24, density 0.19697, nnz 52
- FACTOR 6 (double) SPARSE, size 24x12, density 0.0833333, nnz 24
>>> norm((F+2).toarray() - F.toarray() - 2)
1.85775845048325e-15
>>> F+G+2
Faust size 10x12, density 5.4, nnz_sum 648, 9 factor(s):
- FACTOR 0 (double) SPARSE, size 10x20, density 0.1, nnz 20
- FACTOR 1 (double) SPARSE, size 20x30, density 0.05, nnz 30
- FACTOR 2 (double) SPARSE, size 30x33, density 0.10101, nnz 100
- FACTOR 3 (double) SPARSE, size 33x30, density 0.126263, nnz 125
- FACTOR 4 (double) SPARSE, size 30x31, density 0.11828, nnz 110
- FACTOR 5 (double) SPARSE, size 31x32, density 0.105847, nnz 105
- FACTOR 6 (double) SPARSE, size 32x25, density 0.1225, nnz 98
- FACTOR 7 (double) SPARSE, size 25x24, density 0.06, nnz 36
- FACTOR 8 (double) SPARSE, size 24x12, density 0.0833333, nnz 24
>>> norm((F+G+2).toarray() - F.toarray() - G.toarray() - 2)
7.115828086306871e-15
>>> F+G+F+2+F
Faust size 10x12, density 11.05, nnz_sum 1326, 13 factor(s):
- FACTOR 0 (double) SPARSE, size 10x20, density 0.1, nnz 20
- FACTOR 1 (double) SPARSE, size 20x30, density 0.05, nnz 30
- FACTOR 2 (double) SPARSE, size 30x40, density 0.0333333, nnz 40
- FACTOR 3 (double) SPARSE, size 40x50, density 0.025, nnz 50
- FACTOR 4 (double) SPARSE, size 50x53, density 0.045283, nnz 120
- FACTOR 5 (double) SPARSE, size 53x52, density 0.0634978, nnz 175
- FACTOR 6 (double) SPARSE, size 52x51, density 0.0678733, nnz 180
- FACTOR 7 (double) SPARSE, size 51x55, density 0.0695187, nnz 195
- FACTOR 8 (double) SPARSE, size 55x54, density 0.0717172, nnz 213
- FACTOR 9 (double) SPARSE, size 54x36, density 0.063786, nnz 124
- FACTOR 10 (double) SPARSE, size 36x34, density 0.0743464, nnz 91
- FACTOR 11 (double) SPARSE, size 34x24, density 0.0784314, nnz 64
- FACTOR 12 (double) SPARSE, size 24x12, density 0.0833333, nnz 24
>>> norm((F+G+F+2+F).toarray() - 3*F.toarray() - G.toarray() - 2)
2.7030225852818652e-14

See also

Faust.__sub__()

__div__(s)

Divides F by the scalar s. This method overloads the Python function/operator `/’ (whether s is a float or an integer).

Args:
F: (Faust)

the Faust object.

s: (scalar)

the scalar to divide the Faust object with.

Returns:

the division result as a Faust object.

__floordiv__(s)

F // s.

Warning

this operation is not supported, it raises an exception.

__getitem__(indices)

Indexes or slices a Faust. Returns a Faust representing a submatrix of F.toarray() or a scalar element if that Faust can be reduced to a single element. This function overloads a Python built-in.

Warning

  • This function doesn’t implement F[l1, l2] where l1 and l2 are integer lists, rather use F[l1][:, l2].

  • It is not advised to use this function as an element accessor (e.g. F[0, 0]) because such a use induces to convert the Faust to its dense matrix representation and that is a very expensive computation if used repetitively.

  • Subindexing a Faust which would create an empty Faust will raise an error.

  • ‘Fancy indexing’ must be done with a list not a numpy array.

Args:
F: (Faust)

the Faust object.

indices: (list)

array of length 1 or 2 which elements must be slice, integer or Ellipsis (…) (see examples below). Note that using Ellipsis for more than two indices is forbidden.

Returns:

the Faust object requested or just the corresponding scalar if that Faust has a shape equal to (1, 1).

Raises:

IndexError

Examples:
>>> from pyfaust import rand, seed
>>> import numpy as np
>>> from numpy.random import randint, seed as rseed
>>> seed(42) # just for reproducibility
>>> rseed(42)
>>> F = rand(50, 100)
>>> i1 = randint(1, min(F.shape)-1)
>>> i2 = randint(1, min(F.shape)-1)
>>> # the scalar element located at
>>> # at row i1, column i2 of the F dense matrix
>>> F[i1, i2]
0.0
>>> F[:, i2] # full column i2
Faust size 50x1, density 24.64, nnz_sum 1232, 5 factor(s):
- FACTOR 0 (double) SPARSE, size 50x87, density 0.045977, nnz 200
- FACTOR 1 (double) SPARSE, size 87x63, density 0.0793651, nnz 435
- FACTOR 2 (double) SPARSE, size 63x69, density 0.057971, nnz 252
- FACTOR 3 (double) SPARSE, size 69x60, density 0.0833333, nnz 345
- FACTOR 4 (double) SPARSE, size 60x1, density 0, nnz 0
>>> # from row 2 to 3, each row containing
>>> # only elements from column 1 to 3
>>> F[2:4, 1:4]
Faust size 2x3, density 174.833, nnz_sum 1049, 5 factor(s):
- FACTOR 0 (double) SPARSE, size 2x87, density 0.045977, nnz 8
- FACTOR 1 (double) SPARSE, size 87x63, density 0.0793651, nnz 435
- FACTOR 2 (double) SPARSE, size 63x69, density 0.057971, nnz 252
- FACTOR 3 (double) SPARSE, size 69x60, density 0.0833333, nnz 345
- FACTOR 4 (double) SPARSE, size 60x3, density 0.05, nnz 9
>>> # from row 0 to end row, each row
>>> # containing only elements from column 4 to
>>> # column before the last one.
>>> F[::, 4:-1]
Faust size 50x95, density 0.318947, nnz_sum 1515, 5 factor(s):
- FACTOR 0 (double) SPARSE, size 50x87, density 0.045977, nnz 200
- FACTOR 1 (double) SPARSE, size 87x63, density 0.0793651, nnz 435
- FACTOR 2 (double) SPARSE, size 63x69, density 0.057971, nnz 252
- FACTOR 3 (double) SPARSE, size 69x60, density 0.0833333, nnz 345
- FACTOR 4 (double) SPARSE, size 60x95, density 0.0496491, nnz 283
>>> F[0:i1, ...] # equivalent to F[0:i1, ::]
Faust size 39x100, density 0.381538, nnz_sum 1488, 5 factor(s):
- FACTOR 0 (double) SPARSE, size 39x87, density 0.045977, nnz 156
- FACTOR 1 (double) SPARSE, size 87x63, density 0.0793651, nnz 435
- FACTOR 2 (double) SPARSE, size 63x69, density 0.057971, nnz 252
- FACTOR 3 (double) SPARSE, size 69x60, density 0.0833333, nnz 345
- FACTOR 4 (double) SPARSE, size 60x100, density 0.05, nnz 300
>>> F[2::, :3:] # equivalent to F[2:F.shape[0], 0:3]
Faust size 48x3, density 8.56944, nnz_sum 1234, 5 factor(s):
- FACTOR 0 (double) SPARSE, size 48x87, density 0.045977, nnz 192
- FACTOR 1 (double) SPARSE, size 87x63, density 0.0793651, nnz 435
- FACTOR 2 (double) SPARSE, size 63x69, density 0.057971, nnz 252
- FACTOR 3 (double) SPARSE, size 69x60, density 0.0833333, nnz 345
- FACTOR 4 (double) SPARSE, size 60x3, density 0.0555556, nnz 10
>>> F[0:i2:2, :] # takes every row of even index until the row i2 (excluded)
Faust size 15x100, density 0.928, nnz_sum 1392, 5 factor(s):
- FACTOR 0 (double) SPARSE, size 15x87, density 0.045977, nnz 60
- FACTOR 1 (double) SPARSE, size 87x63, density 0.0793651, nnz 435
- FACTOR 2 (double) SPARSE, size 63x69, density 0.057971, nnz 252
- FACTOR 3 (double) SPARSE, size 69x60, density 0.0833333, nnz 345
- FACTOR 4 (double) SPARSE, size 60x100, density 0.05, nnz 300
>>> F[-1:-3:-1, :] # takes the two last rows in reverse order
Faust size 2x100, density 6.7, nnz_sum 1340, 5 factor(s):
- FACTOR 0 (double) SPARSE, size 2x87, density 0.045977, nnz 8
- FACTOR 1 (double) SPARSE, size 87x63, density 0.0793651, nnz 435
- FACTOR 2 (double) SPARSE, size 63x69, density 0.057971, nnz 252
- FACTOR 3 (double) SPARSE, size 69x60, density 0.0833333, nnz 345
- FACTOR 4 (double) SPARSE, size 60x100, density 0.05, nnz 300
>>> F[i2:0:-2, :] # starts from row i2 and goes backward to take one in two rows until the first one (reversing order of F)
Faust size 15x100, density 0.928, nnz_sum 1392, 5 factor(s):
- FACTOR 0 (double) SPARSE, size 15x87, density 0.045977, nnz 60
- FACTOR 1 (double) SPARSE, size 87x63, density 0.0793651, nnz 435
- FACTOR 2 (double) SPARSE, size 63x69, density 0.057971, nnz 252
- FACTOR 3 (double) SPARSE, size 69x60, density 0.0833333, nnz 345
- FACTOR 4 (double) SPARSE, size 60x100, density 0.05, nnz 300
>>> F[[1, 18, 2], :] # takes in this order the rows 1, 18 and 2
Faust size 3x100, density 4.48, nnz_sum 1344, 5 factor(s):
- FACTOR 0 (double) SPARSE, size 3x87, density 0.045977, nnz 12
- FACTOR 1 (double) SPARSE, size 87x63, density 0.0793651, nnz 435
- FACTOR 2 (double) SPARSE, size 63x69, density 0.057971, nnz 252
- FACTOR 3 (double) SPARSE, size 69x60, density 0.0833333, nnz 345
- FACTOR 4 (double) SPARSE, size 60x100, density 0.05, nnz 300
>>> F[:, [1, 18, 2]] # takes in this order the columns 1, 18 and 2
Faust size 50x3, density 8.27333, nnz_sum 1241, 5 factor(s):
- FACTOR 0 (double) SPARSE, size 50x87, density 0.045977, nnz 200
- FACTOR 1 (double) SPARSE, size 87x63, density 0.0793651, nnz 435
- FACTOR 2 (double) SPARSE, size 63x69, density 0.057971, nnz 252
- FACTOR 3 (double) SPARSE, size 69x60, density 0.0833333, nnz 345
- FACTOR 4 (double) SPARSE, size 60x3, density 0.05, nnz 9
>>> F[[1, 18, 2]][:, [1, 2]] # takes the rows 1, 18 and 2 but keeps only columns 1 and 2 in these rows
Faust size 3x2, density 175, nnz_sum 1050, 5 factor(s):
- FACTOR 0 (double) SPARSE, size 3x87, density 0.045977, nnz 12
- FACTOR 1 (double) SPARSE, size 87x63, density 0.0793651, nnz 435
- FACTOR 2 (double) SPARSE, size 63x69, density 0.057971, nnz 252
- FACTOR 3 (double) SPARSE, size 69x60, density 0.0833333, nnz 345
- FACTOR 4 (double) SPARSE, size 60x2, density 0.05, nnz 6
__iadd__(A)

Inplace operation: F = F + A

__imatmul__(A)

Inplace operation: F = F @ A

__imul__(A)

Inplace operation: F = F * A

__init__(factors=None, filepath=None, **kwargs)

Creates a Faust from a list of factors or alternatively from a file. Other easy ways to create a Faust is to call one of the following functions: pyfaust.rand, pyfaust.dft or pyfaust.wht.

Args:
factors: (list of numpy/scipy array/matrices or a single array/matrix.)

The factors must respect the dimensions needed for the product to be defined <code>(for i in range(0, len(factors)-1): factors[i].shape[1] == factors[i+1].shape[0])</code>. The factors can be sparse or dense matrices (either scipy.sparse.csr_matrix/bsr_matrix or numpy.ndarray with ndim == 2). scipy.sparse.csc_matrix/coo_matrix are supported but converted to csr_matrix on the fly. The Faust will be in the same dtype as the factors (only ‘float32’, ‘float64’/’double’ and ‘complex128’/’complex’ dtype are supported). Passing only an array or a sparse matrix to the constructor is equivalent to passing a list of a single factor.

filepath: (str)

the file from which a Faust is created. The format is Matlab version 5 (.mat extension). The file must have been saved before with Faust.save().()

dev: (str)

‘cpu’ (by default) or ‘gpu’ to instantiate the Faust resp. on CPU or on NVIDIA GPU memory.

kwargs: (dict)

internal purpose arguments.

Note

filepath and factors arguments are mutually exclusive. Either you specify one of them explicitly with the keyword or the first (positional) argument type will determine if it’s a filepath (a str) or a factor list. If both of the keyword arguments are set then filepath will be prioritary.

Examples:
>>> # Example 1 -- Creating a Faust made of CSR matrices and numpy arrays:
>>> from pyfaust import Faust
>>> import numpy as np
>>> from scipy import sparse
>>> factors = []
>>> factors += [sparse.random(100, 100, dtype=np.float64, format='csr', density=0.1)]
>>> factors += [np.random.rand(100, 100).astype(np.float64) ]
>>> factors += [sparse.random(100, 100, dtype=np.float64, format='csr', density=0.1)]
>>> factors += [np.random.rand(100, 100).astype(np.float64) ]
>>> factors += [sparse.random(100, 100, dtype=np.float64, format='csr', density=0.1)]
>>> # define a Faust with those factors
>>> F = Faust(factors)
>>> F.save("F.mat")
>>> # define a Faust from file
>>> H = Faust(filepath="F.mat") # F = Faust("F.mat") does the same
>>> # creating a Faust with only one
>>> # factor
>>> Faust(np.random.rand(10, 10))
Faust size 10x10, density 1, nnz_sum 100, 1 factor(s):
- FACTOR 0 (double) DENSE, size 10x10, density 1, nnz 100
>>> # Example 2 -- Creating a Faust containing one BSR matrix:
>>> from scipy.sparse import bsr_matrix
>>> from pyfaust import Faust
>>> from numpy import allclose
>>> from numpy.random import rand
>>> nzblocks = rand(3, 2, 3) # 3 blocks of size 2x3
>>> # create a scipy BSR matrix
>>> B = bsr_matrix((nzblocks, [0, 1, 2], [0, 1, 2, 3, 3, 3]), shape=(10, 9))
>>> # create the single factor Faust
>>> F = Faust(B)
>>> F
Faust size 10x9, density 0.2, nnz_sum 18, 1 factor(s):
- FACTOR 0 (double) BSR, size 10x9 (blocksize = 2x3), density 0.2, nnz 18 (nnz blocks: 3)
>>> allclose(F.toarray(), B.toarray())
True
>>> # of course it's ok to create a Faust with a BSR and another type of factors
>>> Faust([B, rand(9, 18)])
Faust size 10x18, density 1, nnz_sum 180, 2 factor(s):
- FACTOR 0 (double) BSR, size 10x9 (blocksize = 2x3), density 0.2, nnz 18 (nnz blocks: 3)
- FACTOR 1 (double) DENSE, size 9x18, density 1, nnz 162
__isub__(A)

Inplace operation: F = F - A

__itruediv__(s)

Divides F by the scalar s inplace. This method overloads the Python function/operator `/=’ (whether s is a float or an integer).

Args:
F: (Faust)

the Faust object.

s: (scalar)

the scalar to divide the Faust object with.

Returns:

the division result as a Faust object.

__len__()

Returns the number of factors of F.

Returns:

the number of factors.

Examples:
>>> from pyfaust import rand
>>> F = rand(50, 100)
>>> nf = F.numfactors()
>>> nf
5
>>> nf == len(F)
True
__matmul__(A)

Multiplies F by A which is a dense numpy.matrix/numpy.ndarray or a Faust object. This method overloads a Python function/operator (@).

The primary goal of this function is to implement “fast” multiplication by a Faust, which is the operation performed when A is a dense matrix. In the best case, F @ A is F.rcg() times faster than equivalent F.toarray() @ A.

Other use case is available for this function: - If A is a Faust, no actual multiplication is performed, instead a new Faust is built to implement the multiplication. This Faust verifies that: (F @ A).toarray() == F.toarray() @ A.toarray() (an elementwise non-significant absolute difference between the two members is possible).

Args:
F: (Faust)

the Faust object.

A: a Faust object, a sparse matrix (scipy.sparse.csr_matrix or dia_matrix) or a 2D full matrix (numpy.ndarray, numpy.matrix).

In the latter case, A must be Fortran contiguous (i.e. Column-major order; `order’ argument passed to np.ndararray() must be equal to str ‘F’). .. warning:: if A is a numpy array then A.flags[‘F_CONTIGUOUS’]

must be True (that is, in column-major order) or a on-the-fly conversion will take place (with a little overhead).

Note that performing a Faust-csr_matrix product is often slower than converting first the csr_matrix to a dense representation (toarray()) and then proceed to the Faust-dense_matrix multiplication. In some cases though, it stays quicker: moreover when the Faust is composed of a small number of factors.

Returns:

The result of the multiplication

  • as a numpy.ndarray if A is a ndarray,

  • as a Faust object if A is a Faust.

  • When either F or A is complex, G=F @ A is also complex.

Raises:

ValueError

Examples:
>>> from pyfaust import rand, seed
>>> import numpy as np
>>> F = rand(50, 100)
>>> A = np.random.rand(F.shape[1], 50)
>>> B = F@A # == F*A or F.dot(A)
>>> # is equivalent to B = F.__matmul__(A)
>>> G = rand(F.shape[1], 5)
>>> H = F@G
>>> # H is a Faust because F and G are
The multiplying operand A can’t be a scalar:
>>> from pyfaust import rand
>>> F = rand(5, 50)
>>> F@2
Traceback (most recent call last):
...
ValueError: Scalar operands are not allowed, use '*' instead
__mul__(A)

Multiplies the Faust F by A. This method overloads a Python function/operator (*).

More specifically: - It’s a scalar multiplication if A is a scalar number. - if A is a vector of appropriate size, the function performs a Faust-vector broadcasting.

Args:

F: the Faust object. A: a scalar number or a (vector) numpy.ndarray or a Faust

elementwise multiplication.

Note

to compute the elementwise multiplication F * A column-by-column (hence avoiding calling toarray(), which consumes more memory) enable the environment variable PYFAUST_ELT_WISE_MUL_BY_COL.

Returns: The result of the multiplication as a Faust (either A is a vector or a scalar).

Raises: TypeError

Examples:
>>> from pyfaust import rand
>>> import numpy as np
>>> F = rand(50, 100)
>>> v = np.random.rand(F.shape[1])
>>> B = F*v # Faust vector broadcasting
>>> # is equivalent to B = F.__mul__(v)
>>> F_times_two = F*2

Compute elementwise F * F, gives an array: >>> FFarray = F * F

If the type of A is not supported: >>> import numpy >>> F*’a’ Traceback (most recent call last): … TypeError: must be real number, not str

See also

Faust.__rmul__()

__neg__()

Returns - F (unary operator).

Note

This method overloads the Python unary operator -.

__pos__()

Returns + F (unary operator).

Note

This method overloads the Python unary operator +.

See also

Faust.__add__()

__radd__(lhs_op)

Returns lhs_op+F.

See also

Faust.__add__(),

__repr__()

Returns a str object representing the Faust object. This method overloads a Python function.

Note

Ideally this function is intended to return a valid Python expression but here this is not the case. Only information is displayed.

Args:

F: the Faust object.

Examples:
>>> from pyfaust import rand
>>> F = rand(50, 50)
>>> print(F.__repr__())
Faust size 50x50, density 0.5, nnz_sum 1250, 5 factor(s):
- FACTOR 0 (double) SPARSE, size 50x50, density 0.1, nnz 250
- FACTOR 1 (double) SPARSE, size 50x50, density 0.1, nnz 250
- FACTOR 2 (double) SPARSE, size 50x50, density 0.1, nnz 250
- FACTOR 3 (double) SPARSE, size 50x50, density 0.1, nnz 250
- FACTOR 4 (double) SPARSE, size 50x50, density 0.1, nnz 250
>>> # the same function is called when typing F in a terminal:
>>> F
Faust size 50x50, density 0.5, nnz_sum 1250, 5 factor(s):
- FACTOR 0 (double) SPARSE, size 50x50, density 0.1, nnz 250
- FACTOR 1 (double) SPARSE, size 50x50, density 0.1, nnz 250
- FACTOR 2 (double) SPARSE, size 50x50, density 0.1, nnz 250
- FACTOR 3 (double) SPARSE, size 50x50, density 0.1, nnz 250
- FACTOR 4 (double) SPARSE, size 50x50, density 0.1, nnz 250
__rmatmul__(lhs_op)

Returns lhs_op.__matmul__(F).

Examples:
>>> from pyfaust import rand
>>> import numpy as np
>>> F = rand(50, 100)
>>> A = np.random.rand(50, F.shape[0])
>>> B = A@F # == A*F or pyfaust.dot(A, F)
__rmul__(lhs_op)

lhs_op*F

See also

Faust.__mul__()

__rsub__(lhs_op)

Returns lhs_op-F.

See also

Faust.__sub__()

__str__()

Converts F to its str representation.

Returns:

The str conversion of F.

__sub__(*args)

Subtracts from F one or a sequence of variables. Faust objects, arrays or scalars.

Note

This method overloads the Python function/operator -.

Args:
args: (list[Faust or np.ndarray or scipy.sparse.csr_matrix])

the list of variables to compute the difference with F. These can be Faust objects, arrays (and scipy.sparse.csr_matrix) or scalars. Faust and arrays/matrices must be of compatible sizes.

Returns:

the difference as a Faust object.

Example:
>>> import pyfaust as pf
>>> from numpy.linalg import norm
>>> pf.seed(42) # just for reproducibility
>>> F = pf.rand(10, 12)
>>> F
Faust size 10x12, density 2.04167, nnz_sum 245, 5 factor(s):
- FACTOR 0 (double) SPARSE, size 10x12, density 0.333333, nnz 40
- FACTOR 1 (double) SPARSE, size 12x10, density 0.5, nnz 60
- FACTOR 2 (double) SPARSE, size 10x11, density 0.454545, nnz 50
- FACTOR 3 (double) SPARSE, size 11x10, density 0.5, nnz 55
- FACTOR 4 (double) SPARSE, size 10x12, density 0.333333, nnz 40
>>> G = pf.rand(10, 12)
>>> G
Faust size 10x12, density 2.025, nnz_sum 243, 5 factor(s):
- FACTOR 0 (double) SPARSE, size 10x11, density 0.454545, nnz 50
- FACTOR 1 (double) SPARSE, size 11x10, density 0.5, nnz 55
- FACTOR 2 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 3 (double) SPARSE, size 10x12, density 0.333333, nnz 40
- FACTOR 4 (double) SPARSE, size 12x12, density 0.333333, nnz 48
>>> F-G
Faust size 10x12, density 4.43333, nnz_sum 532, 7 factor(s):
- FACTOR 0 (double) SPARSE, size 10x20, density 0.1, nnz 20
- FACTOR 1 (double) SPARSE, size 20x23, density 0.195652, nnz 90
- FACTOR 2 (double) SPARSE, size 23x20, density 0.25, nnz 115
- FACTOR 3 (double) SPARSE, size 20x21, density 0.238095, nnz 100
- FACTOR 4 (double) SPARSE, size 21x22, density 0.205628, nnz 95
- FACTOR 5 (double) SPARSE, size 22x24, density 0.166667, nnz 88
- FACTOR 6 (double) SPARSE, size 24x12, density 0.0833333, nnz 24
>>> norm((F-G).toarray() - F.toarray() + G.toarray())
3.608651813623577e-15
>>> F-2
Faust size 10x12, density 2.84167, nnz_sum 341, 7 factor(s):
- FACTOR 0 (double) SPARSE, size 10x20, density 0.1, nnz 20
- FACTOR 1 (double) SPARSE, size 20x22, density 0.113636, nnz 50
- FACTOR 2 (double) SPARSE, size 22x20, density 0.159091, nnz 70
- FACTOR 3 (double) SPARSE, size 20x21, density 0.142857, nnz 60
- FACTOR 4 (double) SPARSE, size 21x11, density 0.281385, nnz 65
- FACTOR 5 (double) SPARSE, size 11x24, density 0.19697, nnz 52
- FACTOR 6 (double) SPARSE, size 24x12, density 0.0833333, nnz 24
>>> norm((F-2).toarray() - F.toarray() + 2)
0.0
>>> F-G-2
Faust size 10x12, density 5.4, nnz_sum 648, 9 factor(s):
- FACTOR 0 (double) SPARSE, size 10x20, density 0.1, nnz 20
- FACTOR 1 (double) SPARSE, size 20x30, density 0.05, nnz 30
- FACTOR 2 (double) SPARSE, size 30x33, density 0.10101, nnz 100
- FACTOR 3 (double) SPARSE, size 33x30, density 0.126263, nnz 125
- FACTOR 4 (double) SPARSE, size 30x31, density 0.11828, nnz 110
- FACTOR 5 (double) SPARSE, size 31x32, density 0.105847, nnz 105
- FACTOR 6 (double) SPARSE, size 32x25, density 0.1225, nnz 98
- FACTOR 7 (double) SPARSE, size 25x24, density 0.06, nnz 36
- FACTOR 8 (double) SPARSE, size 24x12, density 0.0833333, nnz 24
>>> norm((F-G-2).toarray() - F.toarray() + G.toarray() + 2)
4.834253232627814e-15
>>> F-G-F-2-F
Faust size 10x12, density 11.05, nnz_sum 1326, 13 factor(s):
- FACTOR 0 (double) SPARSE, size 10x20, density 0.1, nnz 20
- FACTOR 1 (double) SPARSE, size 20x30, density 0.05, nnz 30
- FACTOR 2 (double) SPARSE, size 30x40, density 0.0333333, nnz 40
- FACTOR 3 (double) SPARSE, size 40x50, density 0.025, nnz 50
- FACTOR 4 (double) SPARSE, size 50x53, density 0.045283, nnz 120
- FACTOR 5 (double) SPARSE, size 53x52, density 0.0634978, nnz 175
- FACTOR 6 (double) SPARSE, size 52x51, density 0.0678733, nnz 180
- FACTOR 7 (double) SPARSE, size 51x55, density 0.0695187, nnz 195
- FACTOR 8 (double) SPARSE, size 55x54, density 0.0717172, nnz 213
- FACTOR 9 (double) SPARSE, size 54x36, density 0.063786, nnz 124
- FACTOR 10 (double) SPARSE, size 36x34, density 0.0743464, nnz 91
- FACTOR 11 (double) SPARSE, size 34x24, density 0.0784314, nnz 64
- FACTOR 12 (double) SPARSE, size 24x12, density 0.0833333, nnz 24
>>> norm((F-G-F-2-F).toarray() - F.toarray() + 2*F.toarray() + G.toarray() + 2)
1.4546887564889487e-14

See also

Faust.__add__()

__truediv__(s)

Divides F by the scalar s. This method overloads the Python function/operator `/’ (whether s is a float or an integer).

Args:
F: (Faust)

the Faust object.

s: (scalar)

the scalar to divide the Faust object with.

Returns:

the division result as a Faust object.

__weakref__

list of weak references to the object (if defined)

astype(dtype)

Converts F to the dtype passed as argument in a new Faust.

Args:
dtype: (str)

‘float32’, ‘float64’ or ‘complex’.

Returns:

A Faust copy of F converted to dtype.

Example:
>>> from pyfaust import rand
>>> F = rand(10, 10, dtype='float64')
>>> F.astype('float32')
Faust size 10x10, density 2.5, nnz_sum 250, 5 factor(s):
- FACTOR 0 (float) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 1 (float) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 2 (float) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 3 (float) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 4 (float) SPARSE, size 10x10, density 0.5, nnz 50
>>> F.astype('complex')
Faust size 10x10, density 2.5, nnz_sum 250, 5 factor(s):
- FACTOR 0 (complex) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 1 (complex) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 2 (complex) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 3 (complex) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 4 (complex) SPARSE, size 10x10, density 0.5, nnz 50
average(axis=None, weights=None, sw_returned=False)

Computes the weighted average of F along the specified axis.

Args:
axis: (None or int or tuple of ints)

Axis or axes along which to average the Faust F. The default, axis=None, will average over all of the elements of the input array. If axis is a tuple of ints, averaging is performed on all of the axes specified in the tuple

weights: (np.ndarray)

an array of weights associated with the values in F. Each value in F contributes to the average according to its associated weight. The weights array can either be 1-D (in which case its length must be the size of a along the given axis) or of the same shape as a. If weights=None, then all data in F are assumed to have a weight equal to one. The 1-D calculation is:

avg = sum(F @ weights) / sum(weights)

The only constraint on weights is that sum(weights) must not be 0.

returned: (bool)

True to return the sum of weights in addition to average (as a pair (avg, sum(weights))), False otherwise (default).

Returns:

The Faust average.

Example:
>>> from pyfaust import Faust
>>> import numpy as np
>>> data = np.arange(1, 5).astype('double')
>>> data
array([1., 2., 3., 4.])
>>> F = Faust(data.reshape(1, data.shape[0]))
>>> FA = F.average()
>>> FA
Faust size 1x1, density 9, nnz_sum 9, 3 factor(s):
- FACTOR 0 (double) DENSE, size 1x1, density 1, nnz 1
- FACTOR 1 (double) DENSE, size 1x4, density 1, nnz 4
- FACTOR 2 (double) DENSE, size 4x1, density 1, nnz 4
>>> FA.toarray()
array([[2.5]])
>>> data2 = np.arange(6).reshape((3, 2)).astype('double')
>>> F2 = Faust(data2)
>>> F2.average(axis=1, weights=[1./4, 3./4]).toarray()
array([[0.75],
       [2.75],
       [4.75]])
clone(dev=None)

Clones the Faust (in a new memory space).

Args:
dev: (str)

‘cpu’ to clone on CPU RAM, ‘gpu’ to clone on the GPU device. By default (None), the device is the F.device.

Returns:

The Faust clone.

Example:
>>> from pyfaust import rand, is_gpu_mod_enabled
>>> F = rand(10, 10)
>>> F.clone()
Faust size 10x10, density 2.5, nnz_sum 250, 5 factor(s):
- FACTOR 0 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 1 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 2 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 3 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 4 (double) SPARSE, size 10x10, density 0.5, nnz 50
>>> if is_gpu_mod_enabled(): F.clone(dev='gpu') # only if a NVIDIA compatible GPU is available
>>> #- GPU FACTOR 0 (double) SPARSE size 10 x 10, addr: 0x2c85390, density 0.500000, nnz 50
>>> #- GPU FACTOR 1 (double) SPARSE size 10 x 10, addr: 0x7ff00c0, density 0.500000, nnz 50
>>> #- GPU FACTOR 2 (double) SPARSE size 10 x 10, addr: 0x977f280, density 0.500000, nnz 50
>>> #- GPU FACTOR 3 (double) SPARSE size 10 x 10, addr: 0x9780120, density 0.500000, nnz 50
>>> #- GPU FACTOR 4 (double) SPARSE size 10 x 10, addr: 0x9780fc0, density 0.500000, nnz 50
concatenate(*args, **kwargs)

Concatenates F with len(args) Faust objects, numpy arrays or scipy sparse matrices. The resulting Faust:

C = F.concatenate(G, H, … , axis)

verifies that:

C.toarray() == numpy.concatenate((F.toarray(), G.toarray(), H.toarray(), …), axis)

N.B.: you could have an elementwise non-significant absolute difference between the two members.

Note

it could be wiser to encapsulate a Faust in a

lazylinop.LazyLinearOp for a lazy concatenation.

Args:
F: (Faust)

the Faust to concatenate to.

args: the Fausts or matrices (numpy array or scipy.sparse.csr/csc_matrix)

The objects to be concatenated to F. If args[i] is a matrix it will be Faust-converted on the fly.

axis: (int)

the dimension index (0 or 1) along to concatenate the Faust objects. By default, the axis is 0 (for vertical concatenation).

Returns:

The concatenation result as a Faust.

Raises:

ValueError

Examples:
>>> from pyfaust import rand, seed
>>> F = rand(50, 50)
>>> G = rand(50, 50)
>>> F.concatenate(G) # equivalent to F.concatenate(G, 0)
Faust size 100x50, density 0.52, nnz_sum 2600, 6 factor(s):
- FACTOR 0 (double) SPARSE, size 100x100, density 0.05, nnz 500
- FACTOR 1 (double) SPARSE, size 100x100, density 0.05, nnz 500
- FACTOR 2 (double) SPARSE, size 100x100, density 0.05, nnz 500
- FACTOR 3 (double) SPARSE, size 100x100, density 0.05, nnz 500
- FACTOR 4 (double) SPARSE, size 100x100, density 0.05, nnz 500
- FACTOR 5 (double) SPARSE, size 100x50, density 0.02, nnz 100
>>> F.concatenate(G, axis=1)
Faust size 50x100, density 0.52, nnz_sum 2600, 6 factor(s):
- FACTOR 0 (double) SPARSE, size 50x100, density 0.02, nnz 100
- FACTOR 1 (double) SPARSE, size 100x100, density 0.05, nnz 500
- FACTOR 2 (double) SPARSE, size 100x100, density 0.05, nnz 500
- FACTOR 3 (double) SPARSE, size 100x100, density 0.05, nnz 500
- FACTOR 4 (double) SPARSE, size 100x100, density 0.05, nnz 500
- FACTOR 5 (double) SPARSE, size 100x100, density 0.05, nnz 500
>>> from numpy.random import rand
>>> F.concatenate(rand(34, 50), axis=0) # The random array is auto-converted to a Faust before the vertical concatenation
Faust size 84x50, density 0.773809, nnz_sum 3250, 6 factor(s):
- FACTOR 0 (double) SPARSE, size 84x100, density 0.232143, nnz 1950
- FACTOR 1 (double) SPARSE, size 100x100, density 0.03, nnz 300
- FACTOR 2 (double) SPARSE, size 100x100, density 0.03, nnz 300
- FACTOR 3 (double) SPARSE, size 100x100, density 0.03, nnz 300
- FACTOR 4 (double) SPARSE, size 100x100, density 0.03, nnz 300
- FACTOR 5 (double) SPARSE, size 100x50, density 0.02, nnz 100
>>> from scipy.sparse import rand as sprand
>>> F.concatenate(sprand(50, 24, format='csr'), axis=1) # The sparse random matrix is auto-converted to a Faust before the horizontal concatenation
Faust size 50x74, density 0.422162, nnz_sum 1562, 6 factor(s):
- FACTOR 0 (double) SPARSE, size 50x100, density 0.02, nnz 100
- FACTOR 1 (double) SPARSE, size 100x100, density 0.03, nnz 300
- FACTOR 2 (double) SPARSE, size 100x100, density 0.03, nnz 300
- FACTOR 3 (double) SPARSE, size 100x100, density 0.03, nnz 300
- FACTOR 4 (double) SPARSE, size 100x100, density 0.03, nnz 300
- FACTOR 5 (double) SPARSE, size 100x74, density 0.0354054, nnz 262
>>> F.concatenate(F, G, F, G, rand(34, 50), F, G) # it's allowed to concatenate an arbitrary number of Fausts
Faust size 384x50, density 0.575521, nnz_sum 11050, 6 factor(s):
- FACTOR 0 (double) SPARSE, size 384x400, density 0.0224609, nnz 3450
- FACTOR 1 (double) SPARSE, size 400x400, density 0.01125, nnz 1800
- FACTOR 2 (double) SPARSE, size 400x400, density 0.01125, nnz 1800
- FACTOR 3 (double) SPARSE, size 400x400, density 0.01125, nnz 1800
- FACTOR 4 (double) SPARSE, size 400x400, density 0.01125, nnz 1800
- FACTOR 5 (double) SPARSE, size 400x50, density 0.02, nnz 400
>>> from pyfaust import rand
>>> F = rand(2, 51)
>>> G = rand(2, 25)
>>> F.concatenate(G, 0)
Traceback (most recent call last):
...
ValueError: The dimensions of the two Fausts must agree.
>>> from pyfaust import rand
>>> F = rand(2, 50);
>>> G = rand(2, 50);
>>> F.concatenate(G, axis=5)
Traceback (most recent call last):
...
ValueError: Axis must be 0 or 1.
>>> from pyfaust import rand
>>> F = rand(2, 5);
>>> G = rand(2, 5);
>>> F.concatenate(G, 'a')
Traceback (most recent call last):
...
ValueError: You can't concatenate a Faust with something that is not a Faust, a numpy array or scipy sparse matrix.
conj()

Returns the complex conjugate of F.

Args:

F: the Faust object.

Returns:

a Faust object Fc implementing the conjugate of F.toarray() such that: <code>Fc.toarray() == F.toarray().conj()</code>

Examples:
>>> from pyfaust import rand
>>> F = rand(50, 50, field='complex')
>>> Fc = F.conj()
conjugate()

Returns the complex conjugate of F.

Args:

F: the Faust object.

Returns:

a Faust object Fc implementing the conjugate of F.toarray() such that: <code>Fc.toarray() == F.toarray().conjugate()</code>

Examples:
>>> from pyfaust import rand
>>> F = rand(50, 50, field='complex')
>>> Fc = F.conjugate()
copy(dev='cpu')

Clone alias function (here just to mimic numpy API).

See also

Faust.clone(),

density()

Calculates the density of F such that F.nnz_sum() == F.density() * F.size.

Note

A value of density below one indicates potential memory savings

compared to storing the corresponding dense matrix F.toarray(), as well as potentially faster matrix-vector multiplication when applying F @ x instead of F.toarray() @ x.

Note

A density above one is possible but prevents any saving.

Args:

F: the Faust object.

Returns:

the density value (float).

Examples:
>>> from pyfaust import rand
>>> F = rand(5, 50, density=.5)
>>> dens = F.density()
property device

Returns the device on which the Faust is located (‘cpu’ or ‘gpu’).

Example:
>>> import pyfaust as pf
>>> cpuF = pf.rand(5, 5, dev='cpu')
>>> cpuF.device
'cpu'
>>> if pf.is_gpu_mod_enabled(): gpuF = pf.rand(5, 5, dev='gpu')
>>> gpuF.device if pf.is_gpu_mod_enabled() else None
>>> # it should print 'gpu' if it is available
display()

Displays information about F.

Args:

F: the Faust object.

Examples:
>>> from pyfaust import rand, seed
>>> seed(42) # just for reproducibility
>>> F = rand(50, 100, [1, 2], [50, 100], .5)
>>> F.display()
Faust size 50x100, density 1.3, nnz_sum 6500, 2 factor(s):
- FACTOR 0 (double) SPARSE, size 50x87, density 0.494253, nnz 2150
- FACTOR 1 (double) SPARSE, size 87x100, density 0.5, nnz 4350
>>> F
Faust size 50x100, density 1.3, nnz_sum 6500, 2 factor(s):
- FACTOR 0 (double) SPARSE, size 50x87, density 0.494253, nnz 2150
- FACTOR 1 (double) SPARSE, size 87x100, density 0.5, nnz 4350
>>> print(F)
Faust size 50x100, density 1.3, nnz_sum 6500, 2 factor(s):
- FACTOR 0 (double) SPARSE, size 50x87, density 0.494253, nnz 2150
- FACTOR 1 (double) SPARSE, size 87x100, density 0.5, nnz 4350
dot(A, *args, **kwargs)

Performs equivalent operation of numpy.dot() between the Faust F and A. More specifically:

  • Scalar multiplication if A is a scalar but F*A is preferred.

  • Matrix multiplication if A is a Faust or numpy.ndarray/numpy.matrix but F @ A is preferred.

property dtype

Returns the dtype of the Faust.

This function is intended to be used as a property (see the examples).

Args:

F: the Faust object.

Returns:

the dtype of F, which can be float32, float64 or complex128.

Examples:
>>> from pyfaust import rand
>>> F = rand(5, 10, field='complex')
>>> F.dtype
dtype('complex128')
>>> F = rand(5, 10)
>>> F.dtype
dtype('float64')
>>> G = rand(5, 5, dtype='float32')
>>> G.dtype
dtype('float32')
factors(indices, as_faust=False)

Returns the i-th factor of F or a new Faust composed of F factors whose indices are listed in indices.

Note

Factors are copied in memory.

Args:
F: (Faust)

the Faust object.

indices: (list[int])

the indices of wanted factors.

as_faust: (bool)

True to return a Faust even if a single factor is asked, otherwise (as_faust == False) and a numpy array or a scipy sparse matrix is returned.

Returns:

if indices is a single index and as_faust == False: a copy of the i-th factor. Otherwise a new Faust composed of the factors of F pointed by indices (no copy is made). For a single factor (with as_faust == False), the matrix type is:

  • numpy.ndarray if it is a full storage matrix or,

  • scipy.sparse.csc.matrix_csc if it’s a sparse matrix of a

transposed Faust,

  • scipy.sparse.csr.csr_matrix if it’s a sparse matrix of a

non-transposed Faust.

  • a scipy.sparse.bsr matrix if the factor is a BSR matrix.

Raises:

ValueError.

Example:
>>> from pyfaust import rand
>>> F = rand(5, 10)
>>> f0 = F.factors(0)
>>> G = F.factors(range(3, 5)) # a new Faust composed of the two last factors of F
getH()

Returns the conjugate transpose of F.

Args:

F: the Faust object.

Returns:

a Faust object H implementing the conjugate transpose of F.toarray() such that: H.toarray() == F.toarray().getH()

Examples:
>>> from pyfaust import rand
>>> F = rand(50, 50, field='complex')
>>> H1 = F.getH()
>>> H2 = F.transpose()
>>> H2 = H2.conj()
>>> (H1.toarray() == H2.toarray()).all()
True
get_factor_nonopt(i)

DEPRECATED: use Faust.factors() Returns the i-th factor of F.

Args:
F: (Faust)

the Faust object.

i: (int)

the factor index.

Returns:

a copy of the i-th factor as a dense matrix (of type numpy.ndarray).

Raises:

ValueError.

Examples:
>>> from pyfaust import rand
>>> F = rand(5, 10)
>>> f0 = F.factors(0)
property imag

Returns the imaginary part of F as a Faust.

imshow(name='F')

Displays image of F’s full matrix and its factors.

Args:
F: (Faust)

the Faust object.

name: (str)

the displayed name on the plotted figure.

Examples: >>> from pyfaust import rand >>> import matplotlib.pyplot as plt >>> F = rand(10, 20, density=.5, field=’complex’) >>> F.imshow() >>> plt.show()

See also

Faust.display()

insert(i, new_factor)

Inserts new_factor at index i in a new Faust copy of F.

Note

this is not a true copy, only references to pre-existed factors

are copied.

Args:
i: (int)

The index of insertion.

new_factor: (a numpy array or a scipy matrix)

The factor to insert as the i-th factor of F.

Returns:

a copy of F with new_factor inserted at index i.

Raises:

ValueError: in case new_factor has invalid dimensions.

ValueError: if i is out of range.

TypeError: if i is not an integer.

Example:
>>> import numpy as np
>>> import pyfaust as pf
>>> from scipy.sparse import random
>>> F = pf.rand(10, 10, num_factors=5)
>>> S = random(10, 10, .2)
>>> G = F.insert(2, S)
>>> np.allclose((F.left(1) @ pf.Faust(S) @ F.right(2)).toarray(), G.toarray())
True
static isFaust(obj)

Returns True if obj is a Faust object, False otherwise.

Examples:
>>> from pyfaust import *
>>> Faust.isFaust(2) # isFaust(2) works as well
False
>>> Faust.isFaust(rand(5, 10))
True
isdense()

Returns True if all factors are dense arrays (as np.ndarray-s) False otherwise.

Example:
>>> import pyfaust as pf
>>> pf.seed(42) # just for reproducibility
>>> F = pf.rand(10, 10, fac_type='sparse')
>>> F
Faust size 10x10, density 2.5, nnz_sum 250, 5 factor(s):
- FACTOR 0 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 1 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 2 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 3 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 4 (double) SPARSE, size 10x10, density 0.5, nnz 50
>>> F.isdense()
False
>>> F = pf.rand(10, 10, fac_type='dense')
>>> F
Faust size 10x10, density 2.5, nnz_sum 250, 5 factor(s):
- FACTOR 0 (double) DENSE, size 10x10, density 0.5, nnz 50
- FACTOR 1 (double) DENSE, size 10x10, density 0.5, nnz 50
- FACTOR 2 (double) DENSE, size 10x10, density 0.5, nnz 50
- FACTOR 3 (double) DENSE, size 10x10, density 0.5, nnz 50
- FACTOR 4 (double) DENSE, size 10x10, density 0.5, nnz 50
>>> F.isdense()
True
>>> F = pf.rand(10, 10, fac_type='mixed')
>>> F
Faust size 10x10, density 2.5, nnz_sum 250, 5 factor(s):
- FACTOR 0 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 1 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 2 (double) DENSE, size 10x10, density 0.5, nnz 50
- FACTOR 3 (double) DENSE, size 10x10, density 0.5, nnz 50
- FACTOR 4 (double) DENSE, size 10x10, density 0.5, nnz 50
>>> F.isdense()
False
issparse(csr=True, bsr=False)

Returns True if all F factors are sparse False otherwise.

What a sparse factor is, depends on csr and bsr arguments.

Args:
csr: (bool)

True to consider CSR matrices in F as sparse matrices, False otherwise.

bsr: (bool)

True to consider BSR matrices in F as sparse matrices, False otherwise.

Example:
>>> import pyfaust as pf
>>> F = pf.rand(10, 10, fac_type='sparse')
>>> F
Faust size 10x10, density 2.5, nnz_sum 250, 5 factor(s):
- FACTOR 0 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 1 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 2 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 3 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 4 (double) SPARSE, size 10x10, density 0.5, nnz 50
>>> F.issparse()
True
>>> F.issparse(csr=True)
True
>>> F.issparse(csr=False)
Traceback (most recent call last):
...
ValueError: It doesn't make sense to set csr=False and bsr=False as the function will always return False
>>> F = pf.rand_bsr(10, 10, 2, 2, 2, .1)
>>> F
Faust size 10x10, density 0.24, nnz_sum 24, 2 factor(s):
- FACTOR 0 (double) BSR, size 10x10 (blocksize = 2x2), density 0.12, nnz 12 (nnz blocks: 3)
- FACTOR 1 (double) BSR, size 10x10 (blocksize = 2x2), density 0.12, nnz 12 (nnz blocks: 3)
>>> # default config. recognizes only csr
>>> F.issparse()
False
>>> F.issparse(bsr=True)
True
>>> F = pf.rand(10, 10, fac_type='dense')
>>> F.issparse()
False
left(i, as_faust=False)

Returns the left hand side factors of F from index 0 to i included.

Args:
F: (Faust)

the Faust from which to extract left factors.

i: (int)

the far right index of left factors to extract.

as_faust: (bool)

True to return a Faust even if a single factor is asked (i.e.: F.left(0, as_faust=True) is a Faust, F.left(0) is not).

Returns:

a Faust if the number of factors to be returned is greater than 1 or if as_faust == True, a numpy array or a sparse matrix otherwise.

Examples:
>>> from pyfaust import rand, seed
>>> seed(42) # just for reproducibility
>>> F = rand(5, 10, 5)
>>> LF = F.left(3)
>>> print(F)
Faust size 5x10, density 2.98, nnz_sum 149, 5 factor(s):
- FACTOR 0 (double) SPARSE, size 5x9, density 0.555556, nnz 25
- FACTOR 1 (double) SPARSE, size 9x6, density 0.666667, nnz 36
- FACTOR 2 (double) SPARSE, size 6x7, density 0.714286, nnz 30
- FACTOR 3 (double) SPARSE, size 7x6, density 0.666667, nnz 28
- FACTOR 4 (double) SPARSE, size 6x10, density 0.5, nnz 30
>>> print(LF)
Faust size 5x6, density 3.96667, nnz_sum 119, 4 factor(s):
- FACTOR 0 (double) SPARSE, size 5x9, density 0.555556, nnz 25
- FACTOR 1 (double) SPARSE, size 9x6, density 0.666667, nnz 36
- FACTOR 2 (double) SPARSE, size 6x7, density 0.714286, nnz 30
- FACTOR 3 (double) SPARSE, size 7x6, density 0.666667, nnz 28
static load(filepath)

Loads a Faust from a MAT file.

The format is Matlab format version 5 and the filepath should end with a .mat extension.

The Faust must have been saved before with Faust.save.()

Args:
filepath: (str)

the filepath of the .mat file.

Example:
>>> from pyfaust import rand, Faust
>>> F = rand(5, 10, field='complex')
>>> F.save("F.mat")
>>> G = Faust.load(filepath="F.mat") # equiv. to Faust("F.mat")
load_native()

The format is Matlab format version 5 and the filepath should end with a .mat extension (native C++ version).

The Faust must have been saved before with Faust.save.()

Args:
filepath: (str)

the filepath of the .mat file.

Example:
>>> from pyfaust import rand, Faust
>>> F = rand(5, 10, field='complex')
>>> F.save("F.mat")
>>> G = Faust.load_native(filepath="F.mat") # equiv. to Faust("F.mat")
matvec(x)

This function implements the scipy.sparse.linalg.LinearOperator.matvec function such that scipy.sparse.linalg.aslinearoperator function works on a Faust object.

Example:
>>> import numpy as np
>>> import pyfaust as pf
>>> pf.seed(42) # just for reproducibility
>>> np.random.seed(42)
>>> F = pf.rand(10, 10)
>>> x = np.random.rand(F.shape[1])
>>> F@x
array([30.43795399, 55.16565252, 48.67306554, 34.64963666, 39.76690761,
       47.94326492, 24.18156012, 26.61375659, 43.28975657, 60.90302137])
>>> F.matvec(x)
array([30.43795399, 55.16565252, 48.67306554, 34.64963666, 39.76690761,
       47.94326492, 24.18156012, 26.61375659, 43.28975657, 60.90302137])
property nbytes

Gives the memory size of the Faust in bytes.

Example:
>>> import pyfaust as pf
>>> F = pf.rand(1024, 1024)
>>> F
Faust size 1024x1024, density 0.0244141, nnz_sum 25600, 5 factor(s):
- FACTOR 0 (double) SPARSE, size 1024x1024, density 0.00488281, nnz 5120
- FACTOR 1 (double) SPARSE, size 1024x1024, density 0.00488281, nnz 5120
- FACTOR 2 (double) SPARSE, size 1024x1024, density 0.00488281, nnz 5120
- FACTOR 3 (double) SPARSE, size 1024x1024, density 0.00488281, nnz 5120
- FACTOR 4 (double) SPARSE, size 1024x1024, density 0.00488281, nnz 5120
>>> F.nbytes
327700
>>> F = pf.rand(1024, 1024, fac_type='dense')
>>> F
Faust size 1024x1024, density 0.0244141, nnz_sum 25600, 5 factor(s):
- FACTOR 0 (double) DENSE, size 1024x1024, density 0.00488281, nnz 5120
- FACTOR 1 (double) DENSE, size 1024x1024, density 0.00488281, nnz 5120
- FACTOR 2 (double) DENSE, size 1024x1024, density 0.00488281, nnz 5120
- FACTOR 3 (double) DENSE, size 1024x1024, density 0.00488281, nnz 5120
- FACTOR 4 (double) DENSE, size 1024x1024, density 0.00488281, nnz 5120
>>> F.nbytes
41943040
property ndim

Number of Faust dimensions (always 2).

nnz_sum()

Gives the total number of non-zero elements in the factors of F.

The function sums together the number of non-zero elements of each factor and returns the result. Note that for efficiency the sum is computed at Faust creation time and kept in cache.

Returns:

the number of non-zeros.

Example:
>>> import pyfaust as pf
>>> F = pf.rand(1024, 1024)
>>> F
Faust size 1024x1024, density 0.0244141, nnz_sum 25600, 5 factor(s):
- FACTOR 0 (double) SPARSE, size 1024x1024, density 0.00488281, nnz 5120
- FACTOR 1 (double) SPARSE, size 1024x1024, density 0.00488281, nnz 5120
- FACTOR 2 (double) SPARSE, size 1024x1024, density 0.00488281, nnz 5120
- FACTOR 3 (double) SPARSE, size 1024x1024, density 0.00488281, nnz 5120
- FACTOR 4 (double) SPARSE, size 1024x1024, density 0.00488281, nnz 5120
>>> F.nnz_sum()
25600

See also

Faust.rcg(), Faust.density.()

norm(ord='fro', **kwargs)

Computes the norm of F. Several types of norm are available: 1-norm, 2-norm, inf-norm and Frobenius norm. The norm of F is equal to the numpy.linalg.norm of F.toarray().

Warning

The norm computation time can be expected to be of order n*min(F.shape) with n the time for multipliying F by a vector. Nevertheless, the implementation allows that memory usage remains controlled by avoiding to explicitly compute F.toarray(). Please pay attention to the full_array (and batch_size) arguments for a better understanding.

Args:
F: (Faust)

the Faust object.

ord: (int or str)

the norm order (1, 2, numpy.inf) or “fro” for Frobenius norm (by default the Frobenius norm is computed).

threshold: (float)

power iteration algorithm threshold (default to .001). Used only for norm(2).

max_num_its: (int)

maximum number of iterations for power iteration algorithm (default to 100). Used only for norm(2).

full_array: (bool)

this argument applies only for 1-norm, inf-norm and Frobenius norm. If True the Faust full array is computed before computing the norm otherwise it is not. By default it is set to False. Many configurations exist in which full_array == False can be more efficient but it needs to finetune the batch_size argument.

batch_size: (int)

this argument applies only when the full_array argument is set to False (for the 1-norm, inf-norm and Frobenius norm). It determines the number of Faust columns (resp. rows) that are built in memory in order to compute the Frobenius norm and the 1-norm (resp. the inf-norm). This parameter is primary in the efficiency of the computation and memory consumption. By default, it is set to 1 (which is certainly not the optimal configuration in many cases in matter of computation time but always the best in term of memory cost).

Returns:

the norm (float).

  • If ord == 1,

    the norm is norm(F.toarray(), 1) == max(sum(abs(F.toarray()))).

  • If ord == 2,

    the norm is the maximum singular value of F or approximately norm(F.toarray(), 2) == max(scipy.linalg.svd(F.toarray())[1]). This is the default norm calculated when calling to norm(F).

  • If ord == numpy.inf,

    the norm is norm(F.toarray(), numpy.inf) == max(sum(abs(F.T.toarray())))

  • If ord == ‘fro’,

    the norm is norm(F.toarray(), ‘fro’).

Raises:

ValueError.

Examples:
>>> from pyfaust import rand, seed
>>> import numpy as np
>>> seed(42) # just for reproducibility
>>> F = rand(50, 100, [1, 2], density=.5)
>>> F.norm()
388.2689201639318
>>> F.norm(2)
382.3775910865066
>>> F.norm('fro')
388.2689201639318
>>> F.norm(np.inf)
624.0409076619496
normalize(ord='fro', axis=1)

Normalizes F along the axis dimension using the ord-norm. The function is able to normalize the columns (default axis=1):

NF = F.normalize(ord) is such that for all i in range(0, F.shape[1]) NF.toarray()[:, i] == F.toarray()[:, i]/norm(F.toarray(), ord)

Likewise for normalization of the rows (axis=0):

NF = F.normalize(ord, axis=0) is such that for all i in range(0, F.shape[0]) NF.toarray()[i, :] == F.toarray()[i, :]/norm(F.toarray(), ord)

The variable ord designates one of the Faust.norm() compatible norms.

Args:

ord: the norm order to use (see Faust.norm()). axis: if 1 the columns are normalized, if 0 the rows.

Returns:

the normalized Faust

Example:
>>> from numpy.linalg import norm
>>> import numpy as np
>>> import pyfaust as pf
>>> pf.seed(42) # just for reproducibility
>>> F = pf.rand(10, 10)
>>> F
Faust size 10x10, density 2.5, nnz_sum 250, 5 factor(s):
- FACTOR 0 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 1 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 2 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 3 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 4 (double) SPARSE, size 10x10, density 0.5, nnz 50
>>> # normalize columns according to default fro-norm/2-norm
>>> # then test the second column is properly normalized
>>> nF2 = F.normalize()
>>> nF2
Faust size 10x10, density 2.5, nnz_sum 250, 5 factor(s):
- FACTOR 0 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 1 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 2 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 3 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 4 (double) SPARSE, size 10x10, density 0.5, nnz 50
>>> norm(nF2[:, 1].toarray() - F[:, 1].toarray()/F[:, 1].norm())
1.0385185452638061e-16
>>>
>>> # this time normalize rows using 1-norm and test the third row
>>> nF1 = F.normalize(ord=1, axis=0)
>>> nF1
Faust size 10x10, density 2.5, nnz_sum 250, 5 factor(s):
- FACTOR 0 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 1 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 2 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 3 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 4 (double) SPARSE, size 10x10, density 0.5, nnz 50
>>> norm(nF1[2, :].toarray() - F[2, :].toarray()/F[2, :].norm(ord=1))
2.5438405243138006e-16
>>>
>>> # and the same with inf-norm
>>> nFinf = F.normalize(ord=np.inf, axis=0)
>>> norm(nFinf[2, :].toarray() - F[2, :].toarray()/F[2, :].norm(ord=np.inf))
5.238750013840908e-17

See also

Faust.norm()

numfactors()

Returns the number of factors of F.

Note

using len(F) is more shorter!

Returns:

the number of factors.

Examples: >>> from pyfaust import rand >>> F = rand(100, 100, num_factors=2, density=.5) >>> nf = F.numfactors() >>> nf 2 >>> nf == len(F) True

optimize(transp=False)

Returns a Faust(), optimized with Faust.pruneout(), Faust.optimize_memory(), and Faust.optimize_time().

Args:
transp: (bool)

True in order to optimize the Faust according to its transpose.

Returns:

The optimized Faust.

Note

this function is still experimental, you might use manually Faust.optimize_time(), Faust.optimize_memory() or Faust.pruneout() to be more specific about the optimization to proceed.

Example:

This example shows how Faust.optimize() can diminish the memory size and speed up the Faust.toarray() and Faust-vector multiplication.

>>> import numpy as np
>>> from pyfaust import rand, seed
>>> seed(42) # just for reproducibility
>>> F = rand(1024, 1024, dim_sizes=[1, 1024], num_factors=32, fac_type='mixed')
>>> F.nbytes
13991076
>>> pF = F.optimize()
>>> pF.nbytes
910368
>>> from time import time
>>> t1 = time(); M = F.toarray(); print("F.toarray() time:", time()-t1) 
F.toarray() time: ...
>>> # e.g: F.toarray() time: 0.2779221534729004
>>> t1 = time(); M_ = pF.toarray(); print("pF.toarray() time:", time()-t1) 
pF.toarray() time: ...
>>> # e.g:pF.toarray() time: 0.2017652988433838
>>> np.allclose(M_, M)
True
>>> v = np.random.rand(F.shape[1])
>>> t1 = time(); Fv = F@v;print("F@v time:", time()-t1) 
F@v time: ...
>>> # e.g: F@v time: 0.0016832351684570312
>>> t1 = time(); Fv_ = pF@v; print("pF@v time:", time()-t1) 
pF@v time: ...
>>> # e.g: pF@v time: 0.0002257823944091797
>>> np.allclose(Fv_, Fv)
True
optimize_memory()

Optimizes a Faust by changing the storage format of each factor in order to optimize the memory size.

Returns:

The optimized Faust.

Example:
>>> from pyfaust import rand, seed
>>> seed(42) # just for reproducibility
>>> F = rand(1024, 1024, fac_type='mixed')
>>> F
Faust size 1024x1024, density 0.0244141, nnz_sum 25600, 5 factor(s):
- FACTOR 0 (double) SPARSE, size 1024x1024, density 0.00488281, nnz 5120
- FACTOR 1 (double) SPARSE, size 1024x1024, density 0.00488281, nnz 5120
- FACTOR 2 (double) DENSE, size 1024x1024, density 0.00488281, nnz 5120
- FACTOR 3 (double) SPARSE, size 1024x1024, density 0.00488281, nnz 5120
- FACTOR 4 (double) SPARSE, size 1024x1024, density 0.00488281, nnz 5120
>>> F.nbytes
8650768
>>> pF = F.optimize_memory()
>>> pF
Faust size 1024x1024, density 0.0244141, nnz_sum 25600, 5 factor(s):
- FACTOR 0 (double) SPARSE, size 1024x1024, density 0.00488281, nnz 5120
- FACTOR 1 (double) SPARSE, size 1024x1024, density 0.00488281, nnz 5120
- FACTOR 2 (double) SPARSE, size 1024x1024, density 0.00488281, nnz 5120
- FACTOR 3 (double) SPARSE, size 1024x1024, density 0.00488281, nnz 5120
- FACTOR 4 (double) SPARSE, size 1024x1024, density 0.00488281, nnz 5120
>>> pF.nbytes
327700

See also

Faust.optimize(),

optimize_time(transp=False, inplace=False, nsamples=1, mat=None)

Returns a Faust configured with the quickest Faust-matrix multiplication mode (benchmark ran on the fly).

Note

this function launches a small benchmark on the fly. Basically, the methods available differ by the order used to compute the matrix chain.

The evaluated methods in the benchmark are listed in pyfaust.FaustMulMode. Although depending on the package you installed and the capability of your hardware the methods based on Torch library can be used.

Args:
transp: (bool)

True in order to optimize the Faust according to its transpose.

inplace: (bool)

to optimize the current Faust directly instead of returning a new Faust with the optimization enabled. If True, F is returned otherwise a new Faust object is returned.

nsamples: (int)

the number of Faust-Dense matrix products calculated in order to measure time taken by each method (it could matter to discriminate methods when the performance is similar). By default, only one product is computed to evaluate the method.

mat: (NoneType, np.ndarray, or scipy.sparse.csr_matrix)

Use this argument to run the benchmark on the Faust multiplication by the matrix mat instead of Faust.toarray() (if mat is None). Note that mat must be of the same dtype as F.

Returns:

The optimized Faust.

Example:
>>> from pyfaust import rand, seed
>>> from time import time
>>> import numpy as np
>>> seed(42) # just for reproducibility
>>> F = rand(1024, 1024, dim_sizes=[1, 1024], num_factors=32, fac_type='dense', density=.5)
>>> oF = F.optimize_time()
>>> # possible outout: best method measured in time on operation Faust-toarray is: DYNPROG
>>> t1 = time(); M = F.toarray(); print("F.toarray() time:", time()-t1) 
F.toarray() time:...
>>> # F.toarray() time: 0.2891380786895752
>>> t1 = time(); M_ = oF.toarray(); print("oF.toarray() time:", time()-t1) 
oF.toarray() time:...
>>> # example: oF.toarray() 0.0172119140625
>>> np.allclose(M, M_)
True
pinv()

Computes the (Moore-Penrose) pseudo-inverse of F.toarray().

Warning

this function makes a call to Faust.toarray().()

Returns:

The dense pseudo-inverse matrix.

Example:
>>> from pyfaust import rand
>>> import numpy as np
>>> from numpy.linalg import pinv
>>> F = rand(128, 32)
>>> M = F.toarray()
>>> np.allclose(F.pinv(), pinv(M))
True

See also numpy.linalg.pinv, pyfaust.fact.pinvtj()

power_iteration(threshold=0.001, maxiter=100)

Performs the power iteration algorithm to compute the greatest eigenvalue of the Faust.

For the algorithm to succeed the Faust should be diagonalizable (similar to a digonalizable Faust), ideally, a symmetric positive-definite Faust.

Args:

threshold: the precision required on the eigenvalue. maxiter: the number of iterations above what the algorithm will stop anyway.

Returns:

The greatest eigenvalue approximate.

Examples:
>>> from pyfaust import rand, seed
>>> seed(42) # just for reproducibility
>>> F = rand(50, 50)
>>> F = F@F.H
>>> F.power_iteration()
12653.806783532553
pruneout(thres=0, **kwargs)

Returns a Faust optimized by removing useless zero rows and columns as many times as needed.

Args:
F: (Faust)

the Faust to optimize.

thres: (int)

the threshold in number of nonzeros under what the rows/columns are removed.

Returns:

The optimized Faust.

Example:
>>> from pyfaust import rand, seed
>>> seed(42) # just for reproducibility
>>> F = rand(1024, 1024, dim_sizes=[1, 1024], num_factors=64, fac_type='mixed')
>>> pF = F.pruneout()
>>> F.nbytes
49109760
>>> pF.nbytes
46535972

See also

Faust.optimize()

rcg()

Computes the Relative Complexity Gain.

The RCG is the theoretical gain brought by the Faust representation relatively to its dense matrix equivalent. The higher is the RCG, the more computational savings are made. This gain applies both for storage space and computation time.

Note

F.rcg() == 1/F.density()

Args:

F: the Faust object.

Returns:

the RCG value (float).

Examples:
>>> from pyfaust import rand
>>> F = rand(1024, 1024)
>>> F.rcg()
40.96
>>> F
Faust size 1024x1024, density 0.0244141, nnz_sum 25600, 5 factor(s):
- FACTOR 0 (double) SPARSE, size 1024x1024, density 0.00488281, nnz 5120
- FACTOR 1 (double) SPARSE, size 1024x1024, density 0.00488281, nnz 5120
- FACTOR 2 (double) SPARSE, size 1024x1024, density 0.00488281, nnz 5120
- FACTOR 3 (double) SPARSE, size 1024x1024, density 0.00488281, nnz 5120
- FACTOR 4 (double) SPARSE, size 1024x1024, density 0.00488281, nnz 5120
>>> F.size/F.nnz_sum()
40.96

See also

Faust.density(), Faust.nnz_sum(), Faust.shape.()

property real

Returns the real part of F as a Faust.

replace(i, new_factor)

Replaces the factor of index i by new_factor in a new Faust copy of F.

Note

this is not a true copy, only references to pre-existed factors

are copied.

Args:
i: (int)

The factor index to replace.

new_factor: (a numpy array or a scipy matrix)

The factor replacing the i-th factor of F.

Returns:

a copy of F with the i-th factor replaced by new_factor.

Raises:

ValueError: in case new_factor has invalid dimensions.

ValueError: if i is out of range.

TypeError: if i is not an integer.

Example:
>>> import numpy as np
>>> import pyfaust as pf
>>> from scipy.sparse import random
>>> F = pf.rand(10, 10, num_factors=5)
>>> S = random(10, 10, .2)
>>> G = F.replace(2, S)
>>> np.allclose((F.left(1) @ pf.Faust(S) @ F.right(3)).toarray(), G.toarray())
True
right(i, as_faust=False)

Returns the right hand side factors of F from index i to F.numfactors()-1.

Args:
F: (Faust)

the Faust from which to extract right factors.

i: (int)

the far left index of right factors to extract.

as_faust: (bool)

True to return a Faust even if a single factor is asked (i.e.: F.right(len(F)-1, as_faust=True) is a Faust, F.left(len(F)-1) is not).

Returns:

a Faust if the number of factors to be returned is greater than 1 or if as_faust == True, a numpy array or a sparse matrix otherwise.

Examples:
>>> from pyfaust import rand, seed
>>> seed(42) # just for reproducibility
>>> F = rand(5, 10, 5)
>>> RF = F.right(2)
>>> print(F)
Faust size 5x10, density 2.98, nnz_sum 149, 5 factor(s):
- FACTOR 0 (double) SPARSE, size 5x9, density 0.555556, nnz 25
- FACTOR 1 (double) SPARSE, size 9x6, density 0.666667, nnz 36
- FACTOR 2 (double) SPARSE, size 6x7, density 0.714286, nnz 30
- FACTOR 3 (double) SPARSE, size 7x6, density 0.666667, nnz 28
- FACTOR 4 (double) SPARSE, size 6x10, density 0.5, nnz 30
>>> print(RF)
Faust size 6x10, density 1.46667, nnz_sum 88, 3 factor(s):
- FACTOR 0 (double) SPARSE, size 6x7, density 0.714286, nnz 30
- FACTOR 1 (double) SPARSE, size 7x6, density 0.666667, nnz 28
- FACTOR 2 (double) SPARSE, size 6x10, density 0.5, nnz 30
save(filepath, format='Matlab')

Saves the Faust F into a file.

The file is saved in Matlab format version 5 (.mat extension).

Note

storing F should typically use rcg(F) times less disk space than

storing F.toarray(). See Faust.nbytes() for a precise size.

Args:
F: (Faust)

the Faust object.

filepath: (str)

the path for saving the Faust (should end with .mat if Matlab format is used).

format: (str)

The format to use for writing. By default, it’s “Matlab” to save the Faust in a .mat file (currently only that format is available).

Raises:

ValueError.

Example:
>>> from pyfaust import rand, Faust
>>> F = rand(5, 10, field='complex')
>>> F.save("F.mat")
>>> G = Faust(filepath="F.mat")
property shape

Gives the shape of the Faust F.

This function is intended to be used as a property (see the examples).

The shape is a pair of numbers: the number of rows and the number of columns of F.todense().

Args:

F: the Faust object.

Returns:

the Faust shape tuple, with at first index the number of rows, and at second index the number of columns.

Examples:
>>> from pyfaust import rand
>>> F = rand(50, 50, field='complex')
>>> nrows, ncols = F.shape
>>> nrows = F.shape[0]
>>> ncols = F.shape[1]
property size

Gives the size of the Faust F (that is F.shape[0]*F.shape[1]) .

It’s equivalent to numpy.prod(F.shape)).

This function is intended to be used as a property (see the examples).

Args:

F: the Faust object.

Returns:

The number of elements in the Faust F.

Examples:
>>> from pyfaust import rand
>>> F = rand(50, 50, field='complex')
>>> F.size
2500

See also

Faust.shape()

sum(axis=None, **kwargs)

Sums Faust elements over a given axis.

Args:
axis: (None or int or tuple of ints)

Axis or axes along which the the sum is performed

Returns:

The Faust sum.

Example:
>>> from pyfaust import rand as frand, seed
>>> seed(42) # just for reproducibility
>>> F = frand(10, 10)
>>> F.sum()
Faust size 1x1, density 270, nnz_sum 270, 7 factor(s):
- FACTOR 0 (double) DENSE, size 1x10, density 1, nnz 10
- FACTOR 1 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 2 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 3 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 4 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 5 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 6 (double) DENSE, size 10x1, density 1, nnz 10
>>> F.sum(axis=0).toarray()
array([[ 80.49374154,  45.09382766,  78.8607646 , 136.65920307,
        115.25872767,  62.70720879,  90.48774161,  62.26010951,
          0.        , 158.34355964]])
>>> F.sum(axis=1)
Faust size 10x1, density 26, nnz_sum 260, 6 factor(s):
- FACTOR 0 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 1 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 2 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 3 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 4 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 5 (double) DENSE, size 10x1, density 1, nnz 10
>>> F.sum(axis=1).toarray()
array([[ 60.96835408],
       [112.00753373],
       [ 97.28363377],
       [ 69.45346776],
       [ 80.20182354],
       [ 96.1022012 ],
       [ 48.61063623],
       [ 54.15081295],
       [ 88.20644303],
       [123.1799778 ]])
swap_cols(id1, id2, permutation=False, inplace=False)

Swaps F columns of indices id1 and id2.

Args:
id1: (int)

index of the first column of the swap.

id2: (int)

index of the second column of the swap.

permutation: (bool)

if True then the swap is performed by inserting a permutation matrix to the output Faust. If False, the last matrix in the Faust F sequence is edited to swap the columns.

inplace: (bool)

if True then F is modified instead of returning a new Faust. Otherwise, by default, a new Faust is returned.

Returns:

The column swapped Faust.

Example:
>>> from pyfaust import rand as frand
>>> F = frand(10, 10)
>>> G = F.swap_cols(2, 4)
>>> G
Faust size 10x10, density 2.5, nnz_sum 250, 5 factor(s):
- FACTOR 0 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 1 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 2 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 3 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 4 (double) SPARSE, size 10x10, density 0.5, nnz 50
>>> G[:, 2].toarray() == F[:, 4].toarray()
array([[ True],
       [ True],
       [ True],
       [ True],
       [ True],
       [ True],
       [ True],
       [ True],
       [ True],
       [ True]])
>>> G[:, 4].toarray() == F[:, 2].toarray()
array([[ True],
       [ True],
       [ True],
       [ True],
       [ True],
       [ True],
       [ True],
       [ True],
       [ True],
       [ True]])
>>> H  = F.swap_cols(2, 4, permutation=True)
>>> H
Faust size 10x10, density 2.6, nnz_sum 260, 6 factor(s):
- FACTOR 0 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 1 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 2 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 3 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 4 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 5 (double) SPARSE, size 10x10, density 0.1, nnz 10
swap_rows(id1, id2, permutation=True, inplace=False)

Swaps F rows of indices id1 and id2.

Args:
id1: (int)

index of the first row of the swap.

id2: (int)

index of the second row of the swap.

permutation: (bool)

if True then the swap is performed by inserting a permutation matrix to the output Faust. If False, the last matrix in the Faust F sequence is edited to swap the rows.

inplace: (bool)

if True then F is modified instead of returning a new Faust. Otherwise, by default, a new Faust is returned.

Returns:

The rows swapped Faust.

Example:
>>> from pyfaust import rand as frand
>>> F = frand(10, 10)
>>> G = F.swap_rows(2, 4)
>>> G
Faust size 10x10, density 2.6, nnz_sum 260, 6 factor(s):
- FACTOR 0 (double) SPARSE, size 10x10, density 0.1, nnz 10
- FACTOR 1 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 2 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 3 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 4 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 5 (double) SPARSE, size 10x10, density 0.5, nnz 50
>>> G[2, :].toarray() == F[4, :].toarray()
array([[ True,  True,  True,  True,  True,  True,  True,  True,  True,
         True]])
>>> G[4, :].toarray() == F[2, :].toarray()
array([[ True,  True,  True,  True,  True,  True,  True,  True,  True,
         True]])
>>> H  = F.swap_rows(2, 4, permutation=True)
>>> H
Faust size 10x10, density 2.6, nnz_sum 260, 6 factor(s):
- FACTOR 0 (double) SPARSE, size 10x10, density 0.1, nnz 10
- FACTOR 1 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 2 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 3 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 4 (double) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 5 (double) SPARSE, size 10x10, density 0.5, nnz 50
toarray()

Returns a numpy array for the full matrix implemented by F.

Warning

Using this function is discouraged except for test purposes, as it loses the main potential interests of the Faust structure: compressed memory storage and faster matrix-vector multiplication compared to its equivalent full matrix representation.

Returns:

A numpy ndarray.

Raises:

MemoryError

Warning

running the example below is likely to raise a memory error or freeze your computer for a certain amount of time.

Examples:
>>> from pyfaust import rand
>>> F = rand(10**5, 10**5, 2, 10**5, density=10**-4, fac_type='sparse')
>>> F
Faust size 100000x100000, density 0.00018, nnz_sum 1800000, 2 factor(s):
- FACTOR 0 (double) SPARSE, size 100000x100000, density 9e-05, nnz 900000
- FACTOR 1 (double) SPARSE, size 100000x100000, density 9e-05, nnz 900000
>>> # an attempt to convert F to a dense matrix is most likely to raise a memory error
>>> # the sparse format is the only way to handle such a large Faust
>>> F.toarray()
Traceback (most recent call last):
...
numpy.core._exceptions._ArrayMemoryError: Unable to allocate 74.5 GiB for an array with shape (100000, 100000) and data type float64

See also

Faust.todense()

todense()

Returns a numpy matrix for the full matrix implemented by F.

Warning

Using this function is discouraged except for test purposes, as it loses the main potential interests of the Faust structure: compressed memory storage and faster matrix-vector multiplication compared to its equivalent full matrix representation.

Warning

this function is deprecated in favor to toarray function and will be removed in a future version. The cause behind is that this function returns a numpy.matrix which is deprecated too.

Returns:

A numpy matrix M such that M*x == F*x for any vector x.

Raises:

MemoryError

Warning

running the example below is likely to raise a memory error or freeze your computer for a certain amount of time.

Warning

this function is deprecated and might be deleted in future versions of pyfaust. Please use Faust.toarray() instead.

See also

Faust.toarray()

transpose()

Returns the transpose of F.

Args:

F: the Faust object.

Returns:

a Faust object implementing the transpose of F.toarray(), such that: <code>F.transpose().toarray() == F.toarray().transpose()</code>

Examples:
>>> from pyfaust import rand, seed
>>> F = rand(10, 18)
>>> tF = F.transpose()
>>> tF.shape
(18, 10)
class pyfaust.FaustMulMode

<b/> Enumeration class of all matrix chain multiplication methods available to multiply a Faust to a matrix or to compute Faust.toarray().()

These methods are used by Faust.optimize_time().()

Note

it’s not advisable to use these methods directly. The user should use

Faust.optimize_time() but the access is left open for experimental purpose.

Examples:
>>> from pyfaust import rand as frand, seed as fseed, FaustMulMode
>>> from numpy.random import rand, seed
>>> fseed(42) # just for reproducibility
>>> seed(42)
>>> F = frand(100, 100, 5, [100, 1024])
>>> F.m_faust.set_FM_mul_mode(FaustMulMode.DYNPROG) # method used to compute Faust-matrix product or Faust.toarray()
>>> F @ rand(F.shape[1], 512) # Faust-matrix mul. using method DYNPROG
array([[34.29346113, 33.5135258 , 31.83862847, ..., 32.89901332,
        36.90417709, 32.20140406],
       [45.40983901, 42.52512058, 42.14810308, ..., 44.2648802 ,
        48.4027215 , 40.55809844],
       [27.12859996, 28.26596387, 25.898984  , ..., 27.47460378,
        29.35053152, 25.24167465],
       ...,
       [48.42899773, 47.4001851 , 43.96370573, ..., 47.25218683,
        53.03379773, 46.12690926],
       [39.9583232 , 40.778263  , 36.65671168, ..., 40.69390161,
        43.55280684, 40.37781963],
       [26.92859133, 28.40176389, 24.73304576, ..., 27.72648267,
        28.78612539, 25.82727371]])
>>> F.toarray() # using the same method
array([[1.13340453, 0.94853857, 0.60713635, ..., 1.29155048, 0.98107444,
        0.30254208],
       [1.92753189, 0.45268035, 0.72474175, ..., 0.43439703, 1.21532731,
        0.50115957],
       [1.35028724, 0.30557493, 0.15632569, ..., 0.80602032, 0.56741856,
        0.58193385],
       ...,
       [0.63172715, 1.0883051 , 1.2760964 , ..., 0.45745425, 0.85951258,
        0.25173183],
       [0.85748924, 0.68716077, 0.96293286, ..., 1.09480206, 0.8219215 ,
        0.83602967],
       [0.5461995 , 0.36918089, 0.4556373 , ..., 0.57842966, 0.52784458,
        0.30465166]])
__weakref__

list of weak references to the object (if defined)

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()).

Args:
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.

Example:
>>> 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
pyfaust.bitrev(inds)

Bitreversal permutation.

Args:
inds: (list[int])

the list of indices to bit-reverse.

Returns:

The bit-reversal permutation of inds.

Example:
>>> import numpy as np
>>> from pyfaust import bitrev
>>> bitrev(np.arange(4))
array([0, 2, 1, 3])

See also: https://en.wikipedia.org/wiki/Bit-reversal_permutation.

pyfaust.bitrev_perm(n)

Bitreversal permutation.

Args:
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()

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()).

Args:
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).

Example:
>>> 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
pyfaust.concatenate(_tuple, *args, axis=0, **kwargs)

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

Example:
>>> 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(),

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)\)

Args:
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.

Example:
>>> 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
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.

Args:
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.()

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

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

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)\)

Args:
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.

Example:
>>> 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
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.

Args:
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.

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

Faust identity.

Args:
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
pyfaust.faust_fact(*args, **kwargs)

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

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

Example:
>>> 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()
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.

Example:
>>> 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
pyfaust.implements(numpy_function)

Register an __array_function__ implementation for MyArray objects.

pyfaust.isFaust(obj)

Package alias function of Faust.isFaust.()

See also

Faust.__init__(), Faust.isFaust.(),

pyfaust.is_gpu_mod_enabled()

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

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.

pyfaust.license()

Prints the FAuST license.

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()

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.

Args:

F: The 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.

pyfaust.pinv(F)

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

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.

Args:

num_rows: the Faust number of rows. num_cols: the Faust number of columns. num_factors: if 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_sizes: if 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].

density: the 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_type: the 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).

dtype: the dtype of the Faust (‘float32’, ‘float64’ or ‘complex’). per_row: if 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.

dev: the device on which to create the Faust (‘cpu’ or ‘gpu’). field: (DEPRECATED, use dtype) a str to set the Faust field: ‘real’ or ‘complex’. seed: set 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__, ` :py:func:`pyfaust.rand_bsr()

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

Generates a random Faust composed only of BSR matrices.

Args:
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.

Example:
>>> 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)
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).

Args:

n: order of the butterfly (must be a power of two). diag_opt: if 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.

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.

Example:
>>> 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
pyfaust.toeplitz(c, r=None, dev='cpu', diag_opt=False)

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

Args:
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.

Example:
>>> 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
pyfaust.version()

Returns the FAuST package version.

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.

Example:
>>> 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
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.

Args:
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(),