The Elk Code
fyukawa0.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: fyukawa0
8 ! !INTERFACE:
9 real(8) function fyukawa0(is,l,k)
10 ! !USES:
11 use modmain
12 use 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 ! !DESCRIPTION:
18 ! Calculates the Slater parameters in the unscreened case. See {\it Phys. Rev.
19 ! B} {\bf 52}, 1421 (1995) and {\it Phys. Rev. B} {\bf 80}, 035121 (2009).
20 !
21 ! !REVISION HISTORY:
22 ! Created April 2008 (LN)
23 ! Modified and tested July 2008 (LN and FC)
24 !EOP
25 !BOC
26 implicit none
27 ! arguments
28 integer, intent(in) :: is,l,k
29 ! local variables
30 integer, parameter :: nstart=1
31 integer ias,nr,ir
32 integer ir1,ir2,nr1,nr2
33 real(8) x
34 ! automatic arrays
35 real(8) bint(nrmtmax),cint(nrmtmax),fint(nrmtmax)
36 real(8) blow(nrmtmax),bhigh(nrmtmax)
37 real(8) clow(nrmtmax),chigh(nrmtmax)
38 real(8) ak(nrmtmax),bk(nrmtmax)
39 ias=idxas(1,is)
40 nr=nrmt(is)
41 ak(1:nr)=0.d0
42 bk(1:nr)=0.d0
43 ! calculate unscreened Slater parameters
44 do ir=1,nr
45  bint(ir)=(fdufr(ir,l,ias)*rsp(ir,is))**2
46  x=rsp(ir,is)**k
47  ak(ir)=x
48  bk(ir)=1.d0/(x*rsp(ir,is))
49 end do
50 do ir=nstart,nr
51  nr1=ir-nstart+1
52  nr2=nr-ir+1
53  do ir1=1,nr1
54  ir2=ir1+nstart-1
55  blow(ir1)=bint(ir2)*ak(ir2)
56  end do
57  call fderiv(-1,nr1,rsp(nstart,is),blow,clow)
58  do ir1=1,nr2
59  ir2=ir1+ir-1
60  bhigh(ir1)=bint(ir2)*bk(ir2)
61  end do
62  call fderiv(-1,nr2,rsp(ir,is),bhigh,chigh)
63  cint(ir-nstart+1)=bint(ir)*(bk(ir)*clow(nr1)+ak(ir)*chigh(nr2))
64 end do
65 nr1=nr-nstart+1
66 call fderiv(-1,nr1,rsp(nstart,is),cint,fint)
67 fyukawa0=fint(nr1)
68 end function
69 !EOC
70 
integer, dimension(maxatoms, maxspecies) idxas
Definition: modmain.f90:42
subroutine fderiv(m, n, x, f, g)
Definition: fderiv.f90:10
real(8) function fyukawa0(is, l, k)
Definition: fyukawa0.f90:10
real(8), dimension(:,:,:), allocatable fdufr
Definition: moddftu.f90:58
real(8), dimension(:,:), allocatable rsp
Definition: modmain.f90:135
integer nrmtmax
Definition: modmain.f90:152
integer, dimension(maxspecies) nrmt
Definition: modmain.f90:150