The Elk Code
 
Loading...
Searching...
No Matches
fyukawa.f90
Go to the documentation of this file.
1
2! Copyright (C) 2008 F. Bultmark, F. Cricchio and L. Nordstrom.
3! This file is distributed under the terms of the GNU General Public License.
4! See the file COPYING for license details.
5
6!BOP
7! !ROUTINE: fyukawa
8! !INTERFACE:
9real(8) function fyukawa(is,l,k,lambda)
10! !USES:
11use modmain
12use moddftu
13! !INPUT/OUTPUT PARAMETERS:
14! is : species type (in,integer)
15! l : an angular momentum (in,integer)
16! k : order of Slater parameter (in,integer)
17! lambda : screening length of Yukawa potential (in,real)
18! !DESCRIPTION:
19! Calculates the Slater parameters using a screened Yukawa potential. See
20! {\it Phys. Rev. B} {\bf 52}, 1421 (1995) and {\it Phys. Rev. B} {\bf 80},
21! 035121 (2009).
22!
23! !REVISION HISTORY:
24! Created April 2008 (Lars Nordstrom)
25! Modified and tested July 2008 (LN and FC)
26!EOP
27!BOC
28implicit none
29! arguments
30integer, intent(in) :: is,l,k
31real(8), intent(in) :: lambda
32! local variables
33integer ias,nr,ir
34integer nr1,nr2,ir1,ir2
35real(8) x,t1
36! automatic arrays
37real(8) bint(nrmtmax),cint(nrmtmax),fint(nrmtmax)
38real(8) blow(nrmtmax),bhigh(nrmtmax)
39real(8) clow(nrmtmax),chigh(nrmtmax)
40real(8) a(0:2*l,nrmtmax),b(0:2*l,nrmtmax)
41ias=idxas(1,is)
42nr=nrmt(is)
43! (-1)**k factor takes care of the additional (-1)**k introduced by sbesseli(b)
44t1=lambda*dble((2*k+1)*(-1)**(k+1))
45! zero all quantities
46bint(1:nr)=0.d0
47blow(1:nr)=0.d0
48bhigh(1:nr)=0.d0
49clow(1:nr)=0.d0
50chigh(1:nr)=0.d0
51a(:,1:nr)=0.d0
52b(:,1:nr)=0.d0
53! calculate screened Slater parameters
54do ir=1,nr
55 bint(ir)=(fdufr(ir,l,ias)*rsp(ir,is))**2
56 x=rsp(ir,is)*lambda
57! spherical Bessel and Hankel functions with imaginary arguments
58 call sbesseli(2*l,x,a(:,ir))
59 call shankeli(2*l,x,b(:,ir))
60end do
61do ir=1,nr
62 nr1=ir
63 nr2=nr-ir+1
64! 1st term: r1 < r
65 do ir1=1,nr1
66 ir2=ir1
67 blow(ir1)=bint(ir2)*a(k,ir2)
68 end do
69! integrate 1st term
70 call fderiv(-1,nr1,rsp(1,is),blow,clow)
71! 2nd term : r2 > r
72 do ir1=1,nr2
73 ir2=ir1+ir-1
74 bhigh(ir1)=bint(ir2)*b(k,ir2)
75 end do
76! integrate 2nd term
77 call fderiv(-1,nr2,rsp(ir,is),bhigh,chigh)
78! sum of the two terms
79 cint(ir)=bint(ir)*(b(k,ir)*clow(nr1)+a(k,ir)*chigh(nr2))
80end do
81call fderiv(-1,nr,rsp(1,is),cint,fint)
82fyukawa=t1*fint(nr)
83end function
84!EOC
85
subroutine fderiv(m, n, x, f, g)
Definition fderiv.f90:10
real(8) function fyukawa(is, l, k, lambda)
Definition fyukawa.f90:10
real(8), dimension(:,:,:), allocatable fdufr
Definition moddftu.f90:58
integer, dimension(maxspecies) nrmt
Definition modmain.f90:150
integer, dimension(maxatoms, maxspecies) idxas
Definition modmain.f90:42
real(8), dimension(:,:), allocatable rsp
Definition modmain.f90:135
integer nrmtmax
Definition modmain.f90:152
subroutine sbesseli(lmax, x, jl)
Definition sbesseli.f90:10
subroutine shankeli(lmax, x, hl)
Definition shankeli.f90:10