The Elk Code
 
Loading...
Searching...
No Matches
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
6subroutine genvfxc(tq0,t3hw,gclgq,nm,vchi0,eps0,epsi,vfxc)
7use modmain
8use modtddft
9implicit none
10! arguments
11logical, intent(in) :: tq0,t3hw
12real(8), intent(in) :: gclgq(ngrf)
13integer, intent(in) :: nm
14complex(8), intent(in) :: vchi0(nm,nm,nwrf)
15complex(8), intent(in) :: eps0(nm,nm,nwrf)
16complex(8), intent(in) :: epsi(nm,nm,nwrf)
17complex(8), intent(out) :: vfxc(nm,nm,nwrf)
18! local variables
19integer iw,i,j
20complex(8) z1
21! allocatable arrays
22complex(8), allocatable :: a(:,:)
23! compute v⁻¹⸍² f_xc v⁻¹⸍²
24select case(fxctype(1))
25case(0,1)
26! RPA
27 vfxc(1:nm,1:nm,1:nwrf)=0.d0
28 return
29case(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
39case(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
45case(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
59case default
60 write(*,*)
61 write(*,'("Error(genvfxc): fxctype not defined : ",3I8)') fxctype
62 write(*,*)
63 stop
64end select
65! right multiply by v¹⸍² χ₀ v¹⸍²
66allocate(a(nm,nm))
67do 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)
70end do
71deallocate(a)
72end subroutine
73
subroutine genvfxc(tq0, t3hw, gclgq, nm, vchi0, eps0, epsi, vfxc)
Definition genvfxc.f90:7
subroutine genvfxcg(gclgq, nm, vfxc)
Definition genvfxcg.f90:7
complex(8), parameter zzero
Definition modmain.f90:1238
complex(8), parameter zone
Definition modmain.f90:1238
real(8), parameter fourpi
Definition modmain.f90:1231
complex(8), dimension(:), allocatable wrf
Definition modmain.f90:1167
integer, dimension(3) fxctype
Definition modtddft.f90:12
real(8), dimension(2) fxclrc
Definition modtddft.f90:14