-
3.41.0
|
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... | |
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:
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.
In this documentation, the expression 'full matrix' designates the array Faust.toarray() obtained by the multiplication of the Faust factors.
List of functions that are memory costly:
F[i, j]
/ Faust.__getitem__, but note that slicing is memory efficient through memory views),For more information about FAuST take a look at https://faust.inria.fr.
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.
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. |
Examples
Faust size 10x10, density 1, nnz_sum 100, 1 factor(s):
Faust size 10x9, density 0.2, nnz_sum 18, 1 factor(s):
Faust size 10x18, density 1, nnz_sum 180, 2 factor(s):
def pyfaust.Faust.__add__ | ( | F, | |
* | args, | ||
** | kwargs | ||
) |
Sums F to one or a sequence of variables (Faust objects, arrays or scalars).
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. |
Examples
Faust size 10x12, density 2.04167, nnz_sum 245, 5 factor(s):
Faust size 10x12, density 2.025, nnz_sum 243, 5 factor(s):
Faust size 10x12, density 4.43333, nnz_sum 532, 7 factor(s):
Faust size 10x12, density 2.84167, nnz_sum 341, 7 factor(s):
Faust size 10x12, density 5.4, nnz_sum 648, 9 factor(s):
Faust size 10x12, density 11.05, nnz_sum 1326, 13 factor(s):
def pyfaust.Faust.__floordiv__ | ( | F, | |
s | |||
) |
F // s.
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.
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. |
IndexError |
Examples
Faust size 50x1, density 24.64, nnz_sum 1232, 5 factor(s):
Faust size 2x3, density 174.833, nnz_sum 1049, 5 factor(s):
Faust size 50x95, density 0.318947, nnz_sum 1515, 5 factor(s):
Faust size 39x100, density 0.381538, nnz_sum 1488, 5 factor(s):
Faust size 48x3, density 8.56944, nnz_sum 1234, 5 factor(s):
Faust size 15x100, density 0.928, nnz_sum 1392, 5 factor(s):
Faust size 2x100, density 6.7, nnz_sum 1340, 5 factor(s):
Faust size 15x100, density 0.928, nnz_sum 1392, 5 factor(s):
Faust size 3x100, density 4.48, nnz_sum 1344, 5 factor(s):
Faust size 50x3, density 8.27333, nnz_sum 1241, 5 factor(s):
Faust size 3x2, density 175, nnz_sum 1050, 5 factor(s):
def pyfaust.Faust.__iadd__ | ( | F, | |
A | |||
) |
Inplace F = F + A.
def pyfaust.Faust.__imatmul__ | ( | F, | |
A | |||
) |
Inplace F = F @ A.
def pyfaust.Faust.__imul__ | ( | F, | |
A | |||
) |
Inplace F = F * A.
def pyfaust.Faust.__isub__ | ( | F, | |
A | |||
) |
Inplace F = F - A.
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).
def pyfaust.Faust.__len__ | ( | F | ) |
Returns the number of factors of F.
Examples
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:
(F @ A).toarray() == F.toarray() @ A.toarray()
(an elementwise non-significant absolute difference between the two members is possible).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'). |
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.
ValueError |
def pyfaust.Faust.__mul__ | ( | F, | |
A | |||
) |
Multiplies the Faust F by A.
This method overloads a Python function/operator (*).
More specifically:
F | the Faust object. |
A | a scalar number or a (vector) numpy.ndarray or a Faust elementwise multiplication. |
Returns | The result of the multiplication as a Faust (either A is a vector or a scalar). |
Raises | TypeError |
Examples
Compute elementwise F * F, gives an array:
If the type of A is not supported:
def pyfaust.Faust.__neg__ | ( | F | ) |
Returns F (unary operator).
def pyfaust.Faust.__pos__ | ( | F | ) |
Returns + F (unary operator).
def pyfaust.Faust.__radd__ | ( | F, | |
lhs_op | |||
) |
Returns lhs_op+F.
def pyfaust.Faust.__repr__ | ( | F | ) |
Returns a str object representing the Faust object.
This method overloads a Python function.
F | the Faust object. |
Examples
Faust size 50x50, density 0.5, nnz_sum 1250, 5 factor(s):
Faust size 50x50, density 0.5, nnz_sum 1250, 5 factor(s):
def pyfaust.Faust.__rmatmul__ | ( | F, | |
lhs_op | |||
) |
Returns lhs_op.__matmul__(F).
Examples
def pyfaust.Faust.__rmul__ | ( | F, | |
lhs_op | |||
) |
lhs_op*F
def pyfaust.Faust.__rsub__ | ( | F, | |
lhs_op | |||
) |
Returns lhs_op-F.
def pyfaust.Faust.__str__ | ( | F | ) |
Converts F to its str representation.
def pyfaust.Faust.__sub__ | ( | F, | |
* | args | ||
) |
Subtracts from F one or a sequence of variables.
Faust objects, arrays or scalars.
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. |
Examples
Faust size 10x12, density 2.04167, nnz_sum 245, 5 factor(s):
Faust size 10x12, density 2.025, nnz_sum 243, 5 factor(s):
Faust size 10x12, density 4.43333, nnz_sum 532, 7 factor(s):
Faust size 10x12, density 2.84167, nnz_sum 341, 7 factor(s):
Faust size 10x12, density 5.4, nnz_sum 648, 9 factor(s):
Faust size 10x12, density 11.05, nnz_sum 1326, 13 factor(s):
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).
def pyfaust.Faust.astype | ( | F, | |
dtype | |||
) |
Converts F to the dtype passed as argument in a new Faust.
dtype | (str) 'float32', 'float64' or 'complex'. |
Examples
Faust size 10x10, density 2.5, nnz_sum 250, 5 factor(s):
Faust size 10x10, density 2.5, nnz_sum 250, 5 factor(s):
def pyfaust.Faust.average | ( | F, | |
axis = None , |
|||
weights = None , |
|||
sw_returned = False |
|||
) |
Computes the weighted average of F along the specified axis.
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: |
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). |
Examples
Faust size 1x1, density 9, nnz_sum 9, 3 factor(s):
def pyfaust.Faust.clone | ( | F, | |
dev = None |
|||
) |
Clones the Faust (in a new memory space).
dev | (str) 'cpu' to clone on CPU RAM, 'gpu' to clone on the GPU device. By default (None), the device is the F.device. |
Examples
Faust size 10x10, density 2.5, nnz_sum 250, 5 factor(s):
def pyfaust.Faust.concatenate | ( | F, | |
* | args, | ||
** | kwargs | ||
) |
Concatenates F with len(args) Faust objects, numpy arrays or scipy sparse matrices.
The resulting Faust:
verifies that:
N.B.: you could have an elementwise non-significant absolute difference between the two members.
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). |
ValueError |
Examples
Faust size 100x50, density 0.52, nnz_sum 2600, 6 factor(s):
Faust size 50x100, density 0.52, nnz_sum 2600, 6 factor(s):
Faust size 84x50, density 0.773809, nnz_sum 3250, 6 factor(s):
Faust size 50x74, density 0.422162, nnz_sum 1562, 6 factor(s):
Faust size 384x50, density 0.575521, nnz_sum 11050, 6 factor(s):
def pyfaust.Faust.conj | ( | F | ) |
Returns the complex conjugate of F.
F | the Faust object. |
Fc.toarray() == F.toarray().conj()
Examples
def pyfaust.Faust.conjugate | ( | F | ) |
Returns the complex conjugate of F.
F | the Faust object. |
Fc.toarray() == F.toarray().conjugate()
Examples
def pyfaust.Faust.copy | ( | F, | |
dev = 'cpu' |
|||
) |
Clone alias function (here just to mimic numpy API).
def pyfaust.Faust.density | ( | F | ) |
Calculates the density of F such that F.nnz_sum() == F.density() * F.size
.
F | the Faust object. |
Examples
def pyfaust.Faust.display | ( | F | ) |
Displays information about F.
F | the Faust object. |
Examples
Faust size 50x100, density 1.3, nnz_sum 6500, 2 factor(s):
Faust size 50x100, density 1.3, nnz_sum 6500, 2 factor(s):
Faust size 50x100, density 1.3, nnz_sum 6500, 2 factor(s):
def pyfaust.Faust.dot | ( | F, | |
A, | |||
* | args, | ||
** | kwargs | ||
) |
Performs equivalent operation of numpy.dot() between the Faust F and A.
More specifically:
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.
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. |
ValueError. |
def pyfaust.Faust.getH | ( | F | ) |
Returns the conjugate transpose of F.
F | the Faust object. |
H.toarray() == F.toarray().getH()
Examples
def pyfaust.Faust.imshow | ( | F, | |
name = 'F' |
|||
) |
Displays image of F's full matrix and its factors.
Examples
def pyfaust.Faust.insert | ( | F, | |
i, | |||
new_factor | |||
) |
Inserts new_factor at index i in a new Faust copy of F.
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. |
ValueError | in case new_factor has invalid dimensions. |
ValueError | if i is out of range. |
TypeError | if i is not an integer. |
Examples
def pyfaust.Faust.isdense | ( | F | ) |
Returns True if all factors are dense arrays (as np.ndarray-s) False otherwise.
Examples
Faust size 10x10, density 2.5, nnz_sum 250, 5 factor(s):
Faust size 10x10, density 2.5, nnz_sum 250, 5 factor(s):
Faust size 10x10, density 2.5, nnz_sum 250, 5 factor(s):
|
static |
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.
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
Faust size 10x10, density 2.5, nnz_sum 250, 5 factor(s):
Faust size 10x10, density 0.24, nnz_sum 24, 2 factor(s):
def pyfaust.Faust.left | ( | F, | |
i, | |||
as_faust = False |
|||
) |
Returns the left hand side factors of F from index 0 to i included.
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). |
Examples
Faust size 5x10, density 2.98, nnz_sum 149, 5 factor(s):
Faust size 5x6, density 3.96667, nnz_sum 119, 4 factor(s):
|
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.
filepath | (str) the filepath of the .mat file. |
Examples
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.
filepath | (str) the filepath of the .mat file. |
Examples
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
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.
Examples
Faust size 1024x1024, density 0.0244141, nnz_sum 25600, 5 factor(s):
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().
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). |
norm(F.toarray(), 1) == max(sum(abs(F.toarray())))
.norm(F.toarray(), 2) == max(scipy.linalg.svd(F.toarray())[1])
. This is the default norm calculated when calling to norm(F).norm(F.toarray(), numpy.inf) == max(sum(abs(F.T.toarray())))
If ord == 'fro', the norm is ‘norm(F.toarray(), 'fro’)`.
ValueError. |
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.
ord | the norm order to use (see Faust.norm). |
axis | if 1 the columns are normalized, if 0 the rows. |
Examples
Faust size 10x10, density 2.5, nnz_sum 250, 5 factor(s):
Faust size 10x10, density 2.5, nnz_sum 250, 5 factor(s):
Faust size 10x10, density 2.5, nnz_sum 250, 5 factor(s):
def pyfaust.Faust.numfactors | ( | F | ) |
Returns the number of factors of F.
Examples
def pyfaust.Faust.optimize | ( | F, | |
transp = False |
|||
) |
Returns a Faust, optimized with Faust.pruneout, Faust.optimize_memory, and Faust.optimize_time.
transp | (bool) True in order to optimize the Faust according to its transpose. |
Examples
This example shows how Faust.optimize can diminish the memory size and speed up the Faust.toarray() and Faust-vector multiplication.
F.toarray() time: ...
pF.toarray() time: ...
F@v time: ...
pF@v time: ...
def pyfaust.Faust.optimize_memory | ( | F | ) |
Optimizes a Faust by changing the storage format of each factor in order to optimize the memory size.
Examples
Faust size 1024x1024, density 0.0244141, nnz_sum 25600, 5 factor(s):
Faust size 1024x1024, density 0.0244141, nnz_sum 25600, 5 factor(s):
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).
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.
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. |
Examples
F.toarray() time:...
oF.toarray() time:...
def pyfaust.Faust.pinv | ( | F | ) |
Computes the (Moore-Penrose) pseudo-inverse of F.toarray().
Examples
See also numpy.linalg.pinv, pyfaust.fact.pinvtj
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.
threshold | the precision required on the eigenvalue. |
maxiter | the number of iterations above what the algorithm will stop anyway. |
Examples
def pyfaust.Faust.pruneout | ( | F, | |
thres = 0 , |
|||
** | kwargs | ||
) |
Returns a Faust optimized by removing useless zero rows and columns as many times as needed.
F | (Faust) the Faust to optimize. |
thres | (int) the threshold in number of nonzeros under what the rows/columns are removed. |
Examples
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.
F | the Faust object. |
Examples
Faust size 1024x1024, density 0.0244141, nnz_sum 25600, 5 factor(s):
def pyfaust.Faust.replace | ( | F, | |
i, | |||
new_factor | |||
) |
Replaces the factor of index i by new_factor in a new Faust copy of F.
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. |
ValueError | in case new_factor has invalid dimensions. |
ValueError | if i is out of range. |
TypeError | if i is not an integer. |
Examples
def pyfaust.Faust.right | ( | F, | |
i, | |||
as_faust = False |
|||
) |
Returns the right hand side factors of F from index i to F.numfactors()-1.
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). |
Examples
Faust size 5x10, density 2.98, nnz_sum 149, 5 factor(s):
Faust size 6x10, density 1.46667, nnz_sum 88, 3 factor(s):
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).
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). |
ValueError. |
Examples
def pyfaust.Faust.sum | ( | F, | |
axis = None , |
|||
** | kwargs | ||
) |
Sums Faust elements over a given axis.
axis | (None or int or tuple of ints) Axis or axes along which the the sum is performed |
Examples
Faust size 1x1, density 270, nnz_sum 270, 7 factor(s):
Faust size 10x1, density 26, nnz_sum 260, 6 factor(s):
def pyfaust.Faust.swap_cols | ( | F, | |
id1, | |||
id2, | |||
permutation = False , |
|||
inplace = False |
|||
) |
Swaps F columns of indices id1 and id2.
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. |
Examples
Faust size 10x10, density 2.5, nnz_sum 250, 5 factor(s):
def pyfaust.Faust.swap_rows | ( | F, | |
id1, | |||
id2, | |||
permutation = True , |
|||
inplace = False |
|||
) |
Swaps F rows of indices id1 and id2.
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. |
Examples
Faust size 10x10, density 2.6, nnz_sum 260, 6 factor(s):
Faust size 10x10, density 2.6, nnz_sum 260, 6 factor(s):
def pyfaust.Faust.toarray | ( | F | ) |
Returns a numpy array for the full matrix implemented by F.
MemoryError |
Examples
Faust size 100000x100000, density 0.00018, nnz_sum 1800000, 2 factor(s):
def pyfaust.Faust.todense | ( | F | ) |
Returns a numpy matrix for the full matrix implemented by F.
MemoryError |
def pyfaust.Faust.transpose | ( | F | ) |
Returns the transpose of F.
F | the Faust object. |
F.transpose().toarray() == F.toarray().transpose()
Examples
|
static |
Returns the device on which the Faust is located ('cpu' or 'gpu').
Examples
|
static |
|
static |
Returns the conjugate transpose of F.
This function is intended to be used as a property (see the examples).
H.toarray() == F.toarray().H
Examples
|
static |
Returns the imaginary part of F as a Faust.
|
static |
Gives the memory size of the Faust in bytes.
Examples
Faust size 1024x1024, density 0.0244141, nnz_sum 25600, 5 factor(s):
Faust size 1024x1024, density 0.0244141, nnz_sum 25600, 5 factor(s):
|
static |
Number of Faust dimensions (always 2).
|
static |
Returns the real part of F as a Faust.
|
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().
F | the Faust object. |
Examples
|
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).
F | the Faust object. |
Examples
|
static |
Returns the transpose of F.
This function is intended to be used as a property (see the examples).
F | the Faust object. |
F.T.toarray() == F.toarray().T
Examples