Home > Src > Level_2 > errors_2d.m

errors_2d

PURPOSE ^

ERRORS_2D Computes errors for 2D b.v.p.

SYNOPSIS ^

function [err_inf,err_h1,err_l2]=errors_2d(x,wx,dx,xx,jacx,y,wy,dy,yy,jacy,xy,ww,nov,un,uex,uex_x,uex_y,param)

DESCRIPTION ^

 ERRORS_2D Computes errors for 2D b.v.p.

 Input:
        x = npdx LGL nodes in [-1,1], previously generated by xwlgl
        wx = npdx LGL weigths in [-1,1], previously generated by xwlgl
        dx = LGL 1-st derivative matrix previously generated by derlgl
        xx = 2-indexes array of size (4,ne) 
            xx(1:4,ie)=[x_V1_ie;x_V2_ie;x_V3_ie;x_V4_ie]
            (ne=nex*ney is the global number of spectral elements)
        jacx = array (length(jacx)=ne); jacx(ie)= (x_V2_ie-x_V1_ie)/2
        y = npdx LGL nodes in [-1,1], previously generated by xwlgl
        wy = npdx LGL weigths in [-1,1], previously generated by xwlgl
        dy = LGL 1-st derivative matrix previously generated by derlgl
        yy = 2-indexes array of size (4,ne)
            yy(1:4,ie)=[y_V1_ie;y_V2_ie;y_V3_ie;y_V4_ie]
            (ne=nex*ney is the global number of spectral elements)
        jacy = array (length(jacy)=ne); jacy(ie)= (y_V3_ie-y_V1_ie)/2
         xy = column array with global mesh, length: noe=nov(npdx*npdy,ne)
         ww = column array with global weigths, length: noe=nov(npdx*npdy,ne)
              diag(ww) is the mass matrix associated to SEM discretization
        nov = local -global map, previously generated by cosnov_2d
        un = numerical solution
        uex  = exact solution (uex=@(x,y)[uex(x,y)], with .*, .^, ./)
        uex_x  = exact first x-derivative (uex_x=@(x,y)[uex_x(x,y)], with .*, .^, ./)
        uex_y  = exact first y-derivative (uex_y=@(x,y)[uex_y(x,y)], with .*, .^, ./)
        param(4) = 1: compute errors (L^inf-norm, L2-norm, H1-norm) 
                      on the exact solution
                   2: no  errors are computed
        param(5) = 0: LG quadrature formulas with high precision degree are
                      used to compute norms (exact norms)
                   1: LGL quadrature formulas with npdx,npdy nodes are
                      used to compute norms (discrete norms)
                   (used only if param(4) == 1)
        param(6) = number of nodes for high degree quadrature formula,
                   (used only if param(5) == 0 & param(4) == 1)
        param(7) = 0: absolute errors are computed
                   1: relative errors are computed
                   (used only if param(4) == 1)

 Output: 
 
         err_inf = ||u_ex-u_N|| with respect to discrete maximum-norm
         err_h1 = ||u_ex-u_N|| with respect to discrete H1-norm
         err_l2 = ||u_ex-u_N|| with respect to discrete L2-norm

 Reference: CHQZ2 = C. Canuto, M.Y. Hussaini, A. Quarteroni, T.A. Zang,
                    "Spectral Methods. Fundamentals in Single Domains"
                    Springer Verlag, Berlin Heidelberg New York, 2006.

CROSS-REFERENCE INFORMATION ^

This function calls: This function is called by:

SOURCE CODE ^

0001 function [err_inf,err_h1,err_l2]=errors_2d(x,wx,dx,xx,jacx,y,wy,dy,yy,jacy,...
0002 xy,ww,nov,un,uex,uex_x,uex_y,param)
0003 % ERRORS_2D Computes errors for 2D b.v.p.
0004 %
0005 % Input:
0006 %        x = npdx LGL nodes in [-1,1], previously generated by xwlgl
0007 %        wx = npdx LGL weigths in [-1,1], previously generated by xwlgl
0008 %        dx = LGL 1-st derivative matrix previously generated by derlgl
0009 %        xx = 2-indexes array of size (4,ne)
0010 %            xx(1:4,ie)=[x_V1_ie;x_V2_ie;x_V3_ie;x_V4_ie]
0011 %            (ne=nex*ney is the global number of spectral elements)
0012 %        jacx = array (length(jacx)=ne); jacx(ie)= (x_V2_ie-x_V1_ie)/2
0013 %        y = npdx LGL nodes in [-1,1], previously generated by xwlgl
0014 %        wy = npdx LGL weigths in [-1,1], previously generated by xwlgl
0015 %        dy = LGL 1-st derivative matrix previously generated by derlgl
0016 %        yy = 2-indexes array of size (4,ne)
0017 %            yy(1:4,ie)=[y_V1_ie;y_V2_ie;y_V3_ie;y_V4_ie]
0018 %            (ne=nex*ney is the global number of spectral elements)
0019 %        jacy = array (length(jacy)=ne); jacy(ie)= (y_V3_ie-y_V1_ie)/2
0020 %         xy = column array with global mesh, length: noe=nov(npdx*npdy,ne)
0021 %         ww = column array with global weigths, length: noe=nov(npdx*npdy,ne)
0022 %              diag(ww) is the mass matrix associated to SEM discretization
0023 %        nov = local -global map, previously generated by cosnov_2d
0024 %        un = numerical solution
0025 %        uex  = exact solution (uex=@(x,y)[uex(x,y)], with .*, .^, ./)
0026 %        uex_x  = exact first x-derivative (uex_x=@(x,y)[uex_x(x,y)], with .*, .^, ./)
0027 %        uex_y  = exact first y-derivative (uex_y=@(x,y)[uex_y(x,y)], with .*, .^, ./)
0028 %        param(4) = 1: compute errors (L^inf-norm, L2-norm, H1-norm)
0029 %                      on the exact solution
0030 %                   2: no  errors are computed
0031 %        param(5) = 0: LG quadrature formulas with high precision degree are
0032 %                      used to compute norms (exact norms)
0033 %                   1: LGL quadrature formulas with npdx,npdy nodes are
0034 %                      used to compute norms (discrete norms)
0035 %                   (used only if param(4) == 1)
0036 %        param(6) = number of nodes for high degree quadrature formula,
0037 %                   (used only if param(5) == 0 & param(4) == 1)
0038 %        param(7) = 0: absolute errors are computed
0039 %                   1: relative errors are computed
0040 %                   (used only if param(4) == 1)
0041 %
0042 % Output:
0043 %
0044 %         err_inf = ||u_ex-u_N|| with respect to discrete maximum-norm
0045 %         err_h1 = ||u_ex-u_N|| with respect to discrete H1-norm
0046 %         err_l2 = ||u_ex-u_N|| with respect to discrete L2-norm
0047 %
0048 % Reference: CHQZ2 = C. Canuto, M.Y. Hussaini, A. Quarteroni, T.A. Zang,
0049 %                    "Spectral Methods. Fundamentals in Single Domains"
0050 %                    Springer Verlag, Berlin Heidelberg New York, 2006.
0051 
0052 %   Written by Paola Gervasio
0053 %   $Date: 2007/04/01$
0054 
0055 fdq=param(5); 
0056 if fdq ==0
0057 nq=param(6);
0058 else 
0059 nq=length(wx);
0060 end
0061 errtype=param(7);
0062 
0063 % compute errors on the exact solution
0064 
0065 npdx=length(wx); npdy=length(wy);
0066 [ldnov,ne]=size(nov);
0067 noe=length(un);
0068 
0069 % Evaluate exact solution
0070 
0071 u=zeros(noe,1)+uex(xy(:,1),xy(:,2));
0072 
0073 %  compute difference
0074 
0075 err=u-un;
0076 
0077 % ||u_ex-u_N|| with respect to discrete maximum-norm
0078 
0079 err_inf=norm(err,inf);
0080 if errtype~=0
0081     err_inf=err_inf/norm(u,inf);
0082 end
0083 % ||u_ex-u_N|| with respect to H1 norm
0084 
0085 % fdq=0 LG quadrature formula with nq nodes on each element
0086 % fdq=1 LGL quadrature formula with npdx nodes on each element
0087 
0088 [err_h1]=normah1_2d(fdq,nq,errtype,x,wx,dx,xx,jacx,y,wy,dy,yy,jacy,...
0089 xy,ww,nov,un,u,uex,uex_x,uex_y);
0090 
0091 % ||u_ex-u_N|| with respect to L2 norm
0092 
0093 [err_l2]=normal2_2d(fdq,nq,errtype,x,wx,xx,jacx,y,wy,yy,jacy,...
0094 xy,ww,nov,un,u,uex);
0095 
0096 return

Generated on Fri 21-Sep-2007 10:07:00 by m2html © 2003