DER2LGL Spectral (Legendre Gauss Lobatto) second derivative matrix [d]=der2lgl(x,np) returns the spectral derivative matrix d at the np LGL nodex X (in [-1,1]). np-1 is the polynomial degree used. Legendre Gauss Lobatto (LGL) grid Input: x = array of LGL nodes on [-1,1] (computed by xwlgl) np = number of LGL nodes (=n+1, n=polynomial interpolation degree) Output: d = spectral second derivative matrix: formula (2.3.29), pag. 80 CHQZ2 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 [d] = der2lgl(x,np) 0002 % DER2LGL Spectral (Legendre Gauss Lobatto) second derivative matrix 0003 % 0004 % [d]=der2lgl(x,np) returns the spectral derivative 0005 % matrix d at the np LGL nodex X (in [-1,1]). np-1 is the polynomial 0006 % degree used. Legendre Gauss Lobatto (LGL) grid 0007 % 0008 % 0009 % Input: x = array of LGL nodes on [-1,1] (computed by xwlgl) 0010 % np = number of LGL nodes (=n+1, n=polynomial interpolation degree) 0011 % 0012 % Output: d = spectral second derivative matrix: formula (2.3.29), pag. 80 CHQZ2 0013 % 0014 % Reference: CHQZ2 = C. Canuto, M.Y. Hussaini, A. Quarteroni, T.A. Zang, 0015 % "Spectral Methods. Fundamentals in Single Domains" 0016 % Springer Verlag, Berlin Heidelberg New York, 2006. 0017 0018 % Written by Paola Gervasio 0019 % $Date: 2007/04/01$ 0020 0021 0022 0023 0024 n=np-1; 0025 d=zeros(np,np); 0026 for l =1:np 0027 lnxl = pnleg(x(l),n); 0028 if l~=1 0029 d(1,l)=(-1)^n/lnxl*(n*np*(1+x(l))-4)*0.5/(1+x(l))^2; 0030 end 0031 for j = 2:n 0032 if j ~= l 0033 lnxj = pnleg(x(j),n); 0034 d(j,l) = -2*lnxj/(lnxl*(x(j)-x(l))^2); 0035 else 0036 d(j,j) = pnleg2(x(j),n)/(3*lnxl); 0037 end 0038 end 0039 if l~= np 0040 d(np,l)=1/lnxl*(n*np*(1-x(l))-4)*0.5/(1-x(l))^2; 0041 end 0042 end 0043 d(1,1)=n*(n+1)*(n^2+n-2)/24; 0044 d(np,np)=d(1,1); 0045 0046 return