function [ el2, eh1 ] = errors ( element_area, element_node, node_xy, u, ... element_num, nnodes, quad2_num, node_num, time ) %% ERRORS calculates the L2 and H1-seminorm errors. % % Discussion: % % This routine uses a 13 point quadrature rule in each element, % in order to estimate the values of % % EL2(t) = Sqrt ( Integral ( U(x,y,t) - Uh(x,y,t) )**2 dx dy ) % % EH1(t) = Sqrt ( Integral ( Ux(x,y,t) - Uhx(x,y,t) )**2 + % ( Uy(x,y,t) - Uhy(x,y,t) )**2 dx dy ) % % Here U is the exact solution, and Ux and Uy its spatial derivatives, % as evaluated by a user-supplied routine. % % Uh, Uhx and Uhy are the computed solution and its spatial derivatives, % as specified by the computed finite element solution. % % Licensing: % % This code is distributed under the GNU LGPL license. % % Modified: % % 10 April 2006 % % Author: % % John Burkardt % % Parameters: % % Input, real ELEMENT_AREA(ELEMENT_NUM), the area of elements. % % Input, integer ELEMENT_NODE(NNODES,ELEMENT_NUM); % ELEMENT_NODE(I,J) is the global index of local node I in element J. % % Input, real NODE_XY(2,NODE_NUM), the nodes. % % Input, real U(NODE_NUM), the coefficients of the solution. % % Input, integer ELEMENT_NUM, the number of elements. % % Input, integer NNODES, the number of nodes used to form one element. % % Input, integer QUAD2_NUM, the number of points in the quadrature rule. % This is actually fixed at 13. % % Input, integer NODE_NUM, the number of nodes. % % Input, real TIME, the current time. % % Output, real EL2, the L2 error. % % Output, real EH1, the H1 seminorm error. % % Local Parameters: % % Local, real AR, the weight for a given quadrature point % in a given element. % % Local, real BB, BX, BY, a basis function and its first % derivatives evaluated at a particular quadrature point. % % Local, real EH1, the H1 seminorm error. % % Local, real EL2, the L2 error. % % Local, real UEX, UEXX, UEXY, the exact solution and its first % derivatives evaluated at a particular quadrature point. % % Local, real UH, UHX, UHY, the computed solution and its first % derivatives evaluated at a particular quadrature point. % % Local, real WQE(NQE), stores the quadrature weights. % % Local, real X, Y, the coordinates of a particular % quadrature point. % % Local, real XQE(NQE), YQE(NQE), stores the location % of quadrature points in a given element. % el2 = 0.0; eh1 = 0.0; for element = 1 : element_num [ wqe, xqe, yqe ] = quad_e ( node_xy, element_node, element, ... element_num, nnodes, node_num, quad2_num ); for quad = 1 : quad2_num ar = element_area(element) * wqe(quad); x = xqe(quad); y = yqe(quad); uh = 0.0; dudxh = 0.0; dudyh = 0.0; for in1 = 1 : nnodes i = element_node(in1,element); [ bi, dbidx, dbidy ] = qbf ( x, y, element, in1, node_xy, ... element_node, element_num, nnodes, node_num ); uh = uh + bi * u(i); dudxh = dudxh + dbidx * u(i); dudyh = dudyh + dbidy * u(i); end xy(1,1) = x; xy(2,1) = y; [ u_exact, dudx_exact, dudy_exact ] = exact_u ( 1, xy, time ); el2 = el2 + ar * ( uh - u_exact(1) )^2 ; eh1 = eh1 + ar * ( ( dudxh - dudx_exact(1) )^2 ... + ( dudyh - dudy_exact(1) )^2 ); end end el2 = sqrt ( el2 ); eh1 = sqrt ( eh1 ); fprintf ( 1, ' %12f %12f %12f\n', time, el2, eh1 ); return end