Mercurial > hg > octave-nkf
view doc/interpreter/sparse.txi @ 6210:12b676a0b183 before-graphics-branch
[project @ 2006-12-07 02:37:17 by jwe]
author | jwe |
---|---|
date | Thu, 07 Dec 2006 02:37:17 +0000 |
parents | 55404f3b0da1 |
children | 877b80a8dee7 |
line wrap: on
line source
@c Copyright (C) 2004, 2005 David Bateman @c This is part of the Octave manual. @c For copying conditions, see the file gpl.texi. @ifhtml @set htmltex @end ifhtml @iftex @set htmltex @end iftex @node Sparse Matrices @chapter Sparse Matrices @menu * Basics:: The Creation and Manipulation of Sparse Matrices * Sparse Linear Algebra:: Linear Algebra on Sparse Matrices * Iterative Techniques:: Iterative Techniques applied to Sparse Matrices * Real Life Example:: Real Life Example of the use of Sparse Matrices * Oct-Files:: Using Sparse Matrices in Oct-files * Function Reference:: Documentation from the Specific Sparse Functions @end menu @node Basics, Sparse Linear Algebra, Sparse Matrices, Sparse Matrices @section The Creation and Manipulation of Sparse Matrices The size of mathematical problems that can be treated at any particular time is generally limited by the available computing resources. Both, the speed of the computer and its available memory place limitation on the problem size. There are many classes of mathematical problems which give rise to matrices, where a large number of the elements are zero. In this case it makes sense to have a special matrix type to handle this class of problems where only the non-zero elements of the matrix are stored. Not only does this reduce the amount of memory to store the matrix, but it also means that operations on this type of matrix can take advantage of the a-priori knowledge of the positions of the non-zero elements to accelerate their calculations. A matrix type that stores only the non-zero elements is generally called sparse. It is the purpose of this document to discuss the basics of the storage and creation of sparse matrices and the fundamental operations on them. @menu * Storage:: Storage of Sparse Matrices * Creation:: Creating Sparse Matrices * Information:: Finding out Information about Sparse Matrices * Operators and Functions:: Basic Operators and Functions on Sparse Matrices @end menu @node Storage, Creation, Basics, Basics @subsection Storage of Sparse Matrices It is not strictly speaking necessary for the user to understand how sparse matrices are stored. However, such an understanding will help to get an understanding of the size of sparse matrices. Understanding the storage technique is also necessary for those users wishing to create their own oct-files. There are many different means of storing sparse matrix data. What all of the methods have in common is that they attempt to reduce the complexity and storage given a-priori knowledge of the particular class of problems that will be solved. A good summary of the available techniques for storing sparse matrix is given by Saad @footnote{Youcef Saad "SPARSKIT: A basic toolkit for sparse matrix computation", 1994, @url{ftp://ftp.cs.umn.edu/dept/sparse/SPARSKIT2/DOC/paper.ps}}. With full matrices, knowledge of the point of an element of the matrix within the matrix is implied by its position in the computers memory. However, this is not the case for sparse matrices, and so the positions of the non-zero elements of the matrix must equally be stored. An obvious way to do this is by storing the elements of the matrix as triplets, with two elements being their position in the array (rows and column) and the third being the data itself. This is conceptually easy to grasp, but requires more storage than is strictly needed. The storage technique used within Octave is the compressed column format. In this format the position of each element in a row and the data are stored as previously. However, if we assume that all elements in the same column are stored adjacent in the computers memory, then we only need to store information on the number of non-zero elements in each column, rather than their positions. Thus assuming that the matrix has more non-zero elements than there are columns in the matrix, we win in terms of the amount of memory used. In fact, the column index contains one more element than the number of columns, with the first element always being zero. The advantage of this is a simplification in the code, in that their is no special case for the first or last columns. A short example, demonstrating this in C is. @example for (j = 0; j < nc; j++) for (i = cidx (j); i < cidx(j+1); i++) printf ("non-zero element (%i,%i) is %d\n", ridx(i), j, data(i)); @end example A clear understanding might be had by considering an example of how the above applies to an example matrix. Consider the matrix @example @group 1 2 0 0 0 0 0 3 0 0 0 4 @end group @end example The non-zero elements of this matrix are @example @group (1, 1) @result{} 1 (1, 2) @result{} 2 (2, 4) @result{} 3 (3, 4) @result{} 4 @end group @end example This will be stored as three vectors @var{cidx}, @var{ridx} and @var{data}, representing the column indexing, row indexing and data respectively. The contents of these three vectors for the above matrix will be @example @group @var{cidx} = [0, 1, 2, 2, 4] @var{ridx} = [0, 0, 1, 2] @var{data} = [1, 2, 3, 4] @end group @end example Note that this is the representation of these elements with the first row and column assumed to start at zero, while in Octave itself the row and column indexing starts at one. Thus the number of elements in the @var{i}-th column is given by @code{@var{cidx} (@var{i} + 1) - @var{cidx} (@var{i})}. Although Octave uses a compressed column format, it should be noted that compressed row formats are equally possible. However, in the context of mixed operations between mixed sparse and dense matrices, it makes sense that the elements of the sparse matrices are in the same order as the dense matrices. Octave stores dense matrices in column major ordering, and so sparse matrices are equally stored in this manner. A further constraint on the sparse matrix storage used by Octave is that all elements in the rows are stored in increasing order of their row index, which makes certain operations faster. However, it imposes the need to sort the elements on the creation of sparse matrices. Having dis-ordered elements is potentially an advantage in that it makes operations such as concatenating two sparse matrices together easier and faster, however it adds complexity and speed problems elsewhere. @node Creation, Information, Storage, Basics @subsection Creating Sparse Matrices There are several means to create sparse matrix. @table @asis @item Returned from a function There are many functions that directly return sparse matrices. These include @dfn{speye}, @dfn{sprand}, @dfn{spdiag}, etc. @item Constructed from matrices or vectors The function @dfn{sparse} allows a sparse matrix to be constructed from three vectors representing the row, column and data. Alternatively, the function @dfn{spconvert} uses a three column matrix format to allow easy importation of data from elsewhere. @item Created and then filled The function @dfn{sparse} or @dfn{spalloc} can be used to create an empty matrix that is then filled by the user @item From a user binary program The user can directly create the sparse matrix within an oct-file. @end table 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 @code{speye (@var{n})} or @code{speye (@var{r}, @var{c})}, which creates an @var{n}-by-@var{n} or @var{r}-by-@var{c} sparse identity matrix. Another typical sparse matrix that is often needed is a random distribution of random elements. The functions @dfn{sprand} and @dfn{sprandn} perform this for uniform and normal random distributions of elements. They have exactly the same calling convention, where @code{sprand (@var{r}, @var{c}, @var{d})}, creates an @var{r}-by-@var{c} sparse matrix with a density of filled elements of @var{d}. Other functions of interest that directly creates a sparse matrices, are @dfn{spdiag} or its generalization @dfn{spdiags}, that can take the definition of the diagonals of the matrix and create the sparse matrix that corresponds to this. For example @example s = diag (sparse(randn(1,n)), -1); @end example creates a sparse (@var{n}+1)-by-(@var{n}+1) sparse matrix with a single diagonal defined. 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 @example function x = foo (r, j) idx = randperm (r); x = ([zeros(r-2,1); rand(2,1)])(idx); endfunction ri = []; ci = []; d = []; for j=1:c dtmp = foo (r, j); # Returns vector of length r with (:,j) values idx = find (dtmp != 0.); ri = [ri; idx]; ci = [ci; j*ones(length(idx),1)]; d = [d; dtmp(idx)]; endfor s = sparse (ri, ci, d, r, c); @end example creates an @var{r}-by-@var{c} sparse matrix with a random distribution of 2 elements per row. 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 @dfn{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 @example s = spconvert ([1 2 3 4; 1 3 4 4; 1 2 3 0]') @result{} Compressed Column Sparse (rows=4, cols=4, nnz=3) (1 , 1) -> 1 (2 , 3) -> 2 (3 , 4) -> 3 @end example An example of creating and filling a matrix might be @example 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 @end example 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 @dfn{spalloc} function ignores the @var{nz} argument and does not preassign 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. The above problem can be avoided in oct-files. However, the construction of a sparse matrix from an oct-file is more complex than can be discussed in this brief introduction, and you are referred to section @ref{Oct-Files}, to have a full description of the techniques involved. @node Information, Operators and Functions, Creation, Basics @subsection Finding out Information about Sparse Matrices There are a number of functions that allow information concerning sparse matrices to be obtained. The most basic of these is @dfn{issparse} that identifies whether a particular Octave object is in fact a sparse matrix. Another very basic function is @dfn{nnz} that returns the number of non-zero entries there are in a sparse matrix, while the function @dfn{nzmax} returns the amount of storage allocated to the sparse matrix. Note that Octave tends to crop unused memory at the first opportunity for sparse objects. There are some cases of user created sparse objects where the value returned by @dfn{nzmaz} will not be the same as @dfn{nnz}, but in general they will give the same result. The function @dfn{spstats} returns some basic statistics on the columns of a sparse matrix including the number of elements, the mean and the variance of each column. When solving linear equations involving sparse matrices Octave determines the means to solve the equation based on the type of the matrix as discussed in @ref{Sparse Linear Algebra}. Octave probes the matrix type when the div (/) or ldiv (\) operator is first used with the matrix and then caches the type. However the @dfn{matrix_type} function can be used to determine the type of the sparse matrix prior to use of the div or ldiv operators. For example @example a = tril (sprandn(1024, 1024, 0.02), -1) ... + speye(1024); matrix_type (a); ans = Lower @end example show that Octave correctly determines the matrix type for lower triangular matrices. @dfn{matrix_type} can also be used to force the type of a matrix to be a particular type. For example @example a = matrix_type (tril (sprandn (1024, ... 1024, 0.02), -1) + speye(1024), 'Lower'); @end example This allows the cost of determining the matrix type to be avoided. However, incorrectly defining the matrix type will result in incorrect results from solutions of linear equations, and so it is entirely the responsibility of the user to correctly identify the matrix type There are several graphical means of finding out information about sparse matrices. The first is the @dfn{spy} command, which displays the structure of the non-zero elements of the matrix. @xref{fig:spmatrix}, for an exaple of the use of @dfn{spy}. More advanced graphical information can be obtained with the @dfn{treeplot}, @dfn{etreeplot} and @dfn{gplot} commands. @float Figure,fig:spmatrix @image{spmatrix,8cm} @caption{Structure of simple sparse matrix.} @end float One use of sparse matrices is in graph theory, where the interconnections between nodes is represented as an adjacency matrix. That is, if the i-th node in a graph is connected to the j-th node. Then the ij-th node (and in the case of undirected graphs the ji-th node) of the sparse adjacency matrix is non-zero. If each node is then associated with a set of co-ordinates, then the @dfn{gplot} command can be used to graphically display the interconnections between nodes. As a trivial example of the use of @dfn{gplot}, consider the example @example A = sparse([2,6,1,3,2,4,3,5,4,6,1,5], [1,1,2,2,3,3,4,4,5,5,6,6],1,6,6); xy = [0,4,8,6,4,2;5,0,5,7,5,7]'; gplot(A,xy) @end example which creates an adjacency matrix @code{A} where node 1 is connected to nodes 2 and 6, node 2 with nodes 1 and 3, etc. The co-ordinates of the nodes are given in the n-by-2 matrix @code{xy}. @ifset htmltex @xref{fig:gplot}. @float Figure,fig:gplot @image{gplot,8cm} @caption{Simple use of the @dfn{gplot} command.} @end float @end ifset The dependencies between the nodes of a Cholesky factorization can be calculated in linear time without explicitly needing to calculate the Cholesky factorization by the @code{etree} command. This command returns the elimination tree of the matrix and can be displayed graphically by the command @code{treeplot(etree(A))} if @code{A} is symmetric or @code{treeplot(etree(A+A'))} otherwise. @node Operators and Functions, , Information, Basics @subsection Basic Operators and Functions on Sparse Matrices @menu * Functions:: Sparse Functions * ReturnType:: The Return Types of Operators and Functions * MathConsiderations:: Mathematical Considerations @end menu @node Functions, ReturnType, Operators and Functions, Operators and Functions @subsubsection Sparse Functions An important consideration in the use of the sparse functions of Octave is that many of the internal functions of Octave, such as @dfn{diag}, can not accept sparse matrices as an input. The sparse implementation in Octave therefore uses the @dfn{dispatch} function to overload the normal Octave functions with equivalent functions that work with sparse matrices. However, at any time the sparse matrix specific version of the function can be used by explicitly calling its function name. The table below lists all of the sparse functions of Octave together (with possible future extensions that are currently unimplemented, listed last). Note that in this specific sparse forms of the functions are typically the same as the general versions with a @dfn{sp} prefix. In the table below, and the rest of this article the specific sparse versions of the functions are used. @table @asis @item Generate sparse matrices: @dfn{spalloc}, @dfn{spdiags}, @dfn{speye}, @dfn{sprand}, @dfn{sprandn}, @dfn{sprandsym} @item Sparse matrix conversion: @dfn{full}, @dfn{sparse}, @dfn{spconvert}, @dfn{spfind} @item Manipulate sparse matrices @dfn{issparse}, @dfn{nnz}, @dfn{nonzeros}, @dfn{nzmax}, @dfn{spfun}, @dfn{spones}, @dfn{spy}, @item Graph Theory: @dfn{etree}, @dfn{etreeplot}, @dfn{gplot}, @dfn{treeplot}, (treelayout) @item Sparse matrix reordering: @dfn{ccolamd}, @dfn{colamd}, @dfn{colperm}, @dfn{csymamd}, @dfn{dmperm}, @dfn{symamd}, @dfn{randperm}, (symrcm) @item Linear algebra: @dfn{matrix\_type}, @dfn{spchol}, @dfn{cpcholinv}, @dfn{spchol2inv}, @dfn{spdet}, @dfn{spinv}, @dfn{spkron}, @dfn{splchol}, @dfn{splu}, @dfn{spqr}, (condest, eigs, normest, sprank, svds, spaugment) @item Iterative techniques: @dfn{luinc}, @dfn{pcg}, @dfn{pcr}, (bicg, bicgstab, cholinc, cgs, gmres, lsqr, minres, qmr, symmlq) @item Miscellaneous: @dfn{spparms}, @dfn{symbfact}, @dfn{spstats}, @dfn{spprod}, @dfn{spcumsum}, @dfn{spsum}, @dfn{spsumsq}, @dfn{spmin}, @dfn{spmax}, @dfn{spatan2}, @dfn{spdiag} @end table In addition all of the standard Octave mapper functions (ie. basic math functions that take a single argument) such as @dfn{abs}, etc can accept sparse matrices. The reader is referred to the documentation supplied with these functions within Octave itself for further details. @node ReturnType, MathConsiderations, Functions, Operators and Functions @subsubsection The Return Types of Operators and Functions The two basic reasons to use sparse matrices are to reduce the memory usage and to not have to do calculations on zero elements. The two are closely related in that the computation time on a sparse matrix operator or function is roughly linear with the number of non-zero elements. Therefore, there is a certain density of non-zero elements of a matrix where it no longer makes sense to store it as a sparse matrix, but rather as a full matrix. For this reason operators and functions that have a high probability of returning a full matrix will always return one. For example adding a scalar constant to a sparse matrix will almost always make it a full matrix, and so the example @example speye(3) + 0 @result{} 1 0 0 0 1 0 0 0 1 @end example returns a full matrix as can be seen. Additionally all sparse functions test the amount of memory occupied by the sparse matrix to see if the amount of storage used is larger than the amount used by the full equivalent. Therefore @code{speye (2) * 1} will return a full matrix as the memory used is smaller for the full version than the sparse version. As all of the mixed operators and functions between full and sparse matrices exist, in general this does not cause any problems. However, one area where it does cause a problem is where a sparse matrix is promoted to a full matrix, where subsequent operations would resparsify the matrix. Such cases are rare, but can be artificially created, for example @code{(fliplr(speye(3)) + speye(3)) - speye(3)} gives a full matrix when it should give a sparse one. In general, where such cases occur, they impose only a small memory penalty. There is however one known case where this behavior of Octave's sparse matrices will cause a problem. That is in the handling of the @dfn{diag} function. Whether @dfn{diag} returns a sparse or full matrix depending on the type of its input arguments. So @example a = diag (sparse([1,2,3]), -1); @end example should return a sparse matrix. To ensure this actually happens, the @dfn{sparse} function, and other functions based on it like @dfn{speye}, always returns a sparse matrix, even if the memory used will be larger than its full representation. @node MathConsiderations, , ReturnType, Operators and Functions @subsubsection Mathematical Considerations The attempt has been made to make sparse matrices behave in exactly the same manner as there full counterparts. However, there are certain differences and especially differences with other products sparse implementations. Firstly, the "./" and ".^" operators must be used with care. Consider what the examples @example s = speye (4); a1 = s .^ 2; a2 = s .^ s; a3 = s .^ -2; a4 = s ./ 2; a5 = 2 ./ s; a6 = s ./ s; @end example will give. The first example of @var{s} raised to the power of 2 causes no problems. However @var{s} raised element-wise to itself involves a a large number of terms @code{0 .^ 0} which is 1. There @code{@var{s} .^ @var{s}} is a full matrix. Likewise @code{@var{s} .^ -2} involves terms terms like @code{0 .^ -2} which is infinity, and so @code{@var{s} .^ -2} is equally a full matrix. For the "./" operator @code{@var{s} ./ 2} has no problems, but @code{2 ./ @var{s}} involves a large number of infinity terms as well and is equally a full matrix. The case of @code{@var{s} ./ @var{s}} involves terms like @code{0 ./ 0} which is a @code{NaN} and so this is equally a full matrix with the zero elements of @var{s} filled with @code{NaN} values. The above behavior is consistent with full matrices, but is not consistent with sparse implementations in other products. A particular problem of sparse matrices comes about due to the fact that as the zeros are not stored, the sign-bit of these zeros is equally not stored. In certain cases the sign-bit of zero is important. For example @example a = 0 ./ [-1, 1; 1, -1]; b = 1 ./ a @result{} -Inf Inf Inf -Inf c = 1 ./ sparse (a) @result{} Inf Inf Inf Inf @end example To correct this behavior would mean that zero elements with a negative sign-bit would need to be stored in the matrix to ensure that their sign-bit was respected. This is not done at this time, for reasons of efficient, and so the user is warned that calculations where the sign-bit of zero is important must not be done using sparse matrices. In general any function or operator used on a sparse matrix will result in a sparse matrix with the same or a larger number of non-zero elements than the original matrix. This is particularly true for the important case of sparse matrix factorizations. The usual way to address this is to reorder the matrix, such that its factorization is sparser than the factorization of the original matrix. That is the factorization of @code{L * U = P * S * Q} has sparser terms @code{L} and @code{U} than the equivalent factorization @code{L * U = S}. Several functions are available to reorder depending on the type of the matrix to be factorized. If the matrix is symmetric positive-definite, then @dfn{symamd} or @dfn{csymamd} should be used. Otherwise @dfn{colamd} or @dfn{ccolamd} should be used. For completeness the reordering functions @dfn{colperm} and @dfn{randperm} are also available. @xref{fig:simplematrix}, for an example of the structure of a simple positive definite matrix. @float Figure,fig:simplematrix @image{spmatrix,8cm} @caption{Structure of simple sparse matrix.} @end float The standard Cholesky factorization of this matrix, can be obtained by the same command that would be used for a full matrix. This can be visualized with the command @code{r = chol(A); spy(r);}. @ifset HAVE_CHOLMOD @ifset HAVE_COLAMD @xref{fig:simplechol}. @end ifset @end ifset The original matrix had @ifinfo @ifnothtml 43 @end ifnothtml @end ifinfo @ifset htmltex 598 @end ifset non-zero terms, while this Cholesky factorization has @ifinfo @ifnothtml 71, @end ifnothtml @end ifinfo @ifset htmltex 10200, @end ifset with only half of the symmetric matrix being stored. This is a significant level of fill in, and although not an issue for such a small test case, can represents a large overhead in working with other sparse matrices. The appropriate sparsity preserving permutation of the original matrix is given by @dfn{symamd} and the factorization using this reordering can be visualized using the command @code{q = symamd(A); r = chol(A(q,q)); spy(r)}. This gives @ifinfo @ifnothtml 29 @end ifnothtml @end ifinfo @ifset htmltex 399 @end ifset non-zero terms which is a significant improvement. The Cholesky factorization itself can be used to determine the appropriate sparsity preserving reordering of the matrix during the factorization, In that case this might be obtained with three return arguments as r@code{[r, p, q] = chol(A); spy(r)}. @ifset HAVE_CHOLMOD @ifset HAVE_COLAMD @float Figure,fig:simplechol @image{spchol,8cm} @caption{Structure of the un-permuted Cholesky factorization of the above matrix.} @end float @float Figure,fig:simplecholperm @image{spcholperm,8cm} @caption{Structure of the permuted Cholesky factorization of the above matrix.} @end float @end ifset @end ifset In the case of an asymmetric matrix, the appropriate sparsity preserving permutation is @dfn{colamd} and the factorization using this reordering can be visualized using the command @code{q = colamd(A); [l, u, p] = lu(A(:,q)); spy(l+u)}. Finally, Octave implicitly reorders the matrix when using the div (/) and ldiv (\) operators, and so no the user does not need to explicitly reorder the matrix to maximize performance. @node Sparse Linear Algebra, Iterative Techniques, Basics, Sparse Matrices @section Linear Algebra on Sparse Matrices Octave includes a poly-morphic solver for sparse matrices, where the exact solver used to factorize the matrix, depends on the properties of the sparse matrix itself. Generally, the cost of determining the matrix type is small relative to the cost of factorizing the matrix itself, but in any case the matrix type is cached once it is calculated, so that it is not re-determined each time it is used in a linear equation. The selection tree for how the linear equation is solve is @enumerate 1 @item If the matrix is diagonal, solve directly and goto 8 @item If the matrix is a permuted diagonal, solve directly taking into account the permutations. Goto 8 @item If the matrix is square, banded and if the band density is less than that given by @code{spparms ("bandden")} continue, else goto 4. @enumerate a @item If the matrix is tridiagonal and the right-hand side is not sparse continue, else goto 3b. @enumerate @item If the matrix is hermitian, with a positive real diagonal, attempt Cholesky factorization using @sc{Lapack} xPTSV. @item If the above failed or the matrix is not hermitian with a positive real diagonal use Gaussian elimination with pivoting using @sc{Lapack} xGTSV, and goto 8. @end enumerate @item If the matrix is hermitian with a positive real diagonal, attempt Cholesky factorization using @sc{Lapack} xPBTRF. @item if the above failed or the matrix is not hermitian with a positive real diagonal use Gaussian elimination with pivoting using @sc{Lapack} xGBTRF, and goto 8. @end enumerate @item If the matrix is upper or lower triangular perform a sparse forward or backward substitution, and goto 8 @item If the matrix is a upper triangular matrix with column permutations or lower triangular matrix with row permutations, perform a sparse forward or backward substitution, and goto 8 @item If the matrix is square, hermitian with a real positive diagonal, attempt sparse Cholesky factorization using CHOLMOD. @item If the sparse Cholesky factorization failed or the matrix is not hermitian with a real positive diagonal, and the matrix is square, factorize using UMFPACK. @item If the matrix is not square, or any of the previous solvers flags a singular or near singular matrix, find a minimum norm solution using CXSPARSE@footnote{CHOLMOD, UMFPACK and CXSPARSE are written by Tim Davis and are available at http://www.cise.ufl.edu/research/sparse/}. @end enumerate The band density is defined as the number of non-zero values in the matrix divided by the number of non-zero values in the matrix. The banded matrix solvers can be entirely disabled by using @dfn{spparms} to set @code{bandden} to 1 (i.e. @code{spparms ("bandden", 1)}). The QR solver factorizes the problem with a Dulmage-Mendhelsohn, to seperate the problem into blocks that can be treated as over-determined, multiple well determined blocks, and a final over-determined block. For matrices with blocks of strongly connectted nodes this is a big win as LU decomposition can be used for many blocks. It also significantly improves the chance of finding a solution to over-determined problems rather than just returning a vector of @dfn{NaN}'s. All of the solvers above, can calculate an estimate of the condition number. This can be used to detect numerical stability problems in the solution and force a minimum norm solution to be used. However, for narrow banded, triangular or diagonal matrices, the cost of calculating the condition number is significant, and can in fact exceed the cost of factoring the matrix. Therefore the condition number is not calculated in these case, and octave relies on simplier techniques to detect sinular matrices or the underlying LAPACK code in the case of banded matrices. The user can force the type of the matrix with the @code{matrix_type} function. This overcomes the cost of discovering the type of the matrix. However, it should be noted incorrectly identifying the type of the matrix will lead to unpredictable results, and so @code{matrix_type} should be used with care. @node Iterative Techniques, Real Life Example, Sparse Linear Algebra, Sparse Matrices @section Iterative Techniques applied to sparse matrices There are three functions currently to document here, these being @dfn{luinc}, @dfn{pcg} and @dfn{pcr}. WRITE ME. @node Real Life Example, Oct-Files, Iterative Techniques, Sparse Matrices @section Real Life Example of the use of Sparse Matrices A common application for sparse matrices is in the solution of Finite Element Models. Finite element models allow numerical solution of partial differential equations that do not have closed form solutions, typically because of the complex shape of the domain. In order to motivate this application, we consider the boundary value Laplace equation. This system can model scalar potential fields, such as heat or electrical potential. Given a medium @iftex @tex $\Omega$ @end tex @end iftex @ifinfo Omega @end ifinfo with boundary @iftex @tex $\partial\Omega$ @end tex @end iftex @ifinfo dOmega @end ifinfo . At all points on the @iftex @tex $\partial\Omega$ @end tex @end iftex @ifinfo dOmega @end ifinfo the boundary conditions are known, and we wish to calculate the potential in @iftex @tex $\Omega$ @end tex @end iftex @ifinfo Omega @end ifinfo . Boundary conditions may specify the potential (Dirichlet boundary condition), its normal derivative across the boundary (Neumann boundary condition), or a weighted sum of the potential and its derivative (Cauchy boundary condition). In a thermal model, we want to calculate the temperature in @iftex @tex $\Omega$ @end tex @end iftex @ifinfo Omega @end ifinfo and know the boundary temperature (Dirichlet condition) or heat flux (from which we can calculate the Neumann condition by dividing by the thermal conductivity at the boundary). Similarly, in an electrical model, we want to calculate the voltage in @iftex @tex $\Omega$ @end tex @end iftex @ifinfo Omega @end ifinfo and know the boundary voltage (Dirichlet) or current (Neumann condition after diving by the electrical conductivity). In an electrical model, it is common for much of the boundary to be electrically isolated; this is a Neumann boundary condition with the current equal to zero. The simplest finite element models will divide @iftex @tex $\Omega$ @end tex @end iftex @ifinfo Omega @end ifinfo into simplexes (triangles in 2D, pyramids in 3D). @ifset htmltex We take as an 3D example a cylindrical liquid filled tank with a small non-conductive ball from the EIDORS project@footnote{EIDORS - Electrical Impedance Tomography and Diffuse optical Tomography Reconstruction Software @url{http://eidors3d.sourceforge.net}}. This is model is designed to reflect an application of electrical impedance tomography, where current patterns are applied to such a tank in order to image the internal conductivity distribution. In order to describe the FEM geometry, we have a matrix of vertices @code{nodes} and simplices @code{elems}. @end ifset The following example creates a simple rectangular 2D electrically conductive medium with 10 V and 20 V imposed on opposite sides (Dirichlet boundary conditions). All other edges are electrically isolated. @example node_y= [1;1.2;1.5;1.8;2]*ones(1,11); node_x= ones(5,1)*[1,1.05,1.1,1.2, ... 1.3,1.5,1.7,1.8,1.9,1.95,2]; nodes= [node_x(:), node_y(:)]; [h,w]= size(node_x); elems= []; for idx= 1:w-1 widx= (idx-1)*h; elems= [elems; ... widx+[(1:h-1);(2:h);h+(1:h-1)]'; ... widx+[(2:h);h+(2:h);h+(1:h-1)]' ]; endfor E= size(elems,1); # No. of simplices N= size(nodes,1); # No. of vertices D= size(elems,2); # dimensions+1 @end example This creates a N-by-2 matrix @code{nodes} and a E-by-3 matrix @code{elems} with values, which define finite element triangles: @example nodes(1:7,:)' 1.00 1.00 1.00 1.00 1.00 1.05 1.05 ... 1.00 1.20 1.50 1.80 2.00 1.00 1.20 ... elems(1:7,:)' 1 2 3 4 2 3 4 ... 2 3 4 5 7 8 9 ... 6 7 8 9 6 7 8 ... @end example Using a first order FEM, we approximate the electrical conductivity distribution in @iftex @tex $\Omega$ @end tex @end iftex @ifinfo Omega @end ifinfo as constant on each simplex (represented by the vector @code{conductivity}). Based on the finite element geometry, we first calculate a system (or stiffness) matrix for each simplex (represented as 3-by-3 elements on the diagonal of the element-wise system matrix @code{SE}. Based on @code{SE} and a N-by-DE connectivity matrix @code{C}, representing the connections between simplices and vectices, the global connectivity matrix @code{S} is calculated. @example # Element conductivity conductivity= [1*ones(1,16), ... 2*ones(1,48), 1*ones(1,16)]; # Connectivity matrix C = sparse ((1:D*E), reshape (elems', ... D*E, 1), 1, D*E, N); # Calculate system matrix Siidx = floor ([0:D*E-1]'/D) * D * ... ones(1,D) + ones(D*E,1)*(1:D) ; Sjidx = [1:D*E]'*ones(1,D); Sdata = zeros(D*E,D); dfact = factorial(D-1); for j=1:E a = inv([ones(D,1), ... nodes(elems(j,:), :)]); const = conductivity(j) * 2 / ... dfact / abs(det(a)); Sdata(D*(j-1)+(1:D),:) = const * ... a(2:D,:)' * a(2:D,:); endfor # Element-wise system matrix SE= sparse(Siidx,Sjidx,Sdata); # Global system matrix S= C'* SE *C; @end example The system matrix acts like the conductivity @iftex @tex $S$ @end tex @end iftex @ifinfo @code{S} @end ifinfo in Ohm's law @iftex @tex $SV = I$. @end tex @end iftex @ifinfo @code{S * V = I}. @end ifinfo Based on the Dirichlet and Neumann boundary conditions, we are able to solve for the voltages at each vertex @code{V}. @example # Dirichlet boundary conditions D_nodes=[1:5, 51:55]; D_value=[10*ones(1,5), 20*ones(1,5)]; V= zeros(N,1); V(D_nodes) = D_value; idx = 1:N; # vertices without Dirichlet # boundary condns idx(D_nodes) = []; # Neumann boundary conditions. Note that # N_value must be normalized by the # boundary length and element conductivity N_nodes=[]; N_value=[]; Q = zeros(N,1); Q(N_nodes) = N_value; V(idx) = S(idx,idx) \ ( Q(idx) - ... S(idx,D_nodes) * V(D_nodes)); @end example Finally, in order to display the solution, we show each solved voltage value in the z-axis for each simplex vertex. @ifset htmltex @ifset HAVE_CHOLMOD @ifset HAVE_UMFPACK @ifset HAVE_COLAMD @xref{fig:femmodel}. @end ifset @end ifset @end ifset @end ifset @example elemx = elems(:,[1,2,3,1])'; xelems = reshape (nodes(elemx, 1), 4, E); yelems = reshape (nodes(elemx, 2), 4, E); velems = reshape (V(elemx), 4, E); plot3 (xelems,yelems,velems,'k'); print ('grid.eps'); @end example @ifset htmltex @ifset HAVE_CHOLMOD @ifset HAVE_UMFPACK @ifset HAVE_COLAMD @float Figure,fig:femmodel @image{grid,8cm} @caption{Example finite element model the showing triangular elements. The height of each vertex corresponds to the solution value.} @end float @end ifset @end ifset @end ifset @end ifset @node Oct-Files, Function Reference, Real Life Example, Sparse Matrices @section Using Sparse Matrices in Oct-files An oct-file is a means of writing an Octave function in a compilable language like C++, rather than as a script file. This results in a significant acceleration in the code. It is not the purpose of this section to discuss how to write an oct-file, or discuss what they are. There are already two @footnote{Paul Thomas "Dal Segno al Coda - The octave dynamically linked function cookbook", @url{http://perso.wanadoo.fr/prthomas/intro.html}, and Cristophe Spiel "Del Coda Al Fine - Pushing Octave's Limits", @url{http://octave.sourceforge.net/coda/coda.pdf}} very good references on oct-files themselves. Users who are not familiar with oct-files are urged to read these references to fully understand this chapter. The examples discussed here assume that the oct-file is written entirely in C++. There are three classes of sparse objects that are of interest to the user. @table @asis @item SparseMatrix A double precision sparse matrix class @item SparseComplexMatrix A complex sparse matrix class @item SparseBoolMatrix A boolean sparse matrix class @end table All of these classes inherit from the @code{Sparse<T>} template class, and so all have similar capabilities and usage. The @code{Sparse<T>} class was based on Octave @code{Array<T>} class, and so users familiar with Octave's Array classes will be comfortable with the use of the sparse classes. The sparse classes will not be entirely described in this section, due to their similar with the existing Array classes. However, there are a few differences due the different nature of sparse objects, and these will be described. Firstly, although it is fundamentally possible to have N-dimensional sparse objects, the Octave sparse classes do not allow them at this time. So all operations of the sparse classes must be 2-dimensional. This means that in fact @code{SparseMatrix} is similar to Octave's @code{Matrix} class rather than its @code{NDArray} class. @menu * OctDifferences:: The Differences between the Array and Sparse Classes * OctCreation:: Creating Spare Matrices in Oct-Files * OctUse:: Using Sparse Matrices in Oct-Files @end menu @node OctDifferences, OctCreation, Oct-Files, Oct-Files @subsection The Differences between the Array and Sparse Classes The number of elements in a sparse matrix is considered to be the number of non-zero elements rather than the product of the dimensions. Therefore @example SparseMatrix sm; @dots{} int nel = sm.nelem (); @end example returns the number of non-zero elements. If the user really requires the number of elements in the matrix, including the non-zero elements, they should use @code{numel} rather than @code{nelem}. Note that for very large matrices, where the product of the two dimensions is large that the representation of the an unsigned int, then @code{numel} can overflow. An example is @code{speye(1e6)} which will create a matrix with a million rows and columns, but only a million non-zero elements. Therefore the number of rows by the number of columns in this case is more than two hundred times the maximum value that can be represented by an unsigned int. The use of @code{numel} should therefore be avoided useless it is known it won't overflow. Extreme care must be take with the elem method and the "()" operator, which perform basically the same function. The reason is that if a sparse object is non-const, then Octave will assume that a request for a zero element in a sparse matrix is in fact a request to create this element so it can be filled. Therefore a piece of code like @example SparseMatrix sm; @dots{} for (int j = 0; j < nc; j++) for (int i = 0; i < nr; i++) std::cerr << " (" << i << "," << j << "): " << sm(i,j) << std::endl; @end example is a great way of turning the sparse matrix into a dense one, and a very slow way at that since it reallocates the sparse object at each zero element in the matrix. An easy way of preventing the above from happening is to create a temporary constant version of the sparse matrix. Note that only the container for the sparse matrix will be copied, while the actual representation of the data will be shared between the two versions of the sparse matrix. So this is not a costly operation. For example, the above would become @example SparseMatrix sm; @dots{} const SparseMatrix tmp (sm); for (int j = 0; j < nc; j++) for (int i = 0; i < nr; i++) std::cerr << " (" << i << "," << j << "): " << tmp(i,j) << std::endl; @end example Finally, as the sparse types aren't just represented as a contiguous block of memory, the @code{fortran_vec} method of the @code{Array<T>} is not available. It is however replaced by three separate methods @code{ridx}, @code{cidx} and @code{data}, that access the raw compressed column format that the Octave sparse matrices are stored in. Additionally, these methods can be used in a manner similar to @code{elem}, to allow the matrix to be accessed or filled. However, in that case it is up to the user to respect the sparse matrix compressed column format discussed previous. @node OctCreation, OctUse, OctDifferences, Oct-Files @subsection Creating Spare Matrices in Oct-Files The user has several alternatives in how to create a sparse matrix. They can first create the data as three vectors representing the row and column indexes and the data, and from those create the matrix. Or alternatively, they can create a sparse matrix with the appropriate amount of space and then fill in the values. Both techniques have their advantages and disadvantages. An example of how to create a small sparse matrix with the first technique might be seen the example @example int nz = 4, nr = 3, nc = 4; ColumnVector ridx (nz); ColumnVector cidx (nz); ColumnVector data (nz); ridx(0) = 0; ridx(1) = 0; ridx(2) = 1; ridx(3) = 2; cidx(0) = 0; cidx(1) = 1; cidx(2) = 3; cidx(3) = 3; data(0) = 1; data(1) = 2; data(2) = 3; data(3) = 4; SparseMatrix sm (data, ridx, cidx, nr, nc); @end example which creates the matrix given in section @ref{Storage}. Note that the compressed matrix format is not used at the time of the creation of the matrix itself, however it is used internally. As previously mentioned, the values of the sparse matrix are stored in increasing column-major ordering. Although the data passed by the user does not need to respect this requirement, the pre-sorting the data significantly speeds up the creation of the sparse matrix. The disadvantage of this technique of creating a sparse matrix is that there is a brief time where two copies of the data exists. Therefore for extremely memory constrained problems this might not be the right technique to create the sparse matrix. The alternative is to first create the sparse matrix with the desired number of non-zero elements and then later fill those elements in. The easiest way to do this is @example int nz = 4, nr = 3, nc = 4; SparseMatrix sm (nr, nc, nz); sm(0,0) = 1; sm(0,1) = 2; sm(1,3) = 3; sm(2,3) = 4; @end example That creates the same matrix as previously. Again, although it is not strictly necessary, it is significantly faster if the sparse matrix is created in this manner that the elements are added in column-major ordering. The reason for this is that if the elements are inserted at the end of the current list of known elements then no element in the matrix needs to be moved to allow the new element to be inserted. Only the column indexes need to be updated. There are a few further points to note about this technique of creating a sparse matrix. Firstly, it is not illegal to create a sparse matrix with fewer elements than are actually inserted in the matrix. Therefore @example int nz = 4, nr = 3, nc = 4; SparseMatrix sm (nr, nc, 0); sm(0,0) = 1; sm(0,1) = 2; sm(1,3) = 3; sm(2,3) = 4; @end example is perfectly legal. However it is a very bad idea. The reason is that as each new element is added to the sparse matrix the space allocated to it is increased by reallocating the memory. This is an expensive operation, that will significantly slow this means of creating a sparse matrix. Furthermore, it is not illegal to create a sparse matrix with too much storage, so having @var{nz} above equaling 6 is also legal. The disadvantage is that the matrix occupies more memory than strictly needed. It is not always easy to known the number of non-zero elements prior to filling a matrix. For this reason the additional storage for the sparse matrix can be removed after its creation with the @dfn{maybe_compress} function. Furthermore, the maybe_compress can deallocate the unused storage, but it can equally remove zero elements from the matrix. The removal of zero elements from the matrix is controlled by setting the argument of the @dfn{maybe_compress} function to be 'true'. However, the cost of removing the zeros is high because it implies resorting the elements. Therefore, if possible it is better is the user doesn't add the zeros in the first place. An example of the use of @dfn{maybe_compress} is @example int nz = 6, nr = 3, nc = 4; SparseMatrix sm1 (nr, nc, nz); sm1(0,0) = 1; sm1(0,1) = 2; sm1(1,3) = 3; sm1(2,3) = 4; sm1.maybe_compress (); // No zero elements were added SparseMatrix sm2 (nr, nc, nz); sm2(0,0) = 1; sm2(0,1) = 2; sm(0,2) = 0; sm(1,2) = 0; sm1(1,3) = 3; sm1(2,3) = 4; sm2.maybe_compress (true); // Zero elements were added @end example The use of the @dfn{maybe_compress} function should be avoided if possible, as it will slow the creation of the matrices. A third means of creating a sparse matrix is to work directly with the data in compressed row format. An example of this technique might be @c Note the @verbatim environment is a relatively new addition to texinfo. @c Therefore use the @example environment and replace @, with @@, @c { with @{, etc @example octave_value arg; @dots{} int nz = 6, nr = 3, nc = 4; // Assume we know the max no nz SparseMatrix sm (nr, nc, nz); Matrix m = arg.matrix_value (); int ii = 0; sm.cidx (0) = 0; for (int j = 1; j < nc; j++) @{ for (int i = 0; i < nr; i++) @{ double tmp = foo (m(i,j)); if (tmp != 0.) @{ sm.data(ii) = tmp; sm.ridx(ii) = i; ii++; @} @} sm.cidx(j+1) = ii; @} sm.maybe_compress (); // If don't know a-priori the final no of nz. @end example which is probably the most efficient means of creating the sparse matrix. Finally, it might sometimes arise that the amount of storage initially created is insufficient to completely store the sparse matrix. Therefore, the method @code{change_capacity} exists to reallocate the sparse memory. The above example would then be modified as @example octave_value arg; @dots{} int nz = 6, nr = 3, nc = 4; // Assume we know the max no nz SparseMatrix sm (nr, nc, nz); Matrix m = arg.matrix_value (); int ii = 0; sm.cidx (0) = 0; for (int j = 1; j < nc; j++) @{ for (int i = 0; i < nr; i++) @{ double tmp = foo (m(i,j)); if (tmp != 0.) @{ if (ii == nz) @{ nz += 2; // Add 2 more elements sm.change_capacity (nz); @} sm.data(ii) = tmp; sm.ridx(ii) = i; ii++; @} @} sm.cidx(j+1) = ii; @} sm.maybe_mutate (); // If don't know a-priori the final no of nz. @end example Note that both increasing and decreasing the number of non-zero elements in a sparse matrix is expensive, as it involves memory reallocation. Also as parts of the matrix, though not its entirety, exist as the old and new copy at the same time, additional memory is needed. Therefore if possible this should be avoided. @node OctUse, , OctCreation, Oct-Files @subsection Using Sparse Matrices in Oct-Files Most of the same operators and functions on sparse matrices that are available from the Octave are equally available with oct-files. The basic means of extracting a sparse matrix from an @code{octave_value} and returning them as an @code{octave_value}, can be seen in the following example @example octave_value_list retval; SparseMatrix sm = args(0).sparse_matrix_value (); SparseComplexMatrix scm = args(1).sparse_complex_matrix_value (); SparseBoolMatrix sbm = args(2).sparse_bool_matrix_value (); @dots{} retval(2) = sbm; retval(1) = scm; retval(0) = sm; @end example The conversion to an octave-value is handled by the sparse @code{octave_value} constructors, and so no special care is needed. @node Function Reference, , Oct-Files, Sparse Matrices @section Function Reference @ifset htmltex @subsection Functions by Category @subsubsection Generate sparse matrix @table @asis @item @ref{spdiags} A generalization of the function `spdiag'. @item @ref{speye} Returns a sparse identity matrix. @item @ref{sprand} Generate a random sparse matrix. @item @ref{sprandn} Generate a random sparse matrix. @item @ref{sprandsym} Generate a symmetric random sparse matrix. @end table @subsubsection Sparse matrix conversion @table @asis @item @ref{full} returns a full storage matrix from a sparse one See also: sparse @item @ref{sparse} SPARSE: create a sparse matrix @item @ref{spconvert} This function converts for a simple sparse matrix format easily produced by other programs into Octave's internal sparse format. @item @ref{spfind} SPFIND: a sparse version of the find operator 1. @end table @subsubsection Manipulate sparse matrices @table @asis @item @ref{issparse} Return 1 if the value of the expression EXPR is a sparse matrix. @item @ref{nnz} returns number of non zero elements in SM See also: sparse @item @ref{nonzeros} Returns a vector of the non-zero values of the sparse matrix S @item @ref{nzmax} Returns the amount of storage allocated to the sparse matrix SM. @item @ref{spalloc} Returns an empty sparse matrix of size R-by-C. @item @ref{spfun} Compute `f(X)' for the non-zero values of X This results in a sparse matrix with the same structure as X. @item @ref{spones} Replace the non-zero entries of X with ones. @item @ref{spy} Plot the sparsity pattern of the sparse matrix X @end table @subsubsection Graph Theory @table @asis @item @ref{etree} Returns the elimination tree for the matrix S. @item @ref{etreeplot} Plots the elimination tree of the matrix @var{s} or @code{@var{s}+@var{s}'} if @var{s} in non-symmetric. @item @ref{gplot} Plots a graph defined by @var{A} and @var{xy} in the graph theory sense. @item treelayout @emph{Not implemented} @item @ref{treeplot} Produces a graph of a tree or forest. @end table @subsubsection Sparse matrix reordering @table @asis @item @ref{ccolamd} Constrained column approximate minimum degree permutation. @item @ref{colamd} Column approximate minimum degree permutation. @item @ref{colperm} Returns the column permutations such that the columns of `S (:, P)' are ordered in terms of increase number of non-zero elements. @item @ref{csymamd} For a symmetric positive definite matrix S, returns the permutation vector p such that `S (P, P)' tends to have a sparser Cholesky factor than S. @item @ref{dmperm} Perform a Deulmage-Mendelsohn permutation on the sparse matrix S. @item @ref{symamd} For a symmetric positive definite matrix S, returns the permutation vector p such that `S (P, P)' tends to have a sparser Cholesky factor than S. @item symrcm @emph{Not implemented} @end table @subsubsection Linear algebra @table @asis @item cholinc @emph{Not implemented} @item condest @emph{Not implemented} @item eigs @emph{Not implemented} @item @ref{matrix_type} Identify the matrix type or mark a matrix as a particular type. @item normest @emph{Not implemented} @item @ref{spchol} Compute the Cholesky factor, R, of the symmetric positive definite. @item @ref{spcholinv} Use the Cholesky factorization to compute the inverse of the sparse symmetric positive definite matrix A. @item @ref{spchol2inv} Invert a sparse symmetric, positive definite square matrix from its Cholesky decomposition, U. @item @ref{spdet} Compute the determinant of sparse matrix A using UMFPACK. @item @ref{spinv} Compute the inverse of the square matrix A. @item @ref{spkron} Form the kronecker product of two sparse matrices. @item @ref{splchol} Compute the Cholesky factor, L, of the symmetric positive definite. @item @ref{splu} Compute the LU decomposition of the sparse matrix A, using subroutines from UMFPACK. @item @ref{spqr} Compute the sparse QR factorization of @var{a}, using CSPARSE. @item sprank @emph{Not implemented} @item svds @emph{Not implemented} @end table @subsubsection Iterative techniques @table @asis @item bicg @emph{Not implemented} @item bicgstab @emph{Not implemented} @item cgs @emph{Not implemented} @item gmres @emph{Not implemented} @item @ref{luinc} Produce the incomplete LU factorization of the sparse matrix A. @item lsqr @emph{Not implemented} @item minres @emph{Not implemented} @item pcg Solves the linear system of equations @code{@var{A} * @var{x} = @var{b}} by means of the Preconditioned Conjugate Gradient iterative method. @item pcr Solves the linear system of equations @code{@var{A} * @var{x} = @var{b}} by means of the Preconditioned Conjugate Residual iterative method. @item qmr @emph{Not implemented} @item symmlq @emph{Not implemented} @end table @subsubsection Miscellaneous @table @asis @item spaugment @emph{Not implemented} @item @ref{spparms} Sets or displays the parameters used by the sparse solvers and factorization functions. @item @ref{symbfact} Performs a symbolic factorization analysis on the sparse matrix S. @item @ref{spstats} Return the stats for the non-zero elements of the sparse matrix S COUNT is the number of non-zeros in each column, MEAN is the mean of the non-zeros in each column, and VAR is the variance of the non-zeros in each column @item @ref{spprod} Product of elements along dimension DIM. @item @ref{spcumprod} Cumulative product of elements along dimension DIM. @item @ref{spcumsum} Cumulative sum of elements along dimension DIM. @item @ref{spsum} Sum of elements along dimension DIM. @item @ref{spsumsq} Sum of squares of elements along dimension DIM. @item @ref{spmin} For a vector argument, return the minimum value. @item @ref{spmax} For a vector argument, return the maximum value. @item @ref{spatan2} Compute atan (Y / X) for corresponding sparse matrix elements of Y and X. @item @ref{spdiag} Return a diagonal matrix with the sparse vector V on diagonal K. @end table @subsection Functions Alphabetically @end ifset @menu * ccolamd:: Constrained column approximate minimum degree permutation. * colamd:: Column approximate minimum degree permutation. * colperm:: Returns the column permutations such that the columns of `S (:, P)' are ordered in terms of increase number of non-zero elements. * csymamd:: For a symmetric positive definite matrix S, returns the permutation vector p such that `S (P, P)' tends to have a sparser Cholesky factor than S. * dmperm:: Perfrom a Deulmage-Mendelsohn permutation on the sparse matrix S. * etree:: Returns the elimination tree for the matrix S. * etreeplot:: Plots the elimination tree of the matrix @var{s} or @code{@var{s}+@var{s}'} if @var{s} in non-symmetric. * full:: returns a full storage matrix from a sparse one See also: sparse * gplot:: Plots a graph defined by @var{A} and @var{xy} in the graph theory sense. * issparse:: Return 1 if the value of the expression EXPR is a sparse matrix. * luinc:: Produce the incomplete LU factorization of the sparse A. * matrix_type:: Identify the matrix type or mark a matrix as a particular type. * nnz:: returns number of non zero elements in SM See also: sparse * nonzeros:: Returns a vector of the non-zero values of the sparse matrix S * nzmax:: Returns the amount of storage allocated to the sparse matrix SM. * pcg:: Solves linear system of equations by means of the Preconditioned Conjugate Gradient iterative method. * pcr:: Solves linear system of equations by means of the Preconditioned Conjugate Residual iterative method. * spalloc:: Returns an empty sparse matrix of size R-by-C. * sparse:: SPARSE: create a sparse matrix * spatan2:: Compute atan (Y / X) for corresponding sparse matrix elements of Y and X. * spchol:: Compute the Cholesky factor, R, of the symmetric positive definite. * spcholinv:: Use the Cholesky factorization to compute the inverse of the sparse symmetric positive definite matrix A. * spchol2inv:: Invert a sparse symmetric, positive definite square matrix from its Cholesky decomposition, U. * spconvert:: This function converts for a simple sparse matrix format easily produced by other programs into Octave's internal sparse format. * spcumprod:: Cumulative product of elements along dimension DIM. * spcumsum:: Cumulative sum of elements along dimension DIM. * spdet:: Compute the determinant of sparse matrix A using UMFPACK. * spdiag:: Return a diagonal matrix with the sparse vector V on diagonal K. * spdiags:: A generalization of the function `spdiag'. * speye:: Returns a sparse identity matrix. * spfind:: SPFIND: a sparse version of the find operator 1. * spfun:: Compute `f(X)' for the non-zero values of X This results in a sparse matrix with the same structure as X. * spinv:: Compute the inverse of the square matrix A. * spkron:: Form the kronecker product of two sparse matrices. * splchol:: Compute the Cholesky factor, L, of the symmetric positive definite. * splu:: Compute the LU decomposition of the sparse matrix A, using subroutines from UMFPACK. * spmax:: For a vector argument, return the maximum value. * spmin:: For a vector argument, return the minimum value. * spones:: Replace the non-zero entries of X with ones. * spparms:: Sets or displays the parameters used by the sparse solvers and factorization functions. * spprod:: Product of elements along dimension DIM. * spqr:: Compute the sparse QR factorization of @var{a}, using CSPARSE. * sprand:: Generate a random sparse matrix. * sprandn:: Generate a random sparse matrix. * sprandsym:: Generate a symmetric random sparse matrix. * spstats:: Return the stats for the non-zero elements of the sparse matrix S COUNT is the number of non-zeros in each column, MEAN is the mean of the non-zeros in each column, and VAR is the variance of the non-zeros in each column * spsum:: Sum of elements along dimension DIM. * spsumsq:: Sum of squares of elements along dimension DIM. * spy:: Plot the sparsity pattern of the sparse matrix X * symamd:: For a symmetric positive definite matrix S, returns the permutation vector p such that `S (P, P)' tends to have a sparser Cholesky factor than S. * symbfact:: Performs a symbolic factorization analysis on the sparse matrix S. * treeplot:: Produces a graph of a tree or forest. @end menu @node colamd, ccolamd, , Function Reference @subsubsection colamd @DOCSTRING(colamd) @node ccolamd, colperm, colamd, Function Reference @subsubsection ccolamd @DOCSTRING(ccolamd) @node colperm, csymamd, ccolamd, Function Reference @subsubsection colperm @DOCSTRING(colperm) @node csymamd, dmperm, colperm, Function Reference @subsubsection csymamd @DOCSTRING(csymamd) @node dmperm, etree, csymamd, Function Reference @subsubsection dmperm @DOCSTRING(dmperm) @node etree, etreeplot, dmperm, Function Reference @subsubsection etree @DOCSTRING(etree) @node etreeplot, full, etree, Function Reference @subsubsection etreeplot @DOCSTRING(etreeplot) @node full, gplot, etreeplot, Function Reference @subsubsection full @DOCSTRING(full) @node gplot, issparse, full, Function Reference @subsubsection gplot @DOCSTRING(gplot) @node issparse, luinc, gplot, Function Reference @subsubsection issparse @DOCSTRING(issparse) @node luinc, matrix_type, issparse, Function Reference @subsubsection luinc @DOCSTRING(luinc) @node matrix_type, nnz, luinc, Function Reference @subsubsection matrix_type @DOCSTRING(matrix_type) @node nnz, nonzeros, matrix_type, Function Reference @subsubsection nnz @DOCSTRING(nnz) @node nonzeros, nzmax, nnz, Function Reference @subsubsection nonzeros @DOCSTRING(nonzeros) @node nzmax, pcg, nonzeros, Function Reference @subsubsection nzmax @DOCSTRING(nzmax) @node pcg, pcr, nzmax, Function Reference @subsubsection pcg @DOCSTRING(pcg) @node pcr, spalloc, pcg, Function Reference @subsubsection pcr @DOCSTRING(pcr) @node spalloc, sparse, pcr, Function Reference @subsubsection spalloc @DOCSTRING(spalloc) @node sparse, spatan2, spalloc, Function Reference @subsubsection sparse @DOCSTRING(sparse) @node spatan2, spchol, sparse, Function Reference @subsubsection spatan2 @DOCSTRING(spatan2) @node spchol, spcholinv, spatan2, Function Reference @subsubsection spchol @DOCSTRING(spchol) @node spcholinv, spchol2inv, spchol, Function Reference @subsubsection spcholinv @DOCSTRING(spcholinv) @node spchol2inv, spconvert, spcholinv, Function Reference @subsubsection spchol2inv @DOCSTRING(spchol2inv) @node spconvert, spcumprod, spchol2inv, Function Reference @subsubsection spconvert @DOCSTRING(spconvert) @node spcumprod, spcumsum, spconvert, Function Reference @subsubsection spcumprod @DOCSTRING(spcumprod) @node spcumsum, spdet, spcumprod, Function Reference @subsubsection spcumsum @DOCSTRING(spcumsum) @node spdet, spdiag, spcumsum, Function Reference @subsubsection spdet @DOCSTRING(spdet) @node spdiag, spdiags, spdet, Function Reference @subsubsection spdiag @DOCSTRING(spdiag) @node spdiags, speye, spdiag, Function Reference @subsubsection spdiags @DOCSTRING(spdiags) @node speye, spfind, spdiags, Function Reference @subsubsection speye @DOCSTRING(speye) @node spfind, spfun, speye, Function Reference @subsubsection spfind @DOCSTRING(spfind) @node spfun, spinv, spfind, Function Reference @subsubsection spfun @DOCSTRING(spfun) @node spinv, spkron, spfun, Function Reference @subsubsection spinv @DOCSTRING(spinv) @node spkron, splchol, spinv, Function Reference @subsubsection spkron @DOCSTRING(spkron) @node splchol, splu, spkron, Function Reference @subsubsection splchol @DOCSTRING(splchol) @node splu, spmax, splchol, Function Reference @subsubsection splu @DOCSTRING(splu) @node spmax, spmin, splu, Function Reference @subsubsection spmax @DOCSTRING(spmax) @node spmin, spones, spmax, Function Reference @subsubsection spmin @DOCSTRING(spmin) @node spones, spparms, spmin, Function Reference @subsubsection spones @DOCSTRING(spones) @node spparms, spprod, spones, Function Reference @subsubsection spparms @DOCSTRING(spparms) @node spprod, spqr, spparms, Function Reference @subsubsection spprod @DOCSTRING(spprod) @node spqr, sprand, spprod, Function Reference @subsubsection spqr @DOCSTRING(spqr) @node sprand, sprandn, spqr, Function Reference @subsubsection sprand @DOCSTRING(sprand) @node sprandn, sprandsym, sprand, Function Reference @subsubsection sprandn @DOCSTRING(sprandn) @node sprandsym, spstats, sprandn, Function Reference @subsubsection sprandsym @DOCSTRING(sprandsym) @node spstats, spsum, sprandsym, Function Reference @subsubsection spstats @DOCSTRING(spstats) @node spsum, spsumsq, spstats, Function Reference @subsubsection spsum @DOCSTRING(spsum) @node spsumsq, spy, spsum, Function Reference @subsubsection spsumsq @DOCSTRING(spsumsq) @node spy, symamd, spsumsq, Function Reference @subsubsection spy @DOCSTRING(spy) @node symamd, symbfact, spy, Function Reference @subsubsection symamd @DOCSTRING(symamd) @node symbfact, treeplot, symamd, Function Reference @subsubsection symbfact @DOCSTRING(symbfact) @node treeplot, ,symbfact, Function Reference @subsubsection treeplot @DOCSTRING(treeplot) @c Local Variables: *** @c Mode: texinfo *** @c End: ***