FDMcode
FDMcode
LONG CHEN
We discuss efficient ways of implementing finite difference methods for solving the
Poisson equation on rectangular domains in two and three dimensions. The key is the ma-
trix indexing instead of the traditional linear indexing. With such an indexing system, we
will introduce a matrix-free and a tensor product matrix implementation of finite difference
methods.
Since MATLAB is an interpret language, every line will be complied when it is exe-
cuted. A general guideline for efficient programming in MATLAB is:
avoid large for loops.
A simple modification of the above double loops is to use the vector indexing.
1 i = 2:m-1;
2 j = 2:n-1;
3 Au(i,j) = 4*u(i,j) - u(i-1,j) - u(i+1,j) - u(i,j-1) - u(i,j+1);
To evaluate the right hand side, we can use coordinates (x,y) in the matrix form.
For example, for f (x, y) = 8π 2 sin(2πx) cos(2πy), the h2 scaled right hand side can be
computed as
1 [x,y] = ndgrid(0:h:1,0:h:1);
2 fh2 = hˆ2*8*piˆ2*sin(2*pi*x).*cos(2*pi*y);
Note that .* is used to compute the component-wise product for two matrices. For non-
homogenous boundary conditions, one needs to evaluate boundary values and add to the
right hand side. The evaluation of a function on the whole grid is of complexity O(m × n).
For boundary condition, we can reduce to O(m + n) by restricting to bdidx only.
1
2 LONG CHEN
1 u(bdidx) = sin(2*pi*x(bdidx)).*cos(2*pi*y(bdidx));
The array bdidx for collecting all boundary nodes can be generated as follows
1 isbd = true(size(u));
2 isbd(2:end-1,2:end-1) = false;
3 bdidx = find(isbd(:));
One Jacobi iteration for solving the matrix equation Au = f can be implemented as
1 j = 2:n-1;
2 i = 2:m-1;
3 u(i,j) = (fh2(i,j) + u(i-1,j) + u(i+1,j) + u(i,j-1) + u(i,j+1))/4;
A more efficient iterative methods, Gauss-Seidel (G-S) iteration updates the coordinates
sequentially one at a time. Here is the implementation using for loops.
1 for j = 2:n-1
2 for i = 2:m-1
3 u(i,j) = (fh2(i,j) + u(i-1,j) + u(i+1,j) + u(i,j-1) + u(i,j+1))/4;
4 end
5 end
The ordering does matter in the Gauss-Seidel iteration. The backwards G-S can be imple-
mented by inverse the ordering of i,j indexing.
1 for j = n-1:-1:2
2 for i = m-1:-1:2
3 u(i,j) = (fh2(i,j) + u(i-1,j) + u(i+1,j) + u(i,j-1) + u(i,j+1))/4;
4 end
5 end
Note that for the matrix-free implementation, there is no need to modify the right hand
side for the Dirichlet boundary condition. The boundary values of u is assigned before the
iteration and remains the same since only the interior nodal values are updated during the
iteration. For Neumann boundary conditions, an additional update on boundary nodes is
needed.
The symmetric version Gauss-Seidel will be the combination of one forwards and one
backwards GS iteration and is an SPD operator which can be used in pcg to accelerate the
computation of an approximated solution to the linear system Au = f .
Vectorization of Gauss-Seidel iteration is subtle. If we simply remove the for loops,
it is the Jacobi iteration since the values of u on the right hand side is the old one. To
vectorize G-S, let us first classify the nodes into two category: red nodes and black nodes;
see Fig 1. Black nodes can be identified as mod(i+j,2) == 0. A crucial observation is
that to update red nodes only values of black nodes are needed and vice verse. Then Gauss-
Seidel iteration applied to this red-black ordering can be implemented as Jacobi iterations.
1 [m,n] = size(u);
2 % case 1 (red points): mod(i+j,2) == 0
Red-Black Gauss-Seidel
PROGRAMMING OF FINITE DIFFERENCE METHODS IN MATLAB 3
Red depends
F IGURE 1. only on black,
Red-Black and
Ordering vice-versa.
of vertices
Generalization: multi-color orderings
3 i = 2:2:m-1; j = 2:2:n-1;
4 u(i,j) = (fh2(i,j) + u(i-1,j) + u(i+1,j) + u(i,j-1) + u(i,j+1))/4;
5 i = 3:2:m-1; j = 3:2:n-1;
6 u(i,j) = (fh2(i,j) + u(i-1,j) + u(i+1,j) + u(i,j-1) + u(i,j+1))/4;
7 % case 2 (black points): mod(i+j,2) == 1
8 i = 2:2:m-1; j = 3:2:n-1;
9 u(i,j) = (fh2(i,j) + u(i-1,j) + u(i+1,j) + u(i,j-1) + u(i,j+1))/4;
10 i = 3:2:m-1; j = 2:2:n-1;
11 u(i,j) = (fh2(i,j) + u(i-1,j) + u(i+1,j) + u(i,j-1) + u(i,j+1))/4;
2 1 6 11 2 13 14 15
1.8 1.8
1.6 1.6
2 7 12 10 11 12
1.4 1.4
1.2 1.2
1 3 8 13 1 7 8 9
0.8 0.8
0.6 0.6
4 9 14 4 5 6
0.4 0.4
0.2 0.2
0 5 10 15 0 1 2 3
0 0.5 1 0 0.5 1
8 y =
9 2.0000 2.0000 2.0000
10 1.5000 1.5000 1.5000
11 1.0000 1.0000 1.0000
12 0.5000 0.5000 0.5000
13 0 0 0
In this geometrically consistent system, however, there is an inconsistency of the con-
vention of notation of matrix and Cartesian coordinate. Let us figure out the mapping
between the algebraic index (i,j) and the geometric coordinate (xi , yj ) of a grid point.
In the command
[x,y] = meshgrid(xmin:hx:xmax,ymax:-hy:ymin),
the coordinate of the (i, j)-th grid point is
(xj , yi ) = (xmin + (j − 1)hx , ymin + (n − i + 1)hy ),
which violates the convention of associating index i to xi and j to yj . For a matrix en-
try A(i,j), the 1st index i is the row and the 2nd j is the column while in Cartesian
coordinate, i is usually associated to the x-coordinate and j to the y-coordinate.
The command ndgrid will produce a coordinate consistent matrix in the sense that the
mapping is (i,j) to (xi , yj ) and thus will be called the coordinate consistent indexing.
For example, [x,y] = ndgrid(0:0.5:1,0:0.5:2) will produce two 3 × 5 not 5 × 3
matrices; see Fig. 2(b).
1 >> [x,y] = ndgrid(0:0.5:1,0:0.5:2)
2 x =
3 0 0 0 0 0
4 0.5000 0.5000 0.5000 0.5000 0.5000
5 1.0000 1.0000 1.0000 1.0000 1.0000
6 y =
7 0 0.5000 1.0000 1.5000 2.0000
PROGRAMMING OF FINITE DIFFERENCE METHODS IN MATLAB 5
Then one can easily get the linear indexing of the j-th column of a m × n matrix by
using idxmat(:,j) which is equivalent to sub2ind([m n], 1:m, j*ones(1,m)) but
much easier and intuitive. The price to pay is the extra memory for the full matrix idxmat
which can be minimized using uint32.
For the ndgrid system, to get a geometrically consistent index matrix, we can use the
following command.
idxmat = flipud(transpose(reshape(uint32(1:m*n), n, m))));
For such coordinate consistent system, however, it is recommended to use the subscript
indexing directly.
Similarly we can generate matrices to store the subscripts. For the meshgrid system
1 >> [jj,ii] = meshgrid(1:3,1:5)
2 jj =
3 1 2 3
4 1 2 3
5 1 2 3
6 1 2 3
7 1 2 3
8 ii =
9 1 1 1
10 2 2 2
11 3 3 3
12 4 4 4
13 5 5 5
In the first line, we use size(u) such that it works for both meshgrid and ndgrid system.
PROGRAMMING OF FINITE DIFFERENCE METHODS IN MATLAB 7
The boundary condition can be build into T by changing the entries near the boundary.
Here T corresponds to the homogenous Dirichlet boundary condition. And T is stored as a
sparse matrix to save memory and operations.
For a two dimensional n × n uniform grid, the five point stencil can be decomposed into
(2ui,j − ui−1,j − ui+1,j ) + (2ui,j − ui,j−1 − ui,j+1 )
which can be realized by the left product and right product with the 1-D matrix
Au = u*T + T*u;
For different mesh size or different stencil in x and y-direction, one should generate
specific Tx and Ty and use
1 Au = u*Tx + Ty*u; % meshgrid system
2 Au = Tx*u + u*Ty; % ndgrid system
We can write the matrix as the tensor product of 1-D matrices. Let Am×n and Bp×q be
two matrices. Then the Kronecker (tensor) product of A and B is
a11 B . . . a1n B
A ⊗ B = ... ... ...
am1 B . . . amn B
The matlab command is kron(A,B).
Let Am×m , Bn×n and Xm×n be there matrices. Then it is straightforward to verify the
identities
(1) (AX)(:) = (In ⊗ A) · X(:).
(2) (XB)(:) = (B | ⊗ In ) · X(:).
Here we borrow the notation (:) to change a matrix to a vector by stacking columns from
left to right. Therefore the matrix A for the five point stencil is
A = In ⊗ T + T ⊗ In ,
and the corresponding matlab code is
1 A = kron(speye(nx),Ty) + kron(Tx,speye(ny)); % meshgrid system
2 A = kron(speye(ny),Tx) + kron(Ty,speye(nx)); % ndgrid system
Exercise 3.1. Write out a similar formulae for Neumann boundary condition.
Hint: Change both T and I at boundary indices.
Again in the computation, it is not needed to form A. Instead use the left and right
matrix-product to compute Au if only the matrix-vector product is of interest.
8 LONG CHEN
The tensor product matrix implementation is less obvious since the basic data structure
in MATLAB is matrix not tensor. Denote the stencil matrix in each direction by Ti , i =
1, 2, 3. The first two dimensions can be computed as
1 for k = 1:n3
2 Au(:,:,k) = u(:,:,k)*T2 + T1*u(:,:,k);
3 end
To vectorize the above code, i.e., avoid for loop, one can use reshape which operates
in a column-wise manner. First think about the original data as a long vector by stacking
columns. Then reshape will create the reshaped matrix by transforming consecutive
elements of this long vector into different shape.
We explain the index change by the following example.
1 >> u = reshape(1:3*5*2,3,5,2)
2 u(:,:,1) =
3 1 4 7 10 13
4 2 5 8 11 14
5 3 6 9 12 15
6
7 u(:,:,2) =
8 16 19 22 25 28
9 17 20 23 26 29
PROGRAMMING OF FINITE DIFFERENCE METHODS IN MATLAB 9
10 18 21 24 27 30