The Elk Code
wspline.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 wspline(n,x,wc)
7 implicit none
8 ! arguments
9 integer, intent(in) :: n
10 real(8), intent(in) :: x(n)
11 real(8), intent(out) :: wc(4,3,n)
12 ! local variables
13 integer i,j
14 real(8) f(4),cf(3,4)
15 if (n < 4) then
16  write(*,*)
17  write(*,'("Error(wspline): n < 4 : ",I8)') n
18  write(*,*)
19  stop
20 end if
21 f(1)=1.d0
22 f(2:4)=0.d0
23 call spline(4,x,f,cf)
24 wc(1,1:3,1)=cf(1:3,1)
25 wc(1,1:3,2)=cf(1:3,2)
26 f(1)=0.d0
27 f(2)=1.d0
28 call spline(4,x,f,cf)
29 wc(2,1:3,1)=cf(1:3,1)
30 wc(2,1:3,2)=cf(1:3,2)
31 f(2)=0.d0
32 f(3)=1.d0
33 call spline(4,x,f,cf)
34 wc(3,1:3,1)=cf(1:3,1)
35 wc(3,1:3,2)=cf(1:3,2)
36 f(3)=0.d0
37 f(4)=1.d0
38 call spline(4,x,f,cf)
39 wc(4,1:3,1)=cf(1:3,1)
40 wc(4,1:3,2)=cf(1:3,2)
41 do i=3,n-3
42  j=i-1
43  f(1)=1.d0
44  f(2:4)=0.d0
45  call spline(4,x(j),f,cf)
46  wc(1,1:3,i)=cf(1:3,2)
47  f(1)=0.d0
48  f(2)=1.d0
49  call spline(4,x(j),f,cf)
50  wc(2,1:3,i)=cf(1:3,2)
51  f(2)=0.d0
52  f(3)=1.d0
53  call spline(4,x(j),f,cf)
54  wc(3,1:3,i)=cf(1:3,2)
55  f(3)=0.d0
56  f(4)=1.d0
57  call spline(4,x(j),f,cf)
58  wc(4,1:3,i)=cf(1:3,2)
59 end do
60 j=n-3
61 f(1)=1.d0
62 f(2:4)=0.d0
63 call spline(4,x(j),f,cf)
64 wc(1,1:3,n-2)=cf(1:3,2)
65 wc(1,1:3,n-1)=cf(1:3,3)
66 wc(1,1:3,n)=cf(1:3,4)
67 f(1)=0.d0
68 f(2)=1.d0
69 call spline(4,x(j),f,cf)
70 wc(2,1:3,n-2)=cf(1:3,2)
71 wc(2,1:3,n-1)=cf(1:3,3)
72 wc(2,1:3,n)=cf(1:3,4)
73 f(2)=0.d0
74 f(3)=1.d0
75 call spline(4,x(j),f,cf)
76 wc(3,1:3,n-2)=cf(1:3,2)
77 wc(3,1:3,n-1)=cf(1:3,3)
78 wc(3,1:3,n)=cf(1:3,4)
79 f(3)=0.d0
80 f(4)=1.d0
81 call spline(4,x(j),f,cf)
82 wc(4,1:3,n-2)=cf(1:3,2)
83 wc(4,1:3,n-1)=cf(1:3,3)
84 wc(4,1:3,n)=cf(1:3,4)
85 end subroutine
86 
subroutine spline(n, x, f, cf)
Definition: spline.f90:10
subroutine wspline(n, x, wc)
Definition: wspline.f90:7