The Elk Code
rfinterp.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
3 ! This file is distributed under the terms of the GNU Lesser General Public
4 ! License. See the file COPYING for license details.
5 
6 !BOP
7 ! !ROUTINE: rfinterp
8 ! !INTERFACE:
9 subroutine rfinterp(ni,xi,wci,fi,no,xo,fo)
10 ! !INPUT/OUTPUT PARAMETERS:
11 ! ni : number of input points (in,integer)
12 ! xi : input abscissa array (in,real(ni))
13 ! wci : input spline coefficient weights (in,real(12,ni))
14 ! fi : input data array (in,real(ni))
15 ! no : number of output points (in,integer)
16 ! xo : output abscissa array (in,real(no))
17 ! fo : output interpolated function (out,real(no))
18 ! !DESCRIPTION:
19 ! Given a function defined on a set of input points, this routine uses a
20 ! clamped cubic spline to interpolate the function on a different set of
21 ! points. See routine {\tt spline}.
22 !
23 ! !REVISION HISTORY:
24 ! Created January 2005 (JKD)
25 ! Arguments changed, April 2016 (JKD)
26 !EOP
27 !BOC
28 implicit none
29 ! arguments
30 integer, intent(in) :: ni
31 real(8), intent(in) :: xi(ni),wci(12,ni),fi(ni)
32 integer, intent(in) :: no
33 real(8), intent(in) :: xo(no)
34 real(8), intent(out) :: fo(no)
35 ! local variables
36 integer i,j,k,l
37 real(8) x,dx
38 ! automatic arrays
39 real(8) cf(3,ni)
40 if (ni == 1) then
41  fo(1:no)=fi(1)
42  return
43 end if
44 ! compute the spline coefficients
45 call splinew(ni,wci,fi,cf)
46 ! evaluate spline at output points
47 i=1
48 do l=1,no
49  x=xo(l)
50  if (i >= ni) i=1
51  if (x >= xi(i)) then
52  if (x > xi(i+1)) then
53 ! binary search
54  i=1
55  j=ni
56  do while (j > i+1)
57  k=(i+j)/2
58  if (x < xi(k)) then
59  j=k
60  else
61  i=k
62  end if
63  end do
64  end if
65  end if
66  dx=x-xi(i)
67  fo(l)=fi(i)+dx*(cf(1,i)+dx*(cf(2,i)+dx*cf(3,i)))
68 end do
69 end subroutine
70 !EOC
71 
pure subroutine splinew(n, wc, f, cf)
Definition: splinew.f90:7
subroutine rfinterp(ni, xi, wci, fi, no, xo, fo)
Definition: rfinterp.f90:10