Download this page as a Matlab LiveScript

Table of contents

Bi-linear element quadrature

Compute area(K)=K  dA\textrm{area}\left(K\right)=\int_K^{\;} dA for the quadrilateral element given by the coordinates below, use Gaussian quadrature to integrate!

P = [0 0 
     1 0
     0.5 0.6
     0 1];

Solution:

figure
patch(P(:,1),P(:,2),'w','FaceColor','w'); hold on
plot(P(:,1),P(:,2),'ok','MarkerFaceColor','k','MarkerSize',8,'DisplayName','Nodes')
xlabel('x')
ylabel('y')
np = size(P,1);
xm = mean(P,1);
Tp = P+(xm-P)*0.2;
text(Tp(:,1),Tp(:,2),cellstr(num2str((1:np)')),'FontSize',14)

Using simple polygon area computation we can compute the area of this polygon.

area_geom = polyarea(P(:,1),P(:,2))
area_geom = 0.5500

The area (of any element) is given by the determinant of the Jacobian of the element integrated over the reference element.

area(K)=K  dA=K^det(J)dK^\textrm{area}\left(K\right)=\int_K^{\;} dA=\int_{\hat{K} } \det \left(\mathit{\mathbf{J} }\right)d\hat{K}

We need the Jacobian

J=[xξyξxηyη]=[φξφξ]xc=B^[x1y1xnxn]\mathit{\mathbf{J} }=\left\lbrack \begin{array}{cc} \dfrac{\partial x}{\partial \xi } & \dfrac{\partial y}{\partial \xi }\\ \dfrac{\partial x}{\partial \eta } & \dfrac{\partial y}{\partial \eta } \end{array}\right\rbrack =\left\lbrack \begin{array}{c} \dfrac{\partial \varphi }{\partial \xi }\\ \dfrac{\partial \varphi }{\partial \xi } \end{array}\right\rbrack {\mathit{\mathbf{x} } }_c =\hat{\mathit{\mathbf{B} } } \left\lbrack \begin{array}{cc} x_1 & y_1 \\ \vdots & \vdots \\ x_n & x_n \end{array}\right\rbrack

using the basis functions for the bi-linear element

phi_n = @(xi,eta)[(eta-1.0).*(xi-1.0),-xi.*(eta-1.0),eta.*xi,-eta.*(xi-1.0)];
syms xi eta
phi(xi,eta) = sym(phi_n)
((η1)(ξ1)ξ(η1)ηξη(ξ1)) \left(\begin{array}{cccc} {\left(\eta -1\right)}\,{\left(\xi -1\right)} & -\xi \,{\left(\eta -1\right)} & \eta \,\xi & -\eta \,{\left(\xi -1\right)} \end{array}\right)

The derivatives of the basis functions.

Bhat = [diff(phi,xi)
        diff(phi,eta)]
(η11ηηηξ1ξξ1ξ) \left(\begin{array}{cccc} \eta -1 & 1-\eta & \eta & -\eta \\ \xi -1 & -\xi & \xi & 1-\xi \end{array}\right)
J = Bhat*P
(1η22η5ξ212ξ5) \left(\begin{array}{cc} 1-\dfrac{\eta }{2} & -\dfrac{2\,\eta }{5}\\ -\dfrac{\xi }{2} & 1-\dfrac{2\,\xi }{5} \end{array}\right)
detJ = det(J)
12ξ5η2 1-\dfrac{2\,\xi }{5}-\dfrac{\eta }{2}
area = int( int( detJ ,eta,0,1),xi,0,1 )
1120 \dfrac{11}{20}
vpa(area)
0.55 0.55

Ok, above we used symbolic integration, but we need to use Gauss quadrature to integrate det(J)\det \left(\mathit{\bm{J} }\right), we can see that the function is linear in ξ\xi and η\eta so one-point Gauss is sufficient to integrate it exactly!

Recall

K^  f(ξ,η)dK^=inf(ξi,ηi)wi\int_{\hat{K} }^{\;} f\left(\xi, \eta \right) d\hat{K} =\sum_i^n f\left(\xi_i ,\eta_i \right)w_i where the points are given in the Gauss quadrature table and dK^ d\hat K is the size (area) of the reference element.

we will use ξ=12\xi =\frac{1}{2} and η=12\eta =\frac{1}{2} with wi=1w_i =1

xi = 1/2; eta = 1/2; w = 1;

area = vpa(subs(detJ))
0.550.55

Second order triangle quadrature

The strength of the Gaussian quadrature really shines when it comes to computing the area of curved elements. Remember, Gauss quadrature integrates these elements exactly!

Compute the area of this second order triangle with the nodal coordinates in

P = [0, 0.1
     1 0
     0 1
     0.5 -0.1
     0.5 0.6
     -0.1 0.5]
phi_n = @(xi,eta)[eta.*-3.0-xi.*3.0+eta.^2.*2.0+xi.^2.*2.0+eta.^2.*xi.^2.*1.6e1+1.0,xi.*(xi.*2.0-1.0),eta.*(eta.*2.0-1.0),xi.*(xi+eta.^2.*xi.*4.0-1.0).*-4.0,eta.^2.*xi.^2.*1.6e1,eta.*(eta+eta.*xi.^2.*4.0-1.0).*-4.0];
phi = sym(phi_n);
syms xi eta
Bhat = [diff(phi,xi)
        diff(phi,eta)];

(32ξη2+4ξ34ξ1044ξ(4η2+1)16η2ξ4ξσ2σ2σ1+4η304η1σ1σ144η(4ξ2+1)16ηξ24η)where    σ1=32ηξ2    σ2=32η2ξ \begin{array}{l} \left(\begin{array}{cccccc} 32\,\xi \,\eta^2 +4\,\xi -3 & 4\,\xi -1 & 0 & 4-4\,\xi \,{\left(4\,\eta^2 +1\right)}-16\,\eta^2 \,\xi -4\,\xi & \sigma_2 & -\sigma_2 \\ \sigma_1 +4\,\eta -3 & 0 & 4\,\eta -1 & -\sigma_1 & \sigma_1 & 4-4\,\eta \,{\left(4\,\xi^2 +1\right)}-16\,\eta \,\xi^2 -4\,\eta \end{array}\right)\\ \mathrm{}\\ \textrm{where}\\ \mathrm{}\\ \;\;\sigma_1 =32\,\eta \,\xi^2 \\ \mathrm{}\\ \;\;\sigma_2 =32\,\eta^2 \,\xi \end{array}

J = Bhat*P;

(2ξσ1+56η2ξ5+14ξ5+σ15+8η2ξ7102η5+σ25+8ηξ252512η5σ2+88ηξ25+710)where    σ1=2ξ(4η2+1)    σ2=2η(4ξ2+1) \begin{array}{l} \left(\begin{array}{cc} 2\,\xi -\sigma_1 +\dfrac{56\,\eta^2 \,\xi }{5}+1 & \dfrac{4\,\xi }{5}+\dfrac{\sigma_1 }{5}+8\,\eta^2 \,\xi -\dfrac{7}{10}\\ \dfrac{2\,\eta }{5}+\dfrac{\sigma_2 }{5}+\dfrac{8\,\eta \,\xi^2 }{5}-\dfrac{2}{5} & \dfrac{12\,\eta }{5}-\sigma_2 +\dfrac{88\,\eta \,\xi^2 }{5}+\dfrac{7}{10} \end{array}\right)\\ \mathrm{}\\ \textrm{where}\\ \mathrm{}\\ \;\;\sigma_1 =2\,\xi \,{\left(4\,\eta^2 +1\right)}\\ \mathrm{}\\ \;\;\sigma_2 =2\,\eta \,{\left(4\,\xi^2 +1\right)} \end{array}

detJ = det(J)
32η3ξ5+152η2ξ2596ηξ325+296ηξ22524ηξ25+24η25+12ξ25+2150 -\dfrac{32\,\eta^3 \,\xi }{5}+\dfrac{152\,\eta^2 \,\xi }{25}-\dfrac{96\,\eta \,\xi^3 }{25}+\dfrac{296\,\eta \,\xi^2 }{25}-\dfrac{24\,\eta \,\xi }{25}+\dfrac{24\,\eta }{25}+\dfrac{12\,\xi }{25}+\dfrac{21}{50}

Now we compute the integral of the above expression using a second order Gauss rule since the highest polynomial term is 3.

gp = [0.5, 0
       0.5, 0.5
       0, 0.5];
gw = 1/3*[1;1;1];

The area of the reference triangle is dK^=1/2d \hat K = 1/2, then we use the sum: I=i=1n=3detJ(ξi,ηi)widK^I = \sum_{i=1}^{n=3}\mathrm{det} \bm J (\xi_i,\eta_i) w_i d \hat K

dk_hat = 1/2;
I = 0;
for i = 1:3
    xi = gp(i,1);
    eta = gp(i,2);
    w = gw(i);
    I = I + subs(detJ)*w*dk_hat
end
203300 \dfrac{203}{300}
vpa(I)
0.67666666666666666666666666666667 0.67666666666666666666666666666667

Bi-quadratic element quadrature

with nodes in

P = [0 0
     1 0
     1 1
     0 1
     0.5 -0.1
     0.9 0.5
     0.5 1.1
     0.1 0.5
     0.5 0.5];

Use Gaussian quadrature rules to integrate the area.

Use the code below to get the Gauss points and weights for the quadrilateral element:

order = 2; %Select the Gauss order here
[GP,GW] = GaussQuadrilateral(order);

We see that the area is supposed to be 1 due to symmetry.

phi = @(xi,eta)[(eta.*-3.0+eta.^2.*2.0+1.0).*(xi.*-3.0+xi.^2.*2.0+1.0),xi.*(xi.*2.0-1.0).*(eta.*-3.0+eta.^2.*2.0+1.0),eta.*xi.*(eta.*2.0-1.0).*(xi.*2.0-1.0),eta.*(eta.*2.0-1.0).*(xi.*-3.0+xi.^2.*2.0+1.0),xi.*(xi-1.0).*(eta.*-3.0+eta.^2.*2.0+1.0).*-4.0,eta.*xi.*(xi.*2.0-1.0).*(eta-1.0).*-4.0,eta.*xi.*(eta.*2.0-1.0).*(xi-1.0).*-4.0,eta.*(eta-1.0).*(xi.*-3.0+xi.^2.*2.0+1.0).*-4.0,eta.*xi.*(eta-1.0).*(xi-1.0).*1.6e1];
phi = sym(phi);
 
syms xi eta
Bhat = [diff(phi,xi)
        diff(phi,eta)];
J = Bhat*P;
detJ = det(J)
48η2ξ22548η2ξ25+36η22548ηξ225+48ηξ2536η254ξ225+4ξ25+2925 \dfrac{48\,\eta^2 \,\xi^2 }{25}-\dfrac{48\,\eta^2 \,\xi }{25}+\dfrac{36\,\eta^2 }{25}-\dfrac{48\,\eta \,\xi^2 }{25}+\dfrac{48\,\eta \,\xi }{25}-\dfrac{36\,\eta }{25}-\dfrac{4\,\xi^2 }{25}+\dfrac{4\,\xi }{25}+\dfrac{29}{25}

Here we evaluate the expression numerically:

detJ = matlabFunction(detJ)
detJ = 
    @(eta,xi)eta.*(-3.6e+1./2.5e+1)+xi.*(4.0./2.5e+1)+eta.*xi.*(4.8e+1./2.5e+1)-eta.*xi.^2.*(4.8e+1./2.5e+1)-eta.^2.*xi.*(4.8e+1./2.5e+1)+eta.^2.*(3.6e+1./2.5e+1)-xi.^2.*(4.0./2.5e+1)+eta.^2.*xi.^2.*(4.8e+1./2.5e+1)+2.9e+1./2.5e+1
area = detJ(GP(:,2),GP(:,1)).'*GW
area = 1.0000

function [GP,GW] = GaussQuadrilateral(order)
%GaussQuadrilateral Gauss integration scheme for quadrilateral 2D element
%   [x,y] in [0,1]

[gp,gw] = gauss(order,0,1);

[GPx,GPy] = meshgrid(gp,gp);
[Gw1,Gw2] = meshgrid(gw,gw);
GW = [Gw1(:), Gw2(:)];
GW = prod(GW,2);
GP = [GPx(:), GPy(:)];

end


function [x, w, A] = gauss(n, a, b)
%------------------------------------------------------------------------------
% gauss.m
%------------------------------------------------------------------------------
%
% Purpose:
%
% Generates abscissas and weights on I = [ a, b ] (for Gaussian quadrature).
%
%
% Syntax:
%
% [x, w, A] = gauss(n, a, b);
%
%
% Input:
%
% n    integer    Number of quadrature points.
% a    real       Left endpoint of interval.
% b    real       Right endpoint of interval.
%
%
% Output:
%
% x    real       Quadrature points.
% w    real       Weigths.
% A    real       Area of domain.
%------------------------------------------------------------------------------
% 3-term recurrence coefficients:
n = 1:(n - 1);
beta = 1 ./ sqrt(4 - 1 ./ (n .* n));
% Jacobi matrix:
J = diag(beta, 1) + diag(beta, -1); 
% Eigenvalue decomposition:
%
% e-values are used for abscissas, whereas e-vectors determine weights.
%
[V, D] = eig(J);
x = diag(D);
% Size of domain:
A = b - a;
% Change of interval:
%
% Initally we have I0 = [ -1, 1 ]; now one gets I = [ a, b ] instead.
%
% The abscissas are Legendre points.
%
if ~(a == -1 && b == 1)
  x = 0.5 * A * x + 0.5 * (b + a);
end
% Weigths:
w = V(1, :) .* V(1, :);
end