SPEED
MAKE_LGL_NW.f90 File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine make_lgl_nw (nb_pnt, xq, wq, dd)
 Makes Gauss-Legendre-Lobatto nodes, weigths and spectral derivatives.
 

Function/Subroutine Documentation

◆ make_lgl_nw()

subroutine make_lgl_nw ( integer*4, intent(inout)  nb_pnt,
real*8, dimension(nb_pnt)  xq,
real*8, dimension(nb_pnt)  wq,
real*8, dimension(nb_pnt,nb_pnt)  dd 
)

Makes Gauss-Legendre-Lobatto nodes, weigths and spectral derivatives.

Author
Ilario Mazzieri
Date
September, 2013
Version
1.0
Parameters
[in]nb_pntpolynomial degree
[out]xqLGL nodes
[out]wqLGL weights
[out]ddmatrix for spectral derivatives

Definition at line 28 of file MAKE_LGL_NW.f90.

29
30 implicit real*8 (a-h,o-z)
31
32 parameter(nstep = 1000, acc = 1.d-15)
33
34
35 integer*4, intent(inout) :: nb_pnt
36 real*8, dimension(nb_pnt) :: xq, wq
37 real*8, dimension(nb_pnt) :: dd(nb_pnt,nb_pnt)
38
39 n = nb_pnt-1
40 xq(1) = -1.d0
41 xq(nb_pnt) = 1.d0
42 n2=idint(0.5d0*nb_pnt)
43 dx=2.d0/dfloat(nstep)
44
45 x = -1.d0
46 iroot = 1
47
48 call get_legendre_value(p2,a1,p1,p1der,n,x)
49
50 do while (iroot .lt. n2)
51 x = x + dx
52 call get_legendre_value(p2,a2,p1,p1der,n,x)
53
54 if (dabs(a2) .le. acc) then
55 iroot = iroot + 1
56 xq(iroot) = a2
57 else
58 aa = a1 * a2
59 if (aa .lt. 0.d0) then
60 iroot = iroot + 1
61 xq(iroot) = find_root_pol(x-dx,x,acc,n)
62 endif
63 endif
64
65 a1 = a2
66 enddo
67
68 n_filt = 2*n2
69 if (n_filt .ne. nb_pnt) then ! nb_pnt odd
70 xq(n2+1) = 0.d0
71 do i = 1, n2-1
72 xq(nb_pnt-i) = -xq(i+1)
73 enddo
74 else ! nb_pnt even
75 do i = 1, n2-1
76 xq(nb_pnt-i) = -xq(i+1)
77 enddo
78 endif
79
80
81 xn = dfloat(n)
82 acost = 2.d0/(xn*(xn+1.d0))
83
84 do j = 1,nb_pnt
85 call get_legendre_value(p2,p2der,p1,p1der,n,xq(j))
86 den = p2*p2
87 wq(j) = acost/den
88
89 do i = 1, nb_pnt
90 if (i .ne. j) then
91 call get_legendre_value(pnum,p2der,p1,p1der,n,xq(i))
92 den = p2 * (xq(i)-xq(j))
93 dd(i,j) = pnum / den
94 endif
95 enddo
96
97 enddo
98
99 do j = 2, nb_pnt-1
100 dd(j,j) = 0.d0
101 enddo
102
103 dd(1,1) = -0.25d0 * xn * (xn+1.d0)
104 dd(nb_pnt,nb_pnt) = 0.25d0 * xn * (xn+1.d0)
105
106
107 return
108
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 check_expl(), check_sism(), compute_energy_error(), compute_sdof_input(), get_highest_node(), get_max_values(), get_node_depth_and_vs(), get_node_depth_from_alluvial(), get_node_depth_from_cmplx(), get_node_depth_from_simple(), make_extint_forces(), make_spx_grid_loc(), setup_dg_elem(), and write_output().

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