4
1
mirror of https://github.com/pfloos/quack synced 2024-06-19 11:42:20 +02:00
quack/src/QuAcK/Bethe_Salpeter.f90

92 lines
3.0 KiB
Fortran
Raw Normal View History

2020-01-23 21:22:41 +01:00
subroutine Bethe_Salpeter(TDA,singlet_manifold,triplet_manifold,eta, &
2020-04-22 00:39:52 +02:00
nBas,nC,nO,nV,nR,nS,ERI,eW,eGW,OmRPA,XpY,XmY,rho,EcRPA,EcBSE)
2020-01-14 21:27:34 +01:00
! Compute the Bethe-Salpeter excitation energies
implicit none
include 'parameters.h'
! Input variables
logical,intent(in) :: TDA
logical,intent(in) :: singlet_manifold
logical,intent(in) :: triplet_manifold
2020-01-23 21:22:41 +01:00
double precision,intent(in) :: eta
2020-01-14 21:27:34 +01:00
integer,intent(in) :: nBas,nC,nO,nV,nR,nS
double precision,intent(in) :: eW(nBas)
double precision,intent(in) :: eGW(nBas)
double precision,intent(in) :: ERI(nBas,nBas,nBas,nBas)
2020-04-22 00:39:52 +02:00
double precision :: OmRPA(nS,nspin)
2020-01-14 21:27:34 +01:00
double precision :: XpY(nS,nS,nspin)
double precision :: XmY(nS,nS,nspin)
double precision :: rho(nBas,nBas,nS,nspin)
! Local variables
integer :: ispin
2020-04-22 00:39:52 +02:00
double precision,allocatable :: OmBSE(:,:)
2020-01-14 21:27:34 +01:00
! Output variables
double precision,intent(out) :: EcRPA(nspin)
double precision,intent(out) :: EcBSE(nspin)
2020-04-22 00:39:52 +02:00
! Memory allocation
allocate(OmBSE(nS,nspin))
!-------------------
2020-01-14 21:27:34 +01:00
! Singlet manifold
2020-04-22 00:39:52 +02:00
!-------------------
2020-01-14 21:27:34 +01:00
if(singlet_manifold) then
ispin = 1
EcBSE(ispin) = 0d0
2020-03-13 09:18:18 +01:00
call linear_response(ispin,.true.,TDA,.false.,eta,nBas,nC,nO,nV,nR,nS,1d0,eW,ERI, &
2020-04-22 00:39:52 +02:00
rho(:,:,:,ispin),EcRPA(ispin),OmRPA(:,ispin),XpY(:,:,ispin),XmY(:,:,ispin))
2020-03-13 09:18:18 +01:00
call excitation_density(nBas,nC,nO,nR,nS,ERI,XpY(:,:,ispin),rho(:,:,:,ispin))
2020-01-14 21:27:34 +01:00
2020-04-22 00:39:52 +02:00
OmBSE(:,ispin) = OmRPA(:,ispin)
2020-01-23 21:22:41 +01:00
call linear_response(ispin,.true.,TDA,.true.,eta,nBas,nC,nO,nV,nR,nS,1d0,eGW,ERI, &
2020-04-22 00:39:52 +02:00
rho(:,:,:,ispin),EcBSE(ispin),OmBSE(:,ispin),XpY(:,:,ispin),XmY(:,:,ispin))
call print_excitation('BSE ',ispin,nS,OmBSE(:,ispin))
! Compute dynamic correction for BSE via perturbation theory
2020-04-22 09:55:58 +02:00
call Bethe_Salpeter_dynamic_perturbation(TDA,eta,nBas,nC,nO,nV,nR,nS,OmRPA(:,ispin),OmBSE(:,ispin), &
XpY(:,:,ispin),XmY(:,:,ispin),rho(:,:,:,ispin))
2020-01-14 21:27:34 +01:00
end if
2020-04-22 00:39:52 +02:00
!-------------------
2020-01-14 21:27:34 +01:00
! Triplet manifold
2020-04-22 00:39:52 +02:00
!-------------------
2020-01-14 21:27:34 +01:00
if(triplet_manifold) then
ispin = 2
EcBSE(ispin) = 0d0
2020-01-23 21:22:41 +01:00
call linear_response(ispin,.true.,TDA,.false.,eta,nBas,nC,nO,nV,nR,nS,1d0,eW,ERI, &
2020-04-22 00:39:52 +02:00
rho(:,:,:,ispin),EcRPA(ispin),OmRPA(:,ispin),XpY(:,:,ispin),XmY(:,:,ispin))
2020-01-14 21:27:34 +01:00
call excitation_density(nBas,nC,nO,nR,nS,ERI,XpY(:,:,ispin),rho(:,:,:,ispin))
2020-04-22 00:39:52 +02:00
OmBSE(:,ispin) = OmRPA(:,ispin)
2020-01-23 21:22:41 +01:00
call linear_response(ispin,.true.,TDA,.true.,eta,nBas,nC,nO,nV,nR,nS,1d0,eGW,ERI, &
2020-04-22 00:39:52 +02:00
rho(:,:,:,ispin),EcBSE(ispin),OmBSE(:,ispin),XpY(:,:,ispin),XmY(:,:,ispin))
call print_excitation('BSE ',ispin,nS,OmBSE(:,ispin))
! Compute dynamic correction for BSE via perturbation theory
2020-04-22 09:55:58 +02:00
call Bethe_Salpeter_dynamic_perturbation(TDA,eta,nBas,nC,nO,nV,nR,nS,OmRPA(:,ispin),OmBSE(:,ispin), &
XpY(:,:,ispin),XmY(:,:,ispin),rho(:,:,:,ispin))
2020-04-22 00:39:52 +02:00
2020-01-14 21:27:34 +01:00
end if
end subroutine Bethe_Salpeter