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

430 lines
14 KiB
Fortran
Raw Normal View History

2023-11-13 17:39:30 +01:00
subroutine qsUGW(dotest,maxSCF,thresh,max_diis,doACFDT,exchange_kernel,doXBS,BSE,TDA_W,TDA,dBSE,dTDA,spin_conserved,spin_flip, &
2023-11-09 17:38:30 +01:00
eta,regularize,nNuc,ZNuc,rNuc,ENuc,nBas,nC,nO,nV,nR,nS,EUHF,S,X,T,V,Hc,ERI_AO,ERI_aaaa,ERI_aabb,ERI_bbbb, &
dipole_int_AO,dipole_int_aa,dipole_int_bb,PHF,cHF,eHF)
2020-10-21 12:09:18 +02:00
! Perform a quasiparticle self-consistent GW calculation
implicit none
include 'parameters.h'
! Input variables
2023-11-13 17:39:30 +01:00
logical,intent(in) :: dotest
2020-10-21 12:09:18 +02:00
integer,intent(in) :: maxSCF
integer,intent(in) :: max_diis
double precision,intent(in) :: thresh
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-10-21 12:09:18 +02:00
integer,intent(in) :: nNuc
double precision,intent(in) :: ZNuc(nNuc)
double precision,intent(in) :: rNuc(nNuc,ncart)
double precision,intent(in) :: ENuc
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) :: EUHF
double precision,intent(in) :: eHF(nBas,nspin)
double precision,intent(in) :: cHF(nBas,nBas,nspin)
double precision,intent(in) :: PHF(nBas,nBas,nspin)
double precision,intent(in) :: S(nBas,nBas)
double precision,intent(in) :: T(nBas,nBas)
double precision,intent(in) :: V(nBas,nBas)
double precision,intent(in) :: Hc(nBas,nBas)
double precision,intent(in) :: X(nBas,nBas)
double precision,intent(in) :: ERI_AO(nBas,nBas,nBas,nBas)
double precision,intent(inout):: ERI_aaaa(nBas,nBas,nBas,nBas)
double precision,intent(inout):: ERI_aabb(nBas,nBas,nBas,nBas)
double precision,intent(inout):: ERI_bbbb(nBas,nBas,nBas,nBas)
2020-10-21 12:58:37 +02:00
double precision,intent(in) :: dipole_int_AO(nBas,nBas,ncart)
2023-11-08 15:44:06 +01:00
double precision,intent(inout):: dipole_int_aa(nBas,nBas,ncart)
double precision,intent(inout):: dipole_int_bb(nBas,nBas,ncart)
2020-10-21 12:09:18 +02:00
! Local variables
logical :: dRPA
2020-10-21 12:09:18 +02:00
integer :: nSCF
integer :: nBasSq
integer :: ispin
2023-11-08 15:44:06 +01:00
integer :: ixyz
2020-10-21 12:09:18 +02:00
integer :: is
integer :: n_diis
integer :: nSa,nSb,nSt
2020-10-21 12:09:18 +02:00
double precision :: dipole(ncart)
2020-10-21 20:53:36 +02:00
double precision :: ET(nspin)
double precision :: EV(nspin)
double precision :: EJ(nsp)
2023-11-08 22:58:55 +01:00
double precision :: EK(nspin)
2023-11-27 23:25:10 +01:00
double precision :: EcRPA(nspin)
2021-03-01 14:55:29 +01:00
double precision :: EcGM(nspin)
2020-10-21 20:53:36 +02:00
double precision :: EqsGW
2020-10-21 12:09:18 +02:00
double precision :: EcBSE(nspin)
double precision :: Conv
double precision :: rcond(nspin)
double precision,external :: trace_matrix
2023-11-08 22:58:55 +01:00
double precision,allocatable :: err_diis(:,:,:)
2020-10-21 12:09:18 +02:00
double precision,allocatable :: F_diis(:,:,:)
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-10-21 12:09:18 +02:00
double precision,allocatable :: c(:,:,:)
double precision,allocatable :: cp(:,:,:)
double precision,allocatable :: eGW(:,:)
double precision,allocatable :: P(:,:,:)
double precision,allocatable :: F(:,:,:)
double precision,allocatable :: Fp(:,:,:)
double precision,allocatable :: J(:,:,:)
double precision,allocatable :: K(:,:,:)
double precision,allocatable :: SigC(:,:,:)
double precision,allocatable :: SigCp(:,:,:)
double precision,allocatable :: Z(:,:)
2023-11-08 22:58:55 +01:00
double precision,allocatable :: err(:,:,:)
2020-10-21 12:09:18 +02:00
! Hello world
write(*,*)
2023-11-09 17:22:01 +01:00
write(*,*)'*********************************'
write(*,*)'* Unrestricted qsGW Calculation *'
write(*,*)'*********************************'
2020-10-21 12:09:18 +02:00
write(*,*)
! Warning
write(*,*) '!! ERIs in MO basis will be overwritten in qsUGW !!'
write(*,*)
! Stuff
nBasSq = nBas*nBas
dRPA = .true.
2020-10-21 12:09:18 +02:00
! 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
! Memory allocation
nSa = nS(1)
nSb = nS(2)
nSt = nSa + nSb
2020-10-21 12:09:18 +02:00
allocate(Aph(nSt,nSt),Bph(nSt,nSt),eGW(nBas,nspin),c(nBas,nBas,nspin),cp(nBas,nBas,nspin),P(nBas,nBas,nspin), &
F(nBas,nBas,nspin),Fp(nBas,nBas,nspin),J(nBas,nBas,nspin),K(nBas,nBas,nspin), &
SigC(nBas,nBas,nspin),SigCp(nBas,nBas,nspin),Z(nBas,nspin),Om(nSt),XpY(nSt,nSt),XmY(nSt,nSt), &
rho(nBas,nBas,nSt,nspin),err(nBas,nBas,nspin),err_diis(nBasSq,max_diis,nspin),F_diis(nBasSq,max_diis,nspin))
2020-10-21 12:09:18 +02:00
! Initialization
2023-11-08 22:58:55 +01:00
nSCF = -1
n_diis = 0
ispin = 1
Conv = 1d0
P(:,:,:) = PHF(:,:,:)
eGW(:,:) = eHF(:,:)
c(:,:,:) = cHF(:,:,:)
F_diis(:,:,:) = 0d0
err_diis(:,:,:) = 0d0
rcond(:) = 0d0
2020-10-21 12:09:18 +02:00
!------------------------------------------------------------------------
! Main loop
!------------------------------------------------------------------------
do while(Conv > thresh .and. nSCF < maxSCF)
! Increment
nSCF = nSCF + 1
2023-10-26 09:35:48 +02:00
! Buid Hartree matrix
2020-10-21 12:09:18 +02:00
do is=1,nspin
2023-10-26 09:35:48 +02:00
call Hartree_matrix_AO_basis(nBas,P(:,:,is),ERI_AO(:,:,:,:),J(:,:,is))
2020-10-21 12:09:18 +02:00
end do
! Compute exchange part of the self-energy
do is=1,nspin
call exchange_matrix_AO_basis(nBas,P(:,:,is),ERI_AO(:,:,:,:),K(:,:,is))
end do
!--------------------------------------------------
! AO to MO transformation of two-electron integrals
!--------------------------------------------------
2023-11-08 15:44:06 +01:00
do ixyz=1,ncart
2023-11-10 17:22:51 +01:00
call AOtoMO(nBas,c(:,:,1),dipole_int_AO(:,:,ixyz),dipole_int_aa(:,:,ixyz))
call AOtoMO(nBas,c(:,:,2),dipole_int_AO(:,:,ixyz),dipole_int_bb(:,:,ixyz))
2023-11-08 15:44:06 +01:00
end do
2020-10-21 12:09:18 +02:00
! 4-index transform for (aa|aa) block
2023-11-22 17:02:46 +01:00
call AOtoMO_ERI_UHF(1,1,nBas,c,ERI_AO,ERI_aaaa)
2020-10-21 12:09:18 +02:00
! 4-index transform for (aa|bb) block
2023-11-22 17:02:46 +01:00
call AOtoMO_ERI_UHF(1,2,nBas,c,ERI_AO,ERI_aabb)
2020-10-21 12:09:18 +02:00
! 4-index transform for (bb|bb) block
2023-11-22 17:02:46 +01:00
call AOtoMO_ERI_UHF(2,2,nBas,c,ERI_AO,ERI_bbbb)
2020-10-21 12:09:18 +02:00
! Compute linear response
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-10-21 12:09:18 +02:00
!----------------------!
! Excitation densities !
!----------------------!
call UGW_excitation_density(nBas,nC,nO,nR,nSa,nSb,nSt,ERI_aaaa,ERI_aabb,ERI_bbbb,XpY,rho)
2020-10-21 12:09:18 +02:00
!------------------------------------------------!
! Compute self-energy and renormalization factor !
!------------------------------------------------!
2022-11-30 16:41:19 +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,eGW(:,is),Om,rho(:,:,:,is))
2023-08-01 16:30:41 +02:00
end do
2022-11-30 16:41:19 +01:00
end if
2020-10-21 12:09:18 +02:00
call UGW_self_energy(eta,nBas,nC,nO,nV,nR,nSt,eGW,Om,rho,SigC,Z,EcGM)
2023-08-01 16:30:41 +02:00
2020-10-21 12:09:18 +02:00
! Make correlation self-energy Hermitian and transform it back to AO basis
do is=1,nspin
2023-11-08 22:58:55 +01:00
SigC(:,:,is) = 0.5d0*(SigC(:,:,is) + transpose(SigC(:,:,is)))
2020-10-21 12:09:18 +02:00
end do
do is=1,nspin
2023-11-10 17:22:51 +01:00
call MOtoAO(nBas,S,c(:,:,is),SigC(:,:,is),SigCp(:,:,is))
2020-10-21 12:09:18 +02:00
end do
! Solve the quasi-particle equation
do is=1,nspin
F(:,:,is) = Hc(:,:) + J(:,:,is) + J(:,:,mod(is,2)+1) + K(:,:,is) + SigCp(:,:,is)
end do
! Check convergence
do is=1,nspin
2023-11-08 22:58:55 +01:00
err(:,:,is) = matmul(F(:,:,is),matmul(P(:,:,is),S(:,:))) - matmul(matmul(S(:,:),P(:,:,is)),F(:,:,is))
2020-10-21 12:09:18 +02:00
end do
2023-11-09 22:24:57 +01:00
if(nSCF > 1) Conv = maxval(abs(err))
2020-10-21 12:09:18 +02:00
! DIIS extrapolation
2023-11-09 22:24:57 +01:00
if(max_diis > 1) then
n_diis = min(n_diis+1,max_diis)
do is=1,nspin
if(nO(is) > 1) call DIIS_extrapolation(rcond(is),nBasSq,nBasSq,n_diis,err_diis(:,1:n_diis,is), &
F_diis(:,1:n_diis,is),err(:,:,is),F(:,:,is))
end do
end if
2020-10-21 12:09:18 +02:00
! Transform Fock matrix in orthogonal basis
do is=1,nspin
Fp(:,:,is) = matmul(transpose(X(:,:)),matmul(F(:,:,is),X(:,:)))
end do
! Diagonalize Fock matrix to get eigenvectors and eigenvalues
cp(:,:,:) = Fp(:,:,:)
do is=1,nspin
call diagonalize_matrix(nBas,cp(:,:,is),eGW(:,is))
end do
! Back-transform eigenvectors in non-orthogonal basis
do is=1,nspin
c(:,:,is) = matmul(X(:,:),cp(:,:,is))
end do
2021-03-27 15:03:54 +01:00
! Back-transform self-energy
do is=1,nspin
2023-11-10 17:22:51 +01:00
call AOtoMO(nBas,c(:,:,is),SigCp(:,:,is),SigC(:,:,is))
2021-03-27 15:03:54 +01:00
end do
2020-10-21 12:09:18 +02:00
! Compute density matrix
do is=1,nspin
P(:,:,is) = matmul(c(:,1:nO(is),is),transpose(c(:,1:nO(is),is)))
end do
2020-10-21 20:53:36 +02:00
!------------------------------------------------------------------------
! Compute total energy
!------------------------------------------------------------------------
! Kinetic energy
2020-11-03 21:47:57 +01:00
do is=1,nspin
ET(is) = trace_matrix(nBas,matmul(P(:,:,is),T(:,:)))
2020-10-21 20:53:36 +02:00
end do
! Potential energy
2020-11-03 21:47:57 +01:00
do is=1,nspin
EV(is) = trace_matrix(nBas,matmul(P(:,:,is),V(:,:)))
2020-10-21 20:53:36 +02:00
end do
2023-10-26 09:35:48 +02:00
! Hartree energy
2020-10-21 20:53:36 +02:00
EJ(1) = 0.5d0*trace_matrix(nBas,matmul(P(:,:,1),J(:,:,1)))
2021-03-01 14:55:29 +01:00
EJ(2) = 0.5d0*trace_matrix(nBas,matmul(P(:,:,1),J(:,:,2))) &
+ 0.5d0*trace_matrix(nBas,matmul(P(:,:,2),J(:,:,1)))
2020-10-21 20:53:36 +02:00
EJ(3) = 0.5d0*trace_matrix(nBas,matmul(P(:,:,2),J(:,:,2)))
! Exchange energy
2020-11-03 21:47:57 +01:00
do is=1,nspin
2023-11-08 22:58:55 +01:00
EK(is) = 0.5d0*trace_matrix(nBas,matmul(P(:,:,is),K(:,:,is)))
2020-10-21 20:53:36 +02:00
end do
! Total energy
2023-11-08 22:58:55 +01:00
EqsGW = sum(ET(:)) + sum(EV(:)) + sum(EJ(:)) + sum(EK(:))
2020-10-21 20:53:36 +02:00
!------------------------------------------------------------------------
2020-10-21 12:09:18 +02:00
! Print results
2020-10-21 20:53:36 +02:00
!------------------------------------------------------------------------
2020-10-21 12:09:18 +02:00
2020-10-21 12:58:37 +02:00
call dipole_moment(nBas,P(:,:,1)+P(:,:,2),nNuc,ZNuc,rNuc,dipole_int_AO,dipole)
2023-11-27 23:25:10 +01:00
call print_qsUGW(nBas,nO,nSCF,Conv,thresh,eHF,eGW,c,S,ENuc,ET,EV,EJ,EK,EcGM,EcRPA(ispin),EqsGW,SigCp,Z,dipole)
2020-10-21 12:09:18 +02:00
2023-12-03 18:47:30 +01:00
end do
2020-10-21 12:09:18 +02:00
!------------------------------------------------------------------------
! End main loop
!------------------------------------------------------------------------
! Did it actually converge?
if(nSCF == maxSCF) then
write(*,*)
write(*,*)'!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!'
write(*,*)' Convergence failed '
write(*,*)'!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!'
write(*,*)
stop
2023-12-03 18:47:30 +01:00
end if
2020-10-21 12:09:18 +02:00
! Deallocate memory
2023-11-08 22:58:55 +01:00
deallocate(cp,P,F,Fp,J,K,SigC,SigCp,Z,Om,XpY,XmY,rho,err,err_diis,F_diis)
2020-10-21 12:09:18 +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,c,eGW,eGW,EcBSE)
2020-10-21 12:09:18 +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)
else
EcBSE(2) = 0.0d0
2020-10-21 12:09:18 +02:00
end if
write(*,*)
write(*,*)'-------------------------------------------------------------------------------'
2023-11-27 23:25:10 +01:00
write(*,'(2X,A50,F20.10,A3)') 'Tr@BSE@qsGW@UHF correlation energy (spin-conserved) = ',EcBSE(1),' au'
write(*,'(2X,A50,F20.10,A3)') 'Tr@BSE@qsGW@UHF correlation energy (spin-flip) = ',EcBSE(2),' au'
write(*,'(2X,A50,F20.10,A3)') 'Tr@BSE@qsGW@UHF correlation energy = ',sum(EcBSE),' au'
write(*,'(2X,A50,F20.10,A3)') 'Tr@BSE@qsGW@UHF total energy = ',ENuc + EqsGW + sum(EcBSE),' au'
2020-10-21 12:09:18 +02:00
write(*,*)'-------------------------------------------------------------------------------'
write(*,*)
! Compute the BSE correlation energy via the adiabatic connection
if(doACFDT) then
write(*,*) '--------------------------------------------------------------'
write(*,*) ' Adiabatic connection version of BSE@qsUGW correlation energy '
write(*,*) '--------------------------------------------------------------'
write(*,*)
if(doXBS) then
write(*,*) '*** scaled screening version (XBS) ***'
write(*,*)
end if
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-10-21 12:09:18 +02:00
write(*,*)
write(*,*)'-------------------------------------------------------------------------------'
2023-11-27 23:25:10 +01:00
write(*,'(2X,A50,F20.10,A3)') 'AC@BSE@qsGW@UHF correlation energy (spin-conserved) = ',EcRPA(1),' au'
write(*,'(2X,A50,F20.10,A3)') 'AC@BSE@qsGW@UHF correlation energy (spin-flip) = ',EcRPA(2),' au'
write(*,'(2X,A50,F20.10,A3)') 'AC@BSE@qsGW@UHF correlation energy = ',sum(EcRPA),' au'
write(*,'(2X,A50,F20.10,A3)') 'AC@BSE@qsGW@UHF total energy = ',ENuc + EqsGW + sum(EcRPA),' au'
2020-10-21 12:09:18 +02:00
write(*,*)'-------------------------------------------------------------------------------'
write(*,*)
end if
end if
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','qsGW correlation energy',EcRPA)
call dump_test_value('U','qsGW HOMOa energy',eGW(nO(1),1))
call dump_test_value('U','qsGW LUMOa energy',eGW(nO(1)+1,1))
call dump_test_value('U','qsGW HOMOa energy',eGW(nO(2),2))
call dump_test_value('U','qsGW 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