AMG TEST I: DIFFERENT MESHES
We consider linear finite element discretization of the Poisson equation with homongenous Dirichlet boundary condition on different meshes.
Contents
clear all; close all;
Uniform mesh
[node,elem] = squaremesh([0,1,0,1],0.1); [node,elem] = uniformrefine(node,elem); [node,elem] = uniformrefine(node,elem); showmesh(node,elem); snapnow; [N,itStep,time,err,errHist] = amgtest(node,elem);
Algebraic Multigrid W-cycle Preconditioner with Conjugate Gradient Method nnz/N: 4.90, level: 2, coarse grid 560, nnz/Nc 7.86 #dof: 1521, iter: 10, err = 6.7382e-09, time = 0.125 s Algebraic Multigrid W-cycle Preconditioner with Conjugate Gradient Method nnz/N: 4.95, level: 3, coarse grid 532, nnz/Nc 9.53 #dof: 6241, iter: 12, err = 8.0392e-09, time = 0.286 s Algebraic Multigrid W-cycle Preconditioner with Conjugate Gradient Method nnz/N: 4.97, level: 4, coarse grid 402, nnz/Nc 9.84 #dof: 25281, iter: 14, err = 8.2901e-09, time = 0.93 s Algebraic Multigrid W-cycle Preconditioner with Conjugate Gradient Method nnz/N: 4.99, level: 5, coarse grid 328, nnz/Nc 9.77 #dof: 101761, iter: 16, err = 2.3693e-09, time = 3.73 s Algebraic Multigrid W-cycle Preconditioner with Conjugate Gradient Method nnz/N: 4.99, level: 6, coarse grid 247, nnz/Nc 9.41 #dof: 408321, iter: 17, err = 3.8393e-09, time = 16 s
colname = {'Size','Step','Time','Error'};
disptable(colname, N,[],itStep,[],time,'%4.2f',err,'%0.5e');
Size Step Time Error 1521 10 0.12 6.73819e-09 6241 12 0.29 8.03915e-09 25281 14 0.93 8.29012e-09 101761 16 3.73 2.36934e-09 408321 17 16.03 3.83933e-09
r = showrate(N,time,2); xlabel('N'); ylabel('Time'); title(['Complexity is N^{' num2str(r) '}'] ,'Fontsize', 14);
Circle mesh
close all;
[node,elem] = circlemesh(0,0,1,0.2);
[node,elem] = uniformrefine(node,elem);
[node,elem] = uniformrefine(node,elem);
showmesh(node,elem);
snapnow;
[N,itStep,time,err] = amgtest(node,elem);
- Min quality 0.7571 - Mean quality 0.9696 - Uniformity 4.34%
Algebraic Multigrid W-cycle Preconditioner with Conjugate Gradient Method nnz/N: 6.78, level: 2, coarse grid 271, nnz/Nc 7.82 #dof: 1083, iter: 9, err = 7.4996e-09, time = 0.0406 s Algebraic Multigrid W-cycle Preconditioner with Conjugate Gradient Method nnz/N: 6.89, level: 3, coarse grid 233, nnz/Nc 8.50 #dof: 4453, iter: 12, err = 2.0343e-09, time = 0.152 s Algebraic Multigrid W-cycle Preconditioner with Conjugate Gradient Method nnz/N: 6.95, level: 4, coarse grid 190, nnz/Nc 8.65 #dof: 18057, iter: 13, err = 2.1139e-09, time = 0.498 s Algebraic Multigrid W-cycle Preconditioner with Conjugate Gradient Method nnz/N: 6.97, level: 5, coarse grid 157, nnz/Nc 8.82 #dof: 72721, iter: 14, err = 5.1194e-09, time = 2.16 s Algebraic Multigrid W-cycle Preconditioner with Conjugate Gradient Method nnz/N: 6.99, level: 6, coarse grid 130, nnz/Nc 9.09 #dof: 291873, iter: 15, err = 7.4271e-09, time = 8.67 s
colname = {'Size','Step','Time','Error'};
disptable(colname, N,[],itStep,[],time,'%4.2f',err,'%0.5e');
Size Step Time Error 1083 9 0.04 7.49957e-09 4453 12 0.15 2.03426e-09 18057 13 0.50 2.11393e-09 72721 14 2.16 5.11941e-09 291873 15 8.67 7.42708e-09
r = showrate(N,time,2); xlabel('N'); ylabel('Time'); title(['Complexity is N^{' num2str(r) '}'] ,'Fontsize', 14);
Unstructured mesh
close all; load lakemesh showmesh(node,elem); snapnow; [N,itStep,time,err] = amgtest(node,elem);
Algebraic Multigrid W-cycle Preconditioner with Conjugate Gradient Method nnz/N: 6.08, level: 2, coarse grid 437, nnz/Nc 6.24 #dof: 1770, iter: 9, err = 1.6497e-09, time = 0.0453 s Algebraic Multigrid W-cycle Preconditioner with Conjugate Gradient Method nnz/N: 6.59, level: 3, coarse grid 417, nnz/Nc 6.65 #dof: 7876, iter: 10, err = 6.1181e-09, time = 0.158 s Algebraic Multigrid W-cycle Preconditioner with Conjugate Gradient Method nnz/N: 6.81, level: 4, coarse grid 420, nnz/Nc 7.34 #dof: 33081, iter: 12, err = 5.5876e-09, time = 0.851 s Algebraic Multigrid W-cycle Preconditioner with Conjugate Gradient Method nnz/N: 6.91, level: 5, coarse grid 400, nnz/Nc 7.88 #dof: 135463, iter: 14, err = 2.4224e-09, time = 3.93 s Algebraic Multigrid W-cycle Preconditioner with Conjugate Gradient Method nnz/N: 6.95, level: 6, coarse grid 391, nnz/Nc 7.77 #dof: 548115, iter: 15, err = 4.8749e-09, time = 18.3 s
colname = {'Size','Step','Time','Error'};
disptable(colname, N,[],itStep,[],time,'%4.2f',err,'%0.5e');
Size Step Time Error 1770 9 0.05 1.64974e-09 7876 10 0.16 6.11813e-09 33081 12 0.85 5.58763e-09 135463 14 3.93 2.42238e-09 548115 15 18.28 4.87492e-09
r = showrate(N,time,2); xlabel('N'); ylabel('Time'); title(['Complexity is N^{' num2str(r) '}'],'Fontsize', 14);