10
1
mirror of https://github.com/pfloos/quack synced 2024-12-25 13:53:41 +01:00
QuAcK/src/GW/unrestricted_Bethe_Salpeter.f90

159 lines
5.9 KiB
Fortran
Raw Normal View History

2020-09-22 22:14:08 +02:00
subroutine unrestricted_Bethe_Salpeter(TDA_W,TDA,dBSE,dTDA,evDyn,spin_conserved,spin_flip,eta, &
2020-10-06 14:24:54 +02:00
nBas,nC,nO,nV,nR,nS,S,ERI_aaaa,ERI_aabb,ERI_bbbb, &
dipole_int_aa,dipole_int_bb,cW,eW,eGW,EcBSE)
2020-09-22 22:14:08 +02:00
! Compute the Bethe-Salpeter excitation energies
implicit none
include 'parameters.h'
! Input variables
logical,intent(in) :: TDA_W
logical,intent(in) :: TDA
logical,intent(in) :: dBSE
logical,intent(in) :: dTDA
logical,intent(in) :: evDyn
logical,intent(in) :: spin_conserved
logical,intent(in) :: spin_flip
double precision,intent(in) :: eta
integer,intent(in) :: nBas
integer,intent(in) :: nC(nspin)
integer,intent(in) :: nO(nspin)
integer,intent(in) :: nV(nspin)
integer,intent(in) :: nR(nspin)
2020-09-22 23:08:47 +02:00
integer,intent(in) :: nS(nspin)
2020-10-06 14:24:54 +02:00
double precision,intent(in) :: S(nBas,nBas)
double precision,intent(in) :: cW(nBas,nBas,nspin)
2020-09-22 22:14:08 +02:00
double precision,intent(in) :: eW(nBas,nspin)
double precision,intent(in) :: eGW(nBas,nspin)
double precision,intent(in) :: ERI_aaaa(nBas,nBas,nBas,nBas)
double precision,intent(in) :: ERI_aabb(nBas,nBas,nBas,nBas)
double precision,intent(in) :: ERI_bbbb(nBas,nBas,nBas,nBas)
2020-09-28 22:58:58 +02:00
double precision,intent(in) :: dipole_int_aa(nBas,nBas,ncart)
double precision,intent(in) :: dipole_int_bb(nBas,nBas,ncart)
2020-09-22 22:14:08 +02:00
! Local variables
integer :: ispin
integer :: isp_W
2020-09-23 22:48:54 +02:00
2020-09-22 23:08:47 +02:00
integer :: nS_aa,nS_bb,nS_sc
2020-09-24 14:39:37 +02:00
double precision :: EcRPA
double precision,allocatable :: OmRPA(:)
double precision,allocatable :: XpY_RPA(:,:)
double precision,allocatable :: XmY_RPA(:,:)
double precision,allocatable :: rho_RPA(:,:,:,:)
2020-09-22 23:08:47 +02:00
double precision,allocatable :: OmBSE_sc(:)
double precision,allocatable :: XpY_BSE_sc(:,:)
double precision,allocatable :: XmY_BSE_sc(:,:)
2020-09-22 22:14:08 +02:00
2020-09-23 22:48:54 +02:00
integer :: nS_ab,nS_ba,nS_sf
double precision,allocatable :: OmBSE_sf(:)
double precision,allocatable :: XpY_BSE_sf(:,:)
double precision,allocatable :: XmY_BSE_sf(:,:)
2020-09-22 22:14:08 +02:00
! Output variables
2020-09-22 23:08:47 +02:00
double precision,intent(out) :: EcBSE(nspin)
2020-09-22 22:14:08 +02:00
2020-09-23 22:48:54 +02:00
! Memory allocation
2020-09-22 23:08:47 +02:00
2020-09-23 22:48:54 +02:00
nS_aa = nS(1)
nS_bb = nS(2)
nS_sc = nS_aa + nS_bb
2020-10-07 22:51:30 +02:00
nS_ab = (nO(1) - nC(1))*(nV(2) - nR(2))
nS_ba = (nO(2) - nC(2))*(nV(1) - nR(1))
nS_sf = nS_ab + nS_ba
2020-09-22 23:08:47 +02:00
2020-09-24 14:39:37 +02:00
allocate(OmRPA(nS_sc),XpY_RPA(nS_sc,nS_sc),XmY_RPA(nS_sc,nS_sc),rho_RPA(nBas,nBas,nS_sc,nspin))
!--------------------------!
! Spin-conserved screening !
!--------------------------!
isp_W = 1
EcRPA = 0d0
2020-09-23 22:48:54 +02:00
! Compute spin-conserved RPA screening
2020-09-22 22:14:08 +02:00
2020-09-23 22:48:54 +02:00
call unrestricted_linear_response(isp_W,.true.,TDA_W,.false.,eta,nBas,nC,nO,nV,nR,nS_aa,nS_bb,nS_sc,nS_sc,1d0, &
2020-10-05 23:00:56 +02:00
eW,ERI_aaaa,ERI_aabb,ERI_bbbb,OmRPA,rho_RPA,EcRPA,OmRPA,XpY_RPA,XmY_RPA)
2020-09-24 14:39:37 +02:00
call unrestricted_excitation_density(nBas,nC,nO,nR,nS_aa,nS_bb,nS_sc,ERI_aaaa,ERI_aabb,ERI_bbbb,XpY_RPA,rho_RPA)
2020-09-22 22:14:08 +02:00
2020-09-24 14:39:37 +02:00
!----------------------------!
! Spin-conserved excitations !
!----------------------------!
2020-09-23 14:23:40 +02:00
2020-09-23 22:48:54 +02:00
if(spin_conserved) then
ispin = 1
EcBSE(ispin) = 0d0
allocate(OmBSE_sc(nS_sc),XpY_BSE_sc(nS_sc,nS_sc),XmY_BSE_sc(nS_sc,nS_sc))
2020-09-22 22:14:08 +02:00
2020-09-22 23:08:47 +02:00
! Compute spin-conserved BSE excitation energies
2020-09-22 22:14:08 +02:00
2020-09-23 22:48:54 +02:00
call unrestricted_linear_response(ispin,.true.,TDA,.true.,eta,nBas,nC,nO,nV,nR,nS_aa,nS_bb,nS_sc,nS_sc,1d0, &
2020-10-05 23:00:56 +02:00
eGW,ERI_aaaa,ERI_aabb,ERI_bbbb,OmRPA,rho_RPA,EcBSE(ispin), &
2020-09-22 23:46:33 +02:00
OmBSE_sc,XpY_BSE_sc,XmY_BSE_sc)
2020-10-06 14:24:54 +02:00
call print_excitation('BSE@UGW ',5,nS_sc,OmBSE_sc)
call print_unrestricted_transition_vectors(ispin,nBas,nC,nO,nV,nR,nS,nS_aa,nS_bb,nS_sc,dipole_int_aa,dipole_int_bb, &
cW,S,OmBSE_sc,XpY_BSE_sc,XmY_BSE_sc)
2020-09-22 22:14:08 +02:00
!-------------------------------------------------
! Compute the dynamical screening at the BSE level
!-------------------------------------------------
2020-09-30 00:10:57 +02:00
if(dBSE) &
call unrestricted_Bethe_Salpeter_dynamic_perturbation(ispin,dTDA,eta,nBas,nC,nO,nV,nR,nS,nS_aa,nS_bb,nS_sc,nS_sc, &
2021-02-25 10:55:08 +01:00
eW,eGW,ERI_aaaa,ERI_aabb,ERI_bbbb,dipole_int_aa,dipole_int_bb, &
2020-09-30 00:10:57 +02:00
OmRPA,rho_RPA,OmBSE_sc,XpY_BSE_sc,XmY_BSE_sc)
deallocate(OmBSE_sc,XpY_BSE_sc,XmY_BSE_sc)
2020-09-22 22:14:08 +02:00
end if
!-----------------------!
! Spin-flip excitations !
!-----------------------!
2020-09-23 22:48:54 +02:00
if(spin_flip) then
2020-09-22 22:14:08 +02:00
2020-09-23 22:48:54 +02:00
ispin = 2
EcBSE(ispin) = 0d0
2020-09-22 22:14:08 +02:00
2020-09-23 22:48:54 +02:00
! Memory allocation
allocate(OmBSE_sf(nS_sf),XpY_BSE_sf(nS_sf,nS_sf),XmY_BSE_sf(nS_sf,nS_sf))
2020-09-22 22:14:08 +02:00
2020-10-05 16:58:19 +02:00
! Compute spin-flip BSE excitation energies
2020-09-23 22:48:54 +02:00
call unrestricted_linear_response(ispin,.true.,TDA,.true.,eta,nBas,nC,nO,nV,nR,nS_ab,nS_ba,nS_sf,nS_sc,1d0, &
2020-10-05 23:00:56 +02:00
eGW,ERI_aaaa,ERI_aabb,ERI_bbbb,OmRPA,rho_RPA,EcBSE(ispin), &
2020-09-23 22:48:54 +02:00
OmBSE_sf,XpY_BSE_sf,XmY_BSE_sf)
2020-09-22 22:14:08 +02:00
2020-10-06 14:24:54 +02:00
call print_excitation('BSE@UGW ',6,nS_sf,OmBSE_sf)
call print_unrestricted_transition_vectors(ispin,nBas,nC,nO,nV,nR,nS,nS_ab,nS_ba,nS_sf,dipole_int_aa,dipole_int_bb, &
cW,S,OmBSE_sf,XpY_BSE_sf,XmY_BSE_sf)
2020-09-22 22:14:08 +02:00
!-------------------------------------------------
! Compute the dynamical screening at the BSE level
!-------------------------------------------------
2020-09-30 00:10:57 +02:00
if(dBSE) &
2020-09-30 15:02:48 +02:00
call unrestricted_Bethe_Salpeter_dynamic_perturbation(ispin,dTDA,eta,nBas,nC,nO,nV,nR,nS,nS_ab,nS_ba,nS_sf,nS_sc, &
2021-02-25 10:55:08 +01:00
eW,eGW,ERI_aaaa,ERI_aabb,ERI_bbbb,dipole_int_aa,dipole_int_bb, &
2020-09-30 00:10:57 +02:00
OmRPA,rho_RPA,OmBSE_sf,XpY_BSE_sf,XmY_BSE_sf)
2020-09-22 22:14:08 +02:00
2020-09-30 00:10:57 +02:00
deallocate(OmBSE_sf,XpY_BSE_sf,XmY_BSE_sf)
2020-09-23 22:48:54 +02:00
end if
2020-09-22 22:14:08 +02:00
end subroutine unrestricted_Bethe_Salpeter