The Elk Code
wsplintp.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2019 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 subroutine wsplintp(n,x,wp)
7 implicit none
8 ! arguments
9 integer, intent(in) :: n
10 real(8), intent(in) :: x(n)
11 real(8), intent(out) :: wp(4,n)
12 ! local variables
13 integer i
14 real(8) f(4),t1,t2
15 ! external functions
16 real(8), external :: polynm
17 if (n < 4) then
18  write(*,*)
19  write(*,'("Error(wsplintp): n < 4 : ",I8)') n
20  write(*,*)
21  stop
22 end if
23 wp(:,1)=0.d0
24 f(1)=1.d0
25 f(2:)=0.d0
26 wp(1,2)=polynm(-1,4,x,f,x(2))
27 f(1)=0.d0
28 f(2)=1.d0
29 wp(2,2)=polynm(-1,4,x,f,x(2))
30 f(2)=0.d0
31 f(3)=1.d0
32 wp(3,2)=polynm(-1,4,x,f,x(2))
33 f(3)=0.d0
34 f(4)=1.d0
35 wp(4,2)=polynm(-1,4,x,f,x(2))
36 do i=3,n-1
37  f(1)=1.d0
38  f(2:)=0.d0
39  t1=polynm(-1,4,x(i-2),f,x(i-1))
40  t2=polynm(-1,4,x(i-2),f,x(i))
41  wp(1,i)=t2-t1
42  f(1)=0.d0
43  f(2)=1.d0
44  t1=polynm(-1,4,x(i-2),f,x(i-1))
45  t2=polynm(-1,4,x(i-2),f,x(i))
46  wp(2,i)=t2-t1
47  f(2)=0.d0
48  f(3)=1.d0
49  t1=polynm(-1,4,x(i-2),f,x(i-1))
50  t2=polynm(-1,4,x(i-2),f,x(i))
51  wp(3,i)=t2-t1
52  f(3)=0.d0
53  f(4)=1.d0
54  t1=polynm(-1,4,x(i-2),f,x(i-1))
55  t2=polynm(-1,4,x(i-2),f,x(i))
56  wp(4,i)=t2-t1
57 end do
58 f(1)=1.d0
59 f(2:)=0.d0
60 t1=polynm(-1,4,x(n-3),f,x(n-1))
61 t2=polynm(-1,4,x(n-3),f,x(n))
62 wp(1,n)=t2-t1
63 f(1)=0.d0
64 f(2)=1.d0
65 t1=polynm(-1,4,x(n-3),f,x(n-1))
66 t2=polynm(-1,4,x(n-3),f,x(n))
67 wp(2,n)=t2-t1
68 f(2)=0.d0
69 f(3)=1.d0
70 t1=polynm(-1,4,x(n-3),f,x(n-1))
71 t2=polynm(-1,4,x(n-3),f,x(n))
72 wp(3,n)=t2-t1
73 f(3)=0.d0
74 f(4)=1.d0
75 t1=polynm(-1,4,x(n-3),f,x(n-1))
76 t2=polynm(-1,4,x(n-3),f,x(n))
77 wp(4,n)=t2-t1
78 end subroutine
79 
subroutine wsplintp(n, x, wp)
Definition: wsplintp.f90:7