RATE OF CONVERGENCE OF CUBIC ELEMENT FOR POISSON EQUATION

This example is to show the rate of convergence of cubic finite element approximation of the Poisson equation on the unit square:

$$- \Delta u = f \; \hbox{in } (0,1)^2$$

for the following boundary condition:

  1. Non-empty Dirichlet boundary condition. $u=g_D \hbox{ on }\Gamma_D, \quad \nabla u\cdot n=g_N \hbox{ on }\Gamma_N.$
  2. Pure Neumann boundary condition. $\nabla u\cdot n=g_N \hbox{ on } \partial \Omega$.
  3. Robin boundary condition. $g_R u + \nabla u\cdot n=g_N \hbox{ on }\partial \Omega$

Contents

Setting

[node,elem] = squaremesh([0,1,0,1],0.25);
option.L0 = 2;
option.maxIt = 4;
option.maxN = 1e6;
option.printlevel = 1;
option.elemType = 'P3';

Non-empty Dirichlet boundary condition.

option.plotflag = 1;
pde = sincosdata;
mesh.bdFlag = setboundary(node,elem,'Dirichlet','~(x==0)','Neumann','x==0');
femPoisson(mesh,pde,option);
Multigrid V-cycle Preconditioner with Conjugate Gradient Method
#dof:     2401,  #nnz:    35034, smoothing: (1,1), iter: 18,   err = 3.84e-10,   time = 0.015 s
Multigrid V-cycle Preconditioner with Conjugate Gradient Method
#dof:     9409,  #nnz:   145290, smoothing: (1,1), iter: 18,   err = 4.36e-10,   time = 0.07 s
Multigrid V-cycle Preconditioner with Conjugate Gradient Method
#dof:    37249,  #nnz:   591594, smoothing: (1,1), iter: 18,   err = 4.43e-10,   time = 0.23 s
Multigrid V-cycle Preconditioner with Conjugate Gradient Method
#dof:   148225,  #nnz:  2387370, smoothing: (1,1), iter: 18,   err = 4.41e-10,   time =    1 s

 #Dof        h       ||u-u_h||    ||Du-Du_h||   ||DuI-Du_h|| ||uI-u_h||_{max}

  2401   6.25e-02   1.64347e-06   1.38348e-04   1.26253e-04   3.61495e-06
  9409   3.12e-02   1.01251e-07   1.71699e-05   1.59519e-05   2.26879e-07
 37249   1.56e-02   6.28490e-09   2.13842e-06   2.00358e-06   1.42241e-08
148225   7.81e-03   3.94787e-10   2.66852e-07   2.51090e-07   1.05044e-09

 #Dof    Assemble     Solve      Error      Mesh    

  2401   3.00e-02   1.51e-02   1.00e-02   1.00e-02
  9409   1.50e-01   6.99e-02   2.00e-02   1.00e-02
 37249   3.90e-01   2.32e-01   9.00e-02   3.00e-02
148225   1.24e+00   1.00e+00   3.40e-01   1.40e-01


Pure Neumann boundary condition.

option.plotflag = 0;
pde = sincosNeumanndata;
% pde = sincosdata;
mesh.bdFlag = setboundary(node,elem,'Neumann');
femPoisson(mesh,pde,option);
Multigrid V-cycle Preconditioner with Conjugate Gradient Method
#dof:     2401,  #nnz:    38084, smoothing: (1,1), iter: 20,   err = 3.41e-10,   time = 0.024 s
Multigrid V-cycle Preconditioner with Conjugate Gradient Method
#dof:     9409,  #nnz:   151460, smoothing: (1,1), iter: 20,   err = 7.14e-10,   time = 0.06 s
Multigrid V-cycle Preconditioner with Conjugate Gradient Method
#dof:    37249,  #nnz:   604004, smoothing: (1,1), iter: 21,   err = 3.97e-10,   time = 0.25 s
Multigrid V-cycle Preconditioner with Conjugate Gradient Method
#dof:   148225,  #nnz:  2412260, smoothing: (1,1), iter: 21,   err = 7.45e-10,   time =  1.1 s

 #Dof        h       ||u-u_h||    ||Du-Du_h||   ||DuI-Du_h|| ||uI-u_h||_{max}

  2401   6.25e-02   2.60671e-05   2.15721e-03   2.06194e-03   5.90094e-05
  9409   3.12e-02   1.60744e-06   2.71037e-04   2.58554e-04   3.79705e-06
 37249   1.56e-02   9.99996e-08   3.39784e-05   3.22882e-05   2.39045e-07
148225   7.81e-03   6.23993e-09   4.25404e-06   4.03123e-06   1.49762e-08

 #Dof    Assemble     Solve      Error      Mesh    

  2401   2.00e-02   2.39e-02   1.00e-02   0.00e+00
  9409   5.00e-02   6.04e-02   2.00e-02   0.00e+00
 37249   2.90e-01   2.51e-01   1.20e-01   2.00e-02
148225   1.27e+00   1.15e+00   3.50e-01   5.00e-02


Pure Robin boundary condition.

option.plotflag = 0;
pde = sincosRobindata;
mesh.bdFlag = setboundary(node,elem,'Robin');
femPoisson(mesh,pde,option);
Multigrid V-cycle Preconditioner with Conjugate Gradient Method
#dof:     2401,  #nnz:    38113, smoothing: (1,1), iter: 17,   err = 7.78e-10,   time = 0.02 s
Multigrid V-cycle Preconditioner with Conjugate Gradient Method
#dof:     9409,  #nnz:   151489, smoothing: (1,1), iter: 17,   err = 8.91e-10,   time = 0.054 s
Multigrid V-cycle Preconditioner with Conjugate Gradient Method
#dof:    37249,  #nnz:   604033, smoothing: (1,1), iter: 17,   err = 9.77e-10,   time = 0.23 s
Multigrid V-cycle Preconditioner with Conjugate Gradient Method
#dof:   148225,  #nnz:  2412289, smoothing: (1,1), iter: 17,   err = 9.84e-10,   time =    1 s

 #Dof        h       ||u-u_h||    ||Du-Du_h||   ||DuI-Du_h|| ||uI-u_h||_{max}

  2401   6.25e-02   2.60684e-05   2.15721e-03   2.06171e-03   5.87482e-05
  9409   3.12e-02   1.60746e-06   2.71037e-04   2.58543e-04   3.79065e-06
 37249   1.56e-02   9.99999e-08   3.39784e-05   3.22878e-05   2.38940e-07
148225   7.81e-03   6.23994e-09   4.25404e-06   4.03122e-06   1.49667e-08

 #Dof    Assemble     Solve      Error      Mesh    

  2401   2.00e-02   1.97e-02   1.00e-02   0.00e+00
  9409   6.00e-02   5.44e-02   2.00e-02   1.00e-02
 37249   3.30e-01   2.31e-01   1.10e-01   3.00e-02
148225   1.14e+00   1.02e+00   4.00e-01   5.00e-02


Conclusion

The optimal rate of convergence of the H1-norm (3rd order) and L2-norm (4th order) is observed. The order of |DuI-Duh| is 3rd order and thus no superconvergence exists between nodal interpolate and uh.

MGCG converges uniformly in all cases.