exciting-0.9.150
[exciting.git] / src / rdmdexcdc.f90
blobe6f8c9d3cdf11a489627c2aa96ebf961976598c2
2 ! Copyright (C) 2007-2008 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
6 subroutine rdmdexcdc(ikp,evecsv,dedc)
7 ! calculate the derivative of exchange-correlation energy w.r.t. evecsv
8 use modmain
9 implicit none
10 ! arguments
11 integer, intent(in) :: ikp
12 complex(8), intent(in) :: evecsv(nstsv,nstsv)
13 complex(8), intent(inout) :: dedc(nstsv,nstsv)
14 ! local variables
15 integer ik,jk,iv(3)
16 integer ist1,ist2,ist3,ist4
17 real(8) t1,t2
18 ! allocatable arrays
19 complex(8), allocatable :: vnl(:,:,:,:)
20 ! external functions
21 real(8) r3taxi
22 external r3taxi
23 if (rdmxctype.eq.0) return
24 ! calculate the prefactor
25 if (rdmxctype.eq.1) then
26 t1=1.d0/occmax
27 else if (rdmxctype.eq.2) then
28 if (spinpol) then
29 t1=1.d0
30 else
31 t1=2.d0*(0.25d0)**rdmalpha
32 end if
33 else
34 write(*,*)
35 write(*,'("Error(rdmdexcdc): rdmxctype not defined : ",I8)') rdmxctype
36 write(*,*)
37 stop
38 end if
39 allocate(vnl(nstsv,nstsv,nstsv,nkptnr))
40 ! calculate non-local matrix elements of the type (l-jj-k)
41 call rdmvnlc(ikp,vnl)
42 ! start loop over non-reduced k-points
43 do ik=1,nkptnr
44 ! copy the matrix elements of the type i-jj-i to vnlrdm
45 do ist1=1,nstsv
46 do ist2=1,nstsv
47 vnlrdm(ist1,ikp,ist2,ik)=dble(vnl(ist1,ist1,ist2,ik))
48 end do
49 end do
50 ! find the equivalent reduced k-point
51 iv(:)=ivknr(:,ik)
52 jk=ikmap(iv(1),iv(2),iv(3))
53 do ist1=1,nstsv
54 do ist2=1,nstsv
55 do ist3=1,nstsv
56 do ist4=1,nstsv
57 if (rdmxctype.eq.1) then
58 ! Hartree-Fock functional
59 t2=t1*occsv(ist3,ikp)*occsv(ist4,jk)
60 else if (rdmxctype.eq.2) then
61 ! SDLG functional
62 t2=t1*(occsv(ist3,ikp)*occsv(ist4,jk))**rdmalpha
63 end if
64 dedc(ist2,ist3)=dedc(ist2,ist3)-t2*evecsv(ist2,ist1)* &
65 vnl(ist1,ist3,ist4,ik)
66 end do
67 end do
68 end do
69 end do
70 ! end loop over non-reduced k-points
71 end do
72 deallocate(vnl)
73 return
74 end subroutine