4
1
mirror of https://github.com/pfloos/quack synced 2024-11-09 07:33:55 +01:00
quack/src/GF/BSE2.f90

130 lines
4.2 KiB
Fortran
Raw Normal View History

subroutine BSE2(TDA,dBSE,dTDA,evDyn,singlet,triplet,eta,nBas,nC,nO,nV,nR,nS,ERI,dipole_int,eHF,eGF,EcBSE)
2020-06-01 17:30:02 +02:00
! Compute the second-order Bethe-Salpeter excitation energies
2020-06-01 17:30:02 +02:00
implicit none
include 'parameters.h'
! Input variables
logical,intent(in) :: TDA
2020-06-15 23:04:07 +02:00
logical,intent(in) :: dBSE
logical,intent(in) :: dTDA
2020-06-16 14:02:14 +02:00
logical,intent(in) :: evDyn
logical,intent(in) :: singlet
logical,intent(in) :: triplet
2020-06-01 17:30:02 +02:00
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) :: eHF(nBas)
double precision,intent(in) :: eGF(nBas)
double precision,intent(in) :: ERI(nBas,nBas,nBas,nBas)
double precision,intent(in) :: dipole_int(nBas,nBas,ncart)
2020-06-01 17:30:02 +02:00
! Local variables
integer :: ispin
double precision,allocatable :: OmBSE(:,:)
double precision,allocatable :: XpY(:,:,:)
double precision,allocatable :: XmY(:,:,:)
double precision :: rho
double precision,allocatable :: A_sta(:,:,:)
double precision,allocatable :: B_sta(:,:,:)
2020-06-01 17:30:02 +02:00
! Output variables
double precision,intent(out) :: EcBSE(nspin)
! Memory allocation
allocate(OmBSE(nS,nspin),XpY(nS,nS,nspin),XmY(nS,nS,nspin), &
A_sta(nS,nS,nspin),B_sta(nS,nS,nspin))
2020-06-01 17:30:02 +02:00
!-------------------
! Singlet manifold
!-------------------
if(singlet) then
2020-06-01 17:30:02 +02:00
ispin = 1
EcBSE(ispin) = 0d0
! Compute static kernel
call BSE2_A_matrix_static(ispin,eta,nBas,nC,nO,nV,nR,nS,1d0,ERI,eGF,A_sta(:,:,ispin))
if(.not.TDA) call BSE2_B_matrix_static(ispin,eta,nBas,nC,nO,nV,nR,nS,1d0,ERI,eGF,B_sta(:,:,ispin))
2020-06-01 17:30:02 +02:00
! Compute BSE2 excitation energies
call linear_response_BSE(ispin,.false.,TDA,.true.,eta,nBas,nC,nO,nV,nR,nS,1d0,eGF,ERI,-A_sta(:,:,ispin),-B_sta(:,:,ispin), &
EcBSE(ispin),OmBSE(:,ispin),XpY(:,:,ispin),XmY(:,:,ispin))
2020-06-01 17:30:02 +02:00
call print_excitation('BSE2 ',ispin,nS,OmBSE(:,ispin))
call print_transition_vectors(.true.,nBas,nC,nO,nV,nR,nS,dipole_int,OmBSE(:,ispin),XpY(:,:,ispin),XmY(:,:,ispin))
2021-03-06 23:08:43 +01:00
2020-06-01 17:30:02 +02:00
! Compute dynamic correction for BSE via perturbation theory
2020-06-17 14:32:37 +02:00
if(dBSE) then
if(evDyn) then
call BSE2_dynamic_perturbation_iterative(dTDA,ispin,eta,nBas,nC,nO,nV,nR,nS,ERI,dipole_int,eHF,eGF, &
A_sta(:,:,ispin),B_sta(:,:,ispin),OmBSE(:,ispin),XpY(:,:,ispin),XmY(:,:,ispin))
2020-06-17 14:32:37 +02:00
else
call BSE2_dynamic_perturbation(dTDA,ispin,eta,nBas,nC,nO,nV,nR,nS,ERI,dipole_int,eHF,eGF, &
A_sta(:,:,ispin),B_sta(:,:,ispin),OmBSE(:,ispin),XpY(:,:,ispin),XmY(:,:,ispin))
2020-06-01 17:30:02 +02:00
2020-06-17 14:32:37 +02:00
end if
end if
2020-06-01 17:30:02 +02:00
end if
!-------------------
! Triplet manifold
!-------------------
if(triplet) then
2020-06-01 17:30:02 +02:00
ispin = 2
EcBSE(ispin) = 0d0
! Compute static kernel
call BSE2_A_matrix_static(ispin,eta,nBas,nC,nO,nV,nR,nS,1d0,ERI,eGF,A_sta(:,:,ispin))
if(.not.TDA) call BSE2_B_matrix_static(ispin,eta,nBas,nC,nO,nV,nR,nS,1d0,ERI,eGF,B_sta(:,:,ispin))
2020-06-01 17:30:02 +02:00
! Compute BSE2 excitation energies
2022-02-01 11:42:31 +01:00
call linear_response_BSE(ispin,.false.,TDA,.true.,eta,nBas,nC,nO,nV,nR,nS,1d0,eGF,ERI,-A_sta(:,:,ispin),-B_sta(:,:,ispin), &
EcBSE(ispin),OmBSE(:,ispin),XpY(:,:,ispin),XmY(:,:,ispin))
2020-06-01 17:30:02 +02:00
call print_excitation('BSE2 ',ispin,nS,OmBSE(:,ispin))
call print_transition_vectors(.false.,nBas,nC,nO,nV,nR,nS,dipole_int,OmBSE(:,ispin),XpY(:,:,ispin),XmY(:,:,ispin))
2020-06-01 17:30:02 +02:00
! Compute dynamic correction for BSE via perturbation theory
2020-06-17 14:32:37 +02:00
if(dBSE) then
if(evDyn) then
call BSE2_dynamic_perturbation_iterative(dTDA,ispin,eta,nBas,nC,nO,nV,nR,nS,ERI,dipole_int,eHF,eGF, &
A_sta(:,:,ispin),B_sta(:,:,ispin),OmBSE(:,ispin),XpY(:,:,ispin),XmY(:,:,ispin))
2020-06-17 14:32:37 +02:00
else
call BSE2_dynamic_perturbation(dTDA,ispin,eta,nBas,nC,nO,nV,nR,nS,ERI,dipole_int,eHF,eGF, &
A_sta(:,:,ispin),B_sta(:,:,ispin),OmBSE(:,ispin),XpY(:,:,ispin),XmY(:,:,ispin))
2020-06-17 14:32:37 +02:00
end if
end if
2020-06-01 17:30:02 +02:00
end if
end subroutine BSE2