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

257 lines
8.2 KiB
Fortran
Raw Normal View History

2023-11-20 21:34:08 +01:00
subroutine UG0W0(dotest,doACFDT,exchange_kernel,doXBS,dophBSE,TDA_W,TDA,dBSE,dTDA,spin_conserved,spin_flip, &
linearize,eta,regularize,nBas,nC,nO,nV,nR,nS,ENuc,EUHF,S,ERI_aaaa,ERI_aabb,ERI_bbbb, &
2023-07-30 22:01:44 +02:00
dipole_int_aa,dipole_int_bb,cHF,eHF)
2020-09-21 23:04:26 +02:00
! Perform unrestricted G0W0 calculation
implicit none
include 'parameters.h'
include 'quadrature.h'
! Input variables
2023-11-13 17:39:30 +01:00
logical,intent(in) :: dotest
2020-09-21 23:04:26 +02:00
logical,intent(in) :: doACFDT
logical,intent(in) :: exchange_kernel
logical,intent(in) :: doXBS
2023-11-20 21:34:08 +01:00
logical,intent(in) :: dophBSE
2020-09-21 23:04:26 +02:00
logical,intent(in) :: TDA_W
logical,intent(in) :: TDA
logical,intent(in) :: dBSE
logical,intent(in) :: dTDA
2020-09-22 22:13:51 +02:00
logical,intent(in) :: spin_conserved
logical,intent(in) :: spin_flip
2020-09-21 23:04:26 +02:00
logical,intent(in) :: linearize
double precision,intent(in) :: eta
2021-12-17 11:41:40 +01:00
logical,intent(in) :: regularize
2020-09-21 23:04:26 +02:00
integer,intent(in) :: nBas
integer,intent(in) :: nC(nspin)
integer,intent(in) :: nO(nspin)
integer,intent(in) :: nV(nspin)
integer,intent(in) :: nR(nspin)
integer,intent(in) :: nS(nspin)
double precision,intent(in) :: ENuc
double precision,intent(in) :: EUHF
2020-10-06 14:24:54 +02:00
double precision,intent(in) :: S(nBas,nBas)
2020-09-22 22:13:51 +02:00
double precision,intent(in) :: ERI_aaaa(nBas,nBas,nBas,nBas)
double precision,intent(in) :: ERI_aabb(nBas,nBas,nBas,nBas)
double precision,intent(in) :: ERI_bbbb(nBas,nBas,nBas,nBas)
2020-09-28 22:58:58 +02:00
double precision,intent(in) :: dipole_int_aa(nBas,nBas,ncart)
double precision,intent(in) :: dipole_int_bb(nBas,nBas,ncart)
2020-10-06 14:24:54 +02:00
double precision,intent(in) :: eHF(nBas,nspin)
double precision,intent(in) :: cHF(nBas,nBas,nspin)
2020-09-21 23:04:26 +02:00
! Local variables
logical :: print_W = .true.
2023-10-23 21:36:07 +02:00
logical :: dRPA
2020-09-23 14:23:40 +02:00
integer :: is
2020-09-21 23:04:26 +02:00
integer :: ispin
2020-09-24 14:39:37 +02:00
double precision :: EcRPA
2020-11-04 14:47:45 +01:00
double precision :: EcGM(nspin)
2020-09-22 23:08:47 +02:00
double precision :: EcBSE(nspin)
2020-09-21 23:04:26 +02:00
double precision,allocatable :: SigC(:,:)
double precision,allocatable :: Z(:,:)
integer :: nSa,nSb,nSt
double precision,allocatable :: Aph(:,:)
double precision,allocatable :: Bph(:,:)
2023-08-01 16:30:41 +02:00
double precision,allocatable :: Om(:)
double precision,allocatable :: XpY(:,:)
double precision,allocatable :: XmY(:,:)
double precision,allocatable :: rho(:,:,:,:)
2020-09-21 23:04:26 +02:00
double precision,allocatable :: eGWlin(:,:)
2023-10-24 17:03:41 +02:00
double precision,allocatable :: eGW(:,:)
2020-09-21 23:04:26 +02:00
! Hello world
write(*,*)
2023-11-09 17:22:01 +01:00
write(*,*)'*********************************'
write(*,*)'* Unrestricted G0W0 Calculation *'
write(*,*)'*********************************'
2020-09-21 23:04:26 +02:00
write(*,*)
! Initialization
2020-09-22 15:32:26 +02:00
EcRPA = 0d0
2023-10-23 21:36:07 +02:00
dRPA = .true.
2020-09-21 23:04:26 +02:00
! TDA for W
2020-09-22 22:13:51 +02:00
if(TDA_W) then
write(*,*) 'Tamm-Dancoff approximation for dynamic screening!'
write(*,*)
end if
2020-09-21 23:04:26 +02:00
! TDA
2020-09-22 22:13:51 +02:00
if(TDA) then
write(*,*) 'Tamm-Dancoff approximation activated!'
write(*,*)
end if
2020-09-21 23:04:26 +02:00
! Memory allocation
nSa = nS(1)
nSb = nS(2)
nSt = nSa + nSb
2020-09-21 23:04:26 +02:00
2023-10-24 17:03:41 +02:00
allocate(SigC(nBas,nspin),Z(nBas,nspin),eGWlin(nBas,nspin),eGW(nBas,nspin), &
Aph(nSt,nSt),Bph(nSt,nSt),Om(nSt),XpY(nSt,nSt),XmY(nSt,nSt),rho(nBas,nBas,nSt,nspin))
2020-09-21 23:04:26 +02:00
2020-09-22 15:32:26 +02:00
!-------------------!
! Compute screening !
!-------------------!
2020-09-21 23:04:26 +02:00
2023-10-23 21:36:07 +02:00
! Spin-conserving transitions
2020-09-21 23:04:26 +02:00
2020-09-22 15:32:26 +02:00
ispin = 1
2020-09-21 23:04:26 +02:00
call phULR_A(ispin,dRPA,nBas,nC,nO,nV,nR,nSa,nSb,nSt,1d0,eHF,ERI_aaaa,ERI_aabb,ERI_bbbb,Aph)
if(.not.TDA) call phULR_B(ispin,dRPA,nBas,nC,nO,nV,nR,nSa,nSb,nSt,1d0,ERI_aaaa,ERI_aabb,ERI_bbbb,Bph)
call phULR(TDA_W,nSa,nSb,nSt,Aph,Bph,EcRPA,Om,XpY,XmY)
2023-10-24 17:03:41 +02:00
2023-11-22 10:07:23 +01:00
if(print_W) call print_excitation_energies('phRPA@UHF','spin-conserved',nSt,Om)
2020-09-21 23:04:26 +02:00
2020-09-22 15:32:26 +02:00
!----------------------!
! Excitation densities !
!----------------------!
2020-09-21 23:04:26 +02:00
call UGW_excitation_density(nBas,nC,nO,nR,nSa,nSb,nSt,ERI_aaaa,ERI_aabb,ERI_bbbb,XpY,rho)
2020-09-21 23:04:26 +02:00
2022-01-06 22:10:53 +01:00
!------------------------------------------------!
! Compute self-energy and renormalization factor !
!------------------------------------------------!
2020-09-21 23:04:26 +02:00
2022-01-06 22:10:53 +01:00
if(regularize) then
2023-08-01 16:30:41 +02:00
do is=1,nspin
call GW_regularization(nBas,nC(is),nO(is),nV(is),nR(is),nSt,eHF(:,is),Om,rho(:,:,:,is))
2023-08-01 16:30:41 +02:00
end do
2022-01-06 22:10:53 +01:00
end if
2020-09-21 23:04:26 +02:00
call UGW_self_energy_diag(eta,nBas,nC,nO,nV,nR,nSt,eHF,Om,rho,SigC,Z,EcGM)
2023-08-01 16:30:41 +02:00
2020-09-22 15:32:26 +02:00
!-----------------------------------!
! Solve the quasi-particle equation !
!-----------------------------------!
eGWlin(:,:) = eHF(:,:) + Z(:,:)*SigC(:,:)
2020-09-21 23:04:26 +02:00
if(linearize) then
write(*,*) ' *** Quasiparticle energies obtained by linearization *** '
write(*,*)
eGW(:,:) = eGWlin(:,:)
else
! Find graphical solution of the QP equation
2023-11-27 10:17:14 +01:00
write(*,*) ' *** Quasiparticle energies obtained by root search *** '
2023-10-25 13:19:21 +02:00
write(*,*)
2021-02-15 17:27:06 +01:00
do is=1,nspin
2023-10-25 13:19:21 +02:00
write(*,*)'-----------------------------------------------------'
if(is==1) write(*,*)' Spin-up orbitals '
if(is==2) write(*,*)' Spin-down orbitals '
call UGW_QP_graph(eta,nBas,nC(is),nO(is),nV(is),nR(is),nSt,eHF(:,is), &
2023-10-25 13:19:21 +02:00
Om,rho(:,:,:,is),eGWlin(:,is),eHF(:,is),eGW(:,is),Z(:,is))
2021-02-15 17:27:06 +01:00
end do
2020-09-21 23:04:26 +02:00
end if
2020-09-24 14:39:37 +02:00
! Compute RPA correlation energy
2020-09-21 23:04:26 +02:00
call phULR_A(ispin,dRPA,nBas,nC,nO,nV,nR,nSa,nSb,nSt,1d0,eGW,ERI_aaaa,ERI_aabb,ERI_bbbb,Aph)
if(.not.TDA) call phULR_B(ispin,dRPA,nBas,nC,nO,nV,nR,nSa,nSb,nSt,1d0,ERI_aaaa,ERI_aabb,ERI_bbbb,Bph)
call phULR(TDA_W,nSa,nSb,nSt,Aph,Bph,EcRPA,Om,XpY,XmY)
2020-09-21 23:04:26 +02:00
2020-09-24 14:39:37 +02:00
! Dump results
2020-11-04 14:47:45 +01:00
call print_UG0W0(nBas,nO,eHF,ENuc,EUHF,SigC,Z,eGW,EcRPA,EcGM)
2020-09-21 23:04:26 +02:00
2020-09-22 23:08:47 +02:00
! Free memory
2023-08-01 16:30:41 +02:00
deallocate(Om,XpY,XmY,rho)
2020-09-22 23:08:47 +02:00
2020-09-21 23:04:26 +02:00
! Perform BSE calculation
2023-11-20 21:34:08 +01:00
if(dophBSE) then
2020-09-21 23:04:26 +02:00
2023-07-27 21:59:18 +02:00
call UGW_phBSE(TDA_W,TDA,dBSE,dTDA,spin_conserved,spin_flip,eta,nBas,nC,nO,nV,nR,nS,S, &
ERI_aaaa,ERI_aabb,ERI_bbbb,dipole_int_aa,dipole_int_bb,cHF,eHF,eGW,EcBSE)
2020-09-21 23:04:26 +02:00
2020-10-07 18:26:24 +02:00
if(exchange_kernel) then
EcBSE(1) = 0.5d0*EcBSE(1)
2021-04-02 09:53:23 +02:00
EcBSE(2) = 0.5d0*EcBSE(2)
2020-10-07 18:26:24 +02:00
2021-04-02 09:53:23 +02:00
else
EcBSE(2) = 0.0d0
2020-10-07 18:26:24 +02:00
end if
2020-09-21 23:04:26 +02:00
2020-10-07 18:26:24 +02:00
write(*,*)
write(*,*)'-------------------------------------------------------------------------------'
2023-11-27 23:25:10 +01:00
write(*,'(2X,A50,F20.10,A3)') 'Tr@BSE@G0W0@UHF correlation energy (spin-conserved) = ',EcBSE(1),' au'
write(*,'(2X,A50,F20.10,A3)') 'Tr@BSE@G0W0@UHF correlation energy (spin-flip) = ',EcBSE(2),' au'
write(*,'(2X,A50,F20.10,A3)') 'Tr@BSE@G0W0@UHF correlation energy = ',sum(EcBSE),' au'
write(*,'(2X,A50,F20.10,A3)') 'Tr@BSE@G0W0@UHF total energy = ',ENuc + EUHF + sum(EcBSE),' au'
2020-10-07 18:26:24 +02:00
write(*,*)'-------------------------------------------------------------------------------'
write(*,*)
2020-09-21 23:04:26 +02:00
2020-10-07 18:26:24 +02:00
! Compute the BSE correlation energy via the adiabatic connection
2020-09-21 23:04:26 +02:00
2020-10-07 18:26:24 +02:00
if(doACFDT) then
2020-09-21 23:04:26 +02:00
2021-09-09 15:10:19 +02:00
write(*,*) '------------------------------------------------------------'
write(*,*) 'Adiabatic connection version of BSE@UG0W0 correlation energy'
write(*,*) '------------------------------------------------------------'
2020-10-07 18:26:24 +02:00
write(*,*)
2020-09-21 23:04:26 +02:00
2020-10-07 18:26:24 +02:00
if(doXBS) then
2020-09-21 23:04:26 +02:00
2020-10-07 18:26:24 +02:00
write(*,*) '*** scaled screening version (XBS) ***'
write(*,*)
2020-09-21 23:04:26 +02:00
2020-10-07 18:26:24 +02:00
end if
2020-09-21 23:04:26 +02:00
2023-11-20 21:34:08 +01:00
call UGW_phACFDT(exchange_kernel,doXBS,.true.,TDA_W,TDA,dophBSE,spin_conserved,spin_flip,eta, &
nBas,nC,nO,nV,nR,nS,ERI_aaaa,ERI_aabb,ERI_bbbb,eHF,eGW,EcBSE)
2020-10-07 18:26:24 +02:00
write(*,*)
write(*,*)'-------------------------------------------------------------------------------'
2023-11-27 23:25:10 +01:00
write(*,'(2X,A50,F20.10,A3)') 'AC@BSE@G0W0@UHF correlation energy (spin-conserved) = ',EcBSE(1),' au'
write(*,'(2X,A50,F20.10,A3)') 'AC@BSE@G0W0@UHF correlation energy (spin-flip) = ',EcBSE(2),' au'
write(*,'(2X,A50,F20.10,A3)') 'AC@BSE@G0W0@UHF correlation energy = ',sum(EcBSE),' au'
write(*,'(2X,A50,F20.10,A3)') 'AC@BSE@G0W0@UHF total energy = ',ENuc + EUHF + sum(EcBSE),' au'
2020-10-07 18:26:24 +02:00
write(*,*)'-------------------------------------------------------------------------------'
write(*,*)
end if
2020-09-21 23:04:26 +02:00
2020-09-22 22:13:51 +02:00
end if
2020-09-21 23:04:26 +02:00
2023-11-13 17:39:30 +01:00
! Testing zone
if(dotest) then
2023-11-14 14:31:27 +01:00
call dump_test_value('U','G0W0 correlation energy',EcRPA)
call dump_test_value('U','G0W0 HOMOa energy',eGW(nO(1),1))
call dump_test_value('U','G0W0 LUMOa energy',eGW(nO(1)+1,1))
call dump_test_value('U','G0W0 HOMOa energy',eGW(nO(2),2))
call dump_test_value('U','G0W0 LUMOa energy',eGW(nO(2)+1,2))
2023-11-13 17:39:30 +01:00
end if
2023-07-28 15:00:17 +02:00
end subroutine