Home > Src > Level_0 > xwlgl.m

xwlgl

PURPOSE ^

XWLGL Computes nodes and weights of the Legendre-Gauss-Lobatto quadrature formula.

SYNOPSIS ^

function [x,w] = xwlgl(np,a,b)

DESCRIPTION ^

XWLGL Computes nodes and weights of the Legendre-Gauss-Lobatto quadrature formula.

    [x,w]=xwlgl(np) returns the np weigths and nodes 
    of the corresponding Legendre Gauss-Lobatto quadrature 
    formula in the reference interval [-1,1].

    [x,w]=xwlgl(np,a,b) returns the np weigths and the nodes 
    of the corresponding Legendre Gauss-Lobatto quadrature 
    formula in the  interval [a,b].

 Input: np = number of nodes
        a, b = extrema of the interval

 Output: x(np,1) = LGL nodes  (CHQZ2, (2.3.12), pag. 76)
                   the set is given by (CHQZ2, pag. 92) as well:
               {-1} U {zeros of Jacobi polynomial P_{N-1}^{(1,1)}} U {1}
         w(np,1) = LGL weigths (CHQZ2, (2.3.12), pag. 76)

 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 [x,w] = xwlgl(np,a,b) 
0002 %XWLGL Computes nodes and weights of the Legendre-Gauss-Lobatto quadrature formula.
0003 %
0004 %    [x,w]=xwlgl(np) returns the np weigths and nodes
0005 %    of the corresponding Legendre Gauss-Lobatto quadrature
0006 %    formula in the reference interval [-1,1].
0007 %
0008 %    [x,w]=xwlgl(np,a,b) returns the np weigths and the nodes
0009 %    of the corresponding Legendre Gauss-Lobatto quadrature
0010 %    formula in the  interval [a,b].
0011 %
0012 % Input: np = number of nodes
0013 %        a, b = extrema of the interval
0014 %
0015 % Output: x(np,1) = LGL nodes  (CHQZ2, (2.3.12), pag. 76)
0016 %                   the set is given by (CHQZ2, pag. 92) as well:
0017 %               {-1} U {zeros of Jacobi polynomial P_{N-1}^{(1,1)}} U {1}
0018 %         w(np,1) = LGL weigths (CHQZ2, (2.3.12), pag. 76)
0019 %
0020 % Reference: CHQZ2 = C. Canuto, M.Y. Hussaini, A. Quarteroni, T.A. Zang,
0021 %                    "Spectral Methods. Fundamentals in Single Domains"
0022 %                    Springer Verlag, Berlin Heidelberg New York, 2006.
0023 
0024 %   Written by Paola Gervasio
0025 %   $Date: 2007/04/01$
0026 
0027 
0028 if np <= 1 
0029   es='The number of the quadrature nodes should be greater than 1';
0030   disp(es); x = []; w = [];
0031   return 
0032   
0033 elseif np==2
0034     x=[-1;1];w=[1;1];
0035     
0036 else
0037 x=zeros(np,1);
0038 w=zeros(np,1);
0039 n=np-1;
0040 x(1)=-1;x(np)=1;
0041 x1=jacobi_roots(n-1,1,1);
0042 x(2:n)=x1;
0043 coef=2/(n*np);
0044 w=coef./(pnleg(x,n).^2);
0045 end
0046 
0047 if nargin ==3 
0048    bma=(b-a)*.5;
0049    bpa=(b+a)*.5;
0050    x=bma*x+bpa;
0051    w=w*bma;
0052 end
0053 
0054 return

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