SPEED
MAKE_LGL_NODES.f90 File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine make_lgl_nodes (np, ct)
 Makes Gauss-Legendre-Lobatto nodes.
 

Function/Subroutine Documentation

◆ make_lgl_nodes()

subroutine make_lgl_nodes ( integer*4  np,
real*8, dimension(np)  ct 
)

Makes Gauss-Legendre-Lobatto nodes.

Author
Ilario Mazzieri
Date
September, 2013
Version
1.0
Parameters
[in]npnumber of points
[out]ctLGL nodes

Definition at line 26 of file MAKE_LGL_NODES.f90.

27
28 implicit real*8 (a-h,o-z)
29
30 integer*4 :: np
31
32 real*8, dimension(np) :: ct
33
34 parameter(nstep=1000,acc=1.d-15)
35
36 n=np-1
37 ct(1)=-1.d0
38 ct(np)=1.d0
39 n2=idint(0.5d0*np)
40 dx=2.d0/dfloat(nstep)
41 x=-1.d0
42 iroot=1
43 call get_legendre_value (p2,a1,p1,p1der,n,x)
44 do while (iroot.lt.n2)
45 x=x+dx
46 call get_legendre_value (p2,a2,p1,p1der,n,x)
47 if (dabs(a2).le.acc) then
48 iroot=iroot+1
49 ct(iroot)=a2
50 else
51 aa=a1*a2
52 if (aa.lt.0.d0) then
53 iroot=iroot+1
54 ct(iroot)=find_root_pol(x-dx,x,acc,n)
55 endif
56 endif
57 a1=a2
58 enddo
59
60 n_filt=2*n2
61 if (n_filt.ne.np) then ! np odd
62 ct(n2+1)=0.d0
63 do i=1,n2-1
64 ct(np-i)=-ct(i+1)
65 enddo
66 else ! np even
67 do i=1,n2-1
68 ct(np-i)=-ct(i+1)
69 enddo
70 endif
71
72
73
74 return
75
real *8 function find_root_pol(x1, x2, xacc, n)
Find a root of a polynomial of degree n.
subroutine get_legendre_value(p2, p2der, p1, p1der, n, x)
Computes Legendre polynomial and its derivative on a given point x.

References find_root_pol(), and get_legendre_value().

Referenced by make_loc_matrix_dg().

Here is the call graph for this function:
Here is the caller graph for this function: