The Elk Code
fderiv.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: fderiv
8 ! !INTERFACE:
9 subroutine fderiv(m,n,x,f,g)
10 ! !INPUT/OUTPUT PARAMETERS:
11 ! m : order of derivative (in,integer)
12 ! n : number of points (in,integer)
13 ! x : abscissa array (in,real(n))
14 ! f : function array (in,real(n))
15 ! g : (anti-)derivative of f (out,real(n))
16 ! !DESCRIPTION:
17 ! Given function $f$ defined on a set of points $x_i$ then if $m\ge 0$ this
18 ! routine computes the $m$th derivative of $f$ at each point. If $m=-1$ the
19 ! anti-derivative of $f$ given by
20 ! $$ g(x_i)=\int_{x_1}^{x_i} f(x)\,dx $$
21 ! is calculated. Both derivatives and integrals are computed by first fitting
22 ! the function to a clamped cubic spline.
23 !
24 ! !REVISION HISTORY:
25 ! Created May 2002 (JKD)
26 !EOP
27 !BOC
28 implicit none
29 ! arguments
30 integer, intent(in) :: m,n
31 real(8), intent(in) :: x(n),f(n)
32 real(8), intent(out) :: g(n)
33 ! local variables
34 integer i
35 real(8) sm,dx
36 ! automatic arrays
37 real(8) cf(3,n)
38 if (n < 1) then
39  write(*,*)
40  write(*,'("Error(fderiv): n < 1 : ",I8)') n
41  write(*,*)
42  stop
43 end if
44 ! high accuracy integration/differentiation from spline interpolation
45 call spline(n,x,f,cf)
46 select case(m)
47 case(:-1)
48  g(1)=0.d0
49  sm=0.d0
50  do i=1,n-1
51  dx=x(i+1)-x(i)
52  sm=sm+dx*(f(i) &
53  +dx*(0.5d0*cf(1,i) &
54  +dx*(0.3333333333333333333d0*cf(2,i) &
55  +dx*0.25d0*cf(3,i))))
56  g(i+1)=sm
57  end do
58 case(1)
59  g(1:n)=cf(1,1:n)
60 case(2)
61  g(1:n)=2.d0*cf(2,1:n)
62 case(3)
63  g(1:n)=6.d0*cf(3,1:n)
64 case(4:)
65  g(1:n)=0.d0
66 end select
67 end subroutine
68 !EOC
69 
subroutine spline(n, x, f, cf)
Definition: spline.f90:10
subroutine fderiv(m, n, x, f, g)
Definition: fderiv.f90:10