4
1
mirror of https://github.com/pfloos/quack synced 2024-11-19 04:22:39 +01:00
quack/src/MBPT/UG0W0.f90

247 lines
8.1 KiB
Fortran
Raw Normal View History

2020-09-28 22:58:58 +02:00
subroutine UG0W0(doACFDT,exchange_kernel,doXBS,COHSEX,BSE,TDA_W,TDA,dBSE,dTDA,evDyn,spin_conserved,spin_flip, &
2021-02-15 17:27:06 +01:00
linearize,eta,nBas,nC,nO,nV,nR,nS,ENuc,EUHF,S,ERI,ERI_aaaa,ERI_aabb,ERI_bbbb, &
dipole_int_aa,dipole_int_bb,PHF,cHF,eHF,Vxc,eGW)
2020-09-21 23:04:26 +02:00
! Perform unrestricted G0W0 calculation
implicit none
include 'parameters.h'
include 'quadrature.h'
! Input variables
logical,intent(in) :: doACFDT
logical,intent(in) :: exchange_kernel
logical,intent(in) :: doXBS
logical,intent(in) :: COHSEX
logical,intent(in) :: BSE
logical,intent(in) :: TDA_W
logical,intent(in) :: TDA
logical,intent(in) :: dBSE
logical,intent(in) :: dTDA
logical,intent(in) :: evDyn
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
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)
2021-02-15 17:27:06 +01:00
double precision,intent(in) :: ERI(nBas,nBas,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)
double precision,intent(in) :: PHF(nBas,nBas,nspin)
2021-02-15 17:27:06 +01:00
double precision,intent(in) :: Vxc(nBas,nspin)
2020-09-21 23:04:26 +02:00
! Local variables
logical :: print_W = .true.
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 :: EcAC(nspin)
2021-02-15 17:27:06 +01:00
double precision,allocatable :: SigX(:,:)
2020-09-21 23:04:26 +02:00
double precision,allocatable :: SigC(:,:)
double precision,allocatable :: Z(:,:)
2020-09-22 23:08:47 +02:00
integer :: nS_aa,nS_bb,nS_sc
2020-09-24 14:39:37 +02:00
double precision,allocatable :: OmRPA(:)
double precision,allocatable :: XpY_RPA(:,:)
double precision,allocatable :: XmY_RPA(:,:)
double precision,allocatable :: rho_RPA(:,:,:,:)
2020-09-21 23:04:26 +02:00
double precision,allocatable :: eGWlin(:,:)
! Output variables
double precision :: eGW(nBas,nspin)
! Hello world
write(*,*)
write(*,*)'************************************************'
write(*,*)'| One-shot G0W0 calculation |'
write(*,*)'| *** Unrestricted version *** |'
write(*,*)'************************************************'
write(*,*)
! Initialization
2020-09-22 15:32:26 +02:00
EcRPA = 0d0
2020-09-21 23:04:26 +02:00
! COHSEX approximation
2020-09-22 22:13:51 +02:00
if(COHSEX) then
write(*,*) 'COHSEX approximation activated!'
write(*,*)
end if
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
2020-09-22 23:08:47 +02:00
nS_aa = nS(1)
nS_bb = nS(2)
nS_sc = nS_aa + nS_bb
2020-09-21 23:04:26 +02:00
2021-02-15 17:27:06 +01:00
allocate(SigX(nBas,nspin),SigC(nBas,nspin),Z(nBas,nspin),eGWlin(nBas,nspin), &
OmRPA(nS_sc),XpY_RPA(nS_sc,nS_sc),XmY_RPA(nS_sc,nS_sc),rho_RPA(nBas,nBas,nS_sc,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
2020-09-22 15:32:26 +02:00
! Spin-conserving transition
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
2020-09-23 22:48:54 +02:00
call unrestricted_linear_response(ispin,.true.,TDA_W,.false.,eta,nBas,nC,nO,nV,nR,nS_aa,nS_bb,nS_sc,nS_sc,1d0, &
2020-10-05 23:00:56 +02:00
eHF,ERI_aaaa,ERI_aabb,ERI_bbbb,OmRPA,rho_RPA,EcRPA,OmRPA,XpY_RPA,XmY_RPA)
2020-09-21 23:04:26 +02:00
if(print_W) call print_excitation('RPA@UHF ',5,nS_sc,OmRPA)
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
2020-09-24 14:39:37 +02:00
call unrestricted_excitation_density(nBas,nC,nO,nR,nS_aa,nS_bb,nS_sc,ERI_aaaa,ERI_aabb,ERI_bbbb,XpY_RPA,rho_RPA)
2020-09-21 23:04:26 +02:00
2020-09-22 15:32:26 +02:00
!---------------------!
! Compute self-energy !
!---------------------!
2020-09-21 23:04:26 +02:00
2021-02-15 17:27:06 +01:00
do is=1,nspin
call self_energy_exchange_diag(nBas,cHF(:,:,is),PHF(:,:,is),ERI,SigX(:,is))
end do
2020-11-04 14:47:45 +01:00
call unrestricted_self_energy_correlation_diag(eta,nBas,nC,nO,nV,nR,nS_sc,eHF,OmRPA,rho_RPA,SigC,EcGM)
2020-09-21 23:04:26 +02:00
2020-09-22 15:32:26 +02:00
!--------------------------------!
! Compute renormalization factor !
!--------------------------------!
2020-09-21 23:04:26 +02:00
2020-09-24 14:39:37 +02:00
call unrestricted_renormalization_factor(eta,nBas,nC,nO,nV,nR,nS_sc,eHF,OmRPA,rho_RPA,Z)
2020-09-21 23:04:26 +02:00
2020-09-22 15:32:26 +02:00
!-----------------------------------!
! Solve the quasi-particle equation !
!-----------------------------------!
2021-02-15 17:27:06 +01:00
eGWlin(:,:) = eHF(:,:) + Z(:,:)*(SigX(:,:) + SigC(:,:) - Vxc(:,:))
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
2021-02-15 17:27:06 +01:00
do is=1,nspin
2021-02-15 22:15:30 +01:00
call unrestricted_QP_graph(nBas,nC(is),nO(is),nV(is),nR(is),nS_sc,eta,eHF(:,is),SigX(:,is),Vxc(:,is), &
OmRPA,rho_RPA(:,:,:,is),eGWlin(:,is),eGW(:,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
2020-09-23 22:48:54 +02:00
call unrestricted_linear_response(ispin,.true.,TDA_W,.false.,eta,nBas,nC,nO,nV,nR,nS_aa,nS_bb,nS_sc,nS_sc,1d0, &
2020-10-05 23:00:56 +02:00
eGW,ERI_aaaa,ERI_aabb,ERI_bbbb,OmRPA,rho_RPA,EcRPA,OmRPA,XpY_RPA,XmY_RPA)
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
2020-09-24 14:39:37 +02:00
deallocate(OmRPA,XpY_RPA,XmY_RPA,rho_RPA)
2020-09-22 23:08:47 +02:00
2020-09-21 23:04:26 +02:00
! Perform BSE calculation
2020-09-22 22:13:51 +02:00
if(BSE) then
2020-09-21 23:04:26 +02:00
2020-10-06 14:24:54 +02:00
call unrestricted_Bethe_Salpeter(TDA_W,TDA,dBSE,dTDA,evDyn,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(*,*)'-------------------------------------------------------------------------------'
2020-10-07 22:51:30 +02:00
write(*,'(2X,A50,F20.10)') 'Tr@BSE@UG0W0 correlation energy (spin-conserved) =',EcBSE(1)
write(*,'(2X,A50,F20.10)') 'Tr@BSE@UG0W0 correlation energy (spin-flip) =',EcBSE(2)
write(*,'(2X,A50,F20.10)') 'Tr@BSE@UG0W0 correlation energy =',EcBSE(1) + EcBSE(2)
write(*,'(2X,A50,F20.10)') 'Tr@BSE@UG0W0 total energy =',ENuc + EUHF + EcBSE(1) + EcBSE(2)
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
2020-10-07 22:51:30 +02:00
call unrestricted_ACFDT(exchange_kernel,doXBS,.true.,TDA_W,TDA,BSE,spin_conserved,spin_flip,eta, &
nBas,nC,nO,nV,nR,nS,ERI_aaaa,ERI_aabb,ERI_bbbb,eHF,eGW,EcAC)
2020-10-07 18:26:24 +02:00
write(*,*)
write(*,*)'-------------------------------------------------------------------------------'
2020-10-07 22:51:30 +02:00
write(*,'(2X,A50,F20.10)') 'AC@BSE@UG0W0 correlation energy (spin-conserved) =',EcAC(1)
write(*,'(2X,A50,F20.10)') 'AC@BSE@UG0W0 correlation energy (spin-flip) =',EcAC(2)
write(*,'(2X,A50,F20.10)') 'AC@BSE@UG0W0 correlation energy =',EcAC(1) + EcAC(2)
write(*,'(2X,A50,F20.10)') 'AC@BSE@UG0W0 total energy =',ENuc + EUHF + EcAC(1) + EcAC(2)
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
end subroutine UG0W0