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

324 lines
10 KiB
Fortran
Raw Normal View History

2023-11-13 17:39:30 +01:00
subroutine evUGW(dotest,maxSCF,thresh,max_diis,doACFDT,exchange_kernel,doXBS,BSE,TDA_W,TDA,dBSE,dTDA, &
2024-09-11 17:21:43 +02:00
spin_conserved,spin_flip,linearize,eta,doSRG,nBas,nC,nO,nV,nR,nS,ENuc, &
2023-07-30 22:01:44 +02:00
EUHF,S,ERI_aaaa,ERI_aabb,ERI_bbbb,dipole_int_aa,dipole_int_bb,cHF,eHF)
2020-09-24 22:50:56 +02:00
! Perform self-consistent eigenvalue-only GW calculation
implicit none
include 'parameters.h'
! Input variables
2023-11-13 17:39:30 +01:00
logical,intent(in) :: dotest
2020-09-24 22:50:56 +02:00
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
2023-07-30 22:01:44 +02:00
logical,intent(in) :: linearize
2020-09-24 22:50:56 +02:00
double precision,intent(in) :: eta
2024-09-11 17:21:43 +02:00
logical,intent(in) :: doSRG
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)
double precision,intent(in) :: eHF(nBas,nspin)
double precision,intent(in) :: cHF(nBas,nBas,nspin)
2020-10-06 14:24:54 +02:00
double precision,intent(in) :: S(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 :: dRPA
2020-09-24 22:50:56 +02:00
logical :: linear_mixing
integer :: is
integer :: ispin
integer :: nSCF
integer :: n_diis
double precision :: rcond(nspin)
double precision :: Conv
2023-11-27 23:25:10 +01:00
double precision :: EcRPA(nspin)
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 :: alpha
double precision,allocatable :: error_diis(:,:,:)
double precision,allocatable :: e_diis(:,:,:)
double precision,allocatable :: eGW(:,:)
double precision,allocatable :: eOld(:,:)
double precision,allocatable :: Z(:,:)
integer :: nSa,nSb,nSt
2020-09-24 22:50:56 +02:00
double precision,allocatable :: SigC(:,:)
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-24 22:50:56 +02:00
! Hello world
write(*,*)
2023-11-09 17:22:01 +01:00
write(*,*)'*********************************'
write(*,*)'| Unrestricted evGW Calculation *'
write(*,*)'*********************************'
2020-09-24 22:50:56 +02:00
write(*,*)
! TDA for W
if(TDA_W) then
write(*,*) 'Tamm-Dancoff approximation for dynamic screening!'
write(*,*)
end if
2024-09-11 17:58:36 +02:00
! SRG regularization
2020-09-24 22:50:56 +02:00
2024-09-11 17:58:36 +02:00
if(doSRG) then
write(*,*) '*** SRG regularized qsGW scheme ***'
2020-09-24 22:50:56 +02:00
write(*,*)
2024-09-11 17:58:36 +02:00
2020-09-24 22:50:56 +02:00
end if
! Initialization
2023-11-27 23:25:10 +01:00
EcRPA(:) = 0d0
dRPA = .true.
2020-09-24 22:50:56 +02:00
! Linear mixing
linear_mixing = .false.
alpha = 0.2d0
! Memory allocation
nSa = nS(1)
nSb = nS(2)
nSt = nSa + nSb
2020-09-24 22:50:56 +02:00
allocate(eGW(nBas,nspin),eOld(nBas,nspin),Z(nBas,nspin),SigC(nBas,nspin), &
Aph(nSt,nSt),Bph(nSt,nSt),Om(nSt),XpY(nSt,nSt),XmY(nSt,nSt), &
rho(nBas,nBas,nSt,nspin),error_diis(nBas,max_diis,nspin),e_diis(nBas,max_diis,nspin))
2021-02-15 17:27:06 +01:00
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(:,:) = eHF(:,:)
2020-09-24 22:50:56 +02:00
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
2020-09-24 22:50:56 +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)
2023-11-27 23:25:10 +01:00
call phULR(TDA_W,nSa,nSb,nSt,Aph,Bph,EcRPA(ispin),Om,XpY,XmY)
2020-09-24 22:50:56 +02:00
!----------------------!
! Excitation densities !
!----------------------!
call UGW_excitation_density(nBas,nC,nO,nR,nSa,nSb,nSt,ERI_aaaa,ERI_aabb,ERI_bbbb,XpY,rho)
2020-09-24 22:50:56 +02:00
!------------------------------------------------!
! Compute self-energy and renormalization factor !
!------------------------------------------------!
2024-09-11 17:21:43 +02:00
if(doSRG) then
call UGW_SRG_self_energy_diag(nBas,nC,nO,nV,nR,nSt,eGW,Om,rho,SigC,Z,EcGM)
else
call UGW_self_energy_diag(eta,nBas,nC,nO,nV,nR,nSt,eGW,Om,rho,SigC,Z,EcGM)
end if
2020-09-24 22:50:56 +02:00
2023-08-01 16:30:41 +02:00
2020-09-24 22:50:56 +02:00
!-----------------------------------!
! Solve the quasi-particle equation !
!-----------------------------------!
2023-07-30 22:01:44 +02:00
if(linearize) then
write(*,*) ' *** Quasiparticle energies obtained by linearization *** '
write(*,*)
eGW(:,:) = eHF(:,:) + SigC(:,:)
else
2023-11-27 10:17:14 +01:00
write(*,*) ' *** Quasiparticle energies obtained by root search *** '
2023-07-30 22:01:44 +02:00
write(*,*)
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),eOld(:,is),eOld(:,is),eGW(:,is),Z(:,is))
2023-07-30 22:01:44 +02:00
end do
end if
2020-09-24 22:50:56 +02:00
! Convergence criteria
Conv = maxval(abs(eGW(:,:) - eOld(:,:)))
! Print results
2023-11-27 23:25:10 +01:00
call print_evUGW(nBas,nO,nSCF,Conv,eHF,ENuc,EUHF,SigC,Z,eGW,EcRPA(ispin),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
2023-12-03 18:47:30 +01:00
end do
2020-09-24 22:50:56 +02:00
!------------------------------------------------------------------------
! 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
2023-08-01 16:30:41 +02:00
deallocate(eOld,Z,SigC,Om,XpY,XmY,rho,error_diis,e_diis)
2020-09-24 22:50:56 +02:00
! Perform BSE calculation
if(BSE) then
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,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(*,*)'-------------------------------------------------------------------------------'
2023-11-27 23:25:10 +01:00
write(*,'(2X,A50,F20.10,A3)') 'Tr@BSE@evGW@UHF correlation energy (spin-conserved) =',EcBSE(1),' au'
write(*,'(2X,A50,F20.10,A3)') 'Tr@BSE@evGW@UHF correlation energy (spin-flip) =',EcBSE(2),' au'
write(*,'(2X,A50,F20.10,A3)') 'Tr@BSE@evGW@UHF correlation energy =',sum(EcBSE),' au'
write(*,'(2X,A50,F20.10,A3)') 'Tr@BSE@evGW@UHF total energy =',ENuc + EUHF + sum(EcBSE),' au'
2020-10-07 22:51:30 +02:00
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
2023-07-23 11:16:42 +02:00
call UGW_phACFDT(exchange_kernel,doXBS,.true.,TDA_W,TDA,BSE,spin_conserved,spin_flip, &
2023-11-27 23:25:10 +01:00
eta,nBas,nC,nO,nV,nR,nS,ERI_aaaa,ERI_aabb,ERI_bbbb,eGW,eGW,EcRPA)
2020-09-24 22:50:56 +02:00
2020-10-07 22:51:30 +02:00
write(*,*)
write(*,*)'-------------------------------------------------------------------------------'
2023-11-27 23:25:10 +01:00
write(*,'(2X,A50,F20.10,A3)') 'AC@BSE@evGW@UHF correlation energy (spin-conserved) =',EcRPA(1),' au'
write(*,'(2X,A50,F20.10,A3)') 'AC@BSE@evGW@UHF correlation energy (spin-flip) =',EcRPA(2),' au'
write(*,'(2X,A50,F20.10,A3)') 'AC@BSE@evGW@UHF correlation energy =',sum(EcRPA),' au'
write(*,'(2X,A50,F20.10,A3)') 'AC@BSE@evGW@UHF total energy =',ENuc + EUHF + sum(EcRPA),' au'
2020-10-07 22:51:30 +02:00
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
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','evGW correlation energy',EcRPA)
call dump_test_value('U','evGW HOMOa energy',eGW(nO(1),1))
call dump_test_value('U','evGW LUMOa energy',eGW(nO(1)+1,1))
call dump_test_value('U','evGW HOMOa energy',eGW(nO(2),2))
call dump_test_value('U','evGW 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