RATE OF CONVERGENCE OF MIXED FINITE ELEMENT METHOD (BDM1-P0) FOR POISSON EQUATION
This example is to show the rate of convergence of mixed finite element (RT0-P0) approximation of the Poisson equation on the unit square:
for the following boundary conditions:
- Pure Dirichlet boundary condition. .
- Pure Neumann boundary condition. .
- Mix Dirichlet and Neumann boundary condition.
Written by Ming Wang and improved by Long Chen.
Contents
Setting
[node,elem] = squaremesh([0,1,0,1],0.25);
pde = sincosNeumanndata;
option.L0 = 2;
option.maxIt = 4;
option.printlevel = 1;
option.elemType = 'BDM1';
Non-empty Dirichlet boundary condition.
option.plotflag = 1;
bdFlag = setboundary(node,elem,'Dirichlet');
mfemPoisson(node,elem,pde,bdFlag,option);
Uzawa-type MultiGrid Preconditioned PCG #dof: 2112, #nnz: 14016, V-cycle: 1, iter: 15, err = 8.94e-09, time = 0.05 s Uzawa-type MultiGrid Preconditioned PCG #dof: 8320, #nnz: 55680, V-cycle: 1, iter: 16, err = 3.57e-09, time = 0.32 s Uzawa-type MultiGrid Preconditioned PCG #dof: 33024, #nnz: 221952, V-cycle: 1, iter: 16, err = 6.16e-09, time = 1.5 s Uzawa-type MultiGrid Preconditioned PCG #dof: 131584, #nnz: 886272, V-cycle: 1, iter: 17, err = 3.09e-09, time = 5.3 s Table: Error #Dof h ||u-u_h|| ||u_I-u_h|| ||sigma-sigma_h||||sigma-sigma_h||_{div} 2112 6.25e-02 6.56418e-02 1.19213e-02 8.96176e-02 5.14701e+00 8320 3.12e-02 3.27516e-02 3.03293e-03 2.27048e-02 2.58126e+00 33024 1.56e-02 1.63659e-02 7.61595e-04 5.69906e-03 1.29160e+00 131584 7.81e-03 8.18166e-03 1.90611e-04 1.42668e-03 6.45924e-01 Table: CPU time #Dof Assemble Solve Error Mesh 2112 6.38e-03 5.00e-02 2.00e-02 0.00e+00 8320 2.14e-02 3.20e-01 2.00e-02 0.00e+00 33024 9.15e-02 1.50e+00 1.40e-01 1.00e-02 131584 3.90e-01 5.27e+00 3.70e-01 0.00e+00
Pure Neumann boundary condition.
option.solver = 'tripremixpoisson'; bdFlag = setboundary(node,elem,'Neumann'); mfemPoisson(node,elem,pde,bdFlag,option);
Triangular Preconditioner Preconditioned GMRES #dof: 2112, #nnz: 12999, V-cycle: 1, iter: 26, err = 6.74e-09, time = 0.05 s Triangular Preconditioner Preconditioned GMRES #dof: 8320, #nnz: 53639, V-cycle: 1, iter: 26, err = 6.36e-09, time = 0.19 s Triangular Preconditioner Preconditioned GMRES #dof: 33024, #nnz: 217863, V-cycle: 1, iter: 25, err = 7.62e-09, time = 0.45 s Triangular Preconditioner Preconditioned GMRES #dof: 131584, #nnz: 878087, V-cycle: 1, iter: 25, err = 8.89e-09, time = 1.7 s Table: Error #Dof h ||u-u_h|| ||u_I-u_h|| ||sigma-sigma_h||||sigma-sigma_h||_{div} 2112 6.25e-02 6.59919e-02 1.46218e-02 1.01998e-01 5.14701e+00 8320 3.12e-02 3.27830e-02 3.60835e-03 2.57227e-02 2.58126e+00 33024 1.56e-02 1.63694e-02 8.99012e-04 6.44488e-03 1.29160e+00 131584 7.81e-03 8.18209e-03 2.24561e-04 1.61214e-03 6.45924e-01 Table: CPU time #Dof Assemble Solve Error Mesh 2112 9.63e-03 5.00e-02 1.00e-02 0.00e+00 8320 2.94e-02 1.90e-01 3.00e-02 0.00e+00 33024 1.06e-01 4.50e-01 9.00e-02 1.00e-02 131584 4.34e-01 1.69e+00 3.70e-01 0.00e+00
Mix Dirichlet and Neumann boundary condition.
option.solver = 'uzawapcg'; bdFlag = setboundary(node,elem,'Dirichlet','~(x==0)','Neumann','x==0'); mfemPoisson(node,elem,pde,bdFlag,option);
Uzawa-type MultiGrid Preconditioned PCG #dof: 2112, #nnz: 13760, V-cycle: 1, iter: 15, err = 9.27e-09, time = 0.07 s Uzawa-type MultiGrid Preconditioned PCG #dof: 8320, #nnz: 55168, V-cycle: 1, iter: 16, err = 3.40e-09, time = 0.38 s Uzawa-type MultiGrid Preconditioned PCG #dof: 33024, #nnz: 220928, V-cycle: 1, iter: 16, err = 6.30e-09, time = 1.2 s Uzawa-type MultiGrid Preconditioned PCG #dof: 131584, #nnz: 884224, V-cycle: 1, iter: 17, err = 3.05e-09, time = 4.6 s Table: Error #Dof h ||u-u_h|| ||u_I-u_h|| ||sigma-sigma_h||||sigma-sigma_h||_{div} 2112 6.25e-02 6.57963e-02 1.29671e-02 9.19791e-02 5.14701e+00 8320 3.12e-02 3.27711e-02 3.29173e-03 2.32963e-02 2.58126e+00 33024 1.56e-02 1.63683e-02 8.26130e-04 5.84697e-03 1.29160e+00 131584 7.81e-03 8.18197e-03 2.06734e-04 1.46363e-03 6.45924e-01 Table: CPU time #Dof Assemble Solve Error Mesh 2112 8.09e-03 7.00e-02 1.00e-02 0.00e+00 8320 3.42e-02 3.80e-01 2.00e-02 0.00e+00 33024 9.10e-02 1.22e+00 1.00e-01 1.00e-02 131584 4.30e-01 4.61e+00 3.70e-01 0.00e+00
Conclusion
The optimal rates of convergence for u and sigma are observed, namely, 1st order for L2 norm of u, H(div) norm of sigma, and 2nd order for L2 norm of sigma.
The 2nd order convergent rates between two discrete functions |uI-uh| is known as superconvergence.
The 3rd order convergent rates between two discrete functions |sigmaI-sigmah| for the Pure Neumman boundary condition is known as superconvergence. The 2nd accurate numerical quadrature is required for the integral of rhs to observe such superconvergence.
Triangular preconditioned GMRES and Uzawa preconditioned CG converges uniformly in all cases. Traingular preconditioner is two times faster than PCG although GMRES is used.