The Elk Code
 
Loading...
Searching...
No Matches
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
6subroutine wsplintp(n,x,wp)
7implicit none
8! arguments
9integer, intent(in) :: n
10real(8), intent(in) :: x(n)
11real(8), intent(out) :: wp(4,n)
12! local variables
13integer i
14real(8) f(4),t1,t2
15! external functions
16real(8), external :: polynm
17if (n < 4) then
18 write(*,*)
19 write(*,'("Error(wsplintp): n < 4 : ",I8)') n
20 write(*,*)
21 stop
22end if
23wp(:,1)=0.d0
24f(1)=1.d0
25f(2:)=0.d0
26wp(1,2)=polynm(-1,4,x,f,x(2))
27f(1)=0.d0
28f(2)=1.d0
29wp(2,2)=polynm(-1,4,x,f,x(2))
30f(2)=0.d0
31f(3)=1.d0
32wp(3,2)=polynm(-1,4,x,f,x(2))
33f(3)=0.d0
34f(4)=1.d0
35wp(4,2)=polynm(-1,4,x,f,x(2))
36do 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
57end do
58f(1)=1.d0
59f(2:)=0.d0
60t1=polynm(-1,4,x(n-3),f,x(n-1))
61t2=polynm(-1,4,x(n-3),f,x(n))
62wp(1,n)=t2-t1
63f(1)=0.d0
64f(2)=1.d0
65t1=polynm(-1,4,x(n-3),f,x(n-1))
66t2=polynm(-1,4,x(n-3),f,x(n))
67wp(2,n)=t2-t1
68f(2)=0.d0
69f(3)=1.d0
70t1=polynm(-1,4,x(n-3),f,x(n-1))
71t2=polynm(-1,4,x(n-3),f,x(n))
72wp(3,n)=t2-t1
73f(3)=0.d0
74f(4)=1.d0
75t1=polynm(-1,4,x(n-3),f,x(n-1))
76t2=polynm(-1,4,x(n-3),f,x(n))
77wp(4,n)=t2-t1
78end subroutine
79
pure real(8) function polynm(m, np, xa, ya, x)
Definition polynm.f90:10
subroutine wsplintp(n, x, wp)
Definition wsplintp.f90:7