The Elk Code
 
Loading...
Searching...
No Matches
rdmengyxc.f90
Go to the documentation of this file.
1
2! Copyright (C) 2002-2008 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3! This file is distributed under the terms of the GNU Lesser General Public
4! License. See the file COPYING for license details.
5
6!BOP
7! !ROUTINE: rdmengyxc
8! !INTERFACE:
9subroutine rdmengyxc
10! !USES:
11use modmain
12use modrdm
13! !DESCRIPTION:
14! Calculates RDMFT exchange-correlation energy.
15!
16! !REVISION HISTORY:
17! Created 2008 (Sharma)
18!EOP
19!BOC
20implicit none
21! local variables
22integer ik1,ik2,jk,iv(3)
23integer ist1,ist2
24real(8) t1,t2,t3,t4
25! allocatable arays
26real(8), allocatable :: vcl1221(:,:,:)
27! calculate the prefactor
28if (rdmxctype == 0) then
29 engyx=0.d0
30 return
31else if (rdmxctype == 1) then
32! Hartree-Fock functional
33 t1=0.5d0/occmax
34else if (rdmxctype == 2) then
35! Power functional
36 if (spinpol) then
37 t1=0.5d0
38 else
39 t1=(0.25d0)**rdmalpha
40 end if
41else
42 write(*,*)
43 write(*,'("Error(rdmengyxc): rdmxctype not defined : ",I8)') rdmxctype
44 write(*,*)
45 stop
46end if
47! exchange-correlation energy
48engyx=0.d0
49allocate(vcl1221(nstsv,nstsv,nkpt))
50! start loop over non-reduced k-points
51do ik1=1,nkptnr
52 call getvcl1221(ik1,vcl1221)
53! find the equivalent reduced k-point
54 iv(:)=ivk(:,ik1)
55 jk=ivkik(iv(1),iv(2),iv(3))
56 do ist1=1,nstsv
57! start loop over reduced k-points
58 do ik2=1,nkpt
59 do ist2=1,nstsv
60! Hartree-Fock functional
61 if (rdmxctype == 1) then
62 t2=t1*wkpt(ik2)*occsv(ist2,ik2)*occsv(ist1,jk)
63! Power functional
64 else if (rdmxctype == 2) then
65 t3=occsv(ist2,ik2)*occsv(ist1,jk)
66 t4=sum(abs(vkl(:,ik2)-vkl(:,ik1)))
67 if ((ist2 == ist1).and.(t4 < epslat)) then
68 t2=(0.5d0/occmax)*wkpt(ik2)*t3
69 else
70 t2=t1*wkpt(ik2)*(t3**rdmalpha)
71 end if
72 end if
73 engyx=engyx-t2*vcl1221(ist2,ist1,ik2)
74 end do
75 end do
76 end do
77end do
78deallocate(vcl1221)
79end subroutine
80!EOC
81
subroutine getvcl1221(ikp, vcl1221)
real(8), dimension(:), allocatable wkpt
Definition modmain.f90:475
integer nkptnr
Definition modmain.f90:463
logical spinpol
Definition modmain.f90:228
integer, dimension(:,:), allocatable ivk
Definition modmain.f90:465
real(8) epslat
Definition modmain.f90:24
integer nkpt
Definition modmain.f90:461
real(8) engyx
Definition modmain.f90:972
integer nstsv
Definition modmain.f90:886
real(8), dimension(:,:), allocatable vkl
Definition modmain.f90:471
real(8) occmax
Definition modmain.f90:898
integer, dimension(:,:,:), allocatable ivkik
Definition modmain.f90:467
real(8), dimension(:,:), allocatable occsv
Definition modmain.f90:902
integer rdmxctype
Definition modrdm.f90:21
real(8) rdmalpha
Definition modrdm.f90:29
subroutine rdmengyxc
Definition rdmengyxc.f90:10