Reference

The following provides detailled documentation about types and methods provided by the LazyAlgebra package. This information is also available from the REPL by typing ? followed by the name of a method or a type.

Methods for linear mappings

LazyAlgebra.nrowsFunction
nrows(A)

yields the equivalent number of rows of the linear operator A. Not all operators extend this method.

In the implemented generalization of linear operators, the equivalent number of rows is the number of element of the result of applying the operator be it single- or multi-dimensional.

source
LazyAlgebra.ncolsFunction
ncols(A)

yields the equivalent number of columns of the linear operator A. Not all operators extend this method.

In the implemented generalization of linear operators, the equivalent number of columns is the number of element of an argument of the operator be it single- or multi-dimensional.

source
LazyAlgebra.row_sizeFunction
row_size(A)

yields the dimensions of the result of applying the linear operator A, this is equivalent to output_size(A). Not all operators extend this method.

source
LazyAlgebra.col_sizeFunction
col_size(A)

yields the dimensions of the argument of the linear operator A, this is equivalent to input_size(A). Not all operators extend this method.

source

Sparse operators

Types and compressed storage formats

LazyAlgebra.SparseOperators.SparseOperatorType
SparseOperator{T,M,N}

is the abstract type inherited by sparse operator types. Parameter T is the type of the elements. Parameters M and N are the number of dimensions of the rows and of the columns respectively. Sparse operators are a generalization of sparse matrices in the sense that they implement linear mappings which can be applied to N-dimensonal arguments to produce M-dimensional results (as explained below). See GeneralMatrix for a similar generalization but for dense matrices.

See CompressedSparseOperator for usage of sparse operators implementing compressed storage formats.

source
LazyAlgebra.SparseOperators.CompressedSparseOperatorType
CompressedSparseOperator{F,T,M,N}

is an abstract sub-type of SparseOperator{T,M,N} and is inherited by the concrete types implementing sparse operators with compressed storage in format F.

Format F is specificed as a symbol and can be:

  • :COO for Compressed Sparse Coordinate storage format. This format is not the most efficient, it is mostly used as an intermediate for building a sparse operator in one of the following formats.

  • :CSC for Compressed Sparse Column storage format. This format is very efficient for applying the adjoint of the sparse operator.

  • :CSR for Compressed Sparse Row storage format. This format is very efficient for directly applying the sparse operator.

To construct (or convert to) a sparse operator with compressed storage format F, you can call:

CompressedSparseOperator{F}(args...; kwds...)
CompressedSparseOperator{F,T}(args...; kwds...)
CompressedSparseOperator{F,T,M}(args...; kwds...)
CompressedSparseOperator{F,T,M,N}(args...; kwds...)

where given parameters T, M and N, arguments args... and optional keywords kwds... will be passed to the concrete constructor SparseOperatorCOO, SparseOperatorCSC or SparseOperatorCSR corresponding to the format F.

It is possible to use a compressed sparse operator A as an iterator:

for (Aij,i,j) in A # simple but slow for CSR and CSC
    ...
end

to retrieve the values Aij and respective row i and column j indices for all the entries stored in A. It is however more efficient to access them according to their storage order which depends on the compressed format.

  • If A is in CSC format:

    using LazyAlgebra.SparseMethods
    for j in each_col(A)        # loop over column index
        for k in each_off(A, j) # loop over structural non-zeros in this column
            i   = get_row(A, k) # get row index of entry
            Aij = get_val(A, k) # get value of entry
         end
    end
  • If A is in CSR format:

    using LazyAlgebra.SparseMethods
    for i in each_row(A)        # loop over row index
        for k in each_off(A, i) # loop over structural non-zeros in this row
            j   = get_col(A, k) # get column index of entry
            Aij = get_val(A, k) # get value of entry
         end
    end
  • If A is in COO format:

    using LazyAlgebra.SparseMethods
    for k in each_off(A)
         i   = get_row(A, k) # get row index of entry
         j   = get_col(A, k) # get column index of entry
         Aij = get_val(A, k) # get value of entry
    end

The low-level methods each_row, each_col, each_off, get_row, get_col and get_val are not automatically exported by LazyAlgebra, this is the purpose of the statement using LazyAlgebra.SparseMethods.

source
LazyAlgebra.SparseOperators.SparseOperatorCOOType

Sparse operators in Compressed Sparse Coordinate (COO) format store the significant entries in no particular order, as a vector of values, a vector of linear row indices and a vector of linear column indices. It is even possible to have repeated entries. This format is very useful to build a sparse linear operator. It can be converted to a more efficient format like Compressed Sparse Column (CSC) or Compressed Sparse Row (CSR) for fast application of the sparse linear mapping or of its adjoint.

A sparse operator in COO storage format can be constructed by providing all necessary information:

SparseOperatorCOO(vals, rows, cols, rowsiz, colsiz)

where vals is the vector of values of the sparse entries, rows and cols are integer valued vectors with the linear row and column indices of the sparse entries, rowsiz and colsiz are the sizes of the row and column dimensions. The entries values and respective linear row and column indices of the k-th sparse entry are given by vals[k], rows[k] and cols[k]. For efficiency reasons, sparse operators are currently limited to fast arrays because they can be indexed linearly with no loss of performances. If vals, rows and/or cols are not fast arrays, they will be automatically converted to linearly indexed arrays.

A sparse operator in COO storage format can be directly constructed from a 2-dimensional Julia array A:

SparseOperatorCOO(A, sel = (v,i,j) -> (v != zero(v)))

where optional argument sel is a selector function which is called as sel(v,i,j) with v, i and j the value, the row and the column linear indices for each entries of A and which is assumed to yield true for the entries of A to be selected in the sparse structure and false for the entries of A to discard. The default selector is such that all non-zeros of A are selected.

The element type, say T, for the sparse coefficients can be imposed by rewriting the above examples as:

SparseOperatorCOO{T}(args...)

A sparse operator in COO storage format implementing generalized matrix-vector multiplication can also be directly constructed from a L-dimensional Julia array (with L ≥ 2) A by:

SparseOperatorCOO{T,M}(A[, sel])

with M the number of leading dimensions of A corresponding to the rows of the operator, the trailing N = L - M dimensions being assumed to correspond to the columns of the operator. These dimensions are the size of, respectively, the output and the input arrays when applying the operator. The parameter N may be specified (although it can be automatically determined):

SparseOperatorCOO{T,M,N}(A[, sel])

provided the equality M + N = ndims(A) holds.

A last parameter V can be specified for the type of the vector to store the coefficients of the sparse operator:

SparseOperatorCOO{T,M,N,V}(args...)

provided V implements standard linear indexing. The default is to take V = Vector{T}. As a special case, you can choose a uniform boolean vector from the StructuredArrays package to store the sparse coefficients:

SparseOperatorCOO{T,M,N,UniformVector{Bool}}(args...)

to get a compressed sparse operator in COO format whose values are an immutable uniform vector of true values requiring no storage. This is useful to only store the sparse structure of the operator, that is the indices in COO format of the sparse coefficients not their values.

The SparseOperatorCOO constructor can also be used to convert a sparse operator in another storage format into the COO format. In that case, parameter T may also be specified to convert the type of the sparse coefficients.

source
LazyAlgebra.SparseOperators.SparseOperatorCSCType

Sparse operators in Compressed Sparse Column (CSC) format store the significant entries in a column-wise order, as a vector of values, a vector of corresponding linear row indices and a vector of offsets indicating, for each column, the range of indices in the vectors of values and of row indices. This storage format is very suitable for fast application of the operator, notably its adjoint.

A sparse operator in CSC storage format can be constructed by providing all necessary information:

SparseOperatorCSC(vals, rows, offs, rowsiz, colsiz)

where vals is the vector of values of the sparse entries, rows is an integer valued vector of the linear row indices of the sparse entries, offs is a column-wise table of offsets in these arrays, rowsiz and colsiz are the sizes of the row and column dimensions. The entries values and respective linear row indices of the j-th column are given by vals[k] and rows[k] with k ∈ offs[j]+1:offs[j+1]. The linear column index j is in the range 1:n where n = prod(colsiz) is the equivalent number of columns. For efficiency reasons, sparse operators are currently limited to fast arrays because they can be indexed linearly with no loss of performances. If vals, rows and/or offs are not fast arrays, they will be automatically converted to linearly indexed arrays.

A sparse operator in CSC storage format can be directly constructed from a 2-dimensional Julia array A:

SparseOperatorCSC(A, sel = (v,i,j) -> (v != zero(v)))

where optional argument sel is a selector function which is called as sel(v,i,j) with v, i and j the value, the row and the column linear indices for each entries of A and which is assumed to yield true for the entries of A to be selected in the sparse structure and false for the entries of A to discard. The default selector is such that all non-zeros of A are selected.

The element type, say T, for the sparse coefficients can be imposed by rewriting the above examples as:

SparseOperatorCSC{T}(args...)

A sparse operator in CSC storage format implementing generalized matrix-vector multiplication can also be directly constructed from a L-dimensional Julia array (with L ≥ 2) A by:

SparseOperatorCSC{T,M}(A[, sel])

with M the number of leading dimensions of A corresponding to the rows of the operator, the trailing N = L - M dimensions being assumed to correspond to the columns of the operator. These dimensions are the size of, respectively, the output and the input arrays when applying the operator. The parameter N may be specified (although it can be automatically determined):

SparseOperatorCSC{T,M,N}(A[, sel])

provided the equality M + N = ndims(A) holds.

A last parameter V can be specified for the type of the vector to store the coefficients of the sparse operator:

SparseOperatorCSC{T,M,N,V}(args...)

provided V implements standard linear indexing. The default is to take V = Vector{T}. As a special case, you can choose a uniform boolean vector from the StructuredArrays package to store the sparse coefficients:

SparseOperatorCSC{T,M,N,UniformVector{Bool}}(args...)

to get a compressed sparse operator in CSC format whose values are an immutable uniform vector of true values requiring no storage. This is useful to only store the sparse structure of the operator, that is the indices in CSC format of the sparse coefficients not their values.

The SparseOperatorCSC constructor can also be used to convert a sparse operator in another storage format into the CSC format. In that case, parameter T may also be specified to convert the type of the sparse coefficients.

source
LazyAlgebra.SparseOperators.SparseOperatorCSRType

Sparse operators in Compressed Sparse Row (CSR) format store the significant entries in a row-wise order, as a vector of values, a vector of corresponding linear column indices and a vector of offsets indicating, for each row, the range of indices in the vectors of values and of column indices. This storage format is very suitable for fast application of the operator.

A sparse operator in CSR storage format can be constructed by providing all necessary information:

SparseOperatorCSR(vals, cols, offs, rowsiz, colsiz)

where vals is the vector of values of the sparse entries, cols is an integer valued vector of the linear column indices of the sparse entries, offs is a column-wise table of offsets in these arrays, rowsiz and colsiz are the sizes of the row and column dimensions. The entries values and respective linear column indices of the i-th row are given by vals[k] and cols[k] with k ∈ offs[i]+1:offs[i+1]. The linear row index i is in the range 1:m where m = prod(rowsiz) is the equivalent number of rows. For efficiency reasons, sparse operators are currently limited to fast arrays because they can be indexed linearly with no loss of performances. If vals, cols and/or offs are not fast arrays, they will be automatically converted to linearly indexed arrays.

A sparse operator in CSR storage format can be directly constructed from a 2-dimensional Julia array A:

SparseOperatorCSR(A, sel = (v,i,j) -> (v != zero(v)))

where optional argument sel is a selector function which is called as sel(v,i,j) with v, i and j the value, the row and the column linear indices for each entries of A and which is assumed to yield true for the entries of A to be selected in the sparse structure and false for the entries of A to discard. The default selector is such that all non-zeros of A are selected.

The element type, say T, for the sparse coefficients can be imposed by rewriting the above examples as:

SparseOperatorCSR{T}(args...)

A sparse operator in CSR storage format implementing generalized matrix-vector multiplication can also be directly constructed from a L-dimensional Julia array (with L ≥ 2) A by:

SparseOperatorCSR{T,M}(A[, sel])

with M the number of leading dimensions of A corresponding to the rows of the operator, the trailing N = L - M dimensions being assumed to correspond to the columns of the operator. These dimensions are the size of, respectively, the output and the input arrays when applying the operator. The parameter N may be specified (although it can be automatically determined):

SparseOperatorCSR{T,M,N}(A[, sel])

provided the equality M + N = ndims(A) holds.

A last parameter V can be specified for the type of the vector to store the coefficients of the sparse operator:

SparseOperatorCSR{T,M,N,V}(args...)

provided V implements standard linear indexing. The default is to take V = Vector{T}. As a special case, you can choose a uniform boolean vector from the StructuredArrays package to store the sparse coefficients:

SparseOperatorCSR{T,M,N,UniformVector{Bool}}(args...)

to get a compressed sparse operator in CSR format whose values are an immutable uniform vector of true values requiring no storage. This is useful to only store the sparse structure of the operator, that is the indices in CSR format of the sparse coefficients not their values.

The SparseOperatorCSR constructor can also be used to convert a sparse operator in another storage format into the CSR format. In that case, parameter T may also be specified to convert the type of the sparse coefficients.

source

Methods

LazyAlgebra.SparseOperators.unpack!Function
unpack!(A, S; flatten=false) -> A

unpacks the non-zero coefficients of the sparse operator S into the array A and returns A. Keyword flatten specifies whether to only consider the length of A instead of its dimensions. In any cases, A must have as many elements as length(S) and standard linear indexing.

Just call Array(S) to unpack the coefficients of a sparse operator S without providing the destination array.

source

Low-level interface

These methods are provided by using LazyAlgebra.SparseMethods.

LazyAlgebra.SparseOperators.each_rowFunction
each_row(A)

yields an iterator over the linear row indices of the sparse operator A stored in a Compressed Sparse Row (CSR) format, this includes the adjoint of a sparse operator in Compressed Sparse Column (CSC) format.

source
LazyAlgebra.SparseOperators.each_colFunction
each_col(A)

yields an iterator over the linear column indices of the sparse operator A stored in a Compressed Sparse Column (CSC) format, this includes the adjoint of a sparse operator in Compressed Sparse Row (CSR) format.

source
LazyAlgebra.SparseOperators.each_offFunction

For a sparse operator A stored in a Compressed Sparse Coordinate (COO) format, the call:

each_off(A)

yields an iterator over the indices in the arrays of values and of linear row and column indices for the k-th entry of A.


For a sparse operator A stored in a Compressed Sparse Column (CSC) format, the call:

each_off(A, j)

yields an iterator over the indices in the arrays of values and linear row indices for the j-th column of A.


For a sparse operator A stored in a Compressed Sparse Row (CSR) format, the call:

each_off(A, i)

yields an iterator over the indices in the arrays of values and linear column indices for the i-th row of A.

source
LazyAlgebra.SparseOperators.get_rowFunction
get_row(A, k) -> i

yields the linear row index of the k-th entry of the sparse operator A stored in a Compressed Sparse Column (CSC) or Coordinate (COO) formats (this includes adjoint of sparse operators in CSR format).

source
LazyAlgebra.SparseOperators.get_rowsFunction
get_rows(A)

yields the row indices of the entries of the sparse operator A. The returned array may be shared with A, call copy_rows(A) instead if you want to modify the contents of the returned array with no side effects on A.

source
LazyAlgebra.SparseOperators.get_colFunction
get_col(A, k) -> j

yields the linear column index of the k-th entry of the sparse operator A stored in a Compressed Sparse Row (CSR) or Coordinate (COO) formats (this includes adjoint of sparse operators in CSC format).

source
LazyAlgebra.SparseOperators.get_colsFunction
get_cols(A)

yields the column indices of the entries of the sparse operator A. The returned array may be shared with A, call copy_cols(A) instead if you want to modify the contents of the returned array with no side effects on A.

source
LazyAlgebra.SparseOperators.get_valFunction
get_val(A, k) -> v

yields the value of the k-th entry of the sparse operator A stored in a Compressed Sparse Row (CSR), Compressed Sparse Column (CSC) or Coordinate (COO) format.

Argument may also be the adjoint of a sparse operator:

get_val(A', k) -> conj(get_val(A, k))
source
LazyAlgebra.SparseOperators.get_valsFunction
get_vals(A)

yields the array storing the values of the sparse operator A. The returned array is shared with A, call copy_vals(A) instead if you want to modify the contents of the returned array with no side effects on A.

As a convenience, argument may also be the adjoint of a sparse operator:

get_vals(A') -> get_vals(A)

which yields the unmodified values of A, hence the caller has to take the conjugate of these values. The method get_val(A',k) however takes care of conjugating the values.

source
LazyAlgebra.SparseOperators.set_val!Function
set_val!(A, k, v) -> A

assigns v to the value of the k-th entry of the sparse operator A stored in a Compressed Sparse Row (CSR), Compressed Sparse Column (CSC) or Coordinate (COO) format.

source
LazyAlgebra.SparseOperators.get_offsFunction
get_offs(A)

yields the table of offsets of the sparse operator A. Not all operators extend this method.

Warning

The interpretation of offsets depend on the type of A. For instance, assuming offs = get_offs(A), then the index range of the j-th column of a SparseMatrixCSC is offs[j]:(offs[j+1]-1) while the index range is (offs[j]+1):offs[j+1] for a SparseOperatorCSC. For this reason, it is recommended to call each_off instead or to call get_offs with 2 arguments as shown below.

For a transparent usage of the offsets, the method should be called with 2 arguments:

get_offs(A, i) -> k1, k2

which yields the offsets of the first and last elements in the arrays of values and linear column indices for the i-th row of the sparse operator A stored in a Compressed Sparse Row (CSR) format. If k2 < k1, it means that the i-th row is empty. Calling each_off(A,i) directly yields k1:k2.

get_offs(A, j) -> k1, k2

yields the offsets of the first and last elements in the arrays of values and linear row indices for the j-th column of the sparse operator A stored in a Compressed Sparse Column (CSC) format. If k2 < k1, it means that the j-th column is empty. Calling each_off(A,j) directly yields k1:k2.

source
LazyAlgebra.SparseOperators.copy_rowsFunction
copy_rows(A) -> rows

yields a copy of the linear row indices of entries in sparse operator A. The result is a vector that is not shared by A, the caller may thus modify its contents with no side effects on A.

source
LazyAlgebra.SparseOperators.copy_colsFunction
copy_cols(A) -> cols

yields a copy of the linear column indices of entries in sparse operator A. The result is a vector that is not shared by A, the caller may thus modify its contents with no side effects on A.

source
LazyAlgebra.SparseOperators.copy_valsFunction
copy_vals([T = eltype(A),] A) -> vals

yields a copy of the values of the entries in sparse operator A converted to type T. The result is a vector that is not shared by A, the caller may thus modify its contents with no side effects on A.

source