4
1
mirror of https://github.com/pfloos/quack synced 2024-06-02 11:25:28 +02:00
quack/src/GW/ufG0W0.f90

505 lines
13 KiB
Fortran
Raw Normal View History

2023-11-24 15:31:29 +01:00
subroutine ufG0W0(dotest,TDA_W,nBas,nC,nO,nV,nR,nS,ENuc,ERHF,ERI,eHF)
2021-10-18 21:11:01 +02:00
! Unfold G0W0 equations
implicit none
include 'parameters.h'
! Input variables
2023-11-13 17:39:30 +01:00
logical,intent(in) :: dotest
2023-11-24 15:31:29 +01:00
logical,intent(in) :: TDA_W
2021-10-18 21:11:01 +02:00
integer,intent(in) :: nBas
integer,intent(in) :: nC
integer,intent(in) :: nO
integer,intent(in) :: nV
integer,intent(in) :: nR
integer,intent(in) :: nS
double precision,intent(in) :: ENuc
double precision,intent(in) :: ERHF
double precision,intent(in) :: ERI(nBas,nBas,nBas,nBas)
double precision,intent(in) :: eHF(nBas)
! Local variables
integer :: p
integer :: s
integer :: i,j,k,l
integer :: a,b,c,d
2023-04-20 17:00:34 +02:00
integer :: jb,kc,ia,ja
integer :: klc,kcd,ija,ijb,iab,jab
2021-10-18 21:11:01 +02:00
2023-11-29 17:16:33 +01:00
logical :: print_W = .false.
2023-11-06 11:07:20 +01:00
logical :: dRPA
2023-04-20 17:00:34 +02:00
integer :: ispin
double precision :: EcRPA
2021-10-18 21:11:01 +02:00
integer :: n2h1p,n2p1h,nH
double precision,external :: Kronecker_delta
double precision,allocatable :: H(:,:)
double precision,allocatable :: eGW(:)
double precision,allocatable :: Z(:)
2023-11-06 11:07:20 +01:00
double precision,allocatable :: Aph(:,:)
double precision,allocatable :: Bph(:,:)
double precision,allocatable :: Om(:)
double precision,allocatable :: XpY(:,:)
double precision,allocatable :: XmY(:,:)
double precision,allocatable :: rho(:,:,:)
2021-10-18 21:11:01 +02:00
2022-01-17 08:02:11 +01:00
logical :: verbose = .true.
2023-11-06 11:07:20 +01:00
double precision,parameter :: cutoff1 = 0.01d0
2022-01-17 08:02:11 +01:00
double precision,parameter :: cutoff2 = 0.01d0
2023-11-27 10:49:13 +01:00
double precision :: eF
double precision,parameter :: window = 2d0
2022-01-17 08:02:11 +01:00
2023-11-27 09:56:32 +01:00
double precision :: start_timing,end_timing,timing
2021-10-18 21:11:01 +02:00
! Output variables
! Hello world
write(*,*)
2023-11-24 15:31:29 +01:00
write(*,*)'****************************************'
write(*,*)'* Restricted Upfolded G0W0 Calculation *'
write(*,*)'****************************************'
2021-10-18 21:11:01 +02:00
write(*,*)
! Dimension of the supermatrix
2023-04-20 17:00:34 +02:00
n2h1p = nO*nO*nV
2021-10-18 21:11:01 +02:00
n2p1h = nV*nV*nO
nH = 1 + n2h1p + n2p1h
2023-11-06 11:07:20 +01:00
! Memory allocation
2021-10-18 21:11:01 +02:00
2023-11-24 15:31:29 +01:00
allocate(H(nH,nH),eGW(nH),Z(nH))
2023-11-06 11:07:20 +01:00
! Initialization
dRPA = .true.
EcRPA = 0d0
2021-10-18 21:11:01 +02:00
2023-11-27 10:49:13 +01:00
eF = 0.5d0*(eHF(nO+1) + eHF(nO))
2023-11-27 09:56:32 +01:00
!-------------------------!
! Main loop over orbitals !
!-------------------------!
2021-10-18 21:11:01 +02:00
do p=nO-1,nO
2023-11-24 15:31:29 +01:00
2023-11-27 11:06:33 +01:00
H(:,:) = 0d0
2023-11-27 09:56:32 +01:00
if (TDA_W) then
! TDA for W
write(*,*) 'Tamm-Dancoff approximation actived!'
write(*,*)
!---------------------------!
! Compute GW supermatrix !
!---------------------------!
! !
! | F V2h1p V2p1h | !
! | | !
! H = | V2h1p C2h1p 0 | !
! | | !
! | V2p1h 0 C2p1h | !
! !
!---------------------------!
call wall_time(start_timing)
2023-04-20 17:00:34 +02:00
2023-11-27 09:56:32 +01:00
!---------!
! Block F !
!---------!
H(1,1) = eHF(p)
!-------------!
! Block V2h1p !
!-------------!
ija = 0
do i=nC+1,nO
do j=nC+1,nO
do a=nO+1,nBas-nR
ija = ija + 1
H(1 ,1+ija) = sqrt(2d0)*ERI(p,a,i,j)
H(1+ija,1 ) = sqrt(2d0)*ERI(p,a,i,j)
end do
2021-10-18 21:11:01 +02:00
end do
2023-11-24 15:31:29 +01:00
end do
2023-04-20 17:00:34 +02:00
2023-11-27 09:56:32 +01:00
!-------------!
! Block V2p1h !
!-------------!
iab = 0
do i=nC+1,nO
do a=nO+1,nBas-nR
do b=nO+1,nBas-nR
iab = iab + 1
H(1 ,1+n2h1p+iab) = sqrt(2d0)*ERI(p,i,b,a)
H(1+n2h1p+iab,1 ) = sqrt(2d0)*ERI(p,i,b,a)
end do
2023-04-20 17:00:34 +02:00
end do
2023-11-24 15:31:29 +01:00
end do
2023-11-27 09:56:32 +01:00
!-------------!
! Block C2h1p !
!-------------!
ija = 0
do i=nC+1,nO
do j=nC+1,nO
do a=nO+1,nBas-nR
ija = ija + 1
klc = 0
do k=nC+1,nO
do l=nC+1,nO
do c=nO+1,nBas-nR
klc = klc + 1
H(1+ija,1+klc) &
= ((eHF(i) + eHF(j) - eHF(a))*Kronecker_delta(j,l)*Kronecker_delta(a,c) &
- 2d0*ERI(j,c,a,l))*Kronecker_delta(i,k)
end do
2023-11-24 15:31:29 +01:00
end do
end do
2023-11-27 09:56:32 +01:00
2023-11-24 15:31:29 +01:00
end do
2021-10-18 21:11:01 +02:00
end do
2023-11-24 15:31:29 +01:00
end do
2023-11-27 09:56:32 +01:00
!-------------!
! Block C2p1h !
!-------------!
iab = 0
do i=nC+1,nO
do a=nO+1,nBas-nR
do b=nO+1,nBas-nR
iab = iab + 1
kcd = 0
do k=nC+1,nO
do c=nO+1,nBas-nR
do d=nO+1,nBas-nR
kcd = kcd + 1
H(1+n2h1p+iab,1+n2h1p+kcd) &
= ((eHF(a) + eHF(b) - eHF(i))*Kronecker_delta(i,k)*Kronecker_delta(a,c) &
+ 2d0*ERI(a,k,i,c))*Kronecker_delta(b,d)
end do
2023-11-24 15:31:29 +01:00
end do
end do
2023-11-27 09:56:32 +01:00
2023-11-24 15:31:29 +01:00
end do
end do
end do
2023-11-27 09:56:32 +01:00
call wall_time(end_timing)
timing = end_timing - start_timing
write(*,*)
2023-11-27 13:30:53 +01:00
write(*,'(A65,1X,F9.3,A8)') 'Total CPU time for construction of supermatrix = ',timing,' seconds'
2023-11-27 09:56:32 +01:00
write(*,*)
2023-11-24 15:31:29 +01:00
2023-11-27 09:56:32 +01:00
else
! RPA for W
write(*,*) 'Tamm-Dancoff approximation deactivated!'
write(*,*)
!---------------------------!
! Compute GW supermatrix !
!---------------------------!
! !
! | F W2h1p W2p1h | !
! | | !
! H = | W2h1p D2h1p 0 | !
! | | !
! | W2p1h 0 D2p1h | !
! !
!---------------------------!
! Spin manifold
ispin = 1
2023-11-27 14:29:42 +01:00
! Memory allocation
allocate(Om(nS),Aph(nS,nS),Bph(nS,nS),XpY(nS,nS),XmY(nS,nS),rho(nBas,nBas,nS))
2023-11-27 09:56:32 +01:00
!-------------------!
! Compute screening !
!-------------------!
call phLR_A(ispin,dRPA,nBas,nC,nO,nV,nR,nS,1d0,eHF,ERI,Aph)
if(.not.TDA_W) call phLR_B(ispin,dRPA,nBas,nC,nO,nV,nR,nS,1d0,ERI,Bph)
call phLR(TDA_W,nS,Aph,Bph,EcRPA,Om,XpY,XmY)
2023-11-29 17:16:33 +01:00
if(print_W) call print_excitation_energies('phRPA@RHF','singlet',nS,Om)
2023-11-27 09:56:32 +01:00
!--------------------------!
! Compute spectral weights !
!--------------------------!
call GW_excitation_density(nBas,nC,nO,nR,nS,ERI,XpY,rho)
call wall_time(start_timing)
2023-11-24 15:31:29 +01:00
2023-11-27 09:56:32 +01:00
!---------!
! Block F !
!---------!
H(1,1) = eHF(p)
!-------------!
! Block D2h1p !
!-------------!
ija = 0
do i=nC+1,nO
do ja=1,nS
ija = ija + 1
H(1+ija,1+ija) = eHF(i) - Om(ja)
end do
2023-11-24 15:31:29 +01:00
end do
2023-11-27 09:56:32 +01:00
!-------------!
! Block W2h1p !
!-------------!
ija = 0
do i=nC+1,nO
do ja=1,nS
ija = ija + 1
H(1 ,1+ija) = sqrt(2d0)*rho(p,i,ja)
H(1+ija,1 ) = sqrt(2d0)*rho(p,i,ja)
end do
2023-11-24 15:31:29 +01:00
end do
2023-11-27 09:56:32 +01:00
!-------------!
! Block D2p1h !
!-------------!
iab = 0
do ia=1,nS
do b=nO+1,nBas-nR
iab = iab + 1
H(1+n2h1p+iab,1+n2h1p+iab) = eHF(b) + Om(ia)
end do
2023-11-24 15:31:29 +01:00
end do
2023-11-27 09:56:32 +01:00
!-------------!
! Block W2p1h !
!-------------!
iab = 0
do ia=1,nS
do b=nO+1,nBas-nR
iab = iab + 1
H(1 ,1+n2h1p+iab) = sqrt(2d0)*rho(p,b,ia)
H(1+n2h1p+iab,1 ) = sqrt(2d0)*rho(p,b,ia)
end do
end do
2023-11-27 14:29:42 +01:00
! Memory deallocation
deallocate(Om,Aph,Bph,XpY,XmY,rho)
2023-11-27 09:56:32 +01:00
call wall_time(end_timing)
timing = end_timing - start_timing
write(*,*)
2023-11-27 13:30:53 +01:00
write(*,'(A65,1X,F9.3,A8)') 'Total CPU time for construction of supermatrix = ',timing,' seconds'
2023-11-27 09:56:32 +01:00
write(*,*)
end if
!-------------------------!
! Diagonalize supermatrix !
!-------------------------!
call wall_time(start_timing)
2023-11-24 15:31:29 +01:00
2023-11-27 09:56:32 +01:00
call diagonalize_matrix(nH,H,eGW)
call wall_time(end_timing)
2023-11-24 15:31:29 +01:00
2023-11-27 09:56:32 +01:00
timing = end_timing - start_timing
write(*,*)
2023-11-27 13:30:53 +01:00
write(*,'(A65,1X,F9.3,A8)') 'Total CPU time for diagonalization of supermatrix = ',timing,' seconds'
2023-11-27 09:56:32 +01:00
write(*,*)
2023-11-24 15:31:29 +01:00
2023-11-27 09:56:32 +01:00
!-----------------!
! Compute weights !
!-----------------!
do s=1,nH
Z(s) = H(1,s)**2
2023-11-24 15:31:29 +01:00
end do
2023-11-27 09:56:32 +01:00
!--------------!
! Dump results !
!--------------!
write(*,*)'-------------------------------------------'
write(*,'(1X,A32,I3,A8)')'| G0W0 energies (eV) for orbital',p,' |'
write(*,*)'-------------------------------------------'
write(*,'(1X,A1,1X,A3,1X,A1,1X,A15,1X,A1,1X,A15,1X,A1,1X,A15,1X)') &
'|','#','|','e_QP','|','Z','|'
write(*,*)'-------------------------------------------'
2023-04-20 17:00:34 +02:00
2021-10-18 21:11:01 +02:00
do s=1,nH
2023-11-27 10:49:13 +01:00
if(eGW(s) < eF .and. eGW(s) > eF - window) then
! if(Z(s) > cutoff1) then
2023-11-27 09:56:32 +01:00
write(*,'(1X,A1,1X,I3,1X,A1,1X,F15.6,1X,A1,1X,F15.6,1X,A1,1X)') &
'|',s,'|',eGW(s)*HaToeV,'|',Z(s),'|'
end if
2021-10-18 21:11:01 +02:00
end do
2023-11-27 09:56:32 +01:00
write(*,*)'-------------------------------------------'
write(*,*)
2023-11-24 15:31:29 +01:00
2023-11-27 09:56:32 +01:00
if(verbose) then
2023-11-27 14:29:42 +01:00
if(TDA_W) then
2023-11-27 09:56:32 +01:00
2023-11-27 14:29:42 +01:00
! TDA printing format
2023-11-27 09:56:32 +01:00
2023-11-27 14:29:42 +01:00
do s=1,nH
if(eGW(s) < eF .and. eGW(s) > eF - window) then
write(*,*)'-------------------------------------------------------------'
write(*,'(1X,A7,1X,I3,A6,I3,A1,1X,A7,F12.6,A13,F6.4,1X)') &
'Orbital',p,' and #',s,':','e_QP = ',eGW(s)*HaToeV,' eV and Z = ',Z(s)
write(*,*)'-------------------------------------------------------------'
write(*,'(1X,A20,1X,A20,1X,A15,1X)') &
' Configuration ',' Coefficient ',' Weight '
write(*,*)'-------------------------------------------------------------'
if(p <= nO) &
write(*,'(1X,A7,I3,A16,1X,F15.6,1X,F15.6)') &
' (',p,') ',H(1,s),H(1,s)**2
if(p > nO) &
write(*,'(1X,A16,I3,A7,1X,F15.6,1X,F15.6)') &
' (',p,') ',H(1,s),H(1,s)**2
2023-11-27 09:56:32 +01:00
ija = 0
do i=nC+1,nO
do j=nC+1,nO
do a=nO+1,nBas-nR
ija = ija + 1
2023-11-25 12:36:38 +01:00
2023-11-27 09:56:32 +01:00
if(abs(H(1+ija,s)) > cutoff2) &
write(*,'(1X,A3,I3,A1,I3,A6,I3,A7,1X,F15.6,1X,F15.6)') &
' (',i,',',j,') -> (',a,') ',H(1+ija,s),H(1+ija,s)**2
end do
end do
end do
iab = 0
do i=nC+1,nO
2023-11-25 12:36:38 +01:00
do a=nO+1,nBas-nR
2023-11-27 09:56:32 +01:00
do b=nO+1,nBas-nR
iab = iab + 1
2023-11-25 12:36:38 +01:00
2023-11-27 09:56:32 +01:00
if(abs(H(1+n2h1p+iab,s)) > cutoff2) &
write(*,'(1X,A7,I3,A6,I3,A1,I3,A3,1X,F15.6,1X,F15.6)') &
' (',i,') -> (',a,',',b,') ',H(1+n2h1p+iab,s),H(1+n2h1p+iab,s)**2
end do
2023-11-25 12:36:38 +01:00
end do
2022-01-17 08:02:11 +01:00
end do
2023-11-27 14:29:42 +01:00
write(*,*)'-------------------------------------------------------------'
write(*,*)
end if
end do
2023-11-27 09:56:32 +01:00
2023-11-27 14:29:42 +01:00
else
2023-11-27 09:56:32 +01:00
2023-11-27 14:29:42 +01:00
! non-TDA printing format
2023-11-27 09:56:32 +01:00
2023-11-27 14:29:42 +01:00
do s=1,nH
if(eGW(s) < eF .and. eGW(s) > eF - window) then
write(*,*)'-------------------------------------------------------------'
write(*,'(1X,A7,1X,I3,A6,I3,A1,1X,A7,F12.6,A13,F6.4,1X)') &
'Orbital',p,' and #',s,':','e_QP = ',eGW(s)*HaToeV,' eV and Z = ',Z(s)
write(*,*)'-------------------------------------------------------------'
write(*,'(1X,A20,1X,A20,1X,A15,1X)') &
' Conf. (p,ia) ',' Coefficient ',' Weight '
write(*,*)'-------------------------------------------------------------'
if(p <= nO) &
write(*,'(1X,A7,I3,A16,1X,F15.6,1X,F15.6)') &
' (',p,') ',H(1,s),H(1,s)**2
if(p > nO) &
write(*,'(1X,A16,I3,A7,1X,F15.6,1X,F15.6)') &
' (',p,') ',H(1,s),H(1,s)**2
2023-11-27 09:56:32 +01:00
ija = 0
do i=nC+1,nO
2023-11-27 10:49:13 +01:00
do ja=1,nS
2023-11-27 09:56:32 +01:00
ija = ija + 1
if(abs(H(1+ija,s)) > cutoff2) &
write(*,'(1X,A7,I3,A1,I3,A12,1X,F15.6,1X,F15.6)') &
' (',i,',',ja,') ',H(1+ija,s),H(1+ija,s)**2
end do
end do
iab = 0
do ia=1,nS
2023-11-25 12:36:38 +01:00
do b=nO+1,nBas-nR
iab = iab + 1
2023-11-27 09:56:32 +01:00
if(abs(H(1+n2h1p+iab,s)) > cutoff2) &
write(*,'(1X,A7,I3,A1,I3,A12,1X,F15.6,1X,F15.6)') &
' (',ia,',',b,') ',H(1+n2h1p+iab,s),H(1+n2h1p+iab,s)**2
2023-11-25 12:36:38 +01:00
end do
2022-01-17 08:02:11 +01:00
end do
2023-11-27 14:29:42 +01:00
write(*,*)'-------------------------------------------------------------'
write(*,*)
2023-11-27 09:56:32 +01:00
end if
2023-11-27 14:29:42 +01:00
end do
end if
2023-11-27 09:56:32 +01:00
end if
2023-11-25 12:36:38 +01:00
2023-11-27 09:56:32 +01:00
end do
2021-10-18 21:11:01 +02:00
2023-07-28 15:00:17 +02:00
end subroutine