10
1
mirror of https://github.com/pfloos/quack synced 2024-12-23 12:55:25 +01:00
QuAcK/src/GW/evUGW.f90

291 lines
9.5 KiB
Fortran
Raw Normal View History

2023-07-21 22:54:53 +02:00
subroutine evUGW(maxSCF,thresh,max_diis,doACFDT,exchange_kernel,doXBS,BSE,TDA_W,TDA, &
2023-07-21 10:21:54 +02:00
dBSE,dTDA,spin_conserved,spin_flip,eta,regularize,nBas,nC,nO,nV,nR,nS,ENuc, &
2021-02-15 17:27:06 +01:00
EUHF,S,ERI_AO,ERI_aaaa,ERI_aabb,ERI_bbbb,dipole_int_aa,dipole_int_bb,PHF,cHF,eHF,Vxc,eG0W0)
2020-09-24 22:50:56 +02:00
! Perform self-consistent eigenvalue-only GW calculation
implicit none
include 'parameters.h'
! Input variables
integer,intent(in) :: maxSCF
integer,intent(in) :: max_diis
double precision,intent(in) :: thresh
double precision,intent(in) :: ENuc
2020-10-07 22:51:30 +02:00
double precision,intent(in) :: EUHF
2020-09-24 22:50:56 +02:00
logical,intent(in) :: doACFDT
logical,intent(in) :: exchange_kernel
logical,intent(in) :: doXBS
logical,intent(in) :: BSE
logical,intent(in) :: TDA_W
logical,intent(in) :: TDA
logical,intent(in) :: dBSE
logical,intent(in) :: dTDA
logical,intent(in) :: spin_conserved
logical,intent(in) :: spin_flip
double precision,intent(in) :: eta
2021-12-17 11:41:40 +01:00
logical,intent(in) :: regularize
2020-09-24 22:50:56 +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)
2021-02-15 17:27:06 +01:00
double precision,intent(in) :: PHF(nBas,nBas,nspin)
2020-09-24 22:50:56 +02:00
double precision,intent(in) :: eHF(nBas,nspin)
double precision,intent(in) :: cHF(nBas,nBas,nspin)
2021-02-15 17:27:06 +01:00
double precision,intent(in) :: Vxc(nBas,nspin)
2020-09-24 22:50:56 +02:00
double precision,intent(in) :: eG0W0(nBas,nspin)
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_AO(nBas,nBas,nBas,nBas)
2020-09-24 22:50:56 +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-09-24 22:50:56 +02:00
! Local variables
logical :: linear_mixing
integer :: is
integer :: ispin
integer :: nSCF
integer :: n_diis
double precision :: rcond(nspin)
double precision :: Conv
double precision :: EcRPA
2020-11-04 14:47:45 +01:00
double precision :: EcGM(nspin)
2020-09-24 22:50:56 +02:00
double precision :: EcBSE(nspin)
double precision :: EcAC(nspin)
double precision :: alpha
double precision,allocatable :: error_diis(:,:,:)
double precision,allocatable :: e_diis(:,:,:)
double precision,allocatable :: eGW(:,:)
double precision,allocatable :: eOld(:,:)
double precision,allocatable :: Z(:,:)
integer :: nS_aa,nS_bb,nS_sc
2021-02-15 17:27:06 +01:00
double precision,allocatable :: SigX(:,:)
2020-09-24 22:50:56 +02:00
double precision,allocatable :: SigC(:,:)
double precision,allocatable :: OmRPA(:)
double precision,allocatable :: XpY_RPA(:,:)
double precision,allocatable :: XmY_RPA(:,:)
double precision,allocatable :: rho_RPA(:,:,:,:)
! Hello world
write(*,*)
write(*,*)'************************************************'
write(*,*)'| Self-consistent evGW calculation |'
write(*,*)'************************************************'
write(*,*)
! TDA for W
if(TDA_W) then
write(*,*) 'Tamm-Dancoff approximation for dynamic screening!'
write(*,*)
end if
! TDA
if(TDA) then
write(*,*) 'Tamm-Dancoff approximation activated!'
write(*,*)
end if
! Linear mixing
linear_mixing = .false.
alpha = 0.2d0
! Memory allocation
nS_aa = nS(1)
nS_bb = nS(2)
nS_sc = nS_aa + nS_bb
2021-02-15 17:27:06 +01:00
allocate(eGW(nBas,nspin),eOld(nBas,nspin),Z(nBas,nspin),SigX(nBas,nspin),SigC(nBas,nspin), &
OmRPA(nS_sc),XpY_RPA(nS_sc,nS_sc),XmY_RPA(nS_sc,nS_sc),rho_RPA(nBas,nBas,nS_sc,nspin), &
error_diis(nBas,max_diis,nspin),e_diis(nBas,max_diis,nspin))
! Compute the exchange part of the self-energy
do is=1,nspin
call self_energy_exchange_diag(nBas,cHF(:,:,is),PHF(:,:,is),ERI_AO,SigX(:,is))
end do
2020-09-24 22:50:56 +02:00
! Initialization
nSCF = 0
ispin = 1
n_diis = 0
Conv = 1d0
e_diis(:,:,:) = 0d0
error_diis(:,:,:) = 0d0
eGW(:,:) = eG0W0(:,:)
eOld(:,:) = eGW(:,:)
Z(:,:) = 1d0
2022-01-06 22:36:30 +01:00
rcond(:) = 0d0
2020-09-24 22:50:56 +02:00
!------------------------------------------------------------------------
! Main loop
!------------------------------------------------------------------------
do while(Conv > thresh .and. nSCF <= maxSCF)
! Compute screening
2023-07-21 13:04:29 +02:00
call phULR(ispin,.true.,TDA_W,.false.,eta,nBas,nC,nO,nV,nR,nS_aa,nS_bb,nS_sc,nS_sc,1d0, &
eGW,ERI_aaaa,ERI_aabb,ERI_bbbb,OmRPA,rho_RPA,EcRPA,OmRPA,XpY_RPA,XmY_RPA)
2020-09-24 22:50:56 +02:00
!----------------------!
! Excitation densities !
!----------------------!
call unrestricted_excitation_density(nBas,nC,nO,nR,nS_aa,nS_bb,nS_sc,ERI_aaaa,ERI_aabb,ERI_bbbb,XpY_RPA,rho_RPA)
!------------------------------------------------!
! Compute self-energy and renormalization factor !
!------------------------------------------------!
2022-11-30 16:41:19 +01:00
if(regularize) then
2020-09-24 22:50:56 +02:00
2022-11-30 16:41:19 +01:00
call unrestricted_regularized_self_energy_correlation_diag(eta,nBas,nC,nO,nV,nR,nS_sc,eGW,OmRPA,rho_RPA,SigC,EcGM)
call unrestricted_regularized_renormalization_factor(eta,nBas,nC,nO,nV,nR,nS_sc,eGW,OmRPA,rho_RPA,Z)
2022-11-30 16:41:19 +01:00
else
2022-11-30 16:41:19 +01:00
call unrestricted_self_energy_correlation_diag(eta,nBas,nC,nO,nV,nR,nS_sc,eGW,OmRPA,rho_RPA,SigC,EcGM)
call unrestricted_renormalization_factor(eta,nBas,nC,nO,nV,nR,nS_sc,eGW,OmRPA,rho_RPA,Z)
end if
2020-09-24 22:50:56 +02:00
!-----------------------------------!
! Solve the quasi-particle equation !
!-----------------------------------!
2021-02-15 17:27:06 +01:00
eGW(:,:) = eHF(:,:) + SigX(:,:) + SigC(:,:) - Vxc(:,:)
2020-09-24 22:50:56 +02:00
! Convergence criteria
Conv = maxval(abs(eGW(:,:) - eOld(:,:)))
! Print results
2020-11-04 14:47:45 +01:00
call print_evUGW(nBas,nO,nSCF,Conv,eHF,ENuc,EUHF,SigC,Z,eGW,EcRPA,EcGM)
2020-09-24 22:50:56 +02:00
! Linear mixing or DIIS extrapolation
if(linear_mixing) then
eGW(:,:) = alpha*eGW(:,:) + (1d0 - alpha)*eOld(:,:)
else
n_diis = min(n_diis+1,max_diis)
do is=1,nspin
2020-10-21 12:09:18 +02:00
call DIIS_extrapolation(rcond(ispin),nBas,nBas,n_diis,error_diis(:,1:n_diis,is), &
e_diis(:,1:n_diis,is),eGW(:,is)-eOld(:,is),eGW(:,is))
2020-09-24 22:50:56 +02:00
end do
! Reset DIIS if required
if(minval(rcond(:)) < 1d-15) n_diis = 0
end if
2020-09-24 22:50:56 +02:00
! Save quasiparticles energy for next cycle
eOld(:,:) = eGW(:,:)
! Increment
nSCF = nSCF + 1
enddo
!------------------------------------------------------------------------
! End main loop
!------------------------------------------------------------------------
! Did it actually converge?
if(nSCF == maxSCF+1) then
write(*,*)
write(*,*)'!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!'
write(*,*)' Convergence failed '
write(*,*)'!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!'
write(*,*)
stop
end if
2020-09-24 22:50:56 +02:00
! Deallocate memory
deallocate(eOld,Z,SigC,OmRPA,XpY_RPA,XmY_RPA,rho_RPA,error_diis,e_diis)
! Perform BSE calculation
if(BSE) then
2023-07-21 10:21:54 +02:00
call unrestricted_Bethe_Salpeter(TDA_W,TDA,dBSE,dTDA,spin_conserved,spin_flip,eta,nBas,nC,nO,nV,nR,nS, &
2020-10-06 14:24:54 +02:00
S,ERI_aaaa,ERI_aabb,ERI_bbbb,dipole_int_aa,dipole_int_bb,cHF,eGW,eGW,EcBSE)
2020-09-24 22:50:56 +02:00
2020-10-07 22:51:30 +02:00
if(exchange_kernel) then
2020-09-24 22:50:56 +02:00
2020-10-07 22:51:30 +02:00
EcBSE(1) = 0.5d0*EcBSE(1)
2021-04-02 09:53:23 +02:00
EcBSE(2) = 0.5d0*EcBSE(2)
else
EcBSE(2) = 0.0d0
2020-09-24 22:50:56 +02:00
2020-10-07 22:51:30 +02:00
end if
2020-09-24 22:50:56 +02:00
2020-10-07 22:51:30 +02:00
write(*,*)
write(*,*)'-------------------------------------------------------------------------------'
write(*,'(2X,A50,F20.10)') 'Tr@BSE@evUGW correlation energy (spin-conserved) =',EcBSE(1)
write(*,'(2X,A50,F20.10)') 'Tr@BSE@evUGW correlation energy (spin-flip) =',EcBSE(2)
write(*,'(2X,A50,F20.10)') 'Tr@BSE@evUGW correlation energy =',EcBSE(1) + EcBSE(2)
write(*,'(2X,A50,F20.10)') 'Tr@BSE@evUGW total energy =',ENuc + EUHF + EcBSE(1) + EcBSE(2)
write(*,*)'-------------------------------------------------------------------------------'
write(*,*)
2020-09-24 22:50:56 +02:00
! Compute the BSE correlation energy via the adiabatic connection
2020-10-07 22:51:30 +02:00
if(doACFDT) then
2020-09-24 22:50:56 +02:00
2020-10-07 22:51:30 +02:00
write(*,*) '--------------------------------------------------------------'
write(*,*) ' Adiabatic connection version of BSE@evUGW correlation energy '
write(*,*) '--------------------------------------------------------------'
write(*,*)
2020-09-24 22:50:56 +02:00
2020-10-07 22:51:30 +02:00
if(doXBS) then
2020-09-24 22:50:56 +02:00
2020-10-07 22:51:30 +02:00
write(*,*) '*** scaled screening version (XBS) ***'
write(*,*)
2020-09-24 22:50:56 +02:00
2020-10-07 22:51:30 +02:00
end if
2020-09-24 22:50:56 +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,eGW,eGW,EcAC)
2020-09-24 22:50:56 +02:00
2020-10-07 22:51:30 +02:00
write(*,*)
write(*,*)'-------------------------------------------------------------------------------'
write(*,'(2X,A50,F20.10)') 'AC@BSE@evUGW correlation energy (spin-conserved) =',EcAC(1)
write(*,'(2X,A50,F20.10)') 'AC@BSE@evUGW correlation energy (spin-flip) =',EcAC(2)
write(*,'(2X,A50,F20.10)') 'AC@BSE@evUGW correlation energy =',EcAC(1) + EcAC(2)
write(*,'(2X,A50,F20.10)') 'AC@BSE@evUGW total energy =',ENuc + EUHF + EcAC(1) + EcAC(2)
write(*,*)'-------------------------------------------------------------------------------'
write(*,*)
2020-09-24 22:50:56 +02:00
2020-10-07 22:51:30 +02:00
end if
2020-09-24 22:50:56 +02:00
end if
2020-09-24 22:50:56 +02:00
end subroutine evUGW