NORMAL2_1D Computes L2-norm in 1D [err_h1]=normah1_1d(fdq, nq, errtype, u, uex, ... x, wx, xx, jacx, xy,nov) 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 on the nodes in which the discrete function is knonwn. nq = nodes (in each elemetn) 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|| u = numerical solution, column array of length noe uex = exact solution (uex=@(x)[uex(x)], with .*, .^, ./) x = column array with LGL nodes in [-1,1] wx= column array with LGL weigths in [-1,1] xx = 2-indexes array of size (2,ne): xx(1:2,ie)=[xa_ie;xb_ie] jacx = column array of length = ne, containing jacobians of the maps F_ie:[-1,1]---->[xa_ie,xb_ie] xy = mesh, column array of length noe nov = local to global map. 2-indexes array, size(nov)=[nov,ne] 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.
0001 function [err_l2]=normal2_1d(fdq, nq, errtype, u, uex, ... 0002 x, wx, xx, jacx, xy,nov) 0003 % NORMAL2_1D Computes L2-norm in 1D 0004 % 0005 % [err_h1]=normah1_1d(fdq, nq, errtype, u, uex, ... 0006 % x, wx, xx, jacx, xy,nov) 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 on the 0011 % nodes in which the discrete function is knonwn. 0012 % nq = nodes (in each elemetn) for GL quadrature formulas. Not used if 0013 % fdq=1 0014 % errtype = 0 for absolute error ||u-u_ex|| 0015 % 1 for relative error ||u-u_ex||/||u_ex|| 0016 % u = numerical solution, column array of length noe 0017 % uex = exact solution (uex=@(x)[uex(x)], with .*, .^, ./) 0018 % x = column array with LGL nodes in [-1,1] 0019 % wx= column array with LGL weigths in [-1,1] 0020 % xx = 2-indexes array of size (2,ne): xx(1:2,ie)=[xa_ie;xb_ie] 0021 % jacx = column array of length = ne, containing 0022 % jacobians of the maps F_ie:[-1,1]---->[xa_ie,xb_ie] 0023 % xy = mesh, column array of length noe 0024 % nov = local to global map. 2-indexes array, size(nov)=[nov,ne] 0025 % 0026 % Reference: CHQZ2 = C. Canuto, M.Y. Hussaini, A. Quarteroni, T.A. Zang, 0027 % "Spectral Methods. Fundamentals in Single Domains" 0028 % Springer Verlag, Berlin Heidelberg New York, 2006. 0029 0030 % Written by Paola Gervasio 0031 % $Date: 2007/04/01$ 0032 0033 0034 npdx=length(x); [ldnov,ne]=size(nov); noe=nov(npdx,ne); 0035 num=0; den=0; 0036 err_l2=0; 0037 0038 if fdq==0 0039 0040 % Legendre Gauss nodes and weigths 0041 0042 [xg,wxg] = xwlg(nq,-1,1); 0043 0044 % Evaluation of Lagrange basis polynomials at quadrature nodes xg 0045 % 0046 [phix]= intlag_lgl (x,xg); 0047 0048 % Loop on spectral elements 0049 for ie=1:ne 0050 u_loc=u(nov(1:npdx,ie)); 0051 [norma_loc1,norma_loc2]=normal2_ex_loc(errtype,npdx,uex,... 0052 u_loc,x,jacx(ie),xx(1:2,ie),xg,wxg,nq,phix); 0053 num=num+norma_loc1; 0054 den=den+norma_loc2; 0055 end 0056 0057 0058 elseif fdq==1 0059 for ie=1:ne 0060 u_loc=u(nov(1:npdx,ie)); 0061 xy_loc=xy(nov(1:npdx,ie)); 0062 [norma_loc1,norma_loc2]=normal2_loc(errtype,npdx,uex,... 0063 u_loc,xy_loc, wx,jacx(ie)); 0064 num=num+norma_loc1; 0065 den=den+norma_loc2; 0066 end 0067 end 0068 0069 if errtype==0 0070 err_l2=sqrt(num); 0071 elseif errtype==1 0072 if abs(den)>1.d-14; err_l2=sqrt(num/den); end 0073 end 0074 0075 return 0076 0077 function [norma_loc1,norma_loc2]=normal2_ex_loc(errtype,npdx,uex,... 0078 u, x,jacx,xx,xg,wxg,nq,phix); 0079 0080 % High degree Legendre Gaussian formulas to compute L2-norm error 0081 0082 % mapping quadrature nodes on element ie 0083 0084 xxg=xg*jacx+(xx(2)+xx(1))*.5; 0085 0086 % evaluation of exact solution at quadrature nodes. 0087 0088 U=uex(xxg); 0089 0090 0091 % evaluate numerical solution and its derivative at quadrature nodes. 0092 0093 u_i=phix*u; 0094 % figure(1) 0095 % mesh(xxg,yyg,reshape(U,npx,npy)); 0096 % figure(2) 0097 % mesh(xxg,yyg,uloc_i) 0098 0099 % compute the sum 0100 0101 norma_loc1=sum((U-u_i).^2.*wxg)*jacx; 0102 0103 if errtype==0 0104 norma_loc2=0; 0105 else 0106 norma_loc2=sum(U.^2.*wxg)*jacx; 0107 end 0108 0109 return 0110 0111 function [norma_loc1,norma_loc2]=normal2_loc(errtype,npdx,uex,uloc,... 0112 xy_loc, wx,jacx); 0113 0114 % LGL quadrature formulas on npdx nodes to compute L2-norm error 0115 0116 % evaluation of exact solution at quadrature nodes. 0117 0118 U=uex(xy_loc); 0119 0120 0121 % compute the sum 0122 0123 norma_loc1=sum((U-uloc).^2.*wx)*jacx; 0124 if errtype==0 0125 norma_loc2=0; 0126 else 0127 norma_loc2=sum(U.^2.*wx)*jacx; 0128 end 0129 0130 return 0131