(octave.info)Creating Sparse Matrices


Next: Information Prev: Storage of Sparse Matrices Up: Basics
Enter node , (file) or (file)node

22.1.2 Creating Sparse Matrices
-------------------------------

There are several means to create sparse matrix.

Returned from a function
     There are many functions that directly return sparse matrices.
     These include “speye”, “sprand”, “diag”, etc.

Constructed from matrices or vectors
     The function “sparse” allows a sparse matrix to be constructed from
     three vectors representing the row, column and data.
     Alternatively, the function “spconvert” uses a three column matrix
     format to allow easy importation of data from elsewhere.

Created and then filled
     The function “sparse” or “spalloc” can be used to create an empty
     matrix that is then filled by the user

From a user binary program
     The user can directly create the sparse matrix within an oct-file.

   There are several basic functions to return specific sparse matrices.
For example the sparse identity matrix, is a matrix that is often
needed.  It therefore has its own function to create it as ‘speye (N)’
or ‘speye (R, C)’, which creates an N-by-N or R-by-C sparse identity
matrix.

   Another typical sparse matrix that is often needed is a random
distribution of random elements.  The functions “sprand” and “sprandn”
perform this for uniform and normal random distributions of elements.
They have exactly the same calling convention, where ‘sprand (R, C, D)’,
creates an R-by-C sparse matrix with a density of filled elements of D.

   Other functions of interest that directly create sparse matrices, are
“diag” or its generalization “spdiags”, that can take the definition of
the diagonals of the matrix and create the sparse matrix that
corresponds to this.  For example,

     s = diag (sparse (randn (1,n)), -1);

creates a sparse (N+1)-by-(N+1) sparse matrix with a single diagonal
defined.

 -- : B = spdiags (A)
 -- : [B, D] = spdiags (A)
 -- : B = spdiags (A, D)
 -- : A = spdiags (V, D, A)
 -- : A = spdiags (V, D, M, N)
     A generalization of the function ‘diag’.

     Called with a single input argument, the nonzero diagonals D of A
     are extracted.

     With two arguments the diagonals to extract are given by the vector
     D.

     The other two forms of ‘spdiags’ modify the input matrix by
     replacing the diagonals.  They use the columns of V to replace the
     diagonals represented by the vector D.  If the sparse matrix A is
     defined then the diagonals of this matrix are replaced.  Otherwise
     a matrix of M by N is created with the diagonals given by the
     columns of V.

     Negative values of D represent diagonals below the main diagonal,
     and positive values of D diagonals above the main diagonal.

     For example:

          spdiags (reshape (1:12, 4, 3), [-1 0 1], 5, 4)
             ⇒ 5 10  0  0
                1  6 11  0
                0  2  7 12
                0  0  3  8
                0  0  0  4

     See also: Note: diag.

 -- : S = speye (M, N)
 -- : S = speye (M)
 -- : S = speye (SZ)
     Return a sparse identity matrix of size MxN.

     The implementation is significantly more efficient than ‘sparse
     (eye (M))’ as the full matrix is not constructed.

     Called with a single argument a square matrix of size M-by-M is
     created.  If called with a single vector argument SZ, this argument
     is taken to be the size of the matrix to create.

     See also: Note: sparse, Note: spdiags,
     Note: eye.

 -- : R = spones (S)
     Replace the nonzero entries of S with ones.

     This creates a sparse matrix with the same structure as S.

     See also: Note: sparse, Note: sprand, Note:
     sprandn, Note: sprandsym, *note spfun:
     XREFspfun, Note: spy.

 -- : sprand (M, N, D)
 -- : sprand (M, N, D, RC)
 -- : sprand (S)
     Generate a sparse matrix with uniformly distributed random values.

     The size of the matrix is MxN with a density of values D.  D must
     be between 0 and 1.  Values will be uniformly distributed on the
     interval (0, 1).

     If called with a single matrix argument, a sparse matrix is
     generated with random values wherever the matrix S is nonzero.

     If called with a scalar fourth argument RC, a random sparse matrix
     with reciprocal condition number RC is generated.  If RC is a
     vector, then it specifies the first singular values of the
     generated matrix (‘length (RC) <= min (M, N)’).

     See also: Note: sprandn, *note sprandsym:
     XREFsprandsym, Note: rand.

 -- : sprandn (M, N, D)
 -- : sprandn (M, N, D, RC)
 -- : sprandn (S)
     Generate a sparse matrix with normally distributed random values.

     The size of the matrix is MxN with a density of values D.  D must
     be between 0 and 1.  Values will be normally distributed with a
     mean of 0 and a variance of 1.

     If called with a single matrix argument, a sparse matrix is
     generated with random values wherever the matrix S is nonzero.

     If called with a scalar fourth argument RC, a random sparse matrix
     with reciprocal condition number RC is generated.  If RC is a
     vector, then it specifies the first singular values of the
     generated matrix (‘length (RC) <= min (M, N)’).

     See also: Note: sprand, Note: sprandsym,
     Note: randn.

 -- : sprandsym (N, D)
 -- : sprandsym (S)
     Generate a symmetric random sparse matrix.

     The size of the matrix will be NxN, with a density of values given
     by D.  D must be between 0 and 1 inclusive.  Values will be
     normally distributed with a mean of zero and a variance of 1.

     If called with a single matrix argument, a random sparse matrix is
     generated wherever the matrix S is nonzero in its lower triangular
     part.

     See also: Note: sprand, Note: sprandn,
     Note: spones, Note: sparse.

   The recommended way for the user to create a sparse matrix, is to
create two vectors containing the row and column index of the data and a
third vector of the same size containing the data to be stored.  For
example,

       ri = ci = d = [];
       for j = 1:c
         ri = [ri; randperm(r,n)'];
         ci = [ci; j*ones(n,1)];
         d = [d; rand(n,1)];
       endfor
       s = sparse (ri, ci, d, r, c);

creates an R-by-C sparse matrix with a random distribution of N (<R)
elements per column.  The elements of the vectors do not need to be
sorted in any particular order as Octave will sort them prior to storing
the data.  However, pre-sorting the data will make the creation of the
sparse matrix faster.

   The function “spconvert” takes a three or four column real matrix.
The first two columns represent the row and column index respectively
and the third and four columns, the real and imaginary parts of the
sparse matrix.  The matrix can contain zero elements and the elements
can be sorted in any order.  Adding zero elements is a convenient way to
define the size of the sparse matrix.  For example:

     s = spconvert ([1 2 3 4; 1 3 4 4; 1 2 3 0]')
     ⇒ Compressed Column Sparse (rows=4, cols=4, nnz=3)
           (1 , 1) -> 1
           (2 , 3) -> 2
           (3 , 4) -> 3

   An example of creating and filling a matrix might be

     k = 5;
     nz = r * k;
     s = spalloc (r, c, nz)
     for j = 1:c
       idx = randperm (r);
       s (:, j) = [zeros(r - k, 1); ...
             rand(k, 1)] (idx);
     endfor

   It should be noted, that due to the way that the Octave assignment
functions are written that the assignment will reallocate the memory
used by the sparse matrix at each iteration of the above loop.
Therefore the “spalloc” function ignores the NZ argument and does not
pre-assign the memory for the matrix.  Therefore, it is vitally
important that code using to above structure should be vectorized as
much as possible to minimize the number of assignments and reduce the
number of memory allocations.

 -- : FM = full (SM)
     Return a full storage matrix from a sparse, diagonal, or
     permutation matrix, or a range.

     See also: Note: sparse, Note: issparse.

 -- : S = spalloc (M, N, NZ)
     Create an M-by-N sparse matrix with pre-allocated space for at most
     NZ nonzero elements.

     This is useful for building a matrix incrementally by a sequence of
     indexed assignments.  Subsequent indexed assignments after
     ‘spalloc’ will reuse the pre-allocated memory, provided they are of
     one of the simple forms

        • ‘S(I:J) = X’

        • ‘S(:,I:J) = X’

        • ‘S(K:L,I:J) = X’

     and that the following conditions are met:

        • the assignment does not decrease nnz (S).

        • after the assignment, nnz (S) does not exceed NZ.

        • no index is out of bounds.

     Partial movement of data may still occur, but in general the
     assignment will be more memory and time efficient under these
     circumstances.  In particular, it is possible to efficiently build
     a pre-allocated sparse matrix from a contiguous block of columns.

     The amount of pre-allocated memory for a given matrix may be
     queried using the function ‘nzmax’.

     See also: Note: nzmax, Note: sparse.

 -- : S = sparse (A)
 -- : S = sparse (I, J, SV, M, N)
 -- : S = sparse (I, J, SV)
 -- : S = sparse (M, N)
 -- : S = sparse (I, J, S, M, N, "unique")
 -- : S = sparse (I, J, SV, M, N, NZMAX)
     Create a sparse matrix from a full matrix, or row, column, value
     triplets.

     If A is a full matrix, convert it to a sparse matrix
     representation, removing all zero values in the process.

     Given the integer index vectors I and J, and a 1-by-‘nnz’ vector of
     real or complex values SV, construct the sparse matrix
     ‘S(I(K),J(K)) = SV(K)’ with overall dimensions M and N.  If any of
     SV, I or J are scalars, they are expanded to have a common size.

     If M or N are not specified their values are derived from the
     maximum index in the vectors I and J as given by ‘M = max (I)’, ‘N
     = max (J)’.

     *Note*: if multiple values are specified with the same I, J
     indices, the corresponding value in S will be the sum of the values
     at the repeated location.  See ‘accumarray’ for an example of how
     to produce different behavior, such as taking the minimum instead.

     If the option "unique" is given, and more than one value is
     specified at the same I, J indices, then the last specified value
     will be used.

     ‘sparse (M, N)’ will create an empty MxN sparse matrix and is
     equivalent to ‘sparse ([], [], [], M, N)’

     The argument NZMAX is ignored but accepted for compatibility with
     MATLAB.

     Example 1 (sum at repeated indices):

          I = [1 1 2]; J = [1 1 2]; SV = [3 4 5];
          sparse (I, J, SV, 3, 4)
          ⇒
          Compressed Column Sparse (rows = 3, cols = 4, nnz = 2 [17%])

            (1, 1) ->  7
            (2, 2) ->  5

     Example 2 ("unique" option):

          I = [1 1 2]; J = [1 1 2]; SV = [3 4 5];
          sparse (I, J, SV, 3, 4, "unique")
          ⇒
          Compressed Column Sparse (rows = 3, cols = 4, nnz = 2 [17%])

            (1, 1) ->  4
            (2, 2) ->  5

     See also: Note: full, Note: accumarray,
     Note: spalloc, Note: spdiags, Note:
     speye, Note: spones, *note sprand:
     XREFsprand, Note: sprandn, *note sprandsym:
     XREFsprandsym, Note: spconvert, *note spfun:
     XREFspfun.

 -- : X = spconvert (M)
     Convert a simple sparse matrix format easily generated by other
     programs into Octave’s internal sparse format.

     The input M is either a 3 or 4 column real matrix, containing the
     row, column, real, and imaginary parts of the elements of the
     sparse matrix.  An element with a zero real and imaginary part can
     be used to force a particular matrix size.

     See also: Note: sparse.

   The above problem of memory reallocation can be avoided in oct-files.
However, the construction of a sparse matrix from an oct-file is more
complex than can be discussed here.  Note: External Code Interface,
for a full description of the techniques involved.


automatically generated by info2www version 1.2.2.9