COMPARISON OF CONVERGENCE RATE OF VARIOUS FINITE ELEMENT METHODS
This example is to compare the rate of convergence of various finite element approximation of the Poisson equation on the unit square:

for the Dirichlet boundary condition. Results holds for other boundary conditions.
Contents
Problem
clear all; close all; clc; [node,elem] = squaremesh([0,1,0,1],0.25); pde = sincosdata; bdFlag = setboundary(node,elem,'Dirichlet'); % pde = sincosNeumanndata; % bdFlag = setboundary(node,elem,'Neumann'); option.L0 = 2; option.maxIt = 5; option.printlevel = 1; option.plotflag = 0;
P1 element
err = femPoisson(node,elem,pde,bdFlag,option);
Multigrid V-cycle Preconditioner with Conjugate Gradient Method
#dof: 4225, #nnz: 19593, iter: 9, err = 5.3273e-09, time = 0.14 s
Multigrid V-cycle Preconditioner with Conjugate Gradient Method
#dof: 16641, #nnz: 80137, iter: 10, err = 1.1482e-09, time = 0.17 s
Multigrid V-cycle Preconditioner with Conjugate Gradient Method
#dof: 66049, #nnz: 324105, iter: 10, err = 1.3861e-09, time = 0.62 s
#Dof ||u-u_h|| ||Du-Du_h|| ||DuI-Du_h|| ||uI-u_h||_{max}
289 4.74156e-03 2.17567e-01 7.60779e-03 3.20328e-03
1089 1.19232e-03 1.08979e-01 1.91282e-03 8.02597e-04
4225 2.98519e-04 5.45142e-02 4.78894e-04 2.00761e-04
16641 7.46572e-05 2.72602e-02 1.19767e-04 5.01970e-05
66049 1.86660e-05 1.36305e-02 2.99441e-05 1.25496e-05
The optimal rate of convergence of the H1-norm (1st order), L2-norm (2nd order), and discrte maximum norm is observed. The 2nd order convergent rate between two discrete functions |DuI-Duh| is known as superconvergence.
CR nonconforming P1 element
option.elemType = 'CR';
err = femPoisson(node,elem,pde,bdFlag,option);
Multigrid V-cycle Preconditioner with Conjugate Gradient Method
#dof: 3136, #nnz: 10944, iter: 12, err = 2.8436e-09, time = 0.13 s
Multigrid V-cycle Preconditioner with Conjugate Gradient Method
#dof: 12416, #nnz: 44416, iter: 12, err = 3.0136e-09, time = 0.14 s
Multigrid V-cycle Preconditioner with Conjugate Gradient Method
#dof: 49408, #nnz: 178944, iter: 12, err = 3.0174e-09, time = 0.47 s
Multigrid V-cycle Preconditioner with Conjugate Gradient Method
#dof: 197120, #nnz: 718336, iter: 12, err = 2.8506e-09, time = 2.7 s
#Dof ||u-u_h|| ||Du-Du_h|| ||DuI-Du_h|| ||uI-u_h||_{max}
800 1.20461e-03 1.62309e-01 3.64289e-02 1.55729e-03
3136 3.01947e-04 8.12465e-02 1.81841e-02 3.97659e-04
12416 7.55369e-05 4.06347e-02 9.08831e-03 1.00099e-04
49408 1.88874e-05 2.03188e-02 4.54369e-03 2.50777e-05
197120 4.72202e-06 1.01596e-02 2.27178e-03 6.27347e-06
The optimal rate of convergence of the H1-norm (1st order), L2-norm (2nd order), and discrte maximum norm is observed. But no superconvergence of |DuI-Duh|.
P2 element
option.elemType = 'P2';
err = femPoisson(node,elem,pde,bdFlag,option);
Multigrid V-cycle Preconditioner with Conjugate Gradient Method
#dof: 4225, #nnz: 34721, iter: 9, err = 2.3478e-08, time = 0.15 s
Multigrid V-cycle Preconditioner with Conjugate Gradient Method
#dof: 16641, #nnz: 143137, iter: 11, err = 7.0014e-10, time = 0.29 s
Multigrid V-cycle Preconditioner with Conjugate Gradient Method
#dof: 66049, #nnz: 581153, iter: 12, err = 9.3196e-11, time = 0.89 s
Multigrid V-cycle Preconditioner with Conjugate Gradient Method
#dof: 263169, #nnz: 2341921, iter: 13, err = 1.2514e-11, time = 4.3 s
#Dof ||u-u_h|| ||Du-Du_h|| ||DuI-Du_h|| ||uI-u_h||_{max}
1089 5.72262e-05 8.42020e-03 4.47681e-04 1.36359e-05
4225 7.14446e-06 2.10958e-03 5.75511e-05 8.53377e-07
16641 8.92806e-07 5.27686e-04 7.28791e-06 5.33188e-08
66049 1.11594e-07 1.31940e-04 9.16685e-07 3.33926e-09
263169 1.39490e-08 3.29862e-05 1.14936e-07 2.08791e-10
P3 element
option.elemType = 'P3';
err = femPoisson(node,elem,pde,bdFlag,option);
Multigrid V-cycle Preconditioner with Conjugate Gradient Method
#dof: 2401, #nnz: 34049, iter: 14, err = 5.9261e-08, time = 0.075 s
Multigrid V-cycle Preconditioner with Conjugate Gradient Method
#dof: 9409, #nnz: 143265, iter: 16, err = 5.3279e-09, time = 0.37 s
Multigrid V-cycle Preconditioner with Conjugate Gradient Method
#dof: 37249, #nnz: 587489, iter: 18, err = 4.3333e-10, time = 0.89 s
Multigrid V-cycle Preconditioner with Conjugate Gradient Method
#dof: 148225, #nnz: 2379105, iter: 20, err = 3.6910e-11, time = 3.6 s
#Dof ||u-u_h|| ||Du-Du_h|| ||DuI-Du_h|| ||uI-u_h||_{max}
2401 1.65309e-06 1.39305e-04 1.25209e-04 3.58862e-06
9409 1.01588e-07 1.72328e-05 1.58836e-05 2.27717e-07
37249 6.29771e-09 2.14244e-06 1.99923e-06 1.41694e-08
148225 3.92000e-10 2.67067e-07 2.50731e-07 8.87807e-10