RATE OF CONVERGENCE OF WEAK GALERKIN FINITE ELEMENT METHOD
This example is to show the rate of convergence of CR nonconforming linear finite element approximation of the Poisson equation on the unit square:

for the following boundary condition:
- Non-empty Dirichlet boundary condition.
. - Pure Neumann boundary condition.
. - Robin boundary condition.

Contents
clear all; close all; [node,elem] = squaremesh([0,1,0,1],0.25); option.L0 = 2; option.maxIt = 4; option.printlevel = 1; option.elemType = 'WG'; option.plotflag = 0;
Non-empty Dirichlet boundary condition.
pde = sincosdata; bdFlag = setboundary(node,elem,'Dirichlet','~(x==0)','Neumann','x==0'); err = femPoisson(node,elem,pde,bdFlag,option);
Multigrid V-cycle Preconditioner with Conjugate Gradient Method #dof: 20608, #nnz: 107058, iter: 22, err = 7.5614e-09, time = 0.68 s Multigrid V-cycle Preconditioner with Conjugate Gradient Method #dof: 82176, #nnz: 462050, iter: 22, err = 7.5865e-09, time = 1.4 s
display(' ||u-u_h|| ||Du-Du_h|| ||DuI-Du_h|| ||uI-u_h||_{max}'); format shorte display([err.L2 err.H1 err.uIuhH1 err.uIuhMax]);
||u-u_h|| ||Du-Du_h|| ||DuI-Du_h|| ||uI-u_h||_{max}
ans =
1.9879e-03 1.6246e-01 6.3330e-02 3.2079e-03
4.9724e-04 8.1266e-02 3.1525e-02 8.0320e-04
1.2433e-04 4.0637e-02 1.5745e-02 2.0084e-04
3.1083e-05 2.0319e-02 7.8704e-03 5.0220e-05
Pure Neumann boundary condition.
pde = sincosNeumanndata; % bdFlag = setboundary(node,elem,'Dirichlet','(x==0)','Neumann','~(x==0)'); bdFlag = setboundary(node,elem,'Neumann'); err = femPoisson(node,elem,pde,bdFlag,option);
Multigrid V-cycle Preconditioner with Conjugate Gradient Method #dof: 20608, #nnz: 108249, iter: 22, err = 7.6491e-09, time = 0.57 s Multigrid V-cycle Preconditioner with Conjugate Gradient Method #dof: 82176, #nnz: 464729, iter: 22, err = 8.7841e-09, time = 1.4 s
display(' ||u-u_h|| ||Du-Du_h|| ||DuI-Du_h|| ||uI-u_h||_{max}'); format shorte display([err.L2 err.H1 err.uIuhH1 err.uIuhMax]);
||u-u_h|| ||Du-Du_h|| ||DuI-Du_h|| ||uI-u_h||_{max}
ans =
9.9726e-03 6.4852e-01 2.6172e-01 1.7953e-02
2.4845e-03 3.2489e-01 1.2710e-01 4.3312e-03
8.8151e-04 2.1364e-01 1.1662e-01 9.8067e-02
2.2042e-04 1.0687e-01 5.8311e-02 4.9073e-02
Pure Robin boundary condition.
pdeRobin = sincosRobindata;
bdFlag = setboundary(node,elem,'Robin');
err = femPoisson(node,elem,pdeRobin,bdFlag,option);
Multigrid V-cycle Preconditioner with Conjugate Gradient Method #dof: 20608, #nnz: 108256, iter: 21, err = 6.2675e-09, time = 0.6 s Multigrid V-cycle Preconditioner with Conjugate Gradient Method #dof: 82176, #nnz: 464736, iter: 21, err = 6.9218e-09, time = 1.4 s
display(' ||u-u_h|| ||Du-Du_h|| ||DuI-Du_h|| ||uI-u_h||_{max}'); format shorte display([err.L2 err.H1 err.uIuhH1 err.uIuhMax]);
||u-u_h|| ||Du-Du_h|| ||DuI-Du_h|| ||uI-u_h||_{max}
ans =
9.7330e-03 6.4838e-01 2.6094e-01 1.6932e-02
2.4444e-03 3.2488e-01 1.2708e-01 4.3121e-03
6.1181e-04 1.6253e-01 6.3104e-02 1.0824e-03
1.5300e-04 8.1274e-02 3.1497e-02 2.7081e-04
Conclusion
The optimal rate of convergence of the H1-norm (1st order) and L2-norm (2nd order) is observed. No superconvergence for |DuI-Duh|.
MGCG converges uniformly in all cases.