TABLE OF CONTENTS
ABINIT/der_int [ Functions ]
NAME
der_int
FUNCTION
Given input function f(i) on Teter radial grid, and grid spacing dr(i), compute function derivative df/dr on points from 0 to n. Integrate function f(i) on grid r(i) from r(0) to r(nlast). Note that array dimensions start at 0.
COPYRIGHT
Copyright (C) 1998-2017 ABINIT group (DCA, XG, GMR) This file is distributed under the terms of the GNU General Public License, see ~abinit/COPYING or http://www.gnu.org/copyleft/gpl.txt . For the initials of contributors, see ~abinit/doc/developers/contributors.txt .
INPUTS
f(0 to nlast)=function values on grid r(0 to nlast)=radial grid points dr(0 to nlast)=INVERSE of spacing on grid nlast=radial grid point for upper limit
OUTPUT
df(0 to n)=derivative $ \frac{df}{dr}$ on grid smf= $ \int_{r(0)}^{r(nlast)} f(r) dr $.
PARENTS
psp1lo,psp1nl
CHILDREN
SOURCE
36 #if defined HAVE_CONFIG_H 37 #include "config.h" 38 #endif 39 40 #include "abi_common.h" 41 42 43 subroutine der_int(ff,df,rr,dr,nlast,smf) 44 45 use defs_basis 46 use m_errors 47 use m_profiling_abi 48 49 !This section has been created automatically by the script Abilint (TD). 50 !Do not modify the following lines by hand. 51 #undef ABI_FUNC 52 #define ABI_FUNC 'der_int' 53 !End of the abilint section 54 55 implicit none 56 57 !Arguments ------------------------------------ 58 !nmax sets standard number of grid points ! SHOULD BE REMOVED 59 !scalars 60 integer,parameter :: nmax=2000 61 integer,intent(in) :: nlast 62 real(dp),intent(out) :: smf 63 !no_abirules 64 !Note that dimension here starts at 0 65 real(dp), intent(in) :: dr(0:nmax),ff(0:nmax),rr(0:nmax) 66 real(dp), intent(out) :: df(0:nmax) 67 68 !Local variables------------------------------- 69 !scalars 70 integer :: jj 71 real(dp),parameter :: div12=1.d0/12.d0 72 real(dp) :: hh 73 character(len=500) :: message 74 75 ! ************************************************************************* 76 77 !Check that nlast lie within 0 to nmax 78 if (nlast<0.or.nlast>nmax) then 79 write(message, '(a,i12,a,i12)' )& 80 & ' nlast=',nlast,' lies outside range [0,nmax] with dimension nmax=',nmax 81 MSG_BUG(message) 82 end if 83 84 !Compute derivatives at lower end, near r=0 85 df(0)=-25.d0/12.d0*ff(0)+4.d0*ff(1)-3.d0*ff(2)+4.d0/3.d0*ff(3)& 86 & -1.d0/4.d0*ff(4) 87 df(1)=-1.d0/4.d0*ff(0)-5.d0/6.d0*ff(1)+3.d0/2.d0*ff(2)& 88 & -1.d0/2.d0*ff(3)+1.d0/12.d0*ff(4) 89 90 !Run over range from just past r=0 to near r(n), using central differences 91 do jj=2,nlast-2 92 df(jj)=(ff(jj-2)-8.d0*(ff(jj-1)-ff(jj+1))-ff(jj+2))*div12 93 end do 94 95 !Compute derivative at upper end of range 96 if (nlast < 4) then 97 message = ' der_int: ff does not have enough elements. nlast is too low' 98 MSG_ERROR(message) 99 end if 100 101 df(nlast-1)=-1.d0/12.d0*ff(nlast-4)& 102 & +1.d0/2.d0*ff(nlast-3)& 103 & -3.d0/2.d0*ff(nlast-2)& 104 & +5.d0/6.d0*ff(nlast-1)& 105 & +1.d0/4.d0*ff(nlast) 106 df(nlast)=1.d0/4.d0*ff(nlast-4)& 107 & -4.d0/3.d0*ff(nlast-3)& 108 & +3.d0*ff(nlast-2)& 109 & -4.d0*ff(nlast-1)& 110 & +25.d0/12.d0*ff(nlast) 111 112 !Apply correct normalization over full range 113 do jj=0,nlast 114 df(jj)=df(jj)*dr(jj) 115 end do 116 117 smf=0.d0 118 do jj=0,nlast-1 119 hh=rr(jj+1)-rr(jj) 120 smf=smf+hh*(6.d0*(ff(jj)+ff(jj+1))+hh*(df(jj)-df(jj+1))) 121 end do 122 smf=smf/12.d0 123 124 end subroutine der_int