The Elk Code
genvfxc.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2011 S. Sharma, J. K. Dewhurst 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 subroutine genvfxc(tq0,t3hw,gclgq,nm,vchi0,eps0,epsi,vfxc)
7 use modmain
8 use modtddft
9 implicit none
10 ! arguments
11 logical, intent(in) :: tq0,t3hw
12 real(8), intent(in) :: gclgq(ngrf)
13 integer, intent(in) :: nm
14 complex(8), intent(in) :: vchi0(nm,nm,nwrf)
15 complex(8), intent(in) :: eps0(nm,nm,nwrf)
16 complex(8), intent(in) :: epsi(nm,nm,nwrf)
17 complex(8), intent(out) :: vfxc(nm,nm,nwrf)
18 ! local variables
19 integer iw,i,j
20 complex(8) z1
21 ! allocatable arrays
22 complex(8), allocatable :: a(:,:)
23 ! compute v⁻¹⸍² f_xc v⁻¹⸍²
24 select case(fxctype(1))
25 case(0,1)
26 ! RPA
27  vfxc(1:nm,1:nm,1:nwrf)=0.d0
28  return
29 case(3)
30 ! ALDA
31  if (tq0.and.t3hw) then
32  call genvfxcg(gclgq,nm,vfxc(3,3,1))
33 ! the head and wings are zero
34  vfxc(1:3,1:nm,1:nwrf)=0.d0
35  vfxc(4:,1:3,1:nwrf)=0.d0
36  else
37  call genvfxcg(gclgq,nm,vfxc)
38  end if
39 case(200)
40 ! long-range contribution with dynamic correlations
41  vfxc(1:nm,1:nm,1:nwrf)=0.d0
42  do i=1,nm
43  vfxc(i,i,1:nwrf)=-(fxclrc(1)+fxclrc(2)*dble(wrf(1:nwrf))**2)/fourpi
44  end do
45 case(210,211)
46 ! bootstrap
47  vfxc(1:nm,1:nm,1:nwrf)=0.d0
48  if (tq0.and.t3hw) then
49  z1=(eps0(1,1,1)+eps0(2,2,1)+eps0(3,3,1))/3.d0
50  else
51  z1=eps0(1,1,1)
52  end if
53  z1=-1.d0/(z1-1.d0)
54  do j=1,nm
55  do i=1,nm
56  vfxc(i,j,1:nwrf)=z1*epsi(i,j,1)
57  end do
58  end do
59 case default
60  write(*,*)
61  write(*,'("Error(genvfxc): fxctype not defined : ",3I8)') fxctype
62  write(*,*)
63  stop
64 end select
65 ! right multiply by v¹⸍² χ₀ v¹⸍²
66 allocate(a(nm,nm))
67 do iw=1,nwrf
68  a(1:nm,1:nm)=vfxc(1:nm,1:nm,iw)
69  call zgemm('N','N',nm,nm,nm,zone,a,nm,vchi0(:,:,iw),nm,zzero,vfxc(:,:,iw),nm)
70 end do
71 deallocate(a)
72 end subroutine
73 
complex(8), dimension(:), allocatable wrf
Definition: modmain.f90:1170
complex(8), parameter zone
Definition: modmain.f90:1240
subroutine genvfxcg(gclgq, nm, vfxc)
Definition: genvfxcg.f90:7
real(8), dimension(2) fxclrc
Definition: modtddft.f90:14
complex(8), parameter zzero
Definition: modmain.f90:1240
real(8), parameter fourpi
Definition: modmain.f90:1234
subroutine genvfxc(tq0, t3hw, gclgq, nm, vchi0, eps0, epsi, vfxc)
Definition: genvfxc.f90:7
integer, dimension(3) fxctype
Definition: modtddft.f90:12