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:
9
subroutine
rdmengyxc
10
! !USES:
11
use
modmain
12
use
modrdm
13
! !DESCRIPTION:
14
! Calculates RDMFT exchange-correlation energy.
15
!
16
! !REVISION HISTORY:
17
! Created 2008 (Sharma)
18
!EOP
19
!BOC
20
implicit none
21
! local variables
22
integer
ik1,ik2,jk,iv(3)
23
integer
ist1,ist2
24
real
(8) t1,t2,t3,t4
25
! allocatable arays
26
real
(8),
allocatable
:: vcl1221(:,:,:)
27
! calculate the prefactor
28
if
(
rdmxctype
== 0)
then
29
engyx
=0.d0
30
return
31
else
if
(
rdmxctype
== 1)
then
32
! Hartree-Fock functional
33
t1=0.5d0/
occmax
34
else
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
41
else
42
write
(*,*)
43
write
(*,
'("Error(rdmengyxc): rdmxctype not defined : ",I8)'
)
rdmxctype
44
write
(*,*)
45
stop
46
end if
47
! exchange-correlation energy
48
engyx
=0.d0
49
allocate
(vcl1221(
nstsv
,
nstsv
,
nkpt
))
50
! start loop over non-reduced k-points
51
do
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
77
end do
78
deallocate
(vcl1221)
79
end subroutine
80
!EOC
81
getvcl1221
subroutine getvcl1221(ikp, vcl1221)
Definition
getvcl1221.f90:10
modmain
Definition
modmain.f90:6
modmain::wkpt
real(8), dimension(:), allocatable wkpt
Definition
modmain.f90:475
modmain::nkptnr
integer nkptnr
Definition
modmain.f90:463
modmain::spinpol
logical spinpol
Definition
modmain.f90:228
modmain::ivk
integer, dimension(:,:), allocatable ivk
Definition
modmain.f90:465
modmain::epslat
real(8) epslat
Definition
modmain.f90:24
modmain::nkpt
integer nkpt
Definition
modmain.f90:461
modmain::engyx
real(8) engyx
Definition
modmain.f90:972
modmain::nstsv
integer nstsv
Definition
modmain.f90:886
modmain::vkl
real(8), dimension(:,:), allocatable vkl
Definition
modmain.f90:471
modmain::occmax
real(8) occmax
Definition
modmain.f90:898
modmain::ivkik
integer, dimension(:,:,:), allocatable ivkik
Definition
modmain.f90:467
modmain::occsv
real(8), dimension(:,:), allocatable occsv
Definition
modmain.f90:902
modrdm
Definition
modrdm.f90:6
modrdm::rdmxctype
integer rdmxctype
Definition
modrdm.f90:21
modrdm::rdmalpha
real(8) rdmalpha
Definition
modrdm.f90:29
rdmengyxc
subroutine rdmengyxc
Definition
rdmengyxc.f90:10
rdmengyxc.f90
Generated by
1.9.8