2020-06-15 23:04:07 +02:00
|
|
|
subroutine BSE2_dynamic_perturbation(dTDA,ispin,eta,nBas,nC,nO,nV,nR,nS,ERI,eHF,eGF,OmBSE,XpY,XmY)
|
2020-06-01 17:30:02 +02:00
|
|
|
|
|
|
|
! Compute dynamical effects via perturbation theory for BSE
|
|
|
|
|
|
|
|
implicit none
|
|
|
|
include 'parameters.h'
|
|
|
|
|
|
|
|
! Input variables
|
|
|
|
|
2020-06-15 23:04:07 +02:00
|
|
|
logical,intent(in) :: dTDA
|
2020-06-02 13:20:20 +02:00
|
|
|
integer,intent(in) :: ispin
|
2020-06-01 17:30:02 +02:00
|
|
|
double precision,intent(in) :: eta
|
|
|
|
integer,intent(in) :: nBas
|
|
|
|
integer,intent(in) :: nC
|
|
|
|
integer,intent(in) :: nO
|
|
|
|
integer,intent(in) :: nV
|
|
|
|
integer,intent(in) :: nR
|
|
|
|
integer,intent(in) :: nS
|
|
|
|
|
2020-06-01 22:01:21 +02:00
|
|
|
double precision,intent(in) :: ERI(nBas,nBas,nBas,nBas)
|
2020-06-01 17:30:02 +02:00
|
|
|
double precision,intent(in) :: eHF(nBas)
|
|
|
|
double precision,intent(in) :: eGF(nBas)
|
|
|
|
double precision,intent(in) :: OmBSE(nS)
|
|
|
|
double precision,intent(in) :: XpY(nS,nS)
|
|
|
|
double precision,intent(in) :: XmY(nS,nS)
|
|
|
|
|
|
|
|
! Local variables
|
|
|
|
|
|
|
|
integer :: ia
|
|
|
|
integer,parameter :: maxS = 10
|
|
|
|
double precision :: gapGF
|
|
|
|
|
|
|
|
double precision,allocatable :: OmDyn(:)
|
|
|
|
double precision,allocatable :: ZDyn(:)
|
|
|
|
double precision,allocatable :: X(:)
|
|
|
|
double precision,allocatable :: Y(:)
|
|
|
|
|
2020-06-17 14:32:37 +02:00
|
|
|
double precision,allocatable :: Ap_dyn(:,:)
|
|
|
|
double precision,allocatable :: Am_dyn(:,:)
|
|
|
|
double precision,allocatable :: ZAp_dyn(:,:)
|
|
|
|
double precision,allocatable :: ZAm_dyn(:,:)
|
2020-06-01 17:30:02 +02:00
|
|
|
|
2020-06-01 22:01:21 +02:00
|
|
|
double precision,allocatable :: B_dyn(:,:)
|
2020-06-01 17:30:02 +02:00
|
|
|
|
|
|
|
! Memory allocation
|
|
|
|
|
2020-06-17 14:32:37 +02:00
|
|
|
allocate(OmDyn(nS),ZDyn(nS),X(nS),Y(nS),Ap_dyn(nS,nS),ZAp_dyn(nS,nS))
|
2020-06-01 17:30:02 +02:00
|
|
|
|
2020-06-17 22:10:08 +02:00
|
|
|
if(.not.dTDA) allocate(Am_dyn(nS,nS),ZAm_dyn(nS,nS),B_dyn(nS,nS))
|
2020-06-01 17:30:02 +02:00
|
|
|
|
2020-06-18 08:47:51 +02:00
|
|
|
if(dTDA) then
|
|
|
|
write(*,*)
|
|
|
|
write(*,*) '*** dynamical TDA activated ***'
|
|
|
|
write(*,*)
|
|
|
|
end if
|
|
|
|
|
2020-06-03 12:06:16 +02:00
|
|
|
! Print main components of transition vectors
|
|
|
|
|
|
|
|
call print_transition_vectors(nBas,nC,nO,nV,nR,nS,OmBSE,XpY,XmY)
|
|
|
|
|
2020-06-01 17:30:02 +02:00
|
|
|
gapGF = eGF(nO+1) - eGF(nO)
|
|
|
|
|
|
|
|
write(*,*) '---------------------------------------------------------------------------------------------------'
|
|
|
|
write(*,*) ' First-order dynamical correction to static 2nd-order Bethe-Salpeter excitation energies '
|
|
|
|
write(*,*) '---------------------------------------------------------------------------------------------------'
|
2020-06-03 12:06:16 +02:00
|
|
|
write(*,'(A58,F10.6,A3)') ' BSE neutral excitation must be lower than the GF2 gap = ',gapGF*HaToeV,' eV'
|
|
|
|
write(*,*) '---------------------------------------------------------------------------------------------------'
|
2020-06-01 17:30:02 +02:00
|
|
|
write(*,'(2X,A5,1X,A20,1X,A20,1X,A20,1X,A20)') '#','Static (eV)','Dynamic (eV)','Correction (eV)','Renorm. (eV)'
|
|
|
|
write(*,*) '---------------------------------------------------------------------------------------------------'
|
|
|
|
|
|
|
|
do ia=1,min(nS,maxS)
|
|
|
|
|
|
|
|
X(:) = 0.5d0*(XpY(ia,:) + XmY(ia,:))
|
|
|
|
Y(:) = 0.5d0*(XpY(ia,:) - XmY(ia,:))
|
|
|
|
|
2020-06-01 22:01:21 +02:00
|
|
|
! Resonant part of the BSE correction for dynamical TDA
|
2020-06-01 17:30:02 +02:00
|
|
|
|
2020-06-17 14:32:37 +02:00
|
|
|
call BSE2_A_matrix_dynamic(ispin,eta,nBas,nC,nO,nV,nR,nS,1d0,ERI,eGF,+OmBSE(ia),Ap_dyn,ZAp_dyn)
|
2020-06-01 17:30:02 +02:00
|
|
|
|
2020-06-01 22:01:21 +02:00
|
|
|
if(dTDA) then
|
2020-06-01 17:30:02 +02:00
|
|
|
|
2020-06-17 14:32:37 +02:00
|
|
|
ZDyn(ia) = dot_product(X,matmul(ZAp_dyn,X))
|
|
|
|
OmDyn(ia) = dot_product(X,matmul(Ap_dyn,X))
|
2020-06-01 17:30:02 +02:00
|
|
|
|
|
|
|
else
|
|
|
|
|
2020-06-17 14:32:37 +02:00
|
|
|
! Second part of the resonant and anti-resonant part of the BSE correction (frequency independent)
|
2020-06-01 17:30:02 +02:00
|
|
|
|
2020-06-17 14:32:37 +02:00
|
|
|
call BSE2_A_matrix_dynamic(ispin,eta,nBas,nC,nO,nV,nR,nS,1d0,ERI,eGF,-OmBSE(ia),Am_dyn,ZAm_dyn)
|
2020-06-01 17:30:02 +02:00
|
|
|
|
2020-06-17 22:10:08 +02:00
|
|
|
call BSE2_B_matrix_dynamic(ispin,eta,nBas,nC,nO,nV,nR,nS,1d0,ERI,eGF,B_dyn)
|
2020-06-01 17:30:02 +02:00
|
|
|
|
2020-06-17 14:32:37 +02:00
|
|
|
ZDyn(ia) = dot_product(X,matmul(ZAp_dyn,X)) &
|
2020-06-17 22:10:08 +02:00
|
|
|
+ dot_product(Y,matmul(ZAm_dyn,Y))
|
2020-06-17 14:32:37 +02:00
|
|
|
|
|
|
|
OmDyn(ia) = dot_product(X,matmul(Ap_dyn,X)) &
|
|
|
|
- dot_product(Y,matmul(Am_dyn,Y)) &
|
|
|
|
+ dot_product(X,matmul(B_dyn,Y)) &
|
|
|
|
- dot_product(Y,matmul(B_dyn,X))
|
2020-06-01 17:30:02 +02:00
|
|
|
|
|
|
|
end if
|
|
|
|
|
|
|
|
ZDyn(ia) = 1d0/(1d0 - ZDyn(ia))
|
|
|
|
OmDyn(ia) = ZDyn(ia)*OmDyn(ia)
|
|
|
|
|
|
|
|
write(*,'(2X,I5,5X,F15.6,5X,F15.6,5X,F15.6,5X,F15.6)') &
|
|
|
|
ia,OmBSE(ia)*HaToeV,(OmBSE(ia)+OmDyn(ia))*HaToeV,OmDyn(ia)*HaToeV,ZDyn(ia)
|
|
|
|
|
|
|
|
end do
|
|
|
|
write(*,*) '---------------------------------------------------------------------------------------------------'
|
|
|
|
write(*,*)
|
|
|
|
|
|
|
|
end subroutine BSE2_dynamic_perturbation
|