The Elk Code
Loading...
Searching...
No Matches
genveedu.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2007 F. Bultmark, F. Cricchio and L. Nordstrom.
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: genveedu
8
! !INTERFACE:
9
subroutine
genveedu
(idu,u,j,vee)
10
! !USES:
11
use
modmain
12
use
moddftu
13
! !INPUT/OUTPUT PARAMETERS:
14
! idu : DFT+U entry (in,integer)
15
! u : parameter U (out,real)
16
! j : parameter J (out,real)
17
! vee : Coulomb matrix elements (out,real(-lmaxdm:lmaxdm,-lmaxdm:lmaxdm,
18
! -lmaxdm:lmaxdm,-lmaxdm:lmaxdm))
19
! !DESCRIPTION:
20
! Calculates the Coulomb matrix elements used in DFT+U calculations. See
21
! {\it Phys. Rev. B} {\bf 52}, 5467 (1995).
22
!
23
! !REVISION HISTORY:
24
! Created November 2007 (FC,JKD,FB,LN)
25
! Modified July 2009 (FC)
26
!EOP
27
!BOC
28
implicit none
29
! arguments
30
integer
,
intent(in)
:: idu
31
real
(8),
intent(out)
:: u,j
32
real
(8),
intent(out)
:: vee(-lmaxdm:lmaxdm,-lmaxdm:lmaxdm,-lmaxdm:lmaxdm, &
33
-lmaxdm:lmaxdm)
34
! local variables
35
integer
l,m1,m2,m3,m4,k,q
36
real
(8) sm1,sm2,t1
37
! automatic arrays
38
real
(8) :: f(0:2*lmaxdm)
39
! external functions
40
real
(8),
external
:: gaunt
41
l=
isldu
(2,idu)
42
! calculate Slater integrals
43
call
genfdu
(idu,u,j,f)
44
do
m1=-l,l
45
do
m2=-l,l
46
do
m3=-l,l
47
do
m4=-l,l
48
sm1=0.d0
49
do
k=0,2*l,2
50
sm2=0.d0
51
do
q=-k,k
52
t1=
gaunt
(l,k,l,m1,q,m2)*
gaunt
(l,k,l,m3,-q,m4)
53
if
(mod(q,2) == 0)
then
54
sm2=sm2+t1
55
else
56
sm2=sm2-t1
57
end if
58
end do
59
sm1=sm1+f(k)*sm2/dble(2*k+1)
60
end do
61
vee(m1,m3,m2,m4)=
fourpi
*sm1
62
end do
63
end do
64
end do
65
end do
66
end subroutine
67
!EOC
68
gaunt
real(8) function gaunt(l1, l2, l3, m1, m2, m3)
Definition
gaunt.f90:10
genfdu
subroutine genfdu(idu, u, j, f)
Definition
genfdu.f90:10
genveedu
subroutine genveedu(idu, u, j, vee)
Definition
genveedu.f90:10
moddftu
Definition
moddftu.f90:6
moddftu::isldu
integer, dimension(2, maxdftu) isldu
Definition
moddftu.f90:40
modmain
Definition
modmain.f90:6
modmain::fourpi
real(8), parameter fourpi
Definition
modmain.f90:1231
genveedu.f90
Generated by
1.9.8