-  3.39.22
Public Member Functions | Static Public Member Functions | Protected Attributes | List of all members
matfaust::Faust Class Reference

FAuST Matlab wrapper main class for using multi-layer sparse transforms. More...

Public Member Functions

function Faust (varargin)
 Creates a Faust from a list of factors or alternatively from a file. More...
 
function delete (F)
 Deletes the Faust object F (destructor). More...
 
function uplus (F)
 Returns + F. More...
 
function uminus (F)
 Returns - F. More...
 
function plus (varargin)
 Plus. More...
 
function minus (varargin)
 Minus. More...
 
function mrdivide (F, s)
 / Slash or right Faust divide. More...
 
function mtimes (F, A)
 Multiplies the Faust F by A which is a full matrix, a Faust object or a scalar. More...
 
function full (F)
 The full matrix implemented by F. More...
 
function isreal (F)
 Indicates if F is a real Faust or a complex Faust. More...
 
function real (F)
 Returns the real part of the Faust F. More...
 
function imag (F)
 Returns the imaginary part of the Faust F. More...
 
function single (F)
 Returns a new Faust whose class is single. It allows to convert all the factors to single precision. More...
 
function double (F)
 Returns a new Faust whose class is double. It allows to convert all the factors to double precision. More...
 
function transpose (F)
 The transpose of F. More...
 
function ctranspose (F)
 The conjugate transpose of F. More...
 
function conj (F)
 The complex conjugate of F. More...
 
function pruneout (F, varargin)
 Returns a Faust optimized by removing useless zero rows and columns as many times as needed. More...
 
function optimize_memory (F)
 Optimizes a Faust by changing the storage format of each factor in order to optimize the memory size. More...
 
function optimize (F, varargin)
 Returns a Faust optimized with Faust.pruneout, Faust.optimize_memory and Faust.optimize_time. More...
 
function optimize_time (F, varargin)
 Returns a Faust configured with the quickest Faust-matrix multiplication mode (benchmark ran on the fly). More...
 
function size (F, varargin)
 The size of F. More...
 
function numel (F)
 The number of elements in F. More...
 
function nbytes (F)
 Gives the memory size of the Faust in bytes. More...
 
function end (F, k, n)
 The last index when slicing or indexing a Faust. More...
 
function factors (F, ids, varargin)
 Returns the i-th factor or a new Faust composed of F factors whose indices are listed in indices. More...
 
function left (F, i, varargin)
 Returns the left hand side factors of F from index 1 to i included (in 1-base index). More...
 
function right (F, i, varargin)
 Returns the right hand side factors of F from index i to end (in 1-base index). More...
 
function numfactors (F)
 The number of factors of F. More...
 
function length (F)
 The number of factors of F. More...
 
function replace (F, i, new_factor)
 Replaces the factor of index i by new_factor in a new Faust copy of F. More...
 
function insert (F, i, new_factor)
 Inserts new_factor at index i in a new Faust copy of F. More...
 
function issparse (F, varargin)
 Returns true if F factors are all sparse matrices false otherwise. More...
 
function isdense (F)
 Returns true if F factors are all dense matrices/arrays false otherwise. More...
 
function save (F, filepath)
 Saves the Faust F into a file. More...
 
function subsref (F, S)
 Subscripted reference of a Faust. More...
 
function disp (F)
 Displays information about F. More...
 
function norm (F, varargin)
 The matrix norm of F. More...
 
function power_iteration (F, varargin)
 Performs the power iteration algorithm to compute the greatest eigenvalue of the Faust. More...
 
function normalize (F, varargin)
 Returns the normalized F. More...
 
function nnz_sum (F)
 The total number of nonzeros in the factors of F. More...
 
function density (F)
 The density of F such that nnz_sum(F) == density(F) * numel(F). More...
 
function rcg (F)
 The Relative Complexity Gain of F. More...
 
function cat (varargin)
 Concatenates F with n Faust objects or full/sparse matrices. More...
 
function horzcat (varargin)
 Horizontal concatenation [F,A, B, … ] where F is a Faust and A, B, … are Faust objects or sparse/full matrix. More...
 
function vertcat (varargin)
 Vertical concatenation [F;A;B;…] where F is a Faust and A, B, … are Faust objects or sparse/full matrices. More...
 
function imagesc (F, varargin)
 Displays image of F full matrix and its factors. More...
 
function mldivide (F, B)
 \ Backslash or left full(F) divide. More...
 
function pinv (F)
 Pseudoinverse matrix of full(F). More...
 
function complex (F)
 Converts F to a new complex Faust. More...
 
function clone (F, varargin)
 Clones the Faust (in a new memory space). More...
 
function device (F)
 Returns the Faust device ('cpu' or 'gpu'). More...
 
function class (F)
 Returns the Faust class ('single' or 'double'). More...
 
function get_handle (F)
 Returns the Faust F C++ native object pointer as an integer. More...
 

Static Public Member Functions

static function opt_butterfly (F)
 Optimizes any Faust composed of butterfly factors. More...
 
static function isFaust (obj)
 Returns true if obj is a Faust object, false otherwise. More...
 
static function load (filepath, varargin)
 Loads a Faust from a .mat file. More...
 
static function load_native (filepath)
 Loads a Faust from a .mat file (native C++ version). More...
 

Protected Attributes

Property matrix
 Underlying Faust native object handle. More...
 
Property is_real
 See Faust.isreal. More...
 
Property dev
 See Faust.device. More...
 
Property dtype
 See Faust.class. More...
 

Detailed Description

FAuST Matlab wrapper main class for using multi-layer sparse transforms.

This class provides a Matlab array-like interface for operations with FAuST data structures, which correspond 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) Matlab array.

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 matfaust.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 Matlab builtins have been overloaded to ensure that a Faust is almost handled as a native Matlab matrix.

The main exception is that contrary to a Matlab native array a Faust is immutable. It means that you cannot modify elements of a Faust using the assignment operator ‘=’ like you do with a Matlab matrix (e.g. ‘M(i,j) = 2’). That limitation is the reason why the Matlab built-in ‘SUBSASGN()’ 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 Faust.optimize_time).

Other noticeable limitations are that one cannot:

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

Note
it could be wiser to encapsulate a Faust in a lazylinop.LazyLinearOp for a totally lazy paradigm on all available operations.
Warning
using Faust.full 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 Matlab array Faust.full() obtained by the multiplication of the Faust factors.

List of functions that are memory costly:

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

Constructor & Destructor Documentation

◆ Faust()

function matfaust::Faust::Faust ( varargin  )

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: matfaust.rand(), matfaust.dft() or matfaust.wht().

Usage

    Faust(factors) creates a Faust from a list of factors (1D cell array).

    Faust(filepath) creates a Faust from a previously saved Faust file (filepath is a character array).

Parameters
factors(varargin{1}, first argument) the 1D cell array of factors to initialize the Faust with.
The factors must respect the dimensions needed for the product to be defined (for i=1 to size(factors,2), size(factors{i},2) == size(factors{i+1},1)).
The factors can be sparse, dense or BSR matrices. It might be real matrices in single or double precision but also complex matrices. The created Faust will be consequently a single, double or complex Faust.
To create a BSR matrix you must use the matfaust.create_bsr function.
Passing a single matrix to the constructor instead of a cell array is equivalent to passing a singleton cell array.
filepath(varargin{1}, first argument) the file from which a Faust is created. It must be a character array.
The format is Matlab version 5 (.mat extension).
The file must have been saved before with Faust.save().
'dev',str'gpu or 'cpu' to create the Faust on CPU or GPU (by default on CPU). Examples
>> import matfaust.Faust
>> factors = cell(1,5);
>> is_sparse = false;
>> % alternate sparse and dense factors
>> for i=1:5; if(is_sparse);factors{i} = sprand(100, 100, 0.1); else factors{i} = rand(100, 100); end; is_sparse = ~ is_sparse; end
>> % define a Faust with those factors
>> F = Faust(factors)
F =
Faust size 100x100, density 3.1897, nnz_sum 31897, 5 factor(s):
- FACTOR 0 (double) DENSE, size 100x100, density 1, nnz 10000
- FACTOR 1 (double) SPARSE, size 100x100, density 0.0943, nnz 943
- FACTOR 2 (double) DENSE, size 100x100, density 1, nnz 10000
- FACTOR 3 (double) SPARSE, size 100x100, density 0.0954, nnz 954
- FACTOR 4 (double) DENSE, size 100x100, density 1, nnz 10000
>> save(F, 'F.mat')
>> % define a Faust from file
>> G = Faust('F.mat')
G =
Faust size 100x100, density 3.1897, nnz_sum 31897, 5 factor(s):
- FACTOR 0 (double) DENSE, size 100x100, density 1, nnz 10000
- FACTOR 1 (double) SPARSE, size 100x100, density 0.0943, nnz 943
- FACTOR 2 (double) DENSE, size 100x100, density 1, nnz 10000
- FACTOR 3 (double) SPARSE, size 100x100, density 0.0954, nnz 954
- FACTOR 4 (double) DENSE, size 100x100, density 1, nnz 10000
>> H = Faust(rand(10,10)) % creating a Faust with only one factor
H =
Faust size 10x10, density 1, nnz_sum 100, 1 factor(s):
- FACTOR 0 (double) DENSE, size 10x10, density 1, nnz 100

See also Faust.delete, Faust.save, matfaust.rand, matfaust.dft, matfaust.wht

Member Function Documentation

◆ cat()

function matfaust::Faust::cat ( varargin  )

Concatenates F with n Faust objects or full/sparse matrices.

This function overloads a Matlab built-in function.

The resulting Faust C = cat(DIM,F,G,…) verifies that:

full(C) == full(cat(DIM, full(F), full(G), …))
Note
you could have an elementwise non-significant absolute difference between the two members (not more than eps(1.0)).

Usage

    C=CAT(DIM,F,G) concatenates the Faust F and G, which is a Faust or a matrix, along the dimension DIM. The result is the Faust C.
    CAT(2,F,G) is the same as [F,G].
    CAT(1,F,G) is the same as [F;G].
    CAT(DIM, F, G,…) concatenates an arbitrary number of Fausts and matrices along the DIM-th dimension.

Parameters
DIM(1st arg.) the dimension along which to concatenate. DIM==1 means vertical concatenation, DIM==2 means horizontal concatenation.
F(2nd arg.) the first Faust object.
G,…(3rd to n+3 args) the Faust objects or matrices to concatenate to F. The matrices can be sparse or full.
Return values
Cthe concatenation result as a new Faust.
Note
it could be wiser to encapsulate a Faust in a lazylinop.LazyLinearOp for a lazy concatenation.

Example

>>% in a matlab terminal
>> rng(42); % for reproducibility
>> F = matfaust.rand(50, 50);
>> G = matfaust.rand(50, 50);
>> [F;G] % equiv. to cat(1,F,G)
ans =
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,G] % equiv. to cat(2,F,G)
ans =
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
>> [F;rand(100,50)] % vertical concatenation with auto-conversion of the random matrix to a Faust
ans =
Faust size 150x50, density 0.873333, nnz_sum 6550, 6 factor(s):
- FACTOR 0 (double) SPARSE, size 150x100, density 0.35, nnz 5250
- 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
>> [F;sprand(100,50,.2)] % vertical concatenation with auto-conversion of the random sparse matrix to a Faust
ans =
Faust size 150x50, density 0.329733, nnz_sum 2473, 6 factor(s):
- FACTOR 0 (double) SPARSE, size 150x100, density 0.0782, nnz 1173
- 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
>> [F;G;F;G]; % it's allowed to concatenate an arbitrary number of Fausts
>> [F,G,F,G]; % as long as the dimensions agree
>> [F,G;F,G];

Errors

  • The dimensions of F and G don't agree:
>> F = matfaust.rand(2,51);
>> G = matfaust.rand(2,25);
>> [F;G]
??? Error using mexFaustReal
The dimensions of the two Fausts must agree.
  • The DIM argument is different from 1 and 2:
>> F = matfaust.rand(2,4);
>> G = matfaust.rand(2,4);
>> cat(3,F,G)
??? Error using matfaust.Faust/cat
Wrong first argument: must be an integer between 1 and 2.

See also Faust.vertcat, Faust.horzcat.

◆ class()

function matfaust::Faust::class ( )

Returns the Faust class ('single' or 'double').

Return values
cthe Faust class.

Example:

>> class(matfaust.Faust(rand(10, 10)))
ans =
'double'
>> class(matfaust.Faust(single(rand(10, 10))))
ans =
'single'
>> class(matfaust.Faust(complex(rand(10, 10))))
ans =
'double'

◆ clone()

function matfaust::Faust::clone ( ,
varargin   
)

Clones the Faust (in a new memory space).

Parameters
'dev',str'gpu or 'cpu' to clone the Faust on CPU or GPU (by default on CPU).

Example

>> F = matfaust.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
>> Fc = clone(F)
Fc =
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
>> % and only if your device supports a compatible NVIDIA device
>> Fc_gpu = clone(F, 'dev', 'gpu'); % doctest: +SKIP_IF(~ matfaust.is_gpu_mod_enabled() || ~ matfaust.is_gpu_mod_working())
>> device(gpuF) % doctest: +SKIP_IF(~ matfaust.is_gpu_mod_enabled() || ~ matfaust.is_gpu_mod_working())
ans =
'gpu'
Return values
Fcthe Faust clone.

◆ complex()

function matfaust::Faust::complex ( )

Converts F to a new complex Faust.

This function overloads a Matlab built-in function.

Usage

    cF = COMPLEX(F) for real Faust F returns the complex result cF with real part F and zero matrix as imaginary part of all cF's factors.

Example:

>> dF = matfaust.rand(5, 5, 'field', 'real', 'class', 'double')
dF =
Faust size 5x5, density 5, nnz_sum 125, 5 factor(s):
- FACTOR 0 (double) SPARSE, size 5x5, density 1, nnz 25
- FACTOR 1 (double) SPARSE, size 5x5, density 1, nnz 25
- FACTOR 2 (double) SPARSE, size 5x5, density 1, nnz 25
- FACTOR 3 (double) SPARSE, size 5x5, density 1, nnz 25
- FACTOR 4 (double) SPARSE, size 5x5, density 1, nnz 25
>> cF = complex(dF)
cF =
Faust size 5x5, density 5, nnz_sum 125, 5 factor(s):
- FACTOR 0 (complex) SPARSE, size 5x5, density 1, nnz 25
- FACTOR 1 (complex) SPARSE, size 5x5, density 1, nnz 25
- FACTOR 2 (complex) SPARSE, size 5x5, density 1, nnz 25
- FACTOR 3 (complex) SPARSE, size 5x5, density 1, nnz 25
- FACTOR 4 (complex) SPARSE, size 5x5, density 1, nnz 25

See also Faust.isreal, Faust.conj

◆ conj()

function matfaust::Faust::conj ( )

The complex conjugate of F.

This function overloads a Matlab built-in function.

Usage

    F_conj = conj(F)

Parameters
Fthe Faust object.
Return values
F_conja Faust object implementing the conjugate of full(F), such that:
full(F_conj) == conj(full(F))

Example

>> F = matfaust.rand(5, 5, 'field', 'complex');
>> F_conj = conj(F);
>> all(all(full(F_conj) == conj(full(F))))
ans =
logical
1

See also Faust.transpose, Faust.ctranspose, Faust.complex

◆ ctranspose()

function matfaust::Faust::ctranspose ( )

The conjugate transpose of F.

This function overloads a Matlab built-in function/operator.

Usage

    F_ctrans = ctranspose(F)
    F_ctrans = F'

Parameters
Fthe Faust object.
Return values
F_ctransa Faust object implementing the conjugate transpose of full(F), such that:
full(F_ctrans) == full(F)' == ctranspose(full(F))

Example

>> F = matfaust.rand(5, 10, 'field', 'complex');
>> F_ctrans = F'; % ' == ctranspose
>> size(F_ctrans)
ans =
10 5
>> F_c = transpose(ctranspose(F));
>> all(all(full(F_c) == full(conj(F))))
ans =
logical
1

See also Faust.transpose, Faust.conj, Faust.complex

◆ delete()

function matfaust::Faust::delete ( )

Deletes the Faust object F (destructor).

Usage

    delete(F)

Parameters
Fthe Faust to delete.

Example

>> F = matfaust.rand(5, 10);
>> delete(F) % underlying C++ object wiped out of memory but variable F still exists
>> exist('F')
ans =
1
>> F
??? Invalid or deleted object.
>>
>> F = matfaust.rand(5, 10);
>> clear F % underlying C++ object is removed and variable too
>> exist('F')
ans =
0
>> F
??? Unrecognized function or variable 'F'.

See also Faust.Faust, clear (built-in)

◆ density()

function matfaust::Faust::density ( )

The density of F such that nnz_sum(F) == density(F) * numel(F).

Usage

    dens = density(F)

Note
A value of density below one indicates potential memory savings compared to storing the corresponding dense matrix full(F), as well as potentially faster matrix-vector multiplication when applying F*x instead of full(F)*x.
A density above one is possible but prevents any saving.
Parameters
Fthe Faust object.
Return values
densdensity of F.

Example

>> F = matfaust.rand(10, 10);
>> density(F)
ans =
2.5000

See also Faust.nnz_sum, Faust.rcg, Faust.size, Faust.numel

◆ device()

function matfaust::Faust::device ( )

Returns the Faust device ('cpu' or 'gpu').

Return values
devthe Faust device.

Example:

>> cpuF = matfaust.rand(5, 5, 'dev', 'cpu');
>> device(cpuF)
ans =
'cpu'
>> gpuF = matfaust.rand(5,5, 'dev', 'gpu'); % doctest: +SKIP_IF(~ matfaust.is_gpu_mod_enabled() || ~ matfaust.is_gpu_mod_working())
>> device(gpuF) % doctest: +SKIP_IF(~ matfaust.is_gpu_mod_enabled() || ~ matfaust.is_gpu_mod_working())
ans =
'gpu'

◆ disp()

function matfaust::Faust::disp ( )

Displays information about F.

This function overloads a Matlab built-in function.

Usage

    disp(F)

    F

Parameters
Fthe Faust object.

Example

>> F = matfaust.rand(98, 82, 'num_factors', 2);
>> disp(F) % doctest: +ELLIPSIS
Faust size 98x82, %> %> %>
>> F % doctest: +ELLIPSIS
F =
Faust size 98x82, %> %> %>

See also Faust.nnz_sum, Faust.density, Faust.size, Faust.factors, Faust.numfactors

◆ double()

function matfaust::Faust::double ( )

Returns a new Faust whose class is double. It allows to convert all the factors to double precision.

Usage

    dF = double(F)

Return values
sFthe double class/precision Faust.

Example:

>> sF = matfaust.rand(5, 5, 'field', 'real', 'class', 'single')
sF =
Faust size 5x5, density 5, nnz_sum 125, 5 factor(s):
- FACTOR 0 (float) SPARSE, size 5x5, density 1, nnz 25
- FACTOR 1 (float) SPARSE, size 5x5, density 1, nnz 25
- FACTOR 2 (float) SPARSE, size 5x5, density 1, nnz 25
- FACTOR 3 (float) SPARSE, size 5x5, density 1, nnz 25
- FACTOR 4 (float) SPARSE, size 5x5, density 1, nnz 25
>> dF = double(sF)
dF =
Faust size 5x5, density 5, nnz_sum 125, 5 factor(s):
- FACTOR 0 (double) SPARSE, size 5x5, density 1, nnz 25
- FACTOR 1 (double) SPARSE, size 5x5, density 1, nnz 25
- FACTOR 2 (double) SPARSE, size 5x5, density 1, nnz 25
- FACTOR 3 (double) SPARSE, size 5x5, density 1, nnz 25
- FACTOR 4 (double) SPARSE, size 5x5, density 1, nnz 25

See also Faust.class, Faust.single

◆ end()

function matfaust::Faust::end ( ,
,
 
)

The last index when slicing or indexing a Faust.

size(F, K) == end when indexing the K-th dimension of F \((K \in {1, 2})\).

This function overloads a Matlab built-in function.

Example

% in a matlab terminal
>> rng(42) % for reproducibility
>> F = matfaust.Faust(rand(3, 3));
>> full(F)
ans =
0.3745 0.5987 0.0581
0.9507 0.1560 0.8662
0.7320 0.1560 0.6011
>> full(F(2:end,1))
ans =
0.9507
0.7320
>> full(F(1,1:2:end))
ans =
0.3745 0.0581
>> full(F(1,1:end-1))
ans =
0.3745 0.5987

See also Faust.subsref, Faust.size

◆ factors()

function matfaust::Faust::factors ( ,
ids  ,
varargin   
)

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

Note
Factors are not copied in memory if a subset greater than one is asked or if as_faust == true.

Usage

    factor = factors(F, i) returns the i-th factor of F.
    factor = factors(F, i:j) returns a new Faust formed of the F's factors from the i-th to the j-th included.

Parameters
Fthe Faust object.
idsthe factor indices.
'as_faust',logical(optional): true to return a Faust even if a single factor is asked (length(ids) == 1), otherwise (as_faust == false) and a dense or sparse matrix is returned.
Return values
factorsa matrix copy of the i-th factor if i is a single index or a new Faust composed of i-th to the j-th factors of F. The factors copies keep the storage organization of the source matrix (full or sparse).
Note
Matlab doesn't support float sparse matrices, but matfaust does! Hence when you call Faust.factors on a float sparse Faust to retrieve one factor you'll get a double sparse matrix as a copy of the float sparse matrix.
As well for BSR matrices that aren't not supported by Matlab, the function can't return a bsr matrix so it rather converts it on the fly to a sparse matrix that is finally returned (it applies also for other specialized matrix types).

Example

>> rng(42) % for reproducibility
>> F = matfaust.Faust({sprand(5, 5, .2), sprand(5, 5, .2), sprand(5, 5, .2)});
>> f1 = factors(F, 1)
f1 =
(2,1) 0.9699
(5,1) 0.1818
(1,4) 0.0206
(3,4) 0.8324
(4,5) 0.2123
>> G = factors(F, 2:3) % a new Faust composed of the two last factors of F
G =
Faust size 5x5, density 0.36, nnz_sum 9, 2 factor(s):
- FACTOR 0 (double) SPARSE, size 5x5, density 0.16, nnz 4
- FACTOR 1 (double) SPARSE, size 5x5, density 0.2, nnz 5

See also Faust.numfactors

◆ full()

function matfaust::Faust::full ( )

The full matrix implemented by F.

This function overloads a Matlab built-in function.

Usage

    A = full(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 dense matrix representation.
Return values
Athe Matlab native matrix such that A*x == F*x for any vector x.
Warning
Running the example below is likely to raise a memory error or freeze your computer for a certain amount of time.

Example

% in a matlab terminal
>> F = matfaust.rand(10^5, 10^5, 'num_factors', 2, '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 full matrix is most likely to raise a memory error
>> % the sparse format is the only way to handle such a large Faust
>> full(F)
??? Error using mexFaustReal
std::bad_alloc

◆ get_handle()

function matfaust::Faust::get_handle ( )

Returns the Faust F C++ native object pointer as an integer.

Return values
Hhandle (address as an integer).

Example

% in a matlab terminal
>> F = matfaust.rand(10, 10);
>> get_handle(F) % doctest: +ELLIPSIS
ans =
uint64
%>

◆ horzcat()

function matfaust::Faust::horzcat ( varargin  )

Horizontal concatenation [F,A, B, … ] where F is a Faust and A, B, … are Faust objects or sparse/full matrix.

It's equivalent to cat(2, F, A, B, …).

This function overloads a Matlab built-in function.

Usage

    C=HORZCAT(F,A) concatenates the Faust F and A horizontally. A is a Faust or a sparse/full matrix. The result is the Faust C. HORZCAT(F,A) is the same as [F,G].
    C=HORZCAT(F,A,B,…) concatenates horizontally F to A, B, etc. HORZCAT(F,A, B,…) is the same as [F,A,B,…].

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

See also Faust.vertcat, Faust.cat.

◆ imag()

function matfaust::Faust::imag ( )

Returns the imaginary part of the Faust F.

Example:

% in a matlab terminal
>> F = matfaust.rand(10, 10, 'field', 'complex')
F =
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
>> iF = imag(F)
iF =
Faust size 10x10, density 8, nnz_sum 800, 5 factor(s):
- FACTOR 0 (double) SPARSE, size 10x20, density 0.5, nnz 100
- FACTOR 1 (double) SPARSE, size 20x20, density 0.5, nnz 200
- FACTOR 2 (double) SPARSE, size 20x20, density 0.5, nnz 200
- FACTOR 3 (double) SPARSE, size 20x20, density 0.5, nnz 200
- FACTOR 4 (double) SPARSE, size 20x10, density 0.5, nnz 100
>> norm(full(iF) - imag(full(F))) < 1e-12
ans =
logical
1

See also Faust.real

◆ imagesc()

function matfaust::Faust::imagesc ( ,
varargin   
)

Displays image of F full matrix and its factors.

This function overloads a Matlab built-in function.

Usage

    imagesc(F)
    imagesc(F, 'faust_name') to name the Faust on the plotted figure.

Parameters
faust_name(varargin{1}) The Faust name shown in the figure title (defaultly it's 'F').

Data is scaled to use the full colormap.

Example

% in a matlab terminal (without -nojvm flag)
>> F = matfaust.rand(5, 10);
>> imagesc(F) % imagesc(F, 'my faust') to name the faust differently than 'F'
>> print('test_imsc', '-dpng') % figure saved in test_imsc.png file

See also Faust.disp.

◆ insert()

function matfaust::Faust::insert ( ,
,
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,intthe index of insertion.
new_factor,matrixthe factor to insert as the i-th factor of F.
Return values
Ga copy of F with new_factor inserted at index i.

Example

>> F = matfaust.rand(10, 10, 'num_factors', 5);
>> S = sprand(10, 10, .2);
>> G = insert(F, 3, S);
>> all(all(full(left(F, 2) * matfaust.Faust(S) * right(F, 3)) == full(G)))
ans =
logical
1

See also Faust.factors, Faust.left, Faust.right, Faust.replace

◆ isdense()

function matfaust::Faust::isdense ( )

Returns true if F factors are all dense matrices/arrays false otherwise.

Example

>> F = matfaust.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
>> isdense(F)
ans =
1
>> F = matfaust.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
>> isdense(F)
ans =
0
>> F = matfaust.rand(10, 10, 'fac_type', 'dense') * matfaust.Faust(speye(10))
F =
Faust size 10x10, density 2.6, nnz_sum 260, 6 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
- FACTOR 5 (double) SPARSE, size 10x10, density 0.1, nnz 10
>> isdense(F)
ans =
0

See also Faust.issparse, matfaust.rand, matfaust.rand_bsr

◆ isFaust()

static function matfaust::Faust::isFaust ( obj  )
static

Returns true if obj is a Faust object, false otherwise.

Example

% in a matlab terminal
>> import matfaust.Faust
>> import matfaust.isFaust
>> F = matfaust.rand(5,10);
>> Faust.isFaust(F) % isFaust(F) alone works as well
ans =
logical
1
>> Faust.isFaust(1)
ans =
logical
0

See also Faust.Faust

◆ isreal()

function matfaust::Faust::isreal ( )

Indicates if F is a real Faust or a complex Faust.

This function overloads a Matlab built-in function.

Note
If F is a real Faust, all factors are real. If F is a complex Faust, at least one factor is complex.

Usage

    bool = isreal(F)

Return values
bool1 if F is a real Faust, 0 if it's a complex Faust.

Example:

>> F = matfaust.rand(10, 10, 'field', 'real')
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
>> isreal(F)
ans =
logical
1
>> F = matfaust.rand(10, 10, 'field', 'complex')
F =
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
>> isreal(F)
ans =
logical
0

See also Faust.complex

◆ issparse()

function matfaust::Faust::issparse ( ,
varargin   
)

Returns true if F factors are all sparse matrices false otherwise.

What a sparse factor is, depends on csr and bsr arguments. Defaultly, only a Faust full of CSR matrices is taken as a sparse Faust.

Parameters
csr,booltrue to consider a CSR matrix as a sparse matrix, false otherwise (true by default).
bsr,booltrue to consider a BSR matrix as a sparse matrix, false otherwise (false by default).

Example:

% in a matlab terminal
>> F = matfaust.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
>> issparse(F) % equivalent to issparse(F, 'csr', true)
ans =
1
>> issparse(F, 'csr', false, 'bsr', false)
??? Error using matfaust.Faust/issparse
It doesn't make sense to set csr=false and bsr=false as the function will
always return false
>> F = matfaust.rand_bsr(10, 10, 2, 2, 2, .1)
F =
Faust size 10x10, density 0.6, nnz_sum 60, 5 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)
- FACTOR 2 (double) BSR, size 10x10 (blocksize = 2x2), density 0.12, nnz 12 (nnz blocks: 3)
- FACTOR 3 (double) BSR, size 10x10 (blocksize = 2x2), density 0.12, nnz 12 (nnz blocks: 3)
- FACTOR 4 (double) BSR, size 10x10 (blocksize = 2x2), density 0.12, nnz 12 (nnz blocks: 3)
>> issparse(F) % default config. recognizes only csr
ans =
0
>> issparse(F, 'bsr', true)
ans =
1
>> F = matfaust.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
>> issparse(F)
ans =
0

See also Faust.isdense, matfaust.rand, matfaust.rand_bsr

◆ left()

function matfaust::Faust::left ( ,
,
varargin   
)

Returns the left hand side factors of F from index 1 to i included (in 1-base index).

Parameters
Fthe Faust object.
ithe bound factor index.
'as_faust',logical(optional): true to return a Faust even if a single factor is asked (length(ids) == 1), otherwise (as_faust == false) and a dense or sparse matrix is returned.
Return values
aFaust if the number of factors to be returned is greater than 1, an array or a sparse matrix otherwise.

Example

%in a matlab terminal
>> F = matfaust.rand(8, 5, 'seed', 42)
F =
Faust size 8x5, density 3.675, nnz_sum 147, 5 factor(s):
- FACTOR 0 (double) SPARSE, size 8x7, density 0.714286, nnz 40
- FACTOR 1 (double) SPARSE, size 7x6, density 0.666667, nnz 28
- FACTOR 2 (double) SPARSE, size 6x6, density 0.666667, nnz 24
- FACTOR 3 (double) SPARSE, size 6x5, density 1, nnz 30
- FACTOR 4 (double) SPARSE, size 5x5, density 1, nnz 25
>> LF = left(F, 3)
LF =
Faust size 8x6, density 1.91667, nnz_sum 92, 3 factor(s):
- FACTOR 0 (double) SPARSE, size 8x7, density 0.714286, nnz 40
- FACTOR 1 (double) SPARSE, size 7x6, density 0.666667, nnz 28
- FACTOR 2 (double) SPARSE, size 6x6, density 0.666667, nnz 24

See also Faust.factors, Faust.right

◆ length()

function matfaust::Faust::length ( )

The number of factors of F.

Usage

    A = length(F)

Parameters
Fthe Faust object.
Return values
num_factorsthe number of factors.

Example

>> F = matfaust.rand(5, 10, 'num_factors', 6);
>> nf = length(F)
nf =
6

See also Faust.numfactors

◆ load()

static function matfaust::Faust::load ( filepath  ,
varargin   
)
static

Loads a Faust from a .mat file.

Example

>> import matfaust.Faust
>> F = matfaust.rand(10,10);
>> save(F, 'F.mat')
>> F2 = matfaust.Faust.load('F.mat');
>> F3 = Faust('F.mat');
>> all(all(full(F2) == full(F)))
ans =
logical
1
>> all(all(full(F3) == full(F)))
ans =
logical
1

See also Faust.Faust, Faust.save

◆ load_native()

static function matfaust::Faust::load_native ( filepath  )
static

Loads a Faust from a .mat file (native C++ version).

Example

>> import matfaust.Faust
>> F = matfaust.rand(10,10);
>> save(F, 'F.mat');
>> F2 = matfaust.Faust.load_native('F.mat');
>> F3 = Faust('F.mat');
>> all(all(full(F2) == full(F)))
ans =
logical
1
>> all(all(full(F3) == full(F)))
ans =
logical
1

See also Faust.Faust, Faust.save

◆ minus()

function matfaust::Faust::minus ( varargin  )

Minus.

This function overloads a Matlab built-in function.

Usage

    minus(F,G) or F-G subtracts the Faust G from F, sizes must be compatible.
    minus(F,A) or F-A subtracts a matrix A from F, sizes must be compatible.
    minus(F,s) or F-s subtracts a scalar s from F, such that full(F-s) == full(F)-s.

Parameters
F(first arg.) The Faust object.
G,A,s,…(2nd to n-th args) The variables to subtract from F; Fausts, matrices or scalars.
Return values
Mthe difference as a Faust object.

Example:

>> F = matfaust.rand(10, 12);
>> G = matfaust.rand(10, 12);
>> FmG = F-G; % is a Faust
>> norm(full(FmG) - full(F) + full(G)) < 1e-12
ans =
logical
1
>> Fm2 = F-2; % is a Faust
>> norm(full(Fm2) - full(F) + 2) < 1e-12
ans =
logical
1
>> FmGm2 = F-G-2; % is a Faust
>> norm(full(FmGm2) - full(F) + full(G) + 2) < 1e-12
ans =
logical
1
>> H = F-G-F-2-F; % is a Faust
>> norm(full(H) - full(F) + 2*full(F) + full(G) + 2) < 1e-12
ans =
logical
1

See also Faust.plus

◆ mldivide()

function matfaust::Faust::mldivide ( ,
 
)

\ Backslash or left full(F) divide.

Usage

    X = F\ B is the matrix division of full(F) into B, which is roughly the same as Faust.pinv(F)*B.

Warning
this function makes a call to Faust.full.

See also Faust.pinv, mldivide Matlab built-in.

◆ mrdivide()

function matfaust::Faust::mrdivide ( ,
 
)

/ Slash or right Faust divide.

Usage

    G = F/s is the division of the Faust F by the scalar s, such that full(F)/s == full(F/s)
    G = MRDIVIDE(F,s) is called for the syntax 'F / s' when s is a scalar.

Parameters
FThe Faust object to divide.
sThe scalar to divide F with.
Return values
Gthe division result as a Faust object.

See also Faust.mtimes

◆ mtimes()

function matfaust::Faust::mtimes ( ,
 
)

Multiplies the Faust F by A which is a full matrix, a Faust object or a scalar.

This function overloads a Matlab built-in function.

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

Other use cases are 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:
    full(F*A) == full(F)*full(A)
    Note
    you could have an elementwise non-significant absolute difference between the two members (not more than eps(1.0)).
  • If A is a scalar, F*A is also a Faust such that:
    factors(F*A,1) == factors(F,1)*A
  • Any Faust F (real or complex) can be multiplied by any (real or complex) Faust/array of appropriate dimensions or a scalar.

Usage

    G = mtimes(F, A)
    G = F*A with A and G being full matrices. A could also be a sparse matrix but note that often the Faust-sparse matrix multiplication is slower than performing F*full(A) multiplication. In some cases though, it stays quicker: moreover when the Faust is composed of a small number of factors.
    G = F*s with s a scalar and G a Faust.
    G = s*F with s a scalar and G a Faust.
    G = s*F' with s a scalar multiplying the conjugate-transpose of F.
    G = F'*s with s a scalar multiplying the conjugate-transpose of F.

Parameters
Fthe Faust object.
Aeither a full matrix, a Faust object or a scalar.
Return values
Gthe multiplication result:
  • When A is a full matrix, G = F*A is also a full matrix.
  • When A is a Faust or a scalar, G = F*A is itself a Faust.
  • When either F or A is complex, G=F*A is also complex.

Example

% in a matlab terminal
>> rng(42) % for reproducibility
>> F = matfaust.rand(5, 4, 'seed', 42);
>> A = rand(size(F,2), 5);
>> G = F * A;
>> norm(G - full(F) * A) < 1e-12
ans =
logical
1
>> % this is equivalent to G = mtimes(F, A)
>> G = matfaust.rand(size(F, 2), 15);
>> H = F * G; % H is a Faust because F and G are
>> matfaust.isFaust(H)
ans =
logical
1
>> Fs = F * 2;
>> sF = 2 * F; % sF == Fs, i.e. the Faust F times 2.
>> matfaust.isFaust(Fs) && matfaust.isFaust(sF)
ans =
logical
1

See also Faust.Faust, Faust.rcg, Faust.ctranspose, Faust.complex

◆ nbytes()

function matfaust::Faust::nbytes ( )

Gives the memory size of the Faust in bytes.

Return values
nFaust size in bytes.

Example:

% in a matlab terminal
>> F = matfaust.rand(1024, 1024, 'num_factors', 6, 'seed', 42);
>> nbytes(F)
ans =
393240
>> F = matfaust.rand(1024, 1024, 'num_factors', 6, 'fac_type', 'dense');
>> nbytes(F)
ans =
50331648

◆ nnz_sum()

function matfaust::Faust::nnz_sum ( )

The total number of nonzeros in the factors of F.

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

Usage

    nz = nnz_sum(F)

Parameters
Fthe Faust object.
Return values
nzthe number of non-zeros.

Example:

>> F = matfaust.rand(1024, 1024, 'num_factors', 6)
F =
Faust size 1024x1024, density 0.0292969, nnz_sum 30720, 6 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
- FACTOR 5 (double) SPARSE, size 1024x1024, density 0.00488281, nnz 5120
>> nnz_sum(F)
ans =
30720

See also Faust.rcg, Faust.density.

◆ norm()

function matfaust::Faust::norm ( ,
varargin   
)

The matrix norm of F.

This function overloads a Matlab built-in function.

Several types of norm are available: 1-norm, 2-norm, inf-norm and Frobenius norm.

The norm of F is equal to the norm of full(F).

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 full(F). Please pay attention to the full_array (and batch_size) arguments for a better understanding. Usage

    n = norm(F, 2) the 2-norm or maximum singular value of F: approximately norm(full(F),2) == max(svd(full(F))).

    n = norm(F, 2, 'threshold', 0.001, 'max_num_its', 1000).

    n = norm(F) the same as norm(F, 2).

    n = norm(F, 1) the 1-norm of F: norm(full(F), 1) == max(sum(abs(full(F))))

    n = norm(F, inf) the inf-norm of F: norm(full(F), inf) == max(sum(abs(full(F)')))

    n = norm(F, 'fro') the Frobenius norm of F: norm(full(F), 'fro').

Parameters
Fthe Faust object.
p(optional) the norm order or type. Respectively 1, 2 or inf for the 1-norm, 2-norm and inf-norm or 'fro' for the Frobenius norm (by default the 2-norm is computed).
'threshold',real(optional) power iteration algorithm threshold (default to .001). Used only for norm(2). It's passed in a key-value pair fashion: 'threshold', .001
'max_num_its',int(optional) maximum number of iterations for power iteration algorithm. Used only for norm(2). It's passed in a key-value pair fashion: 'max_num_its', 1000.
'full_array',bool(optional) 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(optional) 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).
Return values
nthe norm (real).

Example

% in a matlab terminal
>> F = matfaust.rand(50, 100, 'num_factors', 2, 'dim_sizes', [50, 100], 'density', .5)
>> norm(F)
ans =
436.4094
>> norm(F,2)
ans =
436.4094
>> norm(F, 'fro')
ans =
442.3946
>> norm(F, inf)
ans =
738.0057

See also Faust.normalize

◆ normalize()

function matfaust::Faust::normalize ( ,
varargin   
)

Returns the normalized F.

Usage

    normalize(F,'norm') normalizes columns with vector 2-norm.
    normalize(F,'norm',p) normalizes columns with vector p-norm (see Faust.norm).
    normalize(F,DIM,'norm', 2) normalizes rows (DIM=1) or columns (DIM=2, default) with vector 2-norm.
    normalize(F,DIM,'norm',p) normalizes rows (DIM=1) or columns (DIM=2, default) with vector p-norm.

Examples:

>> F = matfaust.rand(10, 10, 'seed', 42);
>> % normalize columns according to default 2-norm
>> % then test the second column is properly normalized
>> nF2 = normalize(F); % nF2 is a Faust
>> matfaust.isFaust(nF2)
ans =
logical
1
>> norm(full(nF2(:,2)) - full(F(:,2))/norm(F(:,2))) < 1e-15
ans =
logical
1
>> % this time normalize rows using 1-norm and test the third row
>> nF1 = normalize(F, 1, 'norm', 1);
>> norm(full(nF1(3, :)) - full(F(3,:))/norm(F(3,:), 1)) < 1e-15
ans =
logical
1
>> % and the same except with inf-norm
>> nFinf = normalize(F, 1, 'norm', inf);
>> norm(full(nFinf(3, :)) - full(F(3,:))/norm(F(3,:), inf)) < 1e-15
ans =
logical
1

See also Faust.norm

◆ numel()

function matfaust::Faust::numel ( )

The number of elements in F.

It's equivalent to prod(size(F)).

Usage

    n = numel(F)

Parameters
Fthe Faust object.
Return values
nthe number of elements of the Faust.

Example

>> F = matfaust.rand(5, 10);
>> n = numel(F)
n =
50

See also Faust.size

◆ numfactors()

function matfaust::Faust::numfactors ( )

The number of factors of F.

Note
using length(F) is shorter!

Usage

    A = numfactors(F)

Parameters
Fthe Faust object.
Return values
num_factorsthe number of factors.

Example

>> F = matfaust.rand(5, 10, 'num_factors', 6);
>> nf = numfactors(F)
nf =
6

See also Faust.factors

◆ opt_butterfly()

static function matfaust::Faust::opt_butterfly ( )
static

Optimizes any Faust composed of butterfly factors.

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

Parameters
FThe Faust to optimize. If the factors of F are not set according to a butterfly structure, the result is not defined.
Return values
Gthe optimized Faust.

Example:

% in a matlab terminal
>> F = matfaust.dft(32)
F =
Faust size 32x32, density 0.34375, nnz_sum 352, 6 factor(s):
- FACTOR 0 (complex) SPARSE, size 32x32, density 0.0625, nnz 64
- FACTOR 1 (complex) SPARSE, size 32x32, density 0.0625, nnz 64
- FACTOR 2 (complex) SPARSE, size 32x32, density 0.0625, nnz 64
- FACTOR 3 (complex) SPARSE, size 32x32, density 0.0625, nnz 64
- FACTOR 4 (complex) SPARSE, size 32x32, density 0.0625, nnz 64
- FACTOR 5 (complex) SPARSE, size 32x32, density 0.03125, nnz 32
oF =
Faust size 32x32, density 0.34375, nnz_sum 352, 6 factor(s):
- FACTOR 0 (complex) Butterfly, size 32x32, density 0.0625, nnz 64
- FACTOR 1 (complex) Butterfly, size 32x32, density 0.0625, nnz 64
- FACTOR 2 (complex) Butterfly, size 32x32, density 0.0625, nnz 64
- FACTOR 3 (complex) Butterfly, size 32x32, density 0.0625, nnz 64
- FACTOR 4 (complex) Butterfly, size 32x32, density 0.0625, nnz 64
- FACTOR 5 (complex) PERM, size 32x32, density 0.03125, nnz 32

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

◆ optimize()

function matfaust::Faust::optimize ( ,
varargin   
)

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

Usage
    OF = optimize_time(F) returns an optimized Faust object.
    OF = optimize_time(F, 'transp', true) see Faust.optimize_time

Parameters
'transp',bool(optional) default to false. if true optimize the Faust according to its transpose.
Return values
OFThe optimized Faust.
Note
this function is still experimental, you might use manually Faust.optimize_time, Faust.optimize_memory or Faust.pruneout to be more specific about the optimization to proceed.

Example:

>> F = matfaust.rand(2048, 2048, 'dim_sizes', [1, 1024], 'num_factors', 32, 'fac_type', 'mixed', 'seed', 42);
>> nbytes(F)
ans =
14056612
>> oF = optimize(F);
>> nbytes(oF)
ans =
975904
>> % oF is smaller in memory and should be Faster too
>> tic; full(F); tF = toc;
>> tic; full(oF); toF = toc;
>> toF < tF
ans =
logical
1
>> v = rand(size(F, 1), 1);
>> tic; F*v; tF = toc;
>> tic; oF*v; toF = toc;
>> toF < tF
ans =
logical
1

See also Faust.optimize_time, Faust.optimize_memory, Faust.pruneout, Lives Script about Faust optimizations

◆ optimize_memory()

function matfaust::Faust::optimize_memory ( )

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

Return values
OFThe optimized Faust.

Example:

>> F = matfaust.rand(1024, 1024, 'dim_sizes', [1, 1024], 'num_factors', 32, 'fac_type', 'mixed', 'seed', 42);
>> nbytes(F)
ans =
13991076
>> oF = optimize_memory(F);
>> nbytes(oF)
ans =
972672
>> % oF is smaller in memory

See also Faust.optimize

◆ optimize_time()

function matfaust::Faust::optimize_time ( ,
varargin   
)

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 matfaust.FaustMulMode. Although depending on the package you installed and the capability of your hardware the methods based on Torch library can be used.

Usage
    OF = optimize_time(F) returns a new object with the best product method enabled.
    OF = optimize_time(F, 'inplace', true) modifies directly F instead of creating a new Faust (i.e. OF references the same object as F).
    OF = optimize_time(F, 'transp', true) tries to optimize the Faust transpose instead of the Faust itself.
    OF = optimize_time(F, 'nsamples', 10) benchmarks the product methods by computing 10 products instead of one by default.

Parameters
Fthe Faust object.
'inplace',booldefault to false. If true the current Faust is modified directly.
'transp',booldefault to false. If true optimize the Faust according to its transpose.
'nsamples',intdefault to 1.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 performances are similar). By default, only one product is computed to evaluate the method.
'mat',matrixUse this argument to run the benchmark on the Faust multiplication by the matrix mat instead of Faust.full(). Note that mat must be of the same scalar type as F.
Return values
OFThe optimized Faust.

Example:

>> F = matfaust.rand(2048, 2048, 'dim_sizes', [1, 2048], 'num_factors', 32, 'fac_type', 'dense', 'seed', 42);
>> oF = optimize_time(F); % doctest: +ELLIPSIS
%> >> % oF should be Faster too
>> tF = timeit(@() full(F));
>> toF = timeit(@() full(oF));
>> toF < tF % doctest: +ELLIPSIS
ans =
%>
1

See also Faust.optimize, matfaust.FaustMulMode

◆ pinv()

function matfaust::Faust::pinv ( )

Pseudoinverse matrix of full(F).

This function overloads a Matlab built-in function.

Usage

    X = PINV(F) produces the matrix X of the same dimensions as F' so that F*(F'*X')' == full(F) (or approximately).

Warning
this function makes a call to Faust.full.

Example:

>> F = matfaust.rand(128, 32);
>> M = full(F);
>> norm(pinv(M)-pinv(F))/norm(pinv(M))
ans =
0
>> norm(M-(F*(F'*pinv(F)')))/norm(M) < 1e-12
ans =
logical
1

See also Faust.mldivide, pinv Matlab built-in, matfaust.fact.pinvtj

◆ plus()

function matfaust::Faust::plus ( varargin  )

Plus.

This function overloads a Matlab built-in function.

Usage

    plus(F,G) or F+G adds two Faust together, sizes must be compatible.
    plus(F,A) or F+A adds a Faust and a matrix A, sizes must be compatible.
    plus(F,s) or F+s, with s being a scalar, such that full(F+s) == full(F)+s.

Parameters
F(first arg.) The Faust object.
G,A,s,…(2nd to n-th args) The variables to add to F; Fausts, matrices or scalars.
Return values
Sthe sum as a Faust object.

Examples:

>> F = matfaust.rand(10, 12);
>> G = matfaust.rand(10, 12);
>> FpG = F+G; % is a Faust
>> matfaust.isFaust(FpG)
ans =
logical
1
>> norm(full(FpG) - full(F) - full(G)) < 1e-12
ans =
logical
1
>> Fp2 = F+2; % is a Faust
>> norm(full(Fp2) - full(F) - 2) < 1e-12
ans =
logical
1
>> FpGp2 = F+G+2; % is a Faust
>> norm(full(FpGp2) - full(F) - full(G) - 2) < 1e-12
ans =
logical
1
>> H = F+G+F+2+F; % is a Faust
>> norm(full(H) - 3*full(F) - full(G) - 2) < 1e-12
ans =
logical
1

See also Faust.minus

◆ power_iteration()

function matfaust::Faust::power_iteration ( ,
varargin   
)

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
'threshold',real(optional) the precision required on the eigenvalue. Default value is 1e-3.
'maxiter',integer(optional) the number of iterations above what the algorithm will stop anyway. Default value is 100.
Return values
lambdathe greatest eigenvalue approximate.

Example:

%in a matlab terminal
>> F = matfaust.rand(8, 5, 'seed', 42);
>> l1 = power_iteration(F*F');
>> % the 2-norm of F is its greatest singular value
>> % which is also the square root of the greatest eigenvalue of F*F'
>> abs(sqrt(l1) - norm(full(F), 2)) < 1e-12
ans =
logical
1

◆ pruneout()

function matfaust::Faust::pruneout ( ,
varargin   
)

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

Parameters
Fthe Faust to optimize.
'thres',int(optional) the threshold of number of nonzeros under what the rows/columns are removed.
Return values
GThe optimized Faust.

Example:

>> F = matfaust.rand(1024, 1024, 'dim_sizes', [1, 1024], 'num_factors', 32, 'fac_type', 'mixed', 'seed', 42);
>> pF = pruneout(F);
>> nbytes(pF) < nbytes(F)
ans =
logical
1
>> % pF is smaller in memory and should be Faster too
>> tF = timeit(@() full(F));
>> tpF = timeit(@() full(pF));
>> tpF < tF
ans =
logical
1

See also Faust.optimize, Faust.nbytes

◆ rcg()

function matfaust::Faust::rcg ( )

The Relative Complexity Gain of F.

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
rcg(F) == 1/density(F)

Usage

    gain = rcg(F)

Parameters
Fthe Faust object.
Return values
gainthe RCG value (real).

Example

>> F = matfaust.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
>> rcg(F)
ans =
40.9600
>> numel(F)/nnz_sum(F)
ans =
40.9600

See also Faust.density, Faust.nnz_sum, Faust.size.

◆ real()

function matfaust::Faust::real ( )

Returns the real part of the Faust F.

Example:

% in a matlab terminal
>> F = matfaust.rand(10, 10, 'field', 'complex')
F =
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
>> rF = real(F)
rF =
Faust size 10x10, density 8, nnz_sum 800, 5 factor(s):
- FACTOR 0 (double) SPARSE, size 10x20, density 0.5, nnz 100
- FACTOR 1 (double) SPARSE, size 20x20, density 0.5, nnz 200
- FACTOR 2 (double) SPARSE, size 20x20, density 0.5, nnz 200
- FACTOR 3 (double) SPARSE, size 20x20, density 0.5, nnz 200
- FACTOR 4 (double) SPARSE, size 20x10, density 0.5, nnz 100
>> norm(full(rF) - real(full(F))) < 1e-12
ans =
logical
1

See also Faust.imag

◆ replace()

function matfaust::Faust::replace ( ,
,
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,intthe factor index to replace.
new_factor,matrixthe factor replacing the i-th factor of F.
Return values
Ga copy of F with the i-th factor replaced by new_factor.

Example

>> F = matfaust.rand(10, 10, 'num_factors', 5);
>> S = sprand(10, 10, .2);
>> G = replace(F, 3, S);
>> all(all(full(left(F, 2) * matfaust.Faust(S) * right(F, 4)) == full(G)))
ans =
logical
1

See also Faust.factors, Faust.left, Faust.right, Faust.insert

◆ right()

function matfaust::Faust::right ( ,
,
varargin   
)

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

Parameters
Fthe Faust object.
ithe bound factor index.
'as_faust',logical(optional): true to return a Faust even if a single factor is asked (length(ids) == 1), otherwise (as_faust == false) and a dense or sparse matrix is returned.
Return values
aFaust if the number of factors to be returned is greater than 1, an array or a sparse matrix otherwise.

Example

>> F = matfaust.rand(7,7, 'dim_sizes', [7, 10], 'seed', 42)
F =
Faust size 7x7, density 3.97959, nnz_sum 195, 5 factor(s):
- FACTOR 0 (double) SPARSE, size 7x9, density 0.555556, nnz 35
- FACTOR 1 (double) SPARSE, size 9x8, density 0.625, nnz 45
- FACTOR 2 (double) SPARSE, size 8x8, density 0.625, nnz 40
- FACTOR 3 (double) SPARSE, size 8x7, density 0.714286, nnz 40
- FACTOR 4 (double) SPARSE, size 7x7, density 0.714286, nnz 35
>> RF = right(F, 2)
RF =
Faust size 9x7, density 2.53968, nnz_sum 160, 4 factor(s):
- FACTOR 0 (double) SPARSE, size 9x8, density 0.625, nnz 45
- FACTOR 1 (double) SPARSE, size 8x8, density 0.625, nnz 40
- FACTOR 2 (double) SPARSE, size 8x7, density 0.714286, nnz 40
- FACTOR 3 (double) SPARSE, size 7x7, density 0.714286, nnz 35

See also Faust.factors, Faust.left

◆ save()

function matfaust::Faust::save ( ,
filepath   
)

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 full(F). See Faust.nbytes for a precise size.

Usage

    save(F, filepath)

Parameters
Fthe Faust object.
filepaththe path for saving the Faust. The filename should end with .mat and it must be a character array.

Example

>> import matfaust.Faust
>> F = matfaust.rand(5, 10);
>> save(F, 'F.mat')
>> G = Faust('F.mat');
>> all(all(full(F) == full(G)))
ans =
logical
1

See also Faust.Faust, Faust.rcg, Faust.load, Faust.load_native

◆ single()

function matfaust::Faust::single ( )

Returns a new Faust whose class is single. It allows to convert all the factors to single precision.

Usage

    sF = single(F)

Return values
sFthe single class/precision Faust.

Example:

>> dF = matfaust.rand(5, 5, 'field', 'real', 'class', 'double')
dF =
Faust size 5x5, density 5, nnz_sum 125, 5 factor(s):
- FACTOR 0 (double) SPARSE, size 5x5, density 1, nnz 25
- FACTOR 1 (double) SPARSE, size 5x5, density 1, nnz 25
- FACTOR 2 (double) SPARSE, size 5x5, density 1, nnz 25
- FACTOR 3 (double) SPARSE, size 5x5, density 1, nnz 25
- FACTOR 4 (double) SPARSE, size 5x5, density 1, nnz 25
>> sF = single(dF)
sF =
Faust size 5x5, density 5, nnz_sum 125, 5 factor(s):
- FACTOR 0 (float) SPARSE, size 5x5, density 1, nnz 25
- FACTOR 1 (float) SPARSE, size 5x5, density 1, nnz 25
- FACTOR 2 (float) SPARSE, size 5x5, density 1, nnz 25
- FACTOR 3 (float) SPARSE, size 5x5, density 1, nnz 25
- FACTOR 4 (float) SPARSE, size 5x5, density 1, nnz 25

See also Faust.class, Faust.double

◆ size()

function matfaust::Faust::size ( ,
varargin   
)

The size of F.

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

Usage

    [nrows, ncols] = size(F)
    n = size(F,dim) with n being the size of the dim-th dimension of F.
                 In other words n == nrows if dim == 1, n == ncols if dim == 2.

Parameters
Fthe Faust object.
dim(optional) the index of the dimension to get the size of.

Example

>> F = matfaust.rand(5, 10);
>> [nrows, ncols] = size(F)
nrows =
5
ncols =
10
>> nrows = size(F, 1)
nrows =
5
>> ncols = size(F, 2)
ncols =
10

See also Faust.nnz_sum, Faust.numel

◆ subsref()

function matfaust::Faust::subsref ( ,
 
)

Subscripted reference of a Faust.

The function returns a Faust representing a submatrix of full(F) or a scalar element if that Faust can be reduced to a single element.

This function overloads a Matlab built-in.

Warning
  • It is not advised to use this function as an element accessor (e.g. F(1,1)) 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.

Usage

    G = F(I,J) the Faust representing the sub-matrix of full(F) defined by the subscript vectors I and J (see examples below).

Parameters
Fthe Faust object.
Sthe structure defining the Faust to extract like a submatrix (it's not supposed to be used directly ; usage and examples show how subscript should be used).
Return values
TheFaust object requested or just the corresponding scalar if that Faust has a size equal to [1,1].

Example

>> F = matfaust.rand(50, 100);
>> fF = full(F);
>> i = randi(min(size(F)), 1, 2);
>> i1 = i(1); i2 = i(2);
>> F(i1, i2); % is the scalar element located
>> % at row i1, column i2 of the F dense matrix
>> F(i1,i2) == fF(i1,i2)
ans =
logical
1
>> F(:,i2); % full column i2
>> all(all(full(F(:,i2)) == fF(:,i2)))
ans =
logical
1
>> F(3:4,2:5); % from row 3 to line 4, each row containing only elements from column 2 to 5.
>> all(all(full(F(3:4,2:5)) == fF(3:4,2:5)))
ans =
logical
1
>> F(1:end,5:end-1); % from row 1 to end row, each one containing only elements from column 5 to column before the last one.
>> all(all(full(F(1:end,5:end-1)) == fF(1:end,5:end-1)))
ans =
logical
1
>> F(1:2:end,:); % every row of odd index
>> all(all(full(F(1:2:end,:)) == fF(1:2:end,:)))
ans =
logical
1
>> F(end:-2:1,:); % starts from the last row and goes backward to take one in two rows until the first one (reversing row order of F)
>> all(all(full(F(end:-2:1,:)) == fF(end:-2:1,:)))
ans =
logical
1
>> F([1,18,2],:); % takes in this order the rows 1, 18 and 2.
>> all(all(full(F([1,18,2],:)) == fF([1,18,2],:)))
ans =
logical
1
>> F(:,[1,18,2]); % takes in this order the columns 1, 18 and 2
>> all(all(full(F(:,[1,18,2])) == fF(:,[1,18,2])))
ans =
logical
1
>> F([1,18,2], [1,2]); % takes the rows 1, 18 and 2 but keeps only columns 1 and 2 in these rows.
>> all(all(full(F([1,18,2], [1,2])) == fF([1,18,2], [1,2])))
ans =
logical
1

See also Faust.end.

◆ transpose()

function matfaust::Faust::transpose ( )

The transpose of F.

This function overloads a Matlab built-in function/operator.

Usage

    F_trans = transpose(F)
    F_trans = F.'

Parameters
Fthe Faust object.
Return values
F_transa Faust object implementing the transpose of full(F) such that: full(F_trans) == full(F).' == transpose(full(F)).

Example

>> F = matfaust.rand(8, 2);
>> F_trans = F.';
>> size(F_trans)
ans =
2 8
>> % is equivalent to
>> F_trans = transpose(F);
>> size(F_trans)
ans =
2 8

See also Faust.conj, Faust.ctranspose

◆ uminus()

function matfaust::Faust::uminus ( )

Returns - F.


This function overloads a Matlab built-in function.

Example:

% in a matlab terminal
>> F = matfaust.rand(10, 10);
>> mF = -F;
>> all(all(full(mF) == - full(F)))
ans =
logical
1

See also Faust.sum, Faust.mtimes

◆ uplus()

function matfaust::Faust::uplus ( )

Returns + F.


This function overloads a Matlab built-in function.

Example:

% in a matlab terminal
>> F = matfaust.rand(10, 10);
>> pF = +F;
>> all(all(full(pF) == + full(F)))
ans =
logical
1

See also Faust.plus

◆ vertcat()

function matfaust::Faust::vertcat ( varargin  )

Vertical concatenation [F;A;B;…] where F is a Faust and A, B, … are Faust objects or sparse/full matrices.

It's equivalent to cat(1, F, A, B, …).

This function overloads a Matlab built-in function.

Usage

    C=VERTCAT(F,A) concatenates the Faust F and A vertically. A is a Faust or a sparse/full matrix. The result is the Faust C. VERTCAT(F,A) is the same as [F;G].
    C=VERTCAT(F, A, B,…) concatenates vertically F to A, B etc. VERTCAT(F, A, B,…) is the same as [F;A;B;…].

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

See also Faust.horzcat, Faust.cat.

Member Data Documentation

◆ dev

Property matfaust::Faust::dev
protected

◆ dtype

Property matfaust::Faust::dtype
protected

◆ is_real

Property matfaust::Faust::is_real
protected

◆ matrix

Property matfaust::Faust::matrix
protected

Underlying Faust native object handle.


The documentation for this class was generated from the following file:
matfaust::Faust::length
function length(F)
The number of factors of F.
matfaust::Faust::full
function full(F)
The full matrix implemented by F.
matfaust::Faust::disp
function disp(F)
Displays information about F.
matfaust::Faust::device
function device(F)
Returns the Faust device ('cpu' or 'gpu').
matfaust::Faust::opt_butterfly
static function opt_butterfly(F)
Optimizes any Faust composed of butterfly factors.
matfaust::Faust::isdense
function isdense(F)
Returns true if F factors are all dense matrices/arrays false otherwise.
matfaust::Faust::numel
function numel(F)
The number of elements in F.
matfaust::Faust::cat
function cat(varargin)
Concatenates F with n Faust objects or full/sparse matrices.
matfaust::Faust::optimize_time
function optimize_time(F, varargin)
Returns a Faust configured with the quickest Faust-matrix multiplication mode (benchmark ran on the f...
matfaust::Faust::pinv
function pinv(F)
Pseudoinverse matrix of full(F).
matfaust::Faust::save
function save(F, filepath)
Saves the Faust F into a file.
matfaust::Faust
FAuST Matlab wrapper main class for using multi-layer sparse transforms.
Definition: Faust.m:118
matfaust::Faust::imag
function imag(F)
Returns the imaginary part of the Faust F.
matfaust::Faust::conj
function conj(F)
The complex conjugate of F.
matfaust::Faust::double
function double(F)
Returns a new Faust whose class is double. It allows to convert all the factors to double precision.
matfaust::Faust::right
function right(F, i, varargin)
Returns the right hand side factors of F from index i to end (in 1-base index).
matfaust::Faust::load
static function load(filepath, varargin)
Loads a Faust from a .mat file.
matfaust
The FAuST Matlab Wrapper
Definition: bsl.m:1
matfaust::Faust::get_handle
function get_handle(F)
Returns the Faust F C++ native object pointer as an integer.
matfaust::Faust::numfactors
function numfactors(F)
The number of factors of F.
matfaust::Faust::density
function density(F)
The density of F such that nnz_sum(F) == density(F) * numel(F).
matfaust::Faust::isreal
function isreal(F)
Indicates if F is a real Faust or a complex Faust.
matfaust::Faust::isFaust
static function isFaust(obj)
Returns true if obj is a Faust object, false otherwise.
matfaust::Faust::ctranspose
function ctranspose(F)
The conjugate transpose of F.
matfaust::Faust::load_native
static function load_native(filepath)
Loads a Faust from a .mat file (native C++ version).
matfaust::Faust::power_iteration
function power_iteration(F, varargin)
Performs the power iteration algorithm to compute the greatest eigenvalue of the Faust.
matfaust::Faust::end
function end(F, k, n)
The last index when slicing or indexing a Faust.
matfaust::Faust::replace
function replace(F, i, new_factor)
Replaces the factor of index i by new_factor in a new Faust copy of F.
matfaust::Faust::insert
function insert(F, i, new_factor)
Inserts new_factor at index i in a new Faust copy of F.
matfaust::Faust::pruneout
function pruneout(F, varargin)
Returns a Faust optimized by removing useless zero rows and columns as many times as needed.
matfaust::Faust::clone
function clone(F, varargin)
Clones the Faust (in a new memory space).
matfaust::Faust::mtimes
function mtimes(F, A)
Multiplies the Faust F by A which is a full matrix, a Faust object or a scalar.
matfaust::Faust::nnz_sum
function nnz_sum(F)
The total number of nonzeros in the factors of F.
matfaust::Faust::transpose
function transpose(F)
The transpose of F.
matfaust::rand
function rand(M, N, varargin)
Generates a random Faust.
matfaust::Faust::real
function real(F)
Returns the real part of the Faust F.
matfaust::Faust::size
function size(F, varargin)
The size of F.
matfaust::Faust::complex
function complex(F)
Converts F to a new complex Faust.
matfaust::Faust::normalize
function normalize(F, varargin)
Returns the normalized F.
matfaust::Faust::rcg
function rcg(F)
The Relative Complexity Gain of F.
matfaust::Faust::left
function left(F, i, varargin)
Returns the left hand side factors of F from index 1 to i included (in 1-base index).
matfaust::Faust::matrix
Property matrix
Underlying Faust native object handle.
Definition: Faust.m:123
matfaust::Faust::issparse
function issparse(F, varargin)
Returns true if F factors are all sparse matrices false otherwise.
matfaust::Faust::optimize_memory
function optimize_memory(F)
Optimizes a Faust by changing the storage format of each factor in order to optimize the memory size.
matfaust::Faust::norm
function norm(F, varargin)
The matrix norm of F.
matfaust::Faust::Faust
function Faust(varargin)
Creates a Faust from a list of factors or alternatively from a file.
matfaust::Faust::nbytes
function nbytes(F)
Gives the memory size of the Faust in bytes.
matfaust::Faust::imagesc
function imagesc(F, varargin)
Displays image of F full matrix and its factors.
matfaust::Faust::factors
function factors(F, ids, varargin)
Returns the i-th factor or a new Faust composed of F factors whose indices are listed in indices.
matfaust::Faust::optimize
function optimize(F, varargin)
Returns a Faust optimized with Faust.pruneout, Faust.optimize_memory and Faust.optimize_time.
matfaust::Faust::single
function single(F)
Returns a new Faust whose class is single. It allows to convert all the factors to single precision.