P2 Quadratic Element in 2D

We explain degree of freedoms for quadratic element on triangles. There are two types of dofs: vertex type and edge type. Given a mesh, the required data structure can be constructured by

[elem2dof,edge,bdDof] = dofP2(elem)

Contents

help dofP2
  DOFP2 dof structure for P2 element.
 
   [elem2dof,edge,bdDof] = DOFP2(elem) constructs the dof structure for the
   quadratic element based on a triangle. elem2dof(t,i) is the global index
   of the i-th dof of the t-th element. Each triangle contains 6 dofs which
   are organized according to the order of nodes and edges, i.e.
   elem2dof(t,1:3) is the pointer to dof on nodes and elem2dof(t,4:6) to
   three edges. The i-th edge is opposited to the i-th vertex.
 
   See also dof3P2.
  
   Documentation: <a href="matlab:ifem dofP2doc">P2 Quadratic Element in 2D</a>
 
  Copyright (C) Long Chen. See COPYRIGHT.txt for details. 

Local indexing of DOFs

node = [1,0; 1,1; 0,0];
elem = [1 2 3];
edge = [2 3; 1 3; 1 2];
% elem2dof = 1:6;
figure;
set(gcf,'Units','normal');
set(gcf,'Position',[0,0,0.5,0.3]);
subplot(1,2,1)
showmesh(node,elem);
findnode(node);
findedgedof(node,edge);
subplot(1,2,2)
showmesh(node,elem);
findnode(node);
findedge(node,edge);

The six dofs in a triangle is displayed in the left. The first three are associated to the vertices of the triangle and the second three to the middle points of three edges. The edges are indexed such that the i-th edge is opposite to the i-th vertex for i=1,2,3.

Local bases of P2

The six Lagrange-type bases functions are denoted by $\phi_i, i=1:6$, i.e. $\phi_i(x_j)=\delta _{ij},i,j=1:6$. In barycentric coordinates, they are

$$ \phi_1 = \lambda_1(2\lambda_1 -1),\quad \nabla \phi_1 = \nabla \lambda_1(4\lambda_1-1),$$

$$ \phi_2 = \lambda_2(2\lambda_2 -1),\quad  \nabla \phi_2 = \nabla \lambda_2(4\lambda_2-1),$$

$$ \phi_3 = \lambda_3(2\lambda_3 -1),\quad  \nabla \phi_3 = \nabla \lambda_3(4\lambda_3-1),$$

$$ \phi_4 = 4\lambda_2\lambda_3,\quad  \nabla\phi_4 = 4\left (\lambda_2\nabla \lambda_3 + \lambda_3\nabla \lambda_2\right )\; ,$$

$$ \phi_5 = 4\lambda _3\lambda_1,\quad  \nabla\phi_5= 4\left (\lambda_3\nabla \lambda_1 + \lambda_1\nabla \lambda_3\right )\; ,$$

$$ \phi_6 = 4\lambda _1\lambda_2,\quad  \nabla\phi_6=4\left (\lambda_1\nabla
\lambda_2 + \lambda_2\nabla\lambda_1\right )\; .$$

When transfer to the reference triangle formed by $(0,0),(1,0),(0,1)$, the local bases in x-y coordinate can be obtained by substituting

$$\lambda _1 = x, \quad \lambda _2 = y, \quad \lambda _3 = 1-x-y.$$

Global indexing of DOFs

node = [0,0; 1,0; 1,1; 0,1];
elem = [2,3,1; 4,1,3];
[node,elem] = uniformbisect(node,elem);
figure;
showmesh(node,elem);
findnode(node);
findelem(node,elem,'all','index','FaceColor',[0.5 0.9 0.45]);
[elem2dof,edge,bdDof] = dofP2(elem);
findedgedof(node,edge);

The matrix elem2dof is the local to global index mapping of dofs. The first 3 columns, elem2dof(:,1:3), is the global indices of dofs associated to vertexes and thus elem2dof(:,1:3) = elem. The columns elem2dof(:,4:6) point to indices of dofs on edges. The matrix bdDof contains all dofs on the boundary of the mesh.

display(elem2dof);
display(bdDof);
elem2dof =

           8           6           2          14          15          24
           7           6           4          19          20          23
           5           6           1          11          10          22
           9           6           3          16          18          25
           8           3           6          16          24          17
           7           1           6          11          23          12
           5           2           6          14          22          13
           9           4           6          19          25          21


bdDof =

     1
     2
     3
     4
     5
     7
     8
     9
    10
    12
    13
    15
    17
    18
    20
    21