function mg5ptdemo(p) %MG5PTDEMO graphical demonstration of multigrid convergence % % demonstration of multigrid for 5-point Laplacian: plots % successive iterates of a multigrid v-cycle; uses mesh-size % h=1/2^p, with default value for p. % Douglas N. Arnold, 2009-10-15 if nargin == 0 p=4; end niters = 10; % supply f %f = @(x,y) -2*pi^2*sin(pi*x).*sin(pi*y); f = @(x,y) -190*x.^2.*(1-y)+26.4*exp((x-.25).^2+(y-.25).^2); % set mesh size n = 2^p-1; h = 1/(n+1); x = linspace(0,1,n+2); % x grid points including boundaries y = linspace(0,1,n+2); % y grid points including boundaries [X,Y] = meshgrid(x,y); % 2d arrays of x,y values X = X'; % transpose so that X(i,j),Y(i,j) are Y = Y'; % coordinates of (i,j) grid point, i.e., ((i-1)h,(j-1)h) Iint = 2:n+1; % indices of interior points in x Jint = 2:n+1; % indices of interior points in y Xint = X(Iint,Jint); % interior points Yint = Y(Iint,Jint); fh = f(Xint,Yint); % fh is a mesh function on interior grid points %%%%%% Direct solve %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% % set uh on the boundary to given values (interior values will be ignored) uh = zeros(n+2,n+2); % form 5-point Laplacian matrix (times h^2) I = speye(n); e = ones(n,1); A = spdiags([-e,4*e,-e],[-1,0,1],n,n); J = spdiags([-e,-e],[-1,1],n,n); Lh = kron(I,A) + kron(J,I); fh = f(Xint,Yint); % fh is a mesh function on interior grid points % multiply right hand side by h^2 and adjust for boundary condition gh = h^2*fh; gh(:,1) = gh(:,1) - uh(Iint,1); gh(:,n) = gh(:,n) - uh(Iint,n+2); gh(1,:) = gh(1,:) - uh(1,Jint); gh(n,:) = gh(n,:) - uh(n+2,Jint); % convert the right hand side from a grid function into a vector F = reshape(gh,n*n,1); % Solve the linear system for the grid function uh on the interior, written % as vector U = Lh\F; % convert from vector into grid function use this for uh at the interior % grid points (uh is already given at the boundary grid points) uh(Iint,Jint) = reshape(U,n,n); % plot results: figure(1) clf uhplot=surf(X,Y,uh); set(uhplot,'facecolor','yellow') if p > 5 set(uhplot,'edgecolor','none') end set(gca,'zlim',[-1,1]) l=light('position',[-5,5,10]); lighting gouraud %%%%%% End direct solve %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%% multigrid v-cycle %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% % set up initial iterate (zero grid function) u = zeros(size(fh)); % set up initial iterate (random grid function) %u = 2*rand(size(fh))-1; U = [zeros(1,n+2);zeros(n,1),u,zeros(n,1);zeros(1,n+2)]; figure(2) clf uplot=surf(X,Y,U); set(uplot,'facecolor','yellow') if p > 5 set(uplot,'edgecolor','none') end set(gca,'zlim',[-1,1]) hold on l=light('position',[-5,5,10]); lighting gouraud % print max norm error fprintf('iteration 0. max norm error: %9.5f\n',max(max(abs(uh-U)))) % iterate, updating plot each time for niter = 1:niters pause u = mg5pt(fh,u); clf U = [zeros(1,n+2);zeros(n,1),u,zeros(n,1);zeros(1,n+2)]; uplot=surf(X,Y,U); set(uplot,'facecolor','yellow') if p > 5 set(uplot,'edgecolor','none') end set(gca,'zlim',[-1,1]) l=light('position',[-5,5,10]); lighting gouraud % print max norm error fprintf('iteration %4d. max norm error: %9.5f\n',niter,max(max(abs(uh-U)))) end