The Elk Code
rdmdexcdn.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: rdmdexcdn
8 ! !INTERFACE:
9 subroutine rdmdexcdn(dedn)
10 ! !USES:
11 use modmain
12 use modrdm
13 ! !INPUT/OUTPUT PARAMETERS:
14 ! dedn : energy derivative (inout,real(nstsv,nkpt))
15 ! !DESCRIPTION:
16 ! Calculates the derivative of the exchange-correlation energy w.r.t.
17 ! occupation numbers 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 real(8), intent(inout) :: dedn(nstsv,nkpt)
26 ! local variables
27 integer ik1,ik2,jk,iv(3)
28 integer ist1,ist2
29 ! parameter for calculating the functional derivatives
30 real(8), parameter :: eps=1.d-12
31 real(8) t1,t2,t3,t4
32 ! allocatable arays
33 real(8), allocatable :: vcl1221(:,:,:)
34 if (rdmxctype == 0) return
35 ! calculate the prefactor
36 if (rdmxctype == 1) then
37  t1=1.d0/occmax
38 ! Power functional
39 else if (rdmxctype == 2) then
40  if (spinpol) then
41  t1=rdmalpha
42  else
43  t1=2.d0*rdmalpha*(0.25d0)**rdmalpha
44  end if
45 else
46  write(*,*)
47  write(*,'("Error(rdmdexcdn): rdmxctype not defined : ",I8)') rdmxctype
48  write(*,*)
49  stop
50 end if
51 allocate(vcl1221(nstsv,nstsv,nkpt))
52 ! start loop over non-reduced k-points
53 do ik1=1,nkptnr
54 ! get the Coulomb matrix elements
55  call getvcl1221(ik1,vcl1221)
56 ! find the equivalent reduced k-point
57  iv(:)=ivk(:,ik1)
58  jk=ivkik(iv(1),iv(2),iv(3))
59 ! loop over reduced k-points
60  do ik2=1,nkpt
61  do ist1=1,nstsv
62  do ist2=1,nstsv
63 ! Hartree-Fock functional
64  if (rdmxctype == 1) then
65  t2=t1*occsv(ist1,jk)
66 ! Power functional
67  else if (rdmxctype == 2) then
68  t3=sum(abs(vkl(:,ik2)-vkl(:,ik1)))
69  if ((ist2 == ist1).and.(t3 < epslat)) then
70  t2=(1.d0/occmax)*occsv(ist1,jk)
71  else
72  t3=max(occsv(ist2,ik2),eps)
73  t4=max(occsv(ist1,jk),eps)
74  t2=t1*(t4**rdmalpha)/(t3**(1.d0-rdmalpha))
75  end if
76  end if
77  dedn(ist2,ik2)=dedn(ist2,ik2)+t2*vcl1221(ist2,ist1,ik2)
78  end do
79  end do
80  end do
81 end do
82 deallocate(vcl1221)
83 end subroutine
84 !EOC
85 
logical spinpol
Definition: modmain.f90:228
integer, dimension(:,:,:), allocatable ivkik
Definition: modmain.f90:467
subroutine getvcl1221(ikp, vcl1221)
Definition: getvcl1221.f90:10
subroutine rdmdexcdn(dedn)
Definition: rdmdexcdn.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
integer, dimension(:,:), allocatable ivk
Definition: modmain.f90:465