10
1
mirror of https://github.com/pfloos/quack synced 2024-12-26 14:23:38 +01:00
QuAcK/src/QuAcK/BSE2_dynamic_perturbation.f90

106 lines
3.9 KiB
Fortran
Raw Normal View History

2020-06-01 22:01:21 +02:00
subroutine BSE2_dynamic_perturbation(singlet_manifold,triplet_manifold,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-01 22:01:21 +02:00
logical,intent(in) :: singlet_manifold
logical,intent(in) :: triplet_manifold
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
2020-06-01 22:01:21 +02:00
logical :: dTDA = .true.
2020-06-01 17:30:02 +02:00
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-01 22:01:21 +02:00
double precision,allocatable :: A_dyn(:,:)
double precision,allocatable :: ZA_dyn(:,:)
2020-06-01 17:30:02 +02:00
2020-06-01 22:01:21 +02:00
double precision,allocatable :: B_dyn(:,:)
double precision,allocatable :: ZB_dyn(:,:)
2020-06-01 17:30:02 +02:00
! Memory allocation
2020-06-01 22:01:21 +02:00
allocate(OmDyn(nS),ZDyn(nS),X(nS),Y(nS),A_dyn(nS,nS),ZA_dyn(nS,nS))
2020-06-01 17:30:02 +02:00
2020-06-01 22:01:21 +02:00
if(.not.dTDA) allocate(B_dyn(nS,nS),ZB_dyn(nS,nS))
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(*,*) '---------------------------------------------------------------------------------------------------'
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-01 22:01:21 +02:00
call BSE2_A_matrix_dynamic(singlet_manifold,triplet_manifold,eta,nBas,nC,nO,nV,nR,nS,1d0, &
ERI(:,:,:,:),eHF(:),eGF(:),OmBSE(ia),A_dyn(:,:),ZA_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-01 22:01:21 +02:00
ZDyn(ia) = dot_product(X(:),matmul(ZA_dyn(:,:),X(:)))
OmDyn(ia) = dot_product(X(:),matmul(A_dyn(:,:),X(:)))
2020-06-01 17:30:02 +02:00
else
2020-06-01 22:01:21 +02:00
! Anti-resonant part of the BSE correction
2020-06-01 17:30:02 +02:00
2020-06-01 22:01:21 +02:00
call BSE2_B_matrix_dynamic(singlet_manifold,triplet_manifold,eta,nBas,nC,nO,nV,nR,nS,1d0, &
ERI(:,:,:,:),eHF(:),eGF(:),OmBSE(ia),B_dyn(:,:),ZB_dyn(:,:))
2020-06-01 17:30:02 +02:00
2020-06-01 22:01:21 +02:00
ZDyn(ia) = dot_product(X(:),matmul(ZA_dyn(:,:),X(:))) &
- dot_product(Y(:),matmul(ZA_dyn(:,:),Y(:))) &
+ dot_product(X(:),matmul(ZB_dyn(:,:),Y(:))) &
- dot_product(Y(:),matmul(ZB_dyn(:,:),X(:)))
2020-06-01 17:30:02 +02:00
2020-06-01 22:01:21 +02:00
OmDyn(ia) = dot_product(X(:),matmul(A_dyn(:,:),X(:))) &
- dot_product(Y(:),matmul(A_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)
if(OmBSE(ia) > gapGF) write(*,*) ' !!! BSE2 neutral excitation larger than the GF2 gap !!! '
end do
write(*,*) '---------------------------------------------------------------------------------------------------'
write(*,*)
end subroutine BSE2_dynamic_perturbation