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
11 integer, intent(in
) :: ikp
12 complex(8), intent(in
) :: evecsv(nstsv
,nstsv
)
13 complex(8), intent(inout
) :: dedc(nstsv
,nstsv
)
16 integer ist1
,ist2
,ist3
,ist4
19 complex(8), allocatable
:: vnl(:,:,:,:)
23 if (rdmxctype
.eq
.0) return
24 ! calculate the prefactor
25 if (rdmxctype
.eq
.1) then
27 else if (rdmxctype
.eq
.2) then
31 t1
=2.d0
*(0.25d0)**rdmalpha
35 write(*,'("Error(rdmdexcdc): rdmxctype not defined : ",I8)') rdmxctype
39 allocate(vnl(nstsv
,nstsv
,nstsv
,nkptnr
))
40 ! calculate non-local matrix elements of the type (l-jj-k)
42 ! start loop over non-reduced k-points
44 ! copy the matrix elements of the type i-jj-i to vnlrdm
47 vnlrdm(ist1
,ikp
,ist2
,ik
)=dble(vnl(ist1
,ist1
,ist2
,ik
))
50 ! find the equivalent reduced k-point
52 jk
=ikmap(iv(1),iv(2),iv(3))
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
62 t2
=t1
*(occsv(ist3
,ikp
)*occsv(ist4
,jk
))**rdmalpha
64 dedc(ist2
,ist3
)=dedc(ist2
,ist3
)-t2
*evecsv(ist2
,ist1
)* &
65 vnl(ist1
,ist3
,ist4
,ik
)
70 ! end loop over non-reduced k-points