-  3.41.0
pyfaust.Faust Class Reference

pyfaust main class for using multi-layer sparse transforms. More...

Public Member Functions

def __init__ (F, factors=None, filepath=None, **kwargs)
 Creates a Faust from a list of factors or alternatively from a file. More...
 
def transpose (F)
 Returns the transpose of F. More...
 
def conj (F)
 Returns the complex conjugate of F. More...
 
def conjugate (F)
 Returns the complex conjugate of F. More...
 
def getH (F)
 Returns the conjugate transpose of F. More...
 
def pruneout (F, thres=0, **kwargs)
 Returns a Faust optimized by removing useless zero rows and columns as many times as needed. More...
 
def __repr__ (F)
 Returns a str object representing the Faust object. More...
 
def __str__ (F)
 Converts F to its str representation. More...
 
def display (F)
 Displays information about F. More...
 
def __pos__ (F)
 Returns + F (unary operator). More...
 
def __neg__ (F)
 Returns F (unary operator). More...
 
def __add__ (F, *args, **kwargs)
 Sums F to one or a sequence of variables (Faust objects, arrays or scalars). More...
 
def __radd__ (F, lhs_op)
 Returns lhs_op+F. More...
 
def __sub__ (F, *args)
 Subtracts from F one or a sequence of variables. More...
 
def __rsub__ (F, lhs_op)
 Returns lhs_op-F. More...
 
def __truediv__ (F, s)
 Divides F by the scalar s. More...
 
def __itruediv__ (F, s)
 Divides F by the scalar s inplace. More...
 
def __floordiv__ (F, s)
 F // s. More...
 
def __imatmul__ (F, A)
 Inplace F = F @ A. More...
 
def __imul__ (F, A)
 Inplace F = F * A. More...
 
def __iadd__ (F, A)
 Inplace F = F + A. More...
 
def __isub__ (F, A)
 Inplace F = F - A. More...
 
def __matmul__ (F, A)
 Multiplies F by A which is a dense numpy.matrix/numpy.ndarray or a Faust object. More...
 
def dot (F, A, *args, **kwargs)
 Performs equivalent operation of numpy.dot() between the Faust F and A. More...
 
def matvec (F, x)
 This function implements the scipy.sparse.linalg.LinearOperator.matvec function such that scipy.sparse.linalg.aslinearoperator function works on a Faust object. More...
 
def __mul__ (F, A)
 Multiplies the Faust F by A. More...
 
def __rmul__ (F, lhs_op)
 lhs_op*F More...
 
def __rmatmul__ (F, lhs_op)
 Returns lhs_op.__matmul__(F). More...
 
def concatenate (F, *args, **kwargs)
 Concatenates F with len(args) Faust objects, numpy arrays or scipy sparse matrices. More...
 
def toarray (F)
 Returns a numpy array for the full matrix implemented by F. More...
 
def todense (F)
 Returns a numpy matrix for the full matrix implemented by F. More...
 
def __getitem__ (F, indices)
 Indexes or slices a Faust. More...
 
def nnz_sum (F)
 Gives the total number of non-zero elements in the factors of F. More...
 
def density (F)
 Calculates the density of F such that F.nnz_sum() == F.density() * F.size. More...
 
def rcg (F)
 Computes the Relative Complexity Gain. More...
 
def norm (F, ord='fro', **kwargs)
 Computes the norm of F. More...
 
def power_iteration (self, threshold=1e-3, maxiter=100)
 Performs the power iteration algorithm to compute the greatest eigenvalue of the Faust. More...
 
def normalize (F, ord='fro', axis=1)
 Normalizes F along the axis dimension using the ord-norm. More...
 
def numfactors (F)
 Returns the number of factors of F. More...
 
def __len__ (F)
 Returns the number of factors of F. More...
 
def factors (F, 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. More...
 
def right (F, i, as_faust=False)
 Returns the right hand side factors of F from index i to F.numfactors()-1. More...
 
def left (F, i, as_faust=False)
 Returns the left hand side factors of F from index 0 to i included. More...
 
def replace (F, i, new_factor)
 Replaces the factor of index i by new_factor in a new Faust copy of F. More...
 
def insert (F, i, new_factor)
 Inserts new_factor at index i in a new Faust copy of F. More...
 
def save (F, filepath, format="Matlab")
 Saves the Faust F into a file. More...
 
def load_native (filepath)
 The format is Matlab format version 5 and the filepath should end with a .mat extension (native C++ version). More...
 
def astype (F, dtype)
 Converts F to the dtype passed as argument in a new Faust. More...
 
def imshow (F, name='F')
 Displays image of F's full matrix and its factors. More...
 
def pinv (F)
 Computes the (Moore-Penrose) pseudo-inverse of F.toarray(). More...
 
def issparse (F, csr=True, bsr=False)
 Returns True if all F factors are sparse False otherwise. More...
 
def isdense (F)
 Returns True if all factors are dense arrays (as np.ndarray-s) False otherwise. More...
 
def swap_cols (F, id1, id2, permutation=False, inplace=False)
 Swaps F columns of indices id1 and id2. More...
 
def swap_rows (F, id1, id2, permutation=True, inplace=False)
 Swaps F rows of indices id1 and id2. More...
 
def optimize_memory (F)
 Optimizes a Faust by changing the storage format of each factor in order to optimize the memory size. More...
 
def optimize (F, transp=False)
 Returns a Faust, optimized with Faust.pruneout, Faust.optimize_memory, and Faust.optimize_time. More...
 
def optimize_time (F, transp=False, inplace=False, nsamples=1, mat=None)
 Returns a Faust configured with the quickest Faust-matrix multiplication mode (benchmark ran on the fly). More...
 
def copy (F, dev='cpu')
 Clone alias function (here just to mimic numpy API). More...
 
def clone (F, dev=None)
 Clones the Faust (in a new memory space). More...
 
def sum (F, axis=None, **kwargs)
 Sums Faust elements over a given axis. More...
 
def average (F, axis=None, weights=None, sw_returned=False)
 Computes the weighted average of F along the specified axis. More...
 

Static Public Member Functions

def load (filepath)
 Loads a Faust from a MAT file. More...
 
def isFaust (obj)
 Returns True if obj is a Faust object, False otherwise. More...
 

Properties

 nbytes = property
 Gives the memory size of the Faust in bytes. More...
 
 shape = property
 Gives the shape of the Faust F. More...
 
 ndim = property
 Number of Faust dimensions (always 2). More...
 
 size = property
 Gives the size of the Faust F (that is F.shape[0]*F.shape[1]) . More...
 
 device = property
 Returns the device on which the Faust is located ('cpu' or 'gpu'). More...
 
 T = property
 Returns the transpose of F. More...
 
 H = property
 Returns the conjugate transpose of F. More...
 
 real = property
 Returns the real part of F as a Faust. More...
 
 imag = property
 Returns the imaginary part of F as a Faust. More...
 
 dtype = property
 Returns the dtype of the Faust. More...
 

Detailed Description

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.

See also
pyfaust.Faust.__init__

Constructor & Destructor Documentation

◆ __init__()

def pyfaust.Faust.__init__ (   F,
  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.

Parameters
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 (for i in range(0, len(factors)-1): factors[i].shape[1] == factors[i+1].shape[0]).
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
See also
Faust.save, pyfaust.rand, pyfaust.dft, pyfaust.wht, csr_matrix, bsr_matrix csc_matrix, coo_matrix

Member Function Documentation

◆ __add__()

def pyfaust.Faust.__add__ (   F,
args,
**  kwargs 
)

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

Note
This method overloads the Python function/operator +.
Parameters
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.

Examples

>>> 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
>>> float(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
>>> float(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
>>> float(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
>>> float(norm((F+G+F+2+F).toarray() - 3*F.toarray() - G.toarray() - 2))
2.7030225852818652e-14
See also
Faust.__sub__

◆ __floordiv__()

def pyfaust.Faust.__floordiv__ (   F,
  s 
)

F // s.

Warning
this operation is not supported, it raises an exception.

◆ __getitem__()

def pyfaust.Faust.__getitem__ (   F,
  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.
Parameters
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).
Exceptions
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__()

def pyfaust.Faust.__iadd__ (   F,
  A 
)

Inplace F = F + A.

◆ __imatmul__()

def pyfaust.Faust.__imatmul__ (   F,
  A 
)

Inplace F = F @ A.

◆ __imul__()

def pyfaust.Faust.__imul__ (   F,
  A 
)

Inplace F = F * A.

◆ __isub__()

def pyfaust.Faust.__isub__ (   F,
  A 
)

Inplace F = F - A.

◆ __itruediv__()

def pyfaust.Faust.__itruediv__ (   F,
  s 
)

Divides F by the scalar s inplace.

This method overloads the Python function/operator ‘/=’ (whether s is a float or an integer).

Parameters
F(Faust) the Faust object.
s(scalar) the scalar to divide the Faust object with.
Returns
the division result as a Faust object.
See also
Faust.__mul__, Faust.__truediv__

◆ __len__()

def pyfaust.Faust.__len__ (   F)

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
See also
Faust.factors, Faust.numfactors

◆ __matmul__()

def pyfaust.Faust.__matmul__ (   F,
  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).
Parameters
F(Faust) the Faust object.
Aa 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.


Exceptions
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

See also
Faust.__init__, Faust.rcg, Faust.__mul__, Faust.__rmatmul__, Faust.dot

◆ __mul__()

def pyfaust.Faust.__mul__ (   F,
  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.
Parameters
Fthe Faust object.
Aa 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.
Parameters
ReturnsThe result of the multiplication as a Faust (either A is a vector or a scalar).
RaisesTypeError

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

def pyfaust.Faust.__neg__ (   F)

Returns F (unary operator).

Note
This method overloads the Python unary operator -.
See also
Faust.__sub__, Faust.__mul__

◆ __pos__()

def pyfaust.Faust.__pos__ (   F)

Returns + F (unary operator).

Note
This method overloads the Python unary operator +.
See also
Faust.__add__

◆ __radd__()

def pyfaust.Faust.__radd__ (   F,
  lhs_op 
)

Returns lhs_op+F.

See also
Faust.__add__,

◆ __repr__()

def pyfaust.Faust.__repr__ (   F)

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.
Parameters
Fthe 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
See also
Faust.nnz_sum, Faust.rcg, Faust.shape, Faust.factors, Faust.numfactors, Faust.display

◆ __rmatmul__()

def pyfaust.Faust.__rmatmul__ (   F,
  lhs_op 
)

Returns lhs_op.__matmul__(F).

See also
Faust.__matmul__

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

def pyfaust.Faust.__rmul__ (   F,
  lhs_op 
)

lhs_op*F

See also
Faust.__mul__

◆ __rsub__()

def pyfaust.Faust.__rsub__ (   F,
  lhs_op 
)

Returns lhs_op-F.

See also
Faust.__sub__

◆ __str__()

def pyfaust.Faust.__str__ (   F)

Converts F to its str representation.

Returns
The str conversion of F.

◆ __sub__()

def pyfaust.Faust.__sub__ (   F,
args 
)

Subtracts from F one or a sequence of variables.

Faust objects, arrays or scalars.

Note
This method overloads the Python function/operator -.
Parameters
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.

Examples

>>> 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
>>> float(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
>>> float(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
>>> float(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
>>> float(norm((F-G-F-2-F).toarray() - F.toarray() + 2*F.toarray() + G.toarray() + 2))
1.4546887564889487e-14
See also
Faust.__add__

◆ __truediv__()

def pyfaust.Faust.__truediv__ (   F,
  s 
)

Divides F by the scalar s.

This method overloads the Python function/operator ‘/’ (whether s is a float or an integer).

Parameters
F(Faust) the Faust object.
s(scalar) the scalar to divide the Faust object with.
Returns
the division result as a Faust object.
See also
Faust.__mul__, Faust.__itruediv__

◆ astype()

def pyfaust.Faust.astype (   F,
  dtype 
)

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

Parameters
dtype(str) 'float32', 'float64' or 'complex'.
Returns
A Faust copy of F converted to dtype.

Examples

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

def pyfaust.Faust.average (   F,
  axis = None,
  weights = None,
  sw_returned = False 
)

Computes the weighted average of F along the specified axis.

Parameters
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.

Parameters
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.

Examples

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

def pyfaust.Faust.clone (   F,
  dev = None 
)

Clones the Faust (in a new memory space).

Parameters
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.

Examples

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

def pyfaust.Faust.concatenate (   F,
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.
Parameters
F(Faust) the Faust to concatenate to.
argsthe 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.
Exceptions
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()

def pyfaust.Faust.conj (   F)

Returns the complex conjugate of F.

Parameters
Fthe Faust object.
Returns
a Faust object Fc implementing the conjugate of F.toarray() such that: Fc.toarray() == F.toarray().conj()

Examples

>>> from pyfaust import rand
>>> F = rand(50, 50, field='complex')
>>> Fc = F.conj()
See also
Faust.transpose, Faust.getH, Faust.H,

◆ conjugate()

def pyfaust.Faust.conjugate (   F)

Returns the complex conjugate of F.

Parameters
Fthe Faust object.
Returns
a Faust object Fc implementing the conjugate of F.toarray() such that: Fc.toarray() == F.toarray().conjugate()

Examples

>>> from pyfaust import rand
>>> F = rand(50, 50, field='complex')
>>> Fc = F.conjugate()
See also
Faust.transpose, Faust.getH, Faust.H,

◆ copy()

def pyfaust.Faust.copy (   F,
  dev = 'cpu' 
)

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

See also
Faust.clone,

◆ density()

def pyfaust.Faust.density (   F)

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.
A density above one is possible but prevents any saving.
Parameters
Fthe Faust object.
Returns
the density value (float).

Examples

>>> from pyfaust import rand
>>> F = rand(5, 50, density=.5)
>>> dens = F.density()
See also
Faust.nnz_sum, Faust.rcg, Faust.size, Faust.toarray

◆ display()

def pyfaust.Faust.display (   F)

Displays information about F.

Parameters
Fthe 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
See also
Faust.nnz_sum, Faust.density, Faust.shape, Faust.factors, Faust.numfactors, Faust.__repr__,

◆ dot()

def pyfaust.Faust.dot (   F,
  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.
See also
Faust.__init__, Faust.rcg, Faust.__mul__, Faust.__matmul__, Faust.dot

◆ factors()

def pyfaust.Faust.factors (   F,
  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.
Parameters
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.


Exceptions
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

See also
Faust.numfactors, Faust.transpose, Faust.left, Faust.right

◆ getH()

def pyfaust.Faust.getH (   F)

Returns the conjugate transpose of F.

Parameters
Fthe 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()
>>> bool((H1.toarray() == H2.toarray()).all())
True
See also
Faust.transpose, Faust.conj, Faust.H

◆ imshow()

def pyfaust.Faust.imshow (   F,
  name = 'F' 
)

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

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

def pyfaust.Faust.insert (   F,
  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.
Parameters
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.
Exceptions
ValueErrorin case new_factor has invalid dimensions.
ValueErrorif i is out of range.
TypeErrorif i is not an integer.

Examples

>>> 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
See also
Faust.factors, Faust.left, Faust.right, Faust.replace

◆ isdense()

def pyfaust.Faust.isdense (   F)

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

Examples

>>> 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
See also
Faust.issparse, pyfaust.rand, pyfaust.rand_bsr

◆ isFaust()

def pyfaust.Faust.isFaust (   obj)
static

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

◆ issparse()

def pyfaust.Faust.issparse (   F,
  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.

Parameters
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.

Examples

>>> 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
See also
Faust.isdense, pyfaust.rand, pyfaust.rand_bsr

◆ left()

def pyfaust.Faust.left (   F,
  i,
  as_faust = False 
)

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

Parameters
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
See also
Faust.factors, Faust.right

◆ load()

def pyfaust.Faust.load (   filepath)
static

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.

Parameters
filepath(str) the filepath of the .mat file.

Examples

>>> 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")
See also
Faust.__init__, Faust.save,

◆ load_native()

def pyfaust.Faust.load_native (   filepath)

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.

Parameters
filepath(str) the filepath of the .mat file.

Examples

>>> 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")
See also
Faust.__init__, Faust.save

◆ matvec()

def pyfaust.Faust.matvec (   F,
  x 
)

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

Examples

>>> 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])
See also
Faust.dot, Faust.__matmul__

◆ nnz_sum()

def pyfaust.Faust.nnz_sum (   F)

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.

Examples

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

def pyfaust.Faust.norm (   F,
  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.
Parameters
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’)`.

    Exceptions
    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
    See also
    Faust.normalize

◆ normalize()

def pyfaust.Faust.normalize (   F,
  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.

Parameters
ordthe norm order to use (see Faust.norm).
axisif 1 the columns are normalized, if 0 the rows.
Returns
the normalized Faust

Examples

>>> 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
>>> float(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
>>> float(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)
>>> float(norm(nFinf[2, :].toarray() - F[2, :].toarray()/F[2, :].norm(ord=np.inf)))
5.238750013840908e-17
See also
Faust.norm

◆ numfactors()

def pyfaust.Faust.numfactors (   F)

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
See also
Faust.factors, Faust.__len__

◆ optimize()

def pyfaust.Faust.optimize (   F,
  transp = False 
)

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

Parameters
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.

Examples

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) # doctest:+ELLIPSIS

F.toarray() time: ...

>>> # e.g: F.toarray() time: 0.2779221534729004
>>> t1 = time(); M_ = pF.toarray(); print("pF.toarray() time:", time()-t1) # doctest:+ELLIPSIS

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) # doctest:+ELLIPSIS

F@v time: ...

>>> # e.g: F@v time: 0.0016832351684570312
>>> t1 = time(); Fv_ = pF@v; print("pF@v time:", time()-t1) # doctest:+ELLIPSIS

pF@v time: ...

>>> # e.g: pF@v time: 0.0002257823944091797
>>> np.allclose(Fv_, Fv)
True
See also
Faust.optimize_time, Faust.optimize_memory, Faust.pruneout, Faust.nbytes, Jupyter notebook about Faust, optimizations

◆ optimize_memory()

def pyfaust.Faust.optimize_memory (   F)

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

Returns
The optimized Faust.

Examples

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

def pyfaust.Faust.optimize_time (   F,
  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.

Parameters
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.

Examples

>>> 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) # doctest:+ELLIPSIS

F.toarray() time:...

>>> # F.toarray() time: 0.2891380786895752
>>> t1 = time(); M_ = oF.toarray(); print("oF.toarray() time:", time()-t1) # doctest:+ELLIPSIS

oF.toarray() time:...

>>> # example: oF.toarray() 0.0172119140625
>>> np.allclose(M, M_)
True
See also
Faust.optimize, pyfaust.FaustMulMode

◆ pinv()

def pyfaust.Faust.pinv (   F)

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

Warning
this function makes a call to Faust.toarray().
Returns
The dense pseudo-inverse matrix.

Examples

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

def pyfaust.Faust.power_iteration (   self,
  threshold = 1e-3,
  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.

Parameters
thresholdthe precision required on the eigenvalue.
maxiterthe 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
>>> float(F.power_iteration())
12653.806783532553

◆ pruneout()

def pyfaust.Faust.pruneout (   F,
  thres = 0,
**  kwargs 
)

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

Parameters
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.

Examples

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

def pyfaust.Faust.rcg (   F)

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()
Parameters
Fthe Faust object.
Returns
the RCG value (float).

Examples

>>> from pyfaust import rand
>>> F = rand(1024, 1024)
>>> float(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
>>> float(F.size / F.nnz_sum())
40.96
See also
Faust.density, Faust.nnz_sum, Faust.shape.

◆ replace()

def pyfaust.Faust.replace (   F,
  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.
Parameters
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.
Exceptions
ValueErrorin case new_factor has invalid dimensions.
ValueErrorif i is out of range.
TypeErrorif i is not an integer.

Examples

>>> 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
See also
Faust.factors, Faust.left, Faust.right, Faust.insert

◆ right()

def pyfaust.Faust.right (   F,
  i,
  as_faust = False 
)

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

Parameters
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
See also
Faust.factors, Faust.left, Faust.numfactors,

◆ save()

def pyfaust.Faust.save (   F,
  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.
Parameters
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).
Exceptions
ValueError.

Examples

>>> from pyfaust import rand, Faust
>>> F = rand(5, 10, field='complex')
>>> F.save("F.mat")
>>> G = Faust(filepath="F.mat")
See also
Faust.__init__, Faust.rcg, Faust.load, Faust.load_native

◆ sum()

def pyfaust.Faust.sum (   F,
  axis = None,
**  kwargs 
)

Sums Faust elements over a given axis.

Parameters
axis(None or int or tuple of ints) Axis or axes along which the the sum is performed
Returns
The Faust sum.

Examples

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

def pyfaust.Faust.swap_cols (   F,
  id1,
  id2,
  permutation = False,
  inplace = False 
)

Swaps F columns of indices id1 and id2.

Parameters
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.

Examples

>>> 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
See also
Faust.swap_rows

◆ swap_rows()

def pyfaust.Faust.swap_rows (   F,
  id1,
  id2,
  permutation = True,
  inplace = False 
)

Swaps F rows of indices id1 and id2.

Parameters
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.

Examples

>>> 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
See also
Faust.swap_cols

◆ toarray()

def pyfaust.Faust.toarray (   F)

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

def pyfaust.Faust.todense (   F)

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.
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.
Exceptions
MemoryError
Warning
running the example below is likely to raise a memory error or freeze your computer for a certain amount of time.
this function is deprecated and might be deleted in future versions of pyfaust. Please use Faust.toarray instead.
See also
Faust.toarray

◆ transpose()

def pyfaust.Faust.transpose (   F)

Returns the transpose of F.

Parameters
Fthe Faust object.
Returns
a Faust object implementing the transpose of F.toarray(), such that: F.transpose().toarray() == F.toarray().transpose()

Examples

>>> from pyfaust import rand, seed
>>> F = rand(10, 18)
>>> tF = F.transpose()
>>> tF.shape
(18, 10)
See also
Faust.conj, Faust.getH, Faust.H, Faust.T

Property Documentation

◆ device

pyfaust.Faust.device = property
static

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

Examples

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

◆ dtype

pyfaust.Faust.dtype = property
static

Returns the dtype of the Faust.

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

Parameters
Fthe 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')

◆ H

pyfaust.Faust.H = property
static

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()
>>> bool((H1.toarray() == H2.toarray()).all())
True
See also
Faust.transpose, Faust.conj, Faust.getH

◆ imag

pyfaust.Faust.imag = property
static

Returns the imaginary part of F as a Faust.

◆ nbytes

pyfaust.Faust.nbytes = property
static

Gives the memory size of the Faust in bytes.

Examples

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

◆ ndim

pyfaust.Faust.ndim = property
static

Number of Faust dimensions (always 2).

◆ real

pyfaust.Faust.real = property
static

Returns the real part of F as a Faust.

◆ shape

pyfaust.Faust.shape = property
static

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

Parameters
Fthe 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]
See also
Faust.nnz_sum, Faust.size,

◆ size

pyfaust.Faust.size = property
static

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

Parameters
Fthe Faust object.
Returns
The number of elements in the Faust F.

Examples

>>> from pyfaust import rand
>>> F = rand(50, 50, field='complex')
>>> int(F.size)
2500
See also
Faust.shape

◆ T

pyfaust.Faust.T = property
static

Returns the transpose of F.

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

Parameters
Fthe Faust object.
Returns
a Faust object implementing the transpose of F.toarray(), such that: F.T.toarray() == F.toarray().T

Examples

>>> from pyfaust import rand
>>> F = rand(10, 23)
>>> tF = F.T
>>> tF.shape
(23, 10)
See also
Faust.conj, Faust.getH, Faust.H, Faust.T

The documentation for this class was generated from the following file:
pyfaust.norm
def norm(F, ord='fro', **kwargs)
Returns Faust.norm(F, ord)` ornumpy.linalg.norm(F, ord)`` depending of F type.
Definition: __init__.py:3917
pyfaust.seed
def seed(s)
(Re)Initializes the pyfaust pseudo-random generator.
Definition: __init__.py:5334
pyfaust.rand
def rand(num_rows, num_cols, num_factors=None, dim_sizes=None, density=None, fac_type='sparse', per_row=True, dev='cpu', dtype='float64', field=None, seed=0)
Generates a random Faust.
Definition: __init__.py:5016
pyfaust.pinv
def pinv(F)
A package function alias for the member function Faust.pinv().
Definition: __init__.py:3949
pyfaust.is_gpu_mod_enabled
def is_gpu_mod_enabled()
Returns True if the gpu_mod plug-in has been loaded correctly, False otherwise.
Definition: __init__.py:5251