Solving Poisson PDE with Sparse Matrices
Abstract
In this page, we present the resolution of the Poisson Partial Differential Equation in Scilab with sparse matrices. We show that Scilab 5 can solve in a few seconds sparse linear systems of equations with as many as 250 000 unknowns because Scilab only store nonzero entries. The computations are based on the Scibench module, a toolbox which provides a collection of benchmarks for Scilab.
Contents
Introduction
Sharma and Gobbert analyzed the performance of Scilab for the resolution of sparse linear systems of equations associated with the Poisson equation [1]. In this document, we try to reproduce their experiments.
We consider the Poisson problem with homogeneous Dirichlet boundary conditions and are interested in the numerical solution based on finite differences.
We consider the 2 dimensional Partial Differential Equation:
where the two dimensionnal Laplace operator is
We consider the domain , . The function f is defined by
The solution is
We use a second order finite difference approximation of the Laplace operator based on a grid of N-by-N points.
Sparse matrices can be managed in Scilab since 1989 [6]. In Scilab 5.0, the UMFPACK module was added New Scientific Features in 2008. More details on this topic can be found in [3].
The Scibench external module:
http://atoms.scilab.org/toolboxes/scibench
provides a collection of benchmark scripts to measure the performance of Scilab. For example, it contains benchmarks for the dense matrix-matrix product, the dense backslash operator, the Cholesky decomposition or 2D Lattice Boltzmann simulations.
In order to install this module, we use the statement:
atomsInstall("scibench")
and restart Scilab.
The version 0.6 includes benchmarks for sparse matrices, based on the Poisson equation. The performance presented in this document are measured with Scilab 5.3.2 on Windows XP with a 4GB computer using Intel Xeon E5410 4*2.33 GHz processors.
In this benchmark, we are going to use the "poisson.sce" script provided at:
http://forge.scilab.org/index.php/p/scibench/source/tree/HEAD/demos/poisson.sce
The sparse matrix
The scibench_poissonA(n) statement creates a sparse matrix equation associated with n cells in the X coordinate and n cells in the Y coordinate. Before calling this session, we call the stacksize function in order to let Scilab allocate as much memory as possible. Then we call the PlotSparse function to plot the sparsity pattern of the matrix.
stacksize("max"); A = scibench_poissonA(50); PlotSparse(A)
The previous script produces the following plot.
We emphasize that the previous matrix is a - by - sparse matrix. Even for moderate values of n, this creates huge matrices. Fortunately, only nonzero entries are stored, so that Scilab has no problem to store it, provided that the nonzero entries can be stored in the memory.
-->size(A) ans = 2500. 2500.
The Kronecker operator is used, so that the computation of the matrix is vectorized. To edit the code, we just run the following code.
-->edit scibench_poissonA
At the core of the algorithm, we find the statement
A = I .*. T + T .*. I
where T is a sparse matrix and I is the sparse n-by-n identity matrix. Hence, creating the whole matrix is done with only one statement.
The Poisson Solver
The scibench_poisson function solves the 2D Poisson equation. Its calling sequence is
scibench_poisson ( N , plotgraph , verbose , solver )
where N is the number of cells, plotgraph is a boolean to plot the graphics, verbose is a boolean to print messages and solver is a function which solves the linear equation A*x=b.
The solver argument is designed so that we can customize the linear equation solver which we want to use. For example, if we want to use the sparse backslash operator, all we have to do is to create the mysolverBackslash function as below.
function u=mysolverBackslash(N, b) A = scibench_poissonA(N); u = A\b; endfunction
The following script solves the Poisson equation with N = 50.
scf(); scibench_poisson(50, %t , %t , mysolverBackslash );
The previous script produces the following output, where h is the dimensionnal spacee step and enorminf is the value of the infinite norm of the error.
-->scibench_poisson(50, %t , %t , mysolverBackslash ); N = 50 h = 1.9607843137254902e-002 h^2 = 3.8446751249519417e-004 enorminf = 1.2634082165868810e-003 C = enorminf / h^2 = 3.2861247713424775e+000 wall clock time = 0.15 seconds
The previous script also produces the following plot.
It is then easy to use the pcg function built in Scilab, which uses a pre-conditionned conjugate gradient algorithm. We use the scibench_poissonAu function, which computes the A*u product without actually storing the matrix A.
function u=mysolverPCG(N, b) tol = 0.000001; maxit = 9999; u = zeros(N^2,1); [u,flag,iter,res] = pcg(scibench_poissonAu,b,tol,maxit,[],[],u); endfunction scf(); [timesec,enorminf,h] = scibench_poisson(50, %f , %t , mysolverPCG );
The benchmark
Based on the scibench module, it is easy to compare various sparse linear equation solvers in Scilab. We compared the following functions:
- sparse backslash (which internally use a sparse LU decomposition),
- pcg function,
- umfpack module,
- the TAUCS module.
Both the UMFPACK and TAUCS modules were created by Bruno Pincon [5].
In order to use the UMFPACK module, we created the following mysolverUMF function which solves the A*u=b equation.
function u = mysolverUMF(N,b) A = scibench_poissonA(N); humf = umf_lufact(A); u = umf_lusolve(humf,b) umf_ludel(humf) endfunction
We also created the following wrapper for the TAUCS module.
function u = mysolverTAUCS(N,b) A = scibench_poissonA(N); hchol = taucs_chfact(A); u = taucs_chsolve(hchol,b) taucs_chdel(hchol) endfunction
We were unable to use the TAUCS solver, which makes Scilab unstable in this case. This problem was reported in the following bug report:
http://bugzilla.scilab.org/show_bug.cgi?id=8824
It is straightforward to created increasingly large matrices and to measure the time required to solve the Poisson equation. The script is provided within the demos of the scibench module. The following plot compares the various solvers.
It is obvious that, in this case, the UMFPACK module is much faster.
The fact that we use the pcg function without actually preconditionning the matrix is an obvious limitation of this benchmark. This is why we work on updating the sparse ILU preconditionners in the following Forge project:
http://forge.scilab.org/index.php/p/spilu/
The current work is based on the former Scilin project.
This benchmark shows that we can solve really large systems of equations. The following table displays some performance measures.
Solver |
N |
Matrix Size |
Time (s) |
Backslash |
50 |
2500-by-2500 |
0.15 |
PCG |
50 |
2500-by-2500 |
0.09 |
Backslash |
100 |
10000-by-10000 |
4.32 |
PCG |
100 |
10000-by-10000 |
0.32 |
Backslash |
200 |
40000-by-40000 |
88.89 |
PCG |
200 |
40000-by-40000 |
2.03 |
UMFPACK |
200 |
40000-by-40000 |
0.81 |
PCG |
300 |
90000-by-90000 |
7.27 |
UMFPACK |
300 |
90000-by-90000 |
2.35 |
PCG |
500 |
250000-by-250000 |
44.77 |
For matrices larger than approximately 400, the UMFPACK functions fails to solve the equation because they require more memory than Scilab can provide to it.
--> scibench_poisson(400, %t , %t , mysolverUMF ) !--error 999 umf_lufact: An error occurred: symbolic factorization: not enough memory at line 3 of function mysolverUMF called by : at line 95 of function scibench_poisson called by : scibench_poisson(400, %t , %t , mysolverUMF )
We emphasize that the actual size of the matrix does not matter. What matters is the number of nonzero terms in the matrix. For example, with N=500, the matrix is 250000-by-250000 but has only 1248000 nonzero entries.
-->A = scibench_poissonA ( 500 ); -->size(A) ans = 250000. 250000. -->nnz(A) ans = 1248000.
Conclusion
Scilab 5 can manage sparse matrices and solve partial differential equations such as the Poisson equation for example. Indeed, Scilab provides several sparse linear equation solvers, including a sparse backslash operator, iterative methods (e.g. the pre-conditionned conjugate gradient algorithm) and other direct solvers such as the UMFPACK or TAUCS modules. With these tools we can solve huge systems of linear equations, because Scilab only stores the nonzero entries of these massive matrices. The solvers may fail if the memory required to store the matrix is beyond the capacity of Scilab, but this happens only for huge matrices. In this case, the future Scilab 6 may help to overcome this limitation, by removing the use of the stack which is used by Scilab 5 (see [4] for more details on this topic). Another limitation of Scilab is that there is currently no preconditionner for sparse linear systems of equations. This limitation should be removed once the spilu module is ready for a release.
Bibliography
[1] "A Comparative Evaluation Of Matlab, Octave, Freemat, And Scilab For Research And Teaching", 2010, Neeraj Sharma and Matthias K. Gobbert, http://userpages.umbc.edu/~gobbert/papers/SharmaGobbertTR2010.pdf
[2] "Scibench, A collection of benchmarks", Consortium Scilab - DIGITEO, Michael Baudin, 2011, http://atoms.scilab.org/toolboxes/scibench
[3] "Introduction to Sparse Matrices in Scilab", Consortium Scilab - DIGITEO, Michael Baudin, 2011, http://forge.scilab.org/index.php/p/docscisparse/downloads/
[4] "Programming in Scilab", Consortium Scilab - DIGITEO, Michael Baudin, 2011, http://forge.scilab.org/index.php/p/docprogscilab/downloads/
[5] "Scilab tools for PDE’s : Application to time-reversal", Bruno Pincon and Karim Ramdani, IEEE International Symposium on Computer Aided Control Systems Design, 2004 (PDF) (HTML)
- [6] Basile, un système de CAO pour l'automatique, version 3.0 : guide d'utilisation, François Delebecque, Carlos Klimann, Serge Steer, INRIA, 1989