GNU Octave Manual Version 3 by John W. Eaton, David Bateman, Søren Hauberg Paperback (6"x9"), 568 pages ISBN 095461206X RRP £24.95 ($39.95) |
20.1.2 Creating Sparse Matrices
A sparse matrix can be created in the following ways:
- Returned from a function
-
Many functions return sparse matrices directly. These include
speye
,sprand
,spdiag
, 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 functionspconvert
uses a three column matrix format to allow easy importation of data from elsewhere. - Created and then filled
-
The function
sparse
orspalloc
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 speye (n)
or speye (r,
c)
returns an n-by-n or r-by-c sparse
identity matrix. The function diag
returns a sparse diagonal
matrix when given a sparse input vector. For example,
s = diag (sparse(randn(1,n)), -1);
creates a sparse (n+1)-by-(n+1) sparse matrix with a
single diagonal defined. The function spdiags
described below
is a generalisation of diag
.
- Function File: y = speye (m)
- Function File: y = speye (m, n)
- Function File: y = speye (sz)
- Returns a sparse identity matrix. This 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. Otherwise a matrix of m by n is created. If called with a single vector argument, this argument is taken to be the size of the matrix to create.
- Function File: y = spfun (f,x)
- Compute
f(x)
for the non-zero values of x, This results in a sparse matrix with the same structure as x. The function f can be passed as a string, a function handle or an inline function.
- Function File: y = spones (x)
- Replace the non-zero entries of x with ones. This creates a sparse matrix with the same structure as x.
- Loadable Function: spdiag (v, k)
- Return a diagonal matrix with the sparse vector v on diagonal
k. The second argument is optional. If it is positive, the vector is
placed on the k-th super-diagonal. If it is negative, it is placed
on the -k-th sub-diagonal. The default value of k is 0, and
the vector is placed on the main diagonal. For example,
spdiag ([1, 2, 3], 1) ans = Compressed Column Sparse (rows=4, cols=4, nnz=3) (1 , 2) -> 1 (2 , 3) -> 2 (3 , 4) -> 3
Given a matrix argument, instead of a vector,
spdiag
extracts the k-th diagonal of the sparse matrix.See also diag
- Function File: [b, c] = spdiags (a)
- Function File: b = spdiags (a, c)
- Function File: b = spdiags (v, c, a)
- Function File: b = spdiags (v, c, m, n)
- A generalization of the function
spdiag
. Called with a single input argument, the non-zero diagonals c of A are extracted. With two arguments the diagonals to extract are given by the vector c.The other two forms of
spdiags
modify the input matrix by replacing the diagonals. They use the columns of v to replace the columns represented by the vector c. 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 v.Negative values of c representative diagonals below the main diagonal, and positive values of c 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
- Loadable Function: s = sparse (a)
- Create a sparse matrix from the full matrix a.
- Loadable Function: s = sparse (i, j, sv, m, n, nzmax)
- Create a sparse matrix given integer index vectors i and j,
a 1-by-
nnz
vector of real of complex values sv, overall dimensions m and n of the sparse matrix. The argumentnzmax
is ignored but accepted for compatibility with Matlab.Note: if multiple values are specified with the same i, j indices, the corresponding values in s will be added.
The following are all equivalent:
s = sparse (i, j, s, m, n) s = sparse (i, j, s, m, n, "summation") s = sparse (i, j, s, m, n, "sum")
- Loadable Function: s = sparse (i, j, s, m, n, "unique")
- Same as above, except that if more than two values are specified for the
same i, j indices, the last specified value will be used.
- Loadable Function: s = sparse (i, j, sv)
- Uses
m = max (i)
,n = max (j)
- Loadable Function: s = sparse (m, n)
- Equivalent to
sparse ([], [], [], m, n, 0)
If any of sv, i or j are scalars, they are expanded to have a common size.
See also full
The recommended way to create a sparse matrix is to use the function
sparse
with three vectors containing the row index, column
index and value to be stored. For example, the following commands
create the 3-by-4 sparse matrix shown earlier:
r = 3; c = 4; ri = [1, 1, 2, 3]; ci = [1, 2, 4, 4]; di = [1, 2, 3, 4]; s = sparse (ri, ci, di, r, c);
The indices do not need to be in any particular order as Octave will sort them prior to storing the data.(13)
- Function File: x = spconvert (m)
- This function converts from a simple sparse matrix format easily produced by other programs into Octave's internal sparse format. The input x is either a three or four column real matrix, containing the row index, column index, real and imaginary parts of the elements of the sparse matrix. The matrix can contain zero elements and the elements can be sorted in any order. An element with a zero real and imaginary part can be used to force a particular matrix size.
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
- Function File: s = spalloc (r, c, nz)
- Returns an empty sparse matrix of size r-by-c. This
function is provided only for compatibility reasons and the argument
nz is ignored. Octave resizes the memory required for sparse
matrices on demand whenever they are modified, so
spalloc
does not preassign any memory.It should be noted that this means that code like
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
will reallocate memory at each step. It is therefore vitally important that such code should be vectorized as much as possible, to minimize the number of assignments and reduce the number of memory allocations.(14)
See also sparse, nzmax
- Loadable Function: FM = full (SM)
- returns a full storage matrix from a sparse one
See also sparse
- Function File: sprand (m, n, d)
- Function File: sprand (s)
- Generate a random sparse matrix using uniform variates.
The size of the matrix will be m by n,
with a density of values given by d, which
should be between 0 and 1. Values will be uniformly
distributed between 0 and 1.
The actual density may be slightly lower than d, due to random collisions in the selection of indices. These are unlikely to occur for large, really sparse matrices.
If called with a single matrix argument, a sparse matrix is generated with random values wherever the matrix S is non-zero.
See also sprandn
- Function File: sprandn (m, n, d)
- Function File: sprandn (s)
- Generate a random sparse matrix using normal variates.
The size of the matrix will be m by n,
with a density of values given by d, which
should be between 0 and 1. Values will be normally
distributed with mean of zero and variance 1.
The actual density may be slightly lower than d, due to random collisions in the selection of indices. These are unlikely to occur for large, really sparse matrices.
If called with a single matrix argument, a sparse matrix is generated with random values wherever the matrix S is non-zero.
See also sprand
- Function File: sprandsym (n, d)
- Function File: sprandsym (s)
- Generate a symmetric random sparse matrix. The size of the matrix will be
n by n, with a density of values given by d.
d should be between 0 and 1. Values will be normally
distributed with mean of zero and variance 1.
The actual density may be slightly lower than d, due to random collisions in the selection of indices. These are unlikely to occur for large, really sparse matrices.
If called with a single matrix argument, a sparse matrix is generated with random values wherever the matrix S is non-zero in its lower triangular part.
See also sprand, sprandn
ISBN 095461206X | GNU Octave Manual Version 3 | See the print edition |