DERLGL Spectral (Legendre Gauss Lobatto) derivative matrix [d]=derlgl(x,np) returns the spectral Legendre Gauss Lobatto derivative matrix d at the np LGL nodes x (on [-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 derivative matrix: formula (2.3.28), 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] = derlgl(x,np) 0002 %DERLGL Spectral (Legendre Gauss Lobatto) derivative matrix 0003 % 0004 % [d]=derlgl(x,np) returns the spectral Legendre Gauss Lobatto derivative 0005 % matrix d at the np LGL nodes x (on [-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 derivative matrix: formula (2.3.28), 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 d=zeros(np); 0023 n=np-1; 0024 for j =1:np 0025 lnxj = pnleg(x(j),n); 0026 for i = 1:np 0027 if i ~= j 0028 lnxi = pnleg(x(i),n); 0029 d(i,j) = lnxi/((x(i)-x(j))*lnxj); 0030 end 0031 end 0032 end 0033 d(1,1) = -0.25*n*np; 0034 d(np,np) = -d (1,1); 0035 0036 return