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

139 lines
4.6 KiB
Fortran
Raw Normal View History

subroutine Bethe_Salpeter(TDA_W,TDA,dBSE,dTDA,evDyn,singlet,triplet,eta,nBas,nC,nO,nV,nR,nS,ERI,dipole_int,eW,eGW,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 21:20:01 +02:00
logical,intent(in) :: dBSE
2020-06-14 13:04:16 +02:00
logical,intent(in) :: dTDA
2020-06-14 21:20:01 +02:00
logical,intent(in) :: evDyn
2020-09-24 11:56:06 +02:00
logical,intent(in) :: singlet
logical,intent(in) :: triplet
2020-01-14 21:27:34 +01:00
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)
double precision,intent(in) :: dipole_int(nBas,nBas,ncart)
2020-01-14 21:27:34 +01:00
! Local variables
integer :: ispin
2020-06-14 13:04:16 +02:00
integer :: isp_W
2020-09-24 11:56:06 +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-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(:,:,:)
2020-01-14 21:27:34 +01:00
! Output variables
double precision,intent(out) :: EcBSE(nspin)
2020-04-22 00:39:52 +02:00
! Memory allocation
2020-09-24 11:56:06 +02:00
allocate(OmRPA(nS),XpY_RPA(nS,nS),XmY_RPA(nS,nS),rho_RPA(nBas,nBas,nS))
2020-06-01 11:35:17 +02:00
allocate(OmBSE(nS,nspin),XpY_BSE(nS,nS,nspin),XmY_BSE(nS,nS,nspin))
2020-04-22 00:39:52 +02:00
2020-09-24 11:56:06 +02:00
!---------------------------------
! Compute (singlet) RPA screening
!---------------------------------
isp_W = 1
EcRPA = 0d0
call linear_response(isp_W,.true.,TDA_W,.false.,eta,nBas,nC,nO,nV,nR,nS,1d0,eW,ERI,OmRPA, &
rho_RPA,EcRPA,OmRPA,XpY_RPA,XmY_RPA)
call excitation_density(nBas,nC,nO,nR,nS,ERI,XpY_RPA,rho_RPA)
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
2020-09-24 11:56:06 +02:00
if(singlet) then
2020-01-14 21:27:34 +01:00
ispin = 1
EcBSE(ispin) = 0d0
2020-04-23 23:13:15 +02:00
! Compute BSE excitation energies
2020-01-14 21:27:34 +01:00
2020-09-24 11:56:06 +02:00
call linear_response(ispin,.true.,TDA,.true.,eta,nBas,nC,nO,nV,nR,nS,1d0,eGW,ERI,OmRPA, &
rho_RPA,EcBSE(ispin),OmBSE(:,ispin),XpY_BSE(:,:,ispin),XmY_BSE(:,:,ispin))
call print_excitation('BSE@GW ',ispin,nS,OmBSE(:,ispin))
call print_transition_vectors(.true.,nBas,nC,nO,nV,nR,nS,dipole_int, &
OmBSE(:,ispin),XpY_BSE(:,:,ispin),XmY_BSE(:,:,ispin))
2020-06-01 11:35:17 +02:00
2020-06-14 21:20:01 +02:00
!-------------------------------------------------
2020-06-01 11:35:17 +02:00
! Compute the dynamical screening at the BSE level
2020-06-14 21:20:01 +02:00
!-------------------------------------------------
2020-06-01 11:35:17 +02:00
2020-06-14 21:20:01 +02:00
if(dBSE) then
2020-04-22 00:39:52 +02:00
2020-06-14 21:20:01 +02:00
! Compute dynamic correction for BSE via perturbation theory (iterative or renormalized)
if(evDyn) then
call Bethe_Salpeter_dynamic_perturbation_iterative(dTDA,eta,nBas,nC,nO,nV,nR,nS,eGW,dipole_int,OmRPA,rho_RPA, &
2020-09-24 11:56:06 +02:00
OmBSE(:,ispin),XpY_BSE(:,:,ispin),XmY_BSE(:,:,ispin))
2020-06-01 11:35:17 +02:00
else
2020-06-14 21:20:01 +02:00
2021-02-25 10:55:08 +01:00
call Bethe_Salpeter_dynamic_perturbation(dTDA,eta,nBas,nC,nO,nV,nR,nS,eW,eGW,dipole_int,OmRPA,rho_RPA, &
2020-09-24 11:56:06 +02:00
OmBSE(:,ispin),XpY_BSE(:,:,ispin),XmY_BSE(:,:,ispin))
2020-06-01 11:35:17 +02:00
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
2020-09-24 11:56:06 +02:00
if(triplet) then
2020-01-14 21:27:34 +01:00
ispin = 2
EcBSE(ispin) = 0d0
2020-04-23 23:13:15 +02:00
! Compute BSE excitation energies
2020-01-14 21:27:34 +01:00
2020-09-24 11:56:06 +02:00
call linear_response(ispin,.true.,TDA,.true.,eta,nBas,nC,nO,nV,nR,nS,1d0,eGW,ERI,OmRPA, &
rho_RPA,EcBSE(ispin),OmBSE(:,ispin),XpY_BSE(:,:,ispin),XmY_BSE(:,:,ispin))
call print_excitation('BSE@GW ',ispin,nS,OmBSE(:,ispin))
call print_transition_vectors(.false.,nBas,nC,nO,nV,nR,nS,dipole_int, &
OmBSE(:,ispin),XpY_BSE(:,:,ispin),XmY_BSE(:,:,ispin))
2020-06-01 11:35:17 +02:00
2020-06-14 21:20:01 +02:00
!-------------------------------------------------
2020-06-01 11:35:17 +02:00
! Compute the dynamical screening at the BSE level
2020-06-14 21:20:01 +02:00
!-------------------------------------------------
2020-06-01 11:35:17 +02:00
2020-06-14 21:20:01 +02:00
if(dBSE) then
2020-04-22 00:39:52 +02:00
2020-06-14 21:20:01 +02:00
! Compute dynamic correction for BSE via perturbation theory (iterative or renormalized)
2020-04-22 00:39:52 +02:00
2020-06-14 21:20:01 +02:00
if(evDyn) then
call Bethe_Salpeter_dynamic_perturbation_iterative(dTDA,eta,nBas,nC,nO,nV,nR,nS,eGW,dipole_int,OmRPA,rho_RPA, &
2020-09-24 11:56:06 +02:00
OmBSE(:,ispin),XpY_BSE(:,:,ispin),XmY_BSE(:,:,ispin))
2020-06-01 11:35:17 +02:00
else
2020-06-14 21:20:01 +02:00
2021-02-25 10:55:08 +01:00
call Bethe_Salpeter_dynamic_perturbation(dTDA,eta,nBas,nC,nO,nV,nR,nS,eW,eGW,dipole_int,OmRPA,rho_RPA, &
2020-09-24 11:56:06 +02:00
OmBSE(:,ispin),XpY_BSE(:,:,ispin),XmY_BSE(:,:,ispin))
2020-06-01 11:35:17 +02:00
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