4
1
mirror of https://github.com/pfloos/quack synced 2024-11-06 22:24:03 +01:00
quack/src/QuAcK/Bethe_Salpeter.f90

146 lines
5.3 KiB
Fortran
Raw Normal View History

2020-06-14 13:04:16 +02:00
subroutine Bethe_Salpeter(TDA_W,TDA,dTDA,singlet_manifold,triplet_manifold,eta, &
2020-06-01 17:14:42 +02:00
nBas,nC,nO,nV,nR,nS,ERI,eW,eGW,OmRPA,XpY_RPA,XmY_RPA,rho_RPA,EcRPA,EcBSE)
2020-01-14 21:27:34 +01:00
! Compute the Bethe-Salpeter excitation energies
implicit none
include 'parameters.h'
! Input variables
2020-06-09 21:24:37 +02:00
logical,intent(in) :: TDA_W
2020-01-14 21:27:34 +01:00
logical,intent(in) :: TDA
2020-06-14 13:04:16 +02:00
logical,intent(in) :: dTDA
2020-01-14 21:27:34 +01:00
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-04-23 23:13:15 +02:00
double precision :: XpY_RPA(nS,nS,nspin)
double precision :: XmY_RPA(nS,nS,nspin)
double precision :: rho_RPA(nBas,nBas,nS,nspin)
2020-01-14 21:27:34 +01:00
! Local variables
2020-06-14 15:09:33 +02:00
logical :: evDyn = .true.
2020-06-01 11:35:17 +02:00
logical :: W_BSE = .false.
2020-01-14 21:27:34 +01:00
integer :: ispin
2020-06-14 13:04:16 +02:00
integer :: isp_W
2020-04-22 00:39:52 +02:00
double precision,allocatable :: OmBSE(:,:)
2020-04-23 23:13:15 +02:00
double precision,allocatable :: XpY_BSE(:,:,:)
double precision,allocatable :: XmY_BSE(:,:,:)
double precision,allocatable :: rho_BSE(:,:,:,:)
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
2020-06-01 11:35:17 +02:00
allocate(OmBSE(nS,nspin),XpY_BSE(nS,nS,nspin),XmY_BSE(nS,nS,nspin))
if(W_BSE) allocate(rho_BSE(nBas,nBas,nS,nspin))
2020-04-22 00:39:52 +02:00
!-------------------
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
2020-06-14 13:04:16 +02:00
isp_W = 1
2020-01-14 21:27:34 +01:00
EcBSE(ispin) = 0d0
2020-06-09 21:24:37 +02:00
! Compute (singlet) RPA screening
2020-04-23 23:13:15 +02:00
2020-06-14 13:04:16 +02:00
call linear_response(isp_W,.true.,TDA_W,.false.,eta,nBas,nC,nO,nV,nR,nS,1d0,eW,ERI, &
2020-04-23 23:13:15 +02:00
rho_RPA(:,:,:,ispin),EcRPA(ispin),OmRPA(:,ispin),XpY_RPA(:,:,ispin),XmY_RPA(:,:,ispin))
call excitation_density(nBas,nC,nO,nR,nS,ERI,XpY_RPA(:,:,ispin),rho_RPA(:,:,:,ispin))
! Compute BSE excitation energies
2020-01-14 21:27:34 +01:00
2020-04-22 00:39:52 +02:00
OmBSE(:,ispin) = OmRPA(:,ispin)
2020-04-23 23:13:15 +02:00
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-23 23:13:15 +02:00
rho_RPA(:,:,:,ispin),EcBSE(ispin),OmBSE(:,ispin),XpY_BSE(:,:,ispin),XmY_BSE(:,:,ispin))
2020-04-22 00:39:52 +02:00
call print_excitation('BSE ',ispin,nS,OmBSE(:,ispin))
2020-06-01 11:35:17 +02:00
! Compute the dynamical screening at the BSE level
if(W_BSE) call excitation_density(nBas,nC,nO,nR,nS,ERI,XpY_BSE(:,:,ispin),rho_BSE(:,:,:,ispin))
2020-04-22 00:39:52 +02:00
! Compute dynamic correction for BSE via perturbation theory
2020-05-21 16:13:52 +02:00
if(evDyn) then
2020-01-14 21:27:34 +01:00
2020-06-14 13:04:16 +02:00
call Bethe_Salpeter_dynamic_perturbation_iterative(TDA,dTDA,eta,nBas,nC,nO,nV,nR,nS,eGW(:),OmRPA(:,ispin),OmBSE(:,ispin), &
2020-05-21 16:13:52 +02:00
XpY_BSE(:,:,ispin),XmY_BSE(:,:,ispin),rho_RPA(:,:,:,ispin))
else
2020-06-01 11:35:17 +02:00
if(W_BSE) then
2020-06-14 13:04:16 +02:00
call Bethe_Salpeter_dynamic_perturbation(TDA,dTDA,eta,nBas,nC,nO,nV,nR,nS,eGW(:),OmRPA(:,ispin),OmBSE(:,ispin), &
2020-06-01 11:35:17 +02:00
XpY_BSE(:,:,ispin),XmY_BSE(:,:,ispin),rho_BSE(:,:,:,ispin))
else
2020-06-14 13:04:16 +02:00
call Bethe_Salpeter_dynamic_perturbation(TDA,dTDA,eta,nBas,nC,nO,nV,nR,nS,eGW(:),OmRPA(:,ispin),OmBSE(:,ispin), &
2020-06-01 11:35:17 +02:00
XpY_BSE(:,:,ispin),XmY_BSE(:,:,ispin),rho_RPA(:,:,:,ispin))
end if
2020-05-21 16:13:52 +02:00
end if
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
2020-06-14 13:04:16 +02:00
isp_W = 1
2020-01-14 21:27:34 +01:00
EcBSE(ispin) = 0d0
2020-06-09 21:24:37 +02:00
! Compute (singlet) RPA screening
2020-04-23 23:13:15 +02:00
2020-06-14 13:04:16 +02:00
call linear_response(isp_W,.true.,TDA_W,.false.,eta,nBas,nC,nO,nV,nR,nS,1d0,eW,ERI, &
2020-04-23 23:13:15 +02:00
rho_RPA(:,:,:,ispin),EcRPA(ispin),OmRPA(:,ispin),XpY_RPA(:,:,ispin),XmY_RPA(:,:,ispin))
call excitation_density(nBas,nC,nO,nR,nS,ERI,XpY_RPA(:,:,ispin),rho_RPA(:,:,:,ispin))
! Compute BSE excitation energies
2020-01-14 21:27:34 +01:00
2020-04-22 00:39:52 +02:00
OmBSE(:,ispin) = OmRPA(:,ispin)
2020-04-23 23:13:15 +02:00
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-23 23:13:15 +02:00
rho_RPA(:,:,:,ispin),EcBSE(ispin),OmBSE(:,ispin),XpY_BSE(:,:,ispin),XmY_BSE(:,:,ispin))
2020-04-22 00:39:52 +02:00
call print_excitation('BSE ',ispin,nS,OmBSE(:,ispin))
2020-06-01 11:35:17 +02:00
! Compute the dynamical screening at the BSE level
if(W_BSE) call excitation_density(nBas,nC,nO,nR,nS,ERI,XpY_BSE(:,:,ispin),rho_BSE(:,:,:,ispin))
2020-04-22 00:39:52 +02:00
! Compute dynamic correction for BSE via perturbation theory
2020-05-21 16:13:52 +02:00
if(evDyn) then
2020-06-14 13:04:16 +02:00
call Bethe_Salpeter_dynamic_perturbation_iterative(TDA,dTDA,eta,nBas,nC,nO,nV,nR,nS,eGW(:),OmRPA(:,ispin),OmBSE(:,ispin), &
2020-05-21 16:13:52 +02:00
XpY_BSE(:,:,ispin),XmY_BSE(:,:,ispin),rho_RPA(:,:,:,ispin))
else
2020-04-22 00:39:52 +02:00
2020-06-01 11:35:17 +02:00
if(W_BSE) then
2020-06-14 13:04:16 +02:00
call Bethe_Salpeter_dynamic_perturbation(TDA,dTDA,eta,nBas,nC,nO,nV,nR,nS,eGW(:),OmRPA(:,ispin),OmBSE(:,ispin), &
2020-06-01 11:35:17 +02:00
XpY_BSE(:,:,ispin),XmY_BSE(:,:,ispin),rho_BSE(:,:,:,ispin))
else
2020-06-14 13:04:16 +02:00
call Bethe_Salpeter_dynamic_perturbation(TDA,dTDA,eta,nBas,nC,nO,nV,nR,nS,eGW(:),OmRPA(:,ispin),OmBSE(:,ispin), &
2020-06-01 11:35:17 +02:00
XpY_BSE(:,:,ispin),XmY_BSE(:,:,ispin),rho_RPA(:,:,:,ispin))
end if
2020-05-21 16:13:52 +02:00
end if
2020-01-14 21:27:34 +01:00
end if
end subroutine Bethe_Salpeter