Home > Src > Level_2 > normal2_2d.m

normal2_2d

PURPOSE ^

NORMAL2_2D Computes L2-norm in 2D

SYNOPSIS ^

function [err_l2]=normal2_2d(fdq,nq,errtype,x,wx,xx,jacx,y,wy,yy,jacy,xy,ww,nov,un,u,uex);

DESCRIPTION ^

 NORMAL2_2D   Computes L2-norm in 2D

  [err_l2]=normal2_2d(fdq,nq,errtype,x,wx,xx,jacx,y,wy,yy,jacy,...
                     xy,ww,nov,un,u);

 Input : fdq = 0  uses Legendre Gauss quadrature formulas with nq nodes
               in each element (exactness degree = 2*nq+1)
             = 1  uses Legendre Gauss Lobatto quadrature formulas with npdx/y
                  quadrature nodes in each element.
                Quadrature nodes are the nodes of the mesh.
         nq = nodes (in each element and along each direction) 
              for GL quadrature formulas. Not used if  fdq == 1
         errtype = 0 for absolute error ||u-u_ex||
                   1 for relative error ||u-u_ex||/||u_ex||
         x = column array  with npdx LGL nodes in [-1,1]
         wx= column array  with npdx LGL weigths in [-1,1]
         dx= derivative matrix produced with 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 = npdy LGL nodes in [-1,1], previously generated by xwlgl
         wy = npdy LGL weigths in [-1,1], previously generated by xwlgl
         dy= derivative matrix produced with 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]
         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 to global map. 2-indexes array, size(nov)=[nov,ne]
         un = column array with the numerical solution
         u  = column array with the evaluation of exact solution
         uex  = exact solution (uex=@(x,y)[uex(x,y)], with .*, .^, ./)

 Output: err_l2 = error in 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:

SUBFUNCTIONS ^

SOURCE CODE ^

0001 function [err_l2]=normal2_2d(fdq,nq,errtype,x,wx,xx,jacx,y,wy,yy,jacy,...
0002 xy,ww,nov,un,u,uex);
0003 % NORMAL2_2D   Computes L2-norm in 2D
0004 %
0005 %  [err_l2]=normal2_2d(fdq,nq,errtype,x,wx,xx,jacx,y,wy,yy,jacy,...
0006 %                     xy,ww,nov,un,u);
0007 %
0008 % Input : fdq = 0  uses Legendre Gauss quadrature formulas with nq nodes
0009 %               in each element (exactness degree = 2*nq+1)
0010 %             = 1  uses Legendre Gauss Lobatto quadrature formulas with npdx/y
0011 %                  quadrature nodes in each element.
0012 %                Quadrature nodes are the nodes of the mesh.
0013 %         nq = nodes (in each element and along each direction)
0014 %              for GL quadrature formulas. Not used if  fdq == 1
0015 %         errtype = 0 for absolute error ||u-u_ex||
0016 %                   1 for relative error ||u-u_ex||/||u_ex||
0017 %         x = column array  with npdx LGL nodes in [-1,1]
0018 %         wx= column array  with npdx LGL weigths in [-1,1]
0019 %         dx= derivative matrix produced with derlgl
0020 %         xx = 2-indexes array of size (4,ne)
0021 %            xx(1:4,ie)=[x_V1_ie;x_V2_ie;x_V3_ie;x_V4_ie]
0022 %            (ne=nex*ney is the global number of spectral elements)
0023 %         jacx = array (length(jacx)=ne); jacx(ie)= (x_V2_ie-x_V1_ie)/2
0024 %         y = npdy LGL nodes in [-1,1], previously generated by xwlgl
0025 %         wy = npdy LGL weigths in [-1,1], previously generated by xwlgl
0026 %         dy= derivative matrix produced with derlgl
0027 %         yy = 2-indexes array of size (4,ne):
0028 %            yy(1:4,ie)=[y_V1_ie;y_V2_ie;y_V3_ie;y_V4_ie]
0029 %         jacy = array (length(jacy)=ne); jacy(ie)= (y_V3_ie-y_V1_ie)/2
0030 %         xy = column array with global mesh, length: noe=nov(npdx*npdy,ne)
0031 %         ww = column array with global weigths, length: noe=nov(npdx*npdy,ne)
0032 %              diag(ww) is the mass matrix associated to SEM discretization
0033 %         nov = local to global map. 2-indexes array, size(nov)=[nov,ne]
0034 %         un = column array with the numerical solution
0035 %         u  = column array with the evaluation of exact solution
0036 %         uex  = exact solution (uex=@(x,y)[uex(x,y)], with .*, .^, ./)
0037 %
0038 % Output: err_l2 = error in L2-norm
0039 %
0040 % Reference: CHQZ2 = C. Canuto, M.Y. Hussaini, A. Quarteroni, T.A. Zang,
0041 %                    "Spectral Methods. Fundamentals in Single Domains"
0042 %                    Springer Verlag, Berlin Heidelberg New York, 2006.
0043 
0044 %   Written by Paola Gervasio
0045 %   $Date: 2007/04/01$
0046 
0047 
0048 
0049 npdx=length(x); npdy=length(y);
0050 [ldnov,ne]=size(nov); noe=nov(ldnov,ne);
0051 num=0; den=0;
0052 err_l2=0;
0053 
0054 if fdq==0
0055 
0056 % Legendre Gauss nodes and weigths
0057 
0058 [xg,wxg] = xwlg(nq,-1,1);
0059 [yg,wyg] = xwlg(nq,-1,1);
0060 [wxg1,wyg1]=meshgrid(wxg,wyg); wxyg=wxg1.*wyg1;wxyg=wxyg';wxyg=wxyg(:); 
0061 clear wxg1 wyg1;
0062 
0063 % Evaluation of Lagrange basis polynomials at quadrature nodes xg
0064 %
0065 [phix]= intlag_lgl (x,xg);
0066 [phiy]= intlag_lgl (y,yg);
0067 
0068 % Loop on spectral elements
0069 
0070 for ie=1:ne
0071 un_loc=un(nov(1:ldnov,ie));
0072 [norma_loc1,norma_loc2]=normal2_ex_loc(errtype,nq,uex,un_loc,wxyg,...
0073  jacx(ie),xx(1:4,ie),xg,phix,...
0074  jacy(ie),yy(1:4,ie),yg,phiy);
0075 num=num+norma_loc1;
0076 den=den+norma_loc2;
0077 end
0078 
0079 elseif fdq==1
0080     
0081 num=sum((u-un).^2.*ww);
0082 den=sum(u.^2.*ww);
0083 
0084 end
0085 
0086 
0087 if errtype==0
0088     err_l2=sqrt(num);
0089 elseif errtype==1
0090     if abs(den)>1.d-14; err_l2=sqrt(num/den); end
0091 end
0092 
0093 return
0094 
0095 function [norma_loc1,norma_loc2]=normal2_ex_loc(errtype,nq,uex,...
0096  un,wxyg,jacx,xx,xg,phix,...
0097     jacy,yy,yg,phiy);
0098 
0099 % High degree Legendre Gaussian  formulas to compute H1-norm error
0100 
0101 [nq,npdx]=size(phix);
0102 [nq,npdy]=size(phiy);
0103 
0104 % mapping quadrature nodes on element ie
0105 
0106 xxg=xg*jacx+(xx(2)+xx(1))*.5;
0107 yyg=yg*jacy+(yy(3)+yy(1))*.5;
0108 [xg1,yg1]=meshgrid(xxg,yyg); xg1=xg1'; yg1=yg1'; xyg=[xg1(:),yg1(:)]; 
0109 clear xg1; clear yg1;
0110 
0111 
0112 % evaluation of exact solution at quadrature nodes.
0113 
0114 [U]=uex(xyg(:,1),xyg(:,2));
0115 
0116 
0117 % evaluate numerical solution  at quadrature nodes.
0118 
0119 un=reshape(un,npdx,npdy);
0120 
0121 un_i=phix*(phiy*un')';
0122 
0123 un_i=un_i(:);
0124 
0125 
0126 % compute the sum
0127 
0128 norma_loc1=sum((U-un_i).^2.*wxyg)*jacx*jacy;
0129 
0130 if errtype==0
0131     norma_loc2=0;
0132 else
0133     norma_loc2=sum(U.^2.*wxyg)*jacx*jacy;
0134 end
0135 
0136 return
0137 
0138

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