4
1
mirror of https://github.com/pfloos/quack synced 2024-06-18 03:05:31 +02:00
quack/src/GT/Bethe_Salpeter_Tmatrix.f90

176 lines
6.7 KiB
Fortran
Raw Normal View History

2022-09-09 21:48:50 +02:00
subroutine Bethe_Salpeter_Tmatrix(TDA_T,TDA,dBSE,dTDA,evDyn,singlet,triplet,eta,nBas,nC,nO,nV,nR,nS,nOOab,nVVab,nOOaa,nVVaa, &
Om1ab,X1ab,Y1ab,Om2ab,X2ab,Y2ab,rho1ab,rho2ab,Om1aa,X1aa,Y1aa,Om2aa,X2aa,Y2aa,rho1aa,rho2aa, &
2021-10-16 15:34:34 +02:00
ERI,dipole_int,eT,eGT,EcBSE)
2021-10-15 13:51:04 +02:00
! Compute the Bethe-Salpeter excitation energies with the T-matrix kernel
implicit none
include 'parameters.h'
! Input variables
logical,intent(in) :: TDA_T
logical,intent(in) :: TDA
2021-10-15 22:32:22 +02:00
logical,intent(in) :: dBSE
logical,intent(in) :: dTDA
logical,intent(in) :: evDyn
2021-10-15 13:51:04 +02:00
logical,intent(in) :: singlet
logical,intent(in) :: triplet
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
2022-09-09 21:48:50 +02:00
integer,intent(in) :: nOOab
integer,intent(in) :: nOOaa
integer,intent(in) :: nVVab
integer,intent(in) :: nVVaa
2021-10-16 15:34:34 +02:00
2021-10-15 13:51:04 +02:00
double precision,intent(in) :: eT(nBas)
double precision,intent(in) :: eGT(nBas)
double precision,intent(in) :: ERI(nBas,nBas,nBas,nBas)
double precision,intent(in) :: dipole_int(nBas,nBas,ncart)
2022-09-09 21:48:50 +02:00
double precision,intent(in) :: Om1ab(nVVab)
double precision,intent(in) :: X1ab(nVVab,nVVab)
double precision,intent(in) :: Y1ab(nOOab,nVVab)
double precision,intent(in) :: Om2ab(nOOab)
double precision,intent(in) :: X2ab(nVVab,nOOab)
double precision,intent(in) :: Y2ab(nOOab,nOOab)
double precision,intent(in) :: rho1ab(nBas,nBas,nVVab)
double precision,intent(in) :: rho2ab(nBas,nBas,nOOab)
double precision,intent(in) :: Om1aa(nVVaa)
double precision,intent(in) :: X1aa(nVVaa,nVVaa)
double precision,intent(in) :: Y1aa(nOOaa,nVVaa)
double precision,intent(in) :: Om2aa(nOOaa)
double precision,intent(in) :: X2aa(nVVaa,nOOaa)
double precision,intent(in) :: Y2aa(nOOaa,nOOaa)
double precision,intent(in) :: rho1aa(nBas,nBas,nVVaa)
double precision,intent(in) :: rho2aa(nBas,nBas,nOOaa)
2021-10-16 15:34:34 +02:00
2021-10-15 13:51:04 +02:00
! Local variables
integer :: ispin
integer :: iblock
double precision :: EcRPA(nspin)
2022-09-09 21:48:50 +02:00
double precision,allocatable :: TAab(:,:),TBab(:,:)
double precision,allocatable :: TAaa(:,:),TBaa(:,:)
2021-10-15 13:51:04 +02:00
double precision,allocatable :: OmBSE(:,:)
double precision,allocatable :: XpY_BSE(:,:,:)
double precision,allocatable :: XmY_BSE(:,:,:)
! Output variables
double precision,intent(out) :: EcBSE(nspin)
! Memory allocation
2022-09-09 21:48:50 +02:00
allocate(TAab(nS,nS),TBab(nS,nS),TAaa(nS,nS),TBaa(nS,nS), &
2022-01-07 15:14:00 +01:00
OmBSE(nS,nspin),XpY_BSE(nS,nS,nspin),XmY_BSE(nS,nS,nspin))
2021-10-15 13:51:04 +02:00
2022-01-11 15:23:43 +01:00
!---------------------------------------!
! Compute T-matrix for alpha-beta block !
!---------------------------------------!
2021-10-15 13:51:04 +02:00
ispin = 1
iblock = 3
2022-09-09 21:48:50 +02:00
call linear_response_pp(iblock,TDA_T,nBas,nC,nO,nV,nR,nOOab,nVVab,1d0,eT,ERI, &
Om1ab,X1ab,Y1ab,Om2ab,X2ab,Y2ab,EcRPA(ispin))
2021-10-15 13:51:04 +02:00
2022-09-09 21:48:50 +02:00
call static_Tmatrix_A(eta,nBas,nC,nO,nV,nR,nS,nOOab,nVVab,1d0,Om1ab,rho1ab,Om2ab,rho2ab,TAab)
if(.not.TDA) call static_Tmatrix_B(eta,nBas,nC,nO,nV,nR,nS,nOOab,nVVab,1d0,Om1ab,rho1ab,Om2ab,rho2ab,TBab)
2022-01-02 10:24:30 +01:00
2022-01-11 15:23:43 +01:00
!----------------------------------------!
! Compute T-matrix for alpha-alpha block !
!----------------------------------------!
2021-10-15 13:51:04 +02:00
ispin = 2
iblock = 4
2022-09-09 21:48:50 +02:00
call linear_response_pp(iblock,TDA_T,nBas,nC,nO,nV,nR,nOOaa,nVVaa,1d0,eT,ERI, &
Om1aa,X1aa,Y1aa,Om2aa,X2aa,Y2aa,EcRPA(ispin))
2022-02-15 20:14:01 +01:00
2022-09-09 21:48:50 +02:00
call static_Tmatrix_A(eta,nBas,nC,nO,nV,nR,nS,nOOaa,nVVaa,1d0,Om1aa,rho1aa,Om2aa,rho2aa,TAaa)
if(.not.TDA) call static_Tmatrix_B(eta,nBas,nC,nO,nV,nR,nS,nOOaa,nVVaa,1d0,Om1aa,rho1aa,Om2aa,rho2aa,TBaa)
2022-01-07 15:14:00 +01:00
2022-01-11 15:23:43 +01:00
!------------------!
! Singlet manifold !
!------------------!
2021-10-15 13:51:04 +02:00
if(singlet) then
ispin = 1
EcBSE(ispin) = 0d0
! Compute BSE singlet excitation energies
2022-09-09 21:48:50 +02:00
call linear_response_BSE(ispin,.false.,TDA,.true.,eta,nBas,nC,nO,nV,nR,nS,1d0,eGT,ERI,TAab+TAaa,TBab+TBaa, &
2022-02-15 20:14:01 +01:00
EcBSE(ispin),OmBSE(:,ispin),XpY_BSE(:,:,ispin),XmY_BSE(:,:,ispin))
2021-12-21 13:49:41 +01:00
2021-10-15 13:51:04 +02:00
call print_excitation('BSE@GT ',ispin,nS,OmBSE(:,ispin))
2022-09-09 21:48:50 +02:00
call print_transition_vectors(.true.,nBas,nC,nO,nV,nR,nS,dipole_int,OmBSE(:,ispin),XpY_BSE(:,:,ispin),XmY_BSE(:,:,ispin))
2021-10-15 13:51:04 +02:00
2022-01-17 08:02:11 +01:00
if(dBSE) then
! Compute dynamic correction for BSE via perturbation theory (iterative or renormalized)
if(evDyn) then
print*, ' Iterative dynamical correction for BSE@GT NYI'
! call Bethe_Salpeter_dynamic_perturbation_iterative(dTDA,eta,nBas,nC,nO,nV,nR,nS,eGW,dipole_int,OmRPA,rho_RPA, &
! OmBSE(:,ispin),XpY_BSE(:,:,ispin),XmY_BSE(:,:,ispin))
else
2022-09-09 21:48:50 +02:00
call Bethe_Salpeter_Tmatrix_dynamic_perturbation(ispin,dTDA,eta,nBas,nC,nO,nV,nR,nS,nOOab,nVVab,nOOaa,nVVaa, &
Om1ab,Om2ab,Om1aa,Om2aa,rho1ab,rho2ab,rho1aa,rho2aa,eT,eGT, &
dipole_int,OmBSE(:,ispin),XpY_BSE(:,:,ispin),XmY_BSE(:,:,ispin),TAab,TAaa)
2022-01-17 08:02:11 +01:00
end if
end if
2021-10-15 13:51:04 +02:00
end if
2022-01-11 15:23:43 +01:00
!------------------!
! Triplet manifold !
!------------------!
2021-10-15 13:51:04 +02:00
if(triplet) then
ispin = 2
EcBSE(ispin) = 0d0
! Compute BSE triplet excitation energies
2022-09-09 21:48:50 +02:00
call linear_response_BSE(ispin,.false.,TDA,.true.,eta,nBas,nC,nO,nV,nR,nS,1d0,eGT,ERI,TAaa-TAab,TBaa-TBab, &
EcBSE(ispin),OmBSE(:,ispin),XpY_BSE(:,:,ispin),XmY_BSE(:,:,ispin))
2021-10-15 13:51:04 +02:00
call print_excitation('BSE@GT ',ispin,nS,OmBSE(:,ispin))
2022-09-09 21:48:50 +02:00
call print_transition_vectors(.false.,nBas,nC,nO,nV,nR,nS,dipole_int,OmBSE(:,ispin),XpY_BSE(:,:,ispin),XmY_BSE(:,:,ispin))
2021-10-15 13:51:04 +02:00
2022-01-17 08:02:11 +01:00
if(dBSE) then
! Compute dynamic correction for BSE via perturbation theory (iterative or renormalized)
if(evDyn) then
print*, ' Iterative dynamical correction for BSE@GT NYI'
! call Bethe_Salpeter_dynamic_perturbation_iterative(dTDA,eta,nBas,nC,nO,nV,nR,nS,eGW,dipole_int,OmRPA,rho_RPA, &
! OmBSE(:,ispin),XpY_BSE(:,:,ispin),XmY_BSE(:,:,ispin))
else
2022-09-09 21:48:50 +02:00
call Bethe_Salpeter_Tmatrix_dynamic_perturbation(ispin,dTDA,eta,nBas,nC,nO,nV,nR,nS,nOOab,nVVab,nOOaa,nVVaa, &
Om1ab,Om2ab,Om1aa,Om2aa,rho1ab,rho2ab,rho1aa,rho2aa,eT,eGT, &
dipole_int,OmBSE(:,ispin),XpY_BSE(:,:,ispin),XmY_BSE(:,:,ispin),TAab,TAaa)
2022-01-17 08:02:11 +01:00
end if
2021-10-15 22:32:22 +02:00
end if
2021-10-15 13:51:04 +02:00
end if
end subroutine Bethe_Salpeter_Tmatrix