4
1
mirror of https://github.com/pfloos/quack synced 2024-10-15 20:42:01 +02:00
quack/src/GW/GW_ppBSE.f90

175 lines
5.5 KiB
Fortran
Raw Normal View History

2023-07-04 11:56:00 +02:00
subroutine GW_ppBSE(TDA_W,TDA,singlet,triplet,eta,nBas,nC,nO,nV,nR,nS,ERI,dipole_int,eW,eGW,EcBSE)
2022-08-17 14:32:14 +02:00
! Compute the Bethe-Salpeter excitation energies at the pp level
implicit none
include 'parameters.h'
! Input variables
logical,intent(in) :: TDA_W
logical,intent(in) :: TDA
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
double precision,intent(in) :: eW(nBas)
double precision,intent(in) :: eGW(nBas)
double precision,intent(in) :: ERI(nBas,nBas,nBas,nBas)
double precision,intent(in) :: dipole_int(nBas,nBas,ncart)
! Local variables
integer :: ispin
integer :: isp_W
2023-07-20 12:37:54 +02:00
logical :: dRPA = .false.
logical :: dRPA_W = .true.
2022-08-17 14:32:14 +02:00
integer :: nOO
integer :: nVV
2023-07-21 12:00:08 +02:00
double precision,allocatable :: Aph(:,:)
double precision,allocatable :: Bph(:,:)
2023-07-20 12:37:54 +02:00
2022-08-17 14:32:14 +02:00
double precision :: EcRPA
double precision,allocatable :: OmRPA(:)
double precision,allocatable :: XpY_RPA(:,:)
double precision,allocatable :: XmY_RPA(:,:)
double precision,allocatable :: rho_RPA(:,:,:)
2023-07-20 14:59:02 +02:00
double precision,allocatable :: Bpp(:,:)
double precision,allocatable :: Cpp(:,:)
double precision,allocatable :: Dpp(:,:)
2023-07-20 12:37:54 +02:00
double precision,allocatable :: Om1(:)
2022-08-17 14:32:14 +02:00
double precision,allocatable :: X1(:,:)
double precision,allocatable :: Y1(:,:)
2023-07-20 12:37:54 +02:00
double precision,allocatable :: Om2(:)
2022-08-17 14:32:14 +02:00
double precision,allocatable :: X2(:,:)
double precision,allocatable :: Y2(:,:)
2023-07-21 12:00:08 +02:00
double precision,allocatable :: KB_sta(:,:)
double precision,allocatable :: KC_sta(:,:)
double precision,allocatable :: KD_sta(:,:)
2022-08-17 14:32:14 +02:00
! Output variables
double precision,intent(out) :: EcBSE(nspin)
!---------------------------------
! Compute (singlet) RPA screening
!---------------------------------
isp_W = 1
EcRPA = 0d0
2023-07-21 12:00:08 +02:00
allocate(OmRPA(nS),XpY_RPA(nS,nS),XmY_RPA(nS,nS),rho_RPA(nBas,nBas,nS), &
Aph(nS,nS),Bph(nS,nS))
call phLR_A(isp_W,dRPA_W,nBas,nC,nO,nV,nR,nS,1d0,eW,ERI,Aph)
if(.not.TDA_W) call phLR_B(isp_W,dRPA_W,nBas,nC,nO,nV,nR,nS,1d0,ERI,Bph)
2022-08-17 14:32:14 +02:00
2023-07-21 12:00:08 +02:00
call phLR(TDA_W,nS,Aph,Bph,EcRPA,OmRPA,XpY_RPA,XmY_RPA)
2023-07-04 10:51:03 +02:00
call GW_excitation_density(nBas,nC,nO,nR,nS,ERI,XpY_RPA,rho_RPA)
2022-08-17 14:32:14 +02:00
2023-07-21 12:00:08 +02:00
deallocate(XpY_RPA,XmY_RPA,Aph,Bph)
2023-07-20 12:37:54 +02:00
2022-08-17 14:32:14 +02:00
!-------------------
! Singlet manifold
!-------------------
if(singlet) then
2022-09-09 21:48:50 +02:00
write(*,*) '****************'
write(*,*) '*** Singlets ***'
write(*,*) '****************'
write(*,*)
2022-08-17 14:32:14 +02:00
ispin = 1
2023-07-21 12:00:08 +02:00
EcBSE(ispin) = 0d0
2022-08-17 14:32:14 +02:00
nOO = nO*(nO+1)/2
nVV = nV*(nV+1)/2
2023-07-20 14:59:02 +02:00
allocate(Om1(nVV),X1(nVV,nVV),Y1(nOO,nVV), &
Om2(nOO),X2(nVV,nOO),Y2(nOO,nOO), &
Bpp(nVV,nOO),Cpp(nVV,nVV),Dpp(nOO,nOO), &
2023-07-21 12:00:08 +02:00
KB_sta(nVV,nOO),KC_sta(nVV,nVV),KD_sta(nOO,nOO))
2022-08-17 14:32:14 +02:00
2023-07-20 14:59:02 +02:00
! Compute BSE excitation energies
2023-07-21 12:00:08 +02:00
if(.not.TDA) call GW_ppBSE_static_kernel_B(ispin,eta,nBas,nC,nO,nV,nR,nS,nOO,nVV,1d0,ERI,OmRPA,rho_RPA,KB_sta)
call GW_ppBSE_static_kernel_C(ispin,eta,nBas,nC,nO,nV,nR,nS,nOO,nVV,1d0,ERI,OmRPA,rho_RPA,KC_sta)
call GW_ppBSE_static_kernel_D(ispin,eta,nBas,nC,nO,nV,nR,nS,nOO,nVV,1d0,ERI,OmRPA,rho_RPA,KD_sta)
2022-08-17 14:32:14 +02:00
2023-07-20 14:59:02 +02:00
if(.not.TDA) call ppLR_B(ispin,nBas,nC,nO,nV,nR,nOO,nVV,1d0,ERI,Bpp)
call ppLR_C(ispin,nBas,nC,nO,nV,nR,nVV,1d0,eGW,ERI,Cpp)
call ppLR_D(ispin,nBas,nC,nO,nV,nR,nOO,1d0,eGW,ERI,Dpp)
2022-08-17 14:32:14 +02:00
2023-07-21 12:00:08 +02:00
Bpp(:,:) = Bpp(:,:) + KB_sta(:,:)
Cpp(:,:) = Cpp(:,:) + KC_sta(:,:)
Dpp(:,:) = Dpp(:,:) + KD_sta(:,:)
2023-07-20 14:59:02 +02:00
call ppLR(TDA,nOO,nVV,Bpp,Cpp,Dpp,Om1,X1,Y1,Om2,X2,Y2,EcBSE(ispin))
2022-08-17 14:32:14 +02:00
2023-07-20 12:37:54 +02:00
call print_transition_vectors_pp(.true.,nBas,nC,nO,nV,nR,nOO,nVV,dipole_int,Om1,X1,Y1,Om2,X2,Y2)
2022-08-17 15:52:13 +02:00
2023-07-21 12:00:08 +02:00
deallocate(Om1,X1,Y1,Om2,X2,Y2,Bpp,Cpp,Dpp,KB_sta,KC_sta,KD_sta)
2022-08-17 14:32:14 +02:00
end if
!-------------------
! Triplet manifold
!-------------------
if(triplet) then
2022-09-09 21:48:50 +02:00
write(*,*) '****************'
write(*,*) '*** Triplets ***'
write(*,*) '****************'
write(*,*)
2022-08-17 14:32:14 +02:00
ispin = 2
EcBSE(ispin) = 0d0
nOO = nO*(nO-1)/2
nVV = nV*(nV-1)/2
2023-07-20 14:59:02 +02:00
allocate(Om1(nVV),X1(nVV,nVV),Y1(nOO,nVV), &
Om2(nOO),X2(nVV,nOO),Y2(nOO,nOO), &
Bpp(nVV,nOO),Cpp(nVV,nVV),Dpp(nOO,nOO), &
2023-07-21 12:00:08 +02:00
KB_sta(nVV,nOO),KC_sta(nVV,nVV),KD_sta(nOO,nOO))
2022-08-17 14:32:14 +02:00
2023-07-20 14:59:02 +02:00
! Compute BSE excitation energies
2023-07-21 12:00:08 +02:00
if(.not.TDA) call GW_ppBSE_static_kernel_B(ispin,eta,nBas,nC,nO,nV,nR,nS,nOO,nVV,1d0,ERI,OmRPA,rho_RPA,KB_sta)
call GW_ppBSE_static_kernel_C(ispin,eta,nBas,nC,nO,nV,nR,nS,nOO,nVV,1d0,ERI,OmRPA,rho_RPA,KC_sta)
call GW_ppBSE_static_kernel_D(ispin,eta,nBas,nC,nO,nV,nR,nS,nOO,nVV,1d0,ERI,OmRPA,rho_RPA,KD_sta)
2022-08-17 14:32:14 +02:00
2023-07-20 14:59:02 +02:00
if(.not.TDA) call ppLR_B(ispin,nBas,nC,nO,nV,nR,nOO,nVV,1d0,ERI,Bpp)
call ppLR_C(ispin,nBas,nC,nO,nV,nR,nVV,1d0,eGW,ERI,Cpp)
call ppLR_D(ispin,nBas,nC,nO,nV,nR,nOO,1d0,eGW,ERI,Dpp)
2023-07-21 12:00:08 +02:00
Bpp(:,:) = Bpp(:,:) + KB_sta(:,:)
Cpp(:,:) = Cpp(:,:) + KC_sta(:,:)
Dpp(:,:) = Dpp(:,:) + KD_sta(:,:)
2022-08-17 14:32:14 +02:00
2023-07-20 14:59:02 +02:00
call ppLR(TDA,nOO,nVV,Bpp,Cpp,Dpp,Om1,X1,Y1,Om2,X2,Y2,EcBSE(ispin))
2022-08-17 14:32:14 +02:00
2023-07-20 12:37:54 +02:00
call print_transition_vectors_pp(.false.,nBas,nC,nO,nV,nR,nOO,nVV,dipole_int,Om1,X1,Y1,Om2,X2,Y2)
2022-08-17 15:52:13 +02:00
2023-07-21 12:00:08 +02:00
deallocate(Om1,X1,Y1,Om2,X2,Y2,Bpp,Cpp,Dpp,KB_sta,KC_sta,KD_sta)
2022-08-17 14:32:14 +02:00
end if
2023-07-04 11:56:00 +02:00
end subroutine