The Elk Code
rdmdexcdc.f90
Go to the documentation of this file.
1 
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.
5 
6 !BOP
7 ! !ROUTINE: rdmdexcdc
8 ! !INTERFACE:
9 subroutine rdmdexcdc(dedc)
10 ! !USES:
11 use modmain
12 use modrdm
13 ! !INPUT/OUTPUT PARAMETERS:
14 ! dedc : energy derivative (inout,complex(nstsv,nstsv,nkpt))
15 ! !DESCRIPTION:
16 ! Calculates the derivative of the exchange-correlation energy w.r.t.
17 ! {\tt evecsv} and adds the result to the total.
18 !
19 ! !REVISION HISTORY:
20 ! Created 2008 (Sharma)
21 !EOP
22 !BOC
23 implicit none
24 ! arguments
25 complex(8), intent(inout) :: dedc(nstsv,nstsv,nkpt)
26 ! local variables
27 integer ik1,ik2,jk,iv(3)
28 integer ist1,ist2,ist3,ist4
29 real(8) t1,t2,t3
30 ! allocatable arrays
31 complex(8), allocatable :: vcl1223(:,:,:,:),evecsv(:,:)
32 if (rdmxctype == 0) return
33 ! calculate the prefactor
34 if (rdmxctype == 1) then
35 ! Hartree-Fock functional
36  t1=1.d0/occmax
37 else if (rdmxctype == 2) then
38 ! Power functional
39  if (spinpol) then
40  t1=1.d0
41  else
42  t1=2.d0*(0.25d0)**rdmalpha
43  end if
44 else
45  write(*,*)
46  write(*,'("Error(rdmdexcdc): rdmxctype not defined : ",I8)') rdmxctype
47  write(*,*)
48  stop
49 end if
50 allocate(vcl1223(nstsv,nstsv,nstsv,nkpt))
51 allocate(evecsv(nstsv,nstsv))
52 ! start loop over non-reduced k-points
53 do ik1=1,nkptnr
54 ! get the Coulomb matrix elements
55  call getvcl1223(ik1,vcl1223)
56 ! find the equivalent reduced k-point
57  iv(:)=ivk(:,ik1)
58  jk=ivkik(iv(1),iv(2),iv(3))
59 ! start loop over reduced k-points
60  do ik2=1,nkpt
61 ! get the eigenvectors from file
62  call getevecsv(filext,ik2,vkl(:,ik2),evecsv)
63  do ist4=1,nstsv
64  do ist3=1,nstsv
65  do ist2=1,nstsv
66  do ist1=1,nstsv
67  if (rdmxctype == 1) then
68 ! Hartree-Fock functional
69  t2=t1*occsv(ist3,ik2)*occsv(ist4,jk)
70  else if (rdmxctype == 2) then
71 ! Power functional
72  t3=sum(abs(vkl(:,ik2)-vkl(:,ik1)))
73  if ((ist3 == ist4).and.(t3 < epslat)) then
74  t2=(1.d0/occmax)*occsv(ist4,jk)**2
75  else
76  t2=t1*(occsv(ist3,ik2)*occsv(ist4,jk))**rdmalpha
77  end if
78  end if
79  dedc(ist2,ist3,ik2)=dedc(ist2,ist3,ik2)-t2*evecsv(ist2,ist1)* &
80  vcl1223(ist1,ist3,ist4,ik2)
81  end do
82  end do
83  end do
84  end do
85 ! end loop over reduced k-points
86  end do
87 ! end loop over non-reduced k-points
88 end do
89 deallocate(vcl1223,evecsv)
90 end subroutine
91 !EOC
92 
subroutine getevecsv(fext, ikp, vpl, evecsv)
Definition: getevecsv.f90:7
character(256) filext
Definition: modmain.f90:1301
logical spinpol
Definition: modmain.f90:228
integer, dimension(:,:,:), allocatable ivkik
Definition: modmain.f90:467
subroutine getvcl1223(ikp, vcl1223)
Definition: getvcl1223.f90:10
integer nkptnr
Definition: modmain.f90:463
Definition: modrdm.f90:6
real(8) occmax
Definition: modmain.f90:901
real(8), dimension(:,:), allocatable occsv
Definition: modmain.f90:905
real(8) rdmalpha
Definition: modrdm.f90:29
real(8), dimension(:,:), allocatable vkl
Definition: modmain.f90:471
real(8) epslat
Definition: modmain.f90:24
integer rdmxctype
Definition: modrdm.f90:21
subroutine rdmdexcdc(dedc)
Definition: rdmdexcdc.f90:10
integer, dimension(:,:), allocatable ivk
Definition: modmain.f90:465