-  3.39.12
Namespaces | Classes | Functions
matfaust Namespace Reference

The FAuST Matlab Wrapper More...

Namespaces

 demo
 The matfaust demo namespace.
 
 fact
 The matfaust factorization module.
 
 factparams
 The module for the parametrization of FAuST's algorithms (Palm4MSA and Hierarchical Factorization)
 
 poly
 The matfaust module for polynomial basis as Faust objects.
 
 proj
 This module provides matrix projectors.
 
 tools
 The matfaust tools namespace.
 

Classes

class  Faust
 FAuST Matlab wrapper main class for using multi-layer sparse transforms. More...
 
class  FaustMulMode
 Enumeration class of all matrix chain multiplication methods available to multiply a Faust to a matrix or to compute Faust.full(). These methods are used by Faust.optimize_time(). More...
 

Functions

function anticirc (c, varargin)
 Returns an anticirculant Faust A defined by the vector c (which is the last column of full(A)). More...
 
function bitrev_perm (n)
 Bitreversal permutation. More...
 
function circ (c, varargin)
 Returns a circulant Faust C defined by the vector c (which is the first column of full(C)). More...
 
function create_bsr (M, N, bnnz, bdata, bcolinds, brow_count)
 This function is an helper to create a FAµST BSR matrix (kind of matrix that is not natively supported by Matlab). More...
 
function dct (n, varargin)
 Constructs a Faust implementing the Direct Cosine Transform (Type II) Faust of order n. More...
 
function dft (n, varargin)
 Constructs a Faust F implementing the Discrete Fourier Transform (DFT) of order n. More...
 
function dst (n, varargin)
 Constructs a Faust implementing the Direct Sine Transform (Type II) Faust of order n. More...
 
function enable_gpu_mod (varargin)
 This function loads explicitely the gpu_mod library in memory. More...
 
function eye (varargin)
 Identity Faust. More...
 
function faust_fact (varargin)
 This function is a shorthand for matfaust.fact.hierarchical(). More...
 
function is_gpu_mod_enabled ()
 Returns true if the gpu_mod plug-in has been loaded correctly, false otherwise. More...
 
function is_gpu_mod_working ()
 This function returns true if the gpu_mod is working properly false otherwise. More...
 
function isFaust (obj)
 Package alias of Faust.isFaust. More...
 
function license ()
 Prints the FAuST license. More...
 
function opt_butterfly_faust (F)
 See matfaust.Faust.opt_butterfly alias. More...
 
function rand (M, N, varargin)
 Generates a random Faust. More...
 
function rand_bsr (M, N, BM, BN, varargin)
 Generates a random Faust composed only of BSR matrices. More...
 
function rand_butterfly (n, varargin)
 Constructs a Faust corresponding to the product of log2(n) square factors of size n with butterfly supports and random nonzero coefficients. More...
 
function toeplitz (c, varargin)
 Constructs a toeplitz Faust whose first column is c and first row r. More...
 
function version ()
 Returns the FAuST package version. More...
 
function wht (n, varargin)
 Constructs a Faust implementing the Walsh-Hadamard Transform (WHT) of order n. More...
 

Detailed Description

The FAuST Matlab Wrapper

Function Documentation

◆ anticirc()

function matfaust::anticirc ( ,
varargin   
)

Returns an anticirculant Faust A defined by the vector c (which is the last column of full(A)).

Parameters
cthe vector to define the anticirculant Faust.
'dev',str'gpu' or 'cpu' to create the Faust on CPU or GPU ('cpu' is the default).
'diag_opt',logicalcf. matfaust.circ.

Example

% in a matlab terminal
>> import matfaust.anticirc
>> c = 1:8;
>> A = anticirc(c)

A =

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

  • FACTOR 0 (complex) SPARSE, size 8x8, density 0.25, nnz 16
  • FACTOR 1 (complex) SPARSE, size 8x8, density 0.25, nnz 16
  • FACTOR 2 (complex) SPARSE, size 8x8, density 0.25, nnz 16
  • FACTOR 3 (complex) SPARSE, size 8x8, density 0.25, nnz 16
  • FACTOR 4 (complex) SPARSE, size 8x8, density 0.25, nnz 16
  • FACTOR 5 (complex) SPARSE, size 8x8, density 0.25, nnz 16
>> full_A = full(A);
>> all(full_A(:, end).' - c < 1e-15)
ans =
logical
1
>> real(full_A)
ans =
2.0000 3.0000 4.0000 5.0000 6.0000 7.0000 8.0000 1.0000
3.0000 4.0000 5.0000 6.0000 7.0000 8.0000 1.0000 2.0000
4.0000 5.0000 6.0000 7.0000 8.0000 1.0000 2.0000 3.0000
5.0000 6.0000 7.0000 8.0000 1.0000 2.0000 3.0000 4.0000
6.0000 7.0000 8.0000 1.0000 2.0000 3.0000 4.0000 5.0000
7.0000 8.0000 1.0000 2.0000 3.0000 4.0000 5.0000 6.0000
8.0000 1.0000 2.0000 3.0000 4.0000 5.0000 6.0000 7.0000
1.0000 2.0000 3.0000 4.0000 5.0000 6.0000 7.0000 8.0000
>> % Look at the density of a larger anticirculant Faust
>> % it indicates a speedup of the Faust-matrix/vector product
>> density(anticirc(rand(1, 1024)))
0.0391

See also matfaust.circ, matfaust.toeplitz

◆ bitrev_perm()

function matfaust::bitrev_perm ( )

Bitreversal permutation.

Parameters
nthe size of the permutation, it must be a power of two. P dimensions will be n x n.
Return values
Pa sparse matrix defining the bit-reversal permutation.

See also matfaust.dft

◆ circ()

function matfaust::circ ( ,
varargin   
)

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

Parameters
cthe vector to define the circulant Faust.
'dev',str'gpu' or 'cpu' to create the Faust on CPU or GPU ('cpu' is the default).
'diag_opt',logicalif true then the returned Faust is optimized using matfaust.opt_butterfly_faust (because the DFT is used to implement circ).

Example:

>> import matfaust.circ
>> c = 1:8;
>> C = circ(c)

C =

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

  • FACTOR 0 (complex) SPARSE, size 8x8, density 0.25, nnz 16
  • FACTOR 1 (complex) SPARSE, size 8x8, density 0.25, nnz 16
  • FACTOR 2 (complex) SPARSE, size 8x8, density 0.25, nnz 16
  • FACTOR 3 (complex) SPARSE, size 8x8, density 0.25, nnz 16
  • FACTOR 4 (complex) SPARSE, size 8x8, density 0.25, nnz 16
  • FACTOR 5 (complex) SPARSE, size 8x8, density 0.25, nnz 16
>> full_C = full(C);
>> all(full_C(:, 1).' - c < 1e-15)
ans =
logical
1
>> real(full_C)
ans =
1.0000 8.0000 7.0000 6.0000 5.0000 4.0000 3.0000 2.0000
2.0000 1.0000 8.0000 7.0000 6.0000 5.0000 4.0000 3.0000
3.0000 2.0000 1.0000 8.0000 7.0000 6.0000 5.0000 4.0000
4.0000 3.0000 2.0000 1.0000 8.0000 7.0000 6.0000 5.0000
5.0000 4.0000 3.0000 2.0000 1.0000 8.0000 7.0000 6.0000
6.0000 5.0000 4.0000 3.0000 2.0000 1.0000 8.0000 7.0000
7.0000 6.0000 5.0000 4.0000 3.0000 2.0000 1.0000 8.0000
8.0000 7.0000 6.0000 5.0000 4.0000 3.0000 2.0000 1.0000
>> % Look at the density of a larger circulant Faust
>> % it indicates a speedup of the Faust-matrix/vector product
>> density(circ(rand(1, 1024)))
0.0391

See also matfaust.anticirc, matfaust.toeplitz

◆ create_bsr()

function matfaust::create_bsr ( ,
,
bnnz  ,
bdata  ,
bcolinds  ,
brow_count   
)

This function is an helper to create a FAµST BSR matrix (kind of matrix that is not natively supported by Matlab).

In fact it doesn't directly create a BSR matrix but rather an encoding to create the matrix later through the matfaust.Faust constructor.

Parameters
MThe number of rows of the BSR matrix.
NThe number of columns of the BSR matrix.
bnnzThe number of nonzero blocks of the BSR matrix.
bdataThe horizontal concatenation of the matrix nonzero blocks taken in rowmajor order. Note that all the nonzero blocks have the same size.
bcolindsThe block column indices corresponding (in the same order) to the blocks in bdata. A block column index is not a column index, for example if your matrix contains N columns and each nonzero block contains BN columns, the block column indices must lay in {1, …, N/BN}.
brow_countThe number of nonzero blocks in each block-row of the matrix (a vector of size M/BM, where brow_count(i) is the number of nonzero blocks in the i-th block-row).

Examples

>> bnnz = 3;
>> M = 10;
>> N = 9;
>> nonzero_blocks = [[ 0.6787, 0.7431 0.6555; 0.7577, 0.3922 0.1712] [ 0.7060 0.2769 0.0971; 0.0318 0.0462 0.8235] [0.6948 0.9502 0.4387; 0.3171 0.0344 0.3816]];
>> bsr_mat = matfaust.create_bsr(M, N, bnnz, nonzero_blocks, [2, 3, 1], [2, 0, 1, 0, 0]);
>> F = matfaust.Faust(bsr_mat)
F =
Faust size 10x10, 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)
>> full(F)
ans =
0 0 0 0.6787 0.7431 0.6555 0.7060 0.2769 0.0971
0 0 0 0.7577 0.3922 0.1712 0.0318 0.0462 0.8235
0 0 0 0 0 0 0 0 0
0 0 0 0 0 0 0 0 0
0.6948 0.9502 0.4387 0 0 0 0 0 0
0.3171 0.0344 0.3816 0 0 0 0 0 0
0 0 0 0 0 0 0 0 0
0 0 0 0 0 0 0 0 0
0 0 0 0 0 0 0 0 0
0 0 0 0 0 0 0 0 0

◆ dct()

function matfaust::dct ( ,
varargin   
)

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

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

Parameters
nthe order of the DCT (it must be a power of two).
'dev',str'gpu' or 'cpu' to create the Faust on CPU or GPU ('cpu' is the default).
'normed',booltrue (by default) to normalize the returned Faust as if Faust.normalize() was called, false otherwise.
'class',str'single' or 'double'.
Return values
Dthe DCT Faust.

Example

% in a matlab terminal
>> import matfaust.*
>> F = dct(8);
>> x = ones(8, 1);
>> % apply the DCT to x
>> F*x
ans =
10.2517
0.0000
3.5999
0.0000
2.4054
0.0000
2.0392
0.0000
>> % check the density with a larger DCT Faust of size 1024
>> dct(1024).density()
ans =
0.0610
>> % it is smaller than 1

See also matfaust.dft, matfaust.dst, matfaust.Faust.density

◆ dft()

function matfaust::dft ( ,
varargin   
)

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

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

Usage

    F = dft(n)
    F = dft(n, normed)

Parameters
nthe power of two for a FFT of order n and a factorization in log2(n)+1 factors.
'normed',booltrue (by default) to normalize the returned Faust as if Faust.normalize() was called, false otherwise.
'dev',str'gpu or 'cpu' to create the Faust on CPU or GPU ('cpu' by default).
'diag_opt',boolif true then the returned Faust is optimized using matfaust.opt_butterfly_faust.
Return values
Fthe Faust implementing the FFT transform of dimension n.

Example

% in a matlab terminal
>> import matfaust.*
>> F = dft(1024)

F =

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

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

See also: bitrev_perm, matfaust.wht, matfaust.dct, matfaust.dst, fft, matfaust.rand_butterfly, matfaust.fact.butterfly

◆ dst()

function matfaust::dst ( ,
varargin   
)

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

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

Parameters
nthe order of the DST (it must be a power of two).
'dev',str'gpu' or 'cpu' to create the Faust on CPU or GPU ('cpu' is the default).
'normed',booltrue (by default) to normalize the returned Faust as if Faust.normalize() was called, false otherwise.
'class',str'single' or 'double'.
Return values
Dthe DST Faust.

Example

% in a matlab terminal
>> import matfaust.*
>> F = dst(8);
>> x = 1:8;
>> % apply the DST to x
>> F*x'
ans =
72.0000
-25.7693
0
-2.6938
0
-0.8036
0
-0.2028
>> % check the density with a larger DST Faust of size 1024
>> density(dst(1024))
ans =
0.1581
>> % it is smaller than 1

See also matfaust.dft, matfaust.dct, matfaust.Faust.density

◆ enable_gpu_mod()

function matfaust::enable_gpu_mod ( varargin  )

This function loads explicitely the gpu_mod library in memory.

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

Parameters
libpaththe absolute or relative path where to find the dynamic library (gm) to load. By default, it's none to auto-find the library (if possible).
backendthe GPU backend to use, only 'cuda' is available for now.
silentif True nothing or almost will be displayed on loading, otherwise all messages are visible.

◆ eye()

function matfaust::eye ( varargin  )

Identity Faust.

Usage

    eye(m,n) or eye([m,n]) forms a M-by-N Faust F = Faust(speye(M,N)).
    eye(m) is a short for eye(m,m).
    eye(S, 'complex') with S the size, does the same as above but returns a complex Faust.
    eye(S, 'complex', 'dev', 'gpu') or eye(S, 'dev', 'gpu') same as above but creates the Faust on GPU.

Parameters
'dev',str'gpu or 'cpu' to create the Faust on CPU or GPU (by default on CPU).
'class',str'double' (by default) or 'single' to select the scalar type used for the Faust generated.

Example

% in a matlab terminal
>> matfaust.eye(4)
ans =
Faust size 4x4, density 0.25, nnz_sum 4, 1 factor(s):
- FACTOR 0 (real) SPARSE, size 4x4, density 0.25, nnz 4
identity matrix flag
>> full(matfaust.eye(4))
ans =
1 0 0 0
0 1 0 0
0 0 1 0
0 0 0 1
>> full(matfaust.eye(4,5))
ans =
1 0 0 0 0
0 1 0 0 0
0 0 1 0 0
0 0 0 1 0
>> full(matfaust.eye([5,4]))
ans =
1 0 0 0
0 1 0 0
0 0 1 0
0 0 0 1
0 0 0 0
>> full(matfaust.eye([5,4],'complex'))
ans =
1.0000 + 0.0000i 0.0000 + 0.0000i 0.0000 + 0.0000i 0.0000 + 0.0000i
0.0000 + 0.0000i 1.0000 + 0.0000i 0.0000 + 0.0000i 0.0000 + 0.0000i
0.0000 + 0.0000i 0.0000 + 0.0000i 1.0000 + 0.0000i 0.0000 + 0.0000i
0.0000 + 0.0000i 0.0000 + 0.0000i 0.0000 + 0.0000i 1.0000 + 0.0000i
0.0000 + 0.0000i 0.0000 + 0.0000i 0.0000 + 0.0000i 0.0000 + 0.0000i
>> full(matfaust.eye([4],'complex'))
ans =
1.0000 + 0.0000i 0.0000 + 0.0000i 0.0000 + 0.0000i 0.0000 + 0.0000i
0.0000 + 0.0000i 1.0000 + 0.0000i 0.0000 + 0.0000i 0.0000 + 0.0000i
0.0000 + 0.0000i 0.0000 + 0.0000i 1.0000 + 0.0000i 0.0000 + 0.0000i
0.0000 + 0.0000i 0.0000 + 0.0000i 0.0000 + 0.0000i 1.0000 + 0.0000i

◆ faust_fact()

function matfaust::faust_fact ( varargin  )

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

See also matfaust.fact.hierarchical

◆ is_gpu_mod_enabled()

function matfaust::is_gpu_mod_enabled ( )

Returns true if the gpu_mod plug-in has been loaded correctly, false otherwise.

See also matfaust.is_gpu_mod_working

◆ is_gpu_mod_working()

function matfaust::is_gpu_mod_working ( )

This function returns true if the gpu_mod is working properly false otherwise.

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

◆ isFaust()

function matfaust::isFaust ( obj  )

Package alias of Faust.isFaust.

See also Faust.Faust, Faust.isFaust

◆ license()

function matfaust::license ( )

Prints the FAuST license.

◆ opt_butterfly_faust()

function matfaust::opt_butterfly_faust ( )

◆ rand()

function matfaust::rand ( ,
,
varargin   
)

Generates a random Faust.

Warning
if this function is imported through 'import matfaust.rand' or 'import matfaust.*' the Matlab builtin function rand will be unreachable (matfaust.rand will be called instead). So, generally, it is not advisable to import this function, rather directly call matfaust.rand without any import.

Usage

    rand(M,N) generates a M-by-N Faust object. The numbers of rows and columns of intermediary factors are all randomly chosen between M and N (included). The number of factors is 5. The factors are sparse and real. The nnz per row of factors is 5.

    rand(M, N, 'num_factors', NF, 'dim_sizes', S) with NF and S two integers, generates a M-by-N Faust of NF factors. The factor size is S for both dimensions (except the first factor number of rows and the last factor number of columns which are respectively M and N). The factors are sparse and real.

    rand(M, N, 'num_factors', [N1, N2], 'dim_sizes', S) same as above except that here the number of factors is randomly chosen between N1 and N2 inclusively.

    rand(M, N, 'num_factors', [N1, N2], [S1, S2]) or rand(M, N, 'num_factors', NF, 'dim_sizes', [S1, S2]) same as above except that here the intermediary factor dimension sizes are random; the number of rows and columns are both randomly chosen between S1 and S2 inclusively.

    rand(M, N, 'num_factors', NF, 'dim_sizes', S, 'density', D) or rand(M, N, 'num_factors', [N1, N2], 'dim_sizes', [S1, S2], 'density', D) same as above but specifying D the approximate density of each factor.

    rand(M, N, 'num_factors', NF, 'dim_sizes', S, 'density', D, 'per_row', true) or rand(M, N, 'num_factors', [N1, N2], 'dim_sizes', [S1, S2], 'density',D, 'per_row', true) same as above but specifying D, the density of each factor per row ('per_row', true) or per column ('per_row', false).

    rand(M, N, 'num_factors', NF, 'dim_sizes', S, 'density', D, 'fac_type', 'dense') or rand(M, N, 'num_factors', [N1, N2], 'dim_sizes', [S1, S2], 'density', D, 'fac_type', 'dense') same as above but generating only dense matrices as factors.

    rand(M, N, 'num_factors', NF, 'dim_sizes', S, 'density', D, 'mixed') or rand(M, N, 'num_factors', [N1, N2], 'dim_sizes', [S1, S2], 'density', D, 'fac_type', 'sparse') same as above but generating either sparse or dense matrices as factors.

    rand(M, N, 'num_factors', NF, 'dim_sizes', S, 'density', D, 'fac_type', 'sparse', 'field' , 'complex'), rand(M, N, 'num_factors', [N1, N2], 'dim_sizes', [S1, S2], 'density', D, 'fac_type', 'sparse', 'per_row', false), rand(M, N, 'num_factors', NF, 'dim_sizes', S, 'density', D, 'fac_type', 'dense', 'field' , 'complex') or rand(M, N, 'num_factors', [N1, N2], 'dim_sizes', [S1, S2], D, 'fac_type', 'dense', 'field', 'complex') same as above but generating a complex Faust, that is, matrices defined over the complex field.

Parameters
Mthe Faust number of rows.
Nthe Faust number of columns.
'num_factors',NFif it's an integer it is the number of random factors to set in the Faust. If NF is a vector of 2 integers then the number of factors is set randomly between NF(1) and NF(2) (inclusively). Defaultly, a 5 factors long Faust is generated.
'dim_sizes',Sif it's an integer all Faust factors are square matrices (except maybe the first and last ones, depending on M and N). The size of the intermediary square factors is S**2. If it's a vector of 2 integers then the number of rows and columns are both a random number between size_dims(1) and size_dims(2) (inclusively). Defaultly, dim_sizes is [M, N].
'density',Dthe approximate density of generated factors. D must be a floating point number greater than 0 and lower of equal to 1. The default value is such that each factor gets 5 non-zero elements per row, if per_row is true, or per column otherwise. A density of zero is equivalent to the default case.
'per_row',boolthis argument is to specify the density per row or per column. By default the density is set per row and is such that the Faust's factors will have 5 non-zero elements per row.
fac_type,strthe storage representation of factors. str must be 'sparse', 'dense' or 'mixed'. The latter designates a mix of dense and sparse matrices in the generated Faust (the choice is made according to a uniform distribution). The default value is 'sparse'.
'field',strstr is either 'real' or 'complex' (the Faust field). The default value is 'real'.
'dev',str'gpu or 'cpu' to create the random Faust on CPU or GPU (by default on CPU).
'class',str'double' (by default) or 'single' to select the scalar type used for the Faust generated.
'seed',intseed to initialize the PRNG used to generate the Faust factors. If 0 (default) the seed will be initialized with a random value depending of the time clock.
Return values
Fthe random Faust.

Example 1

>> F = matfaust.rand(10,5)
F =
Faust size 10x5, density 4.32, nnz_sum 216, 5 factor(s):
- FACTOR 0 (real) SPARSE, size 10x9, density 0.555556, nnz 50
- FACTOR 1 (real) SPARSE, size 9x6, density 0.666667, nnz 36
- FACTOR 2 (real) SPARSE, size 6x10, density 0.5, nnz 30
- FACTOR 3 (real) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 4 (real) SPARSE, size 10x5, density 1, nnz 50

Example 2

% in a matlab terminal
>> G = matfaust.rand(10, 10, 'num_factors', 2, 'dim_sizes', 10, 'density', .5, 'mixed', 'complex')
G =
Faust size 10x10, density 1, nnz_sum 100, 2 factor(s):
- FACTOR 0 (complex) SPARSE, size 10x10, density 0.5, nnz 50
- FACTOR 1 (complex) DENSE, size 10x10, density 0.5, nnz 50

Example 3

>> H = matfaust.rand(10, 10, 'num_factors', [2, 5], 'dim_sizes', [10, 20], 'density', .5, 'fac_type', 'dense')
H =
Faust size 5x10, density 4.64, nnz_sum 232, 3 factor(s):
- FACTOR 0 (real) DENSE, size 5x14, density 0.5, nnz 35
- FACTOR 1 (real) DENSE, size 14x17, density 0.470588, nnz 112
- FACTOR 2 (real) DENSE, size 17x10, density 0.5, nnz 85

See also Faust.Faust, matfaust.rand_bsr.

◆ rand_bsr()

function matfaust::rand_bsr ( ,
,
BM  ,
BN  ,
varargin   
)

Generates a random Faust composed only of BSR matrices.

Parameters
Mthe number of rows of the random Faust.
Nthe number of columns of the random Faust.
BMthe nonzero block number of rows (must divide M).
BNthe nonzero block number of columns (must divide N).
'num_factors',NFif it's an integer it will be the number of random factors to set in the Faust. If NF is a vector of 2 integers then the number of factors will be set randomly between NF(1) and NF(2) (inclusively).
'density',Dthe approximate density of generated factors. D must be a floating point number between 0 and 1.
'field',strstr is either 'real' or 'complex' to set the Faust field. The default value is 'real'.
'dev',str'cpu' or 'gpu' to create the random Faust on CPU or GPU (by default on CPU).
'class',str'double' (by default) or 'single' to select the scalar type used for the generated Faust.
Return values
Fthe random Faust.

Example 1

>> matfaust.rand_bsr(10, 10, 2, 2)
ans =
Faust size 10x10, density 0.6, nnz_sum 60, 5 factor(s):
- FACTOR 0 (double) BSR, size 10x10, density 0.12, nnz 12
- FACTOR 1 (double) BSR, size 10x10, density 0.12, nnz 12
- FACTOR 2 (double) BSR, size 10x10, density 0.12, nnz 12
- FACTOR 3 (double) BSR, size 10x10, density 0.12, nnz 12
- FACTOR 4 (double) BSR, size 10x10, density 0.12, nnz 12

Example 2

>> matfaust.rand_bsr(128, 128, 32, 32, 'num_factors', [6, 10], 'field', 'real', 'class', 'double', 'density', .8)
ans =
Faust size 128x128, density 8.125, nnz_sum 133120, 10 factor(s):
- FACTOR 0 (double) BSR, size 128x128, density 0.8125, nnz 13312
- FACTOR 1 (double) BSR, size 128x128, density 0.8125, nnz 13312
- FACTOR 2 (double) BSR, size 128x128, density 0.8125, nnz 13312
- FACTOR 3 (double) BSR, size 128x128, density 0.8125, nnz 13312
- FACTOR 4 (double) BSR, size 128x128, density 0.8125, nnz 13312
- FACTOR 5 (double) BSR, size 128x128, density 0.8125, nnz 13312
- FACTOR 6 (double) BSR, size 128x128, density 0.8125, nnz 13312
- FACTOR 7 (double) BSR, size 128x128, density 0.8125, nnz 13312
- FACTOR 8 (double) BSR, size 128x128, density 0.8125, nnz 13312
- FACTOR 9 (double) BSR, size 128x128, density 0.8125, nnz 13312

See also Faust.Faust, matfaust.rand.

◆ rand_butterfly()

function matfaust::rand_butterfly ( ,
varargin   
)

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

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

Parameters
norder of the butterfly (must be a power of two).
diag_opt,logicalif true then the returned Faust is optimized using matfaust.opt_butterfly_faust.
'dev',str'gpu or 'cpu' to create the Faust on CPU or GPU ('cpu' by default).
'field',strstr is either 'real' or 'complex' (the Faust field). The default value is 'real'.
'class',str'double' (by default) or 'single' to select the scalar type used for the Faust generated.
Return values
Fa random butterfly support Faust.

See also: matfaust.fact.butterfly, matfaust.rand_butterfly, matfaust.dft, matfaust.opt_butterfly_faust

◆ toeplitz()

function matfaust::toeplitz ( ,
varargin   
)

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

Usage

    T = toeplitz(c), T is a symmetric Toeplitz Faust whose the first column is c.
    T = toeplitz(c, r), T is a Toeplitz Faust whose the first column is c and the first row is [c(1), r(2:)]

Parameters
cthe first column of the toeplitz Faust.
r(2nd argument) the first row of the toeplitz Faust. Defaulty r = conj(c). r(1) is ignored, the first row is always [c(1), r(2:)].
'dev',str'gpu' or 'cpu' to create the Faust on CPU or GPU ('cpu' is the default).
'diag_opt',logicalcf. matfaust.circ.
Return values
Tthe toeplitz Faust.

Example

% in a matlab terminal
>> import matfaust.toeplitz
>> c = 1:10;
>> T = toeplitz(c)
T =

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

  • FACTOR 0 (complex) SPARSE, size 10x32, density 0.0625, nnz 20, addr: 0x7f278c6aec20
  • FACTOR 1 (complex) SPARSE, size 32x32, density 0.0625, nnz 64, addr: 0x7f278ded1440
  • FACTOR 2 (complex) SPARSE, size 32x32, density 0.0625, nnz 64, addr: 0x7f278defd560
  • FACTOR 3 (complex) SPARSE, size 32x32, density 0.0625, nnz 64, addr: 0x7f278f687dc0
  • FACTOR 4 (complex) SPARSE, size 32x32, density 0.0625, nnz 64, addr: 0x7f27899a94d0
  • FACTOR 5 (complex) SPARSE, size 32x32, density 0.0625, nnz 64, addr: 0x7f2724d04e30
  • FACTOR 6 (complex) SPARSE, size 32x32, density 0.0625, nnz 64, addr: 0x7f278deb8fb0
  • FACTOR 7 (complex) SPARSE, size 32x32, density 0.0625, nnz 64, addr: 0x7f278f69e8f0
  • FACTOR 8 (complex) SPARSE, size 32x32, density 0.0625, nnz 64, addr: 0x7f278ded1840
  • FACTOR 9 (complex) SPARSE, size 32x10, density 0.0625, nnz 20, addr: 0x7f278c6b2070
>> full_T = full(T);
>> all(full_T(:,1).' - c < 1e-15)
ans =
logical
1
>> all(full_T(1, :) - c < 1e-15)
ans =
logical
1
>> r = 11:20;
>> T2 = toeplitz(c, r)

T2 =

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

  • FACTOR 0 (complex) SPARSE, size 10x32, density 0.0625, nnz 20, addr: 0x7f278dee2640
  • FACTOR 1 (complex) SPARSE, size 32x32, density 0.0625, nnz 64, addr: 0x7f278c694ff0
  • FACTOR 2 (complex) SPARSE, size 32x32, density 0.0625, nnz 64, addr: 0x7f278debd710
  • FACTOR 3 (complex) SPARSE, size 32x32, density 0.0625, nnz 64, addr: 0x7f278deb5e70
  • FACTOR 4 (complex) SPARSE, size 32x32, density 0.0625, nnz 64, addr: 0x7f278c6dc0d0
  • FACTOR 5 (complex) SPARSE, size 32x32, density 0.0625, nnz 64, addr: 0x7f278decd230
  • FACTOR 6 (complex) SPARSE, size 32x32, density 0.0625, nnz 64, addr: 0x7f278deb44a0
  • FACTOR 7 (complex) SPARSE, size 32x32, density 0.0625, nnz 64, addr: 0x7f278dee9c00
  • FACTOR 8 (complex) SPARSE, size 32x32, density 0.0625, nnz 64, addr: 0x7f278deb6090
  • FACTOR 9 (complex) SPARSE, size 32x10, density 0.0625, nnz 20, addr: 0x7f278c682e60
>> full_T2 = full(T2);
>> all(full_T2(1, :) - [c(1) , r(2:end)] < 1e-15)
ans =
logical
1
>> all(full_T2(:, 1).' - c < 1e-14)
ans =
logical
1
>> % Look at the density of a larger toeplitz Faust
>> % it indicates a speedup of the Faust-matrix/vector product
>> density(toeplitz(rand(1, 1024)))
0.0820

See also matfaust.circ, matfaust.anticirc

◆ version()

function matfaust::version ( )

Returns the FAuST package version.

◆ wht()

function matfaust::wht ( ,
varargin   
)

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

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

Usage

    H = wht(n)
    H = wht(n, 'normed', bool)
    H = wht(n, 'normed', bool, 'dev', str) str might be 'cpu' or 'gpu'.

Parameters
norder of the WHT (must be a power of two).
'normed',bool(optional) true (by default) to normalize the returned Faust as if Faust.normalize() was called, false otherwise.
'dev',str(optional) 'cpu' to create a CPU Faust (default choice) and 'gpu' for a GPU Faust.
'dtype',str(optional) 'double' (default choice) or 'float' to select the scalar type of the generated Faust.
Return values
Hthe Faust implementing the Hadamard transform of dimension n.

Example

% in a matlab terminal
>> import matfaust.*
>> H = wht(1024) % is equal to
>> H = normalize(wht(1024, 'normed', false, 'dev', 'cpu'))

H =

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

  • FACTOR 0 (real) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 1 (real) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 2 (real) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 3 (real) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 4 (real) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 5 (real) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 6 (real) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 7 (real) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 8 (real) SPARSE, size 1024x1024, density 0.00195312, nnz 2048
  • FACTOR 9 (real) SPARSE, size 1024x1024, density 0.00195312, nnz 2048

    See also: hadamard, matfaust.dft, matfaust.rand_butterfly, matfaust.fact.butterfly

matfaust::anticirc
function anticirc(c, varargin)
Returns an anticirculant Faust A defined by the vector c (which is the last column of full(A)).
matfaust::dct
function dct(n, varargin)
Constructs a Faust implementing the Direct Cosine Transform (Type II) Faust of order n.
matfaust::dft
function dft(n, varargin)
Constructs a Faust F implementing the Discrete Fourier Transform (DFT) of order n.
matfaust::wht
function wht(n, varargin)
Constructs a Faust implementing the Walsh-Hadamard Transform (WHT) of order n.
matfaust::eye
function eye(varargin)
Identity Faust.
matfaust::Faust
FAuST Matlab wrapper main class for using multi-layer sparse transforms.
Definition: Faust.m:118
matfaust
The FAuST Matlab Wrapper
Definition: bsl.m:1
matfaust::poly::next
function next(F)
Gives the next Faust basis of dimension (n+1) from the Faust F polynomial basis of dimension n.
matfaust::circ
function circ(c, varargin)
Returns a circulant Faust C defined by the vector c (which is the first column of full(C)).
matfaust::rand
function rand(M, N, varargin)
Generates a random Faust.
matfaust::toeplitz
function toeplitz(c, varargin)
Constructs a toeplitz Faust whose first column is c and first row r.
matfaust::dst
function dst(n, varargin)
Constructs a Faust implementing the Direct Sine Transform (Type II) Faust of order n.