4
1
mirror of https://github.com/pfloos/quack synced 2024-10-06 00:06:25 +02:00
quack/src/RPA/phACFDT.f90

182 lines
5.7 KiB
Fortran
Raw Normal View History

2023-07-12 23:16:37 +02:00
subroutine phACFDT(exchange_kernel,doXBS,dRPA,TDA_W,TDA,BSE,singlet,triplet,eta,nBas,nC,nO,nV,nR,nS,ERI,eW,e,EcAC)
2020-01-17 15:31:38 +01:00
! Compute the correlation energy via the adiabatic connection fluctuation dissipation theorem
implicit none
include 'parameters.h'
include 'quadrature.h'
! Input variables
logical,intent(in) :: doXBS
logical,intent(in) :: exchange_kernel
logical,intent(in) :: dRPA
2020-06-09 21:24:37 +02:00
logical,intent(in) :: TDA_W
2020-01-17 15:31:38 +01:00
logical,intent(in) :: TDA
logical,intent(in) :: BSE
2020-09-24 11:56:06 +02:00
logical,intent(in) :: singlet
logical,intent(in) :: triplet
2020-01-17 15:31:38 +01:00
2020-01-23 21:22:41 +01:00
double precision,intent(in) :: eta
2023-07-17 15:17:02 +02:00
integer,intent(in) :: nBas
integer,intent(in) :: nC
integer,intent(in) :: nO
integer,intent(in) :: nV
integer,intent(in) :: nR
integer,intent(in) :: nS
2020-03-13 10:25:20 +01:00
double precision,intent(in) :: eW(nBas)
2020-01-17 15:31:38 +01:00
double precision,intent(in) :: e(nBas)
double precision,intent(in) :: ERI(nBas,nBas,nBas,nBas)
! Local variables
integer :: ispin
2020-06-14 13:04:16 +02:00
integer :: isp_W
2020-01-17 15:31:38 +01:00
integer :: iAC
double precision :: lambda
double precision,allocatable :: Ec(:,:)
2020-09-24 11:56:06 +02:00
double precision :: EcRPA
2023-07-17 15:17:02 +02:00
double precision,allocatable :: KA(:,:)
double precision,allocatable :: KB(:,:)
2020-09-24 11:56:06 +02:00
double precision,allocatable :: OmRPA(:)
double precision,allocatable :: XpY_RPA(:,:)
double precision,allocatable :: XmY_RPA(:,:)
double precision,allocatable :: rho_RPA(:,:,:)
2023-07-17 15:17:02 +02:00
double precision,allocatable :: Om(:)
double precision,allocatable :: XpY(:,:)
double precision,allocatable :: XmY(:,:)
2020-09-24 11:56:06 +02:00
2020-01-17 15:31:38 +01:00
! Output variables
double precision,intent(out) :: EcAC(nspin)
! Memory allocation
allocate(Ec(nAC,nspin))
2023-07-17 15:17:02 +02:00
allocate(KA(nS,nS),KB(nS,nS),OmRPA(nS),XpY_RPA(nS,nS),XmY_RPA(nS,nS),rho_RPA(nBas,nBas,nS))
allocate(Om(nS),XpY(nS,nS),XmY(nS,nS))
2020-01-17 15:31:38 +01:00
! Antisymmetrized kernel version
if(exchange_kernel) then
write(*,*)
write(*,*) '*** Exchange kernel version ***'
write(*,*)
end if
2020-01-28 09:00:11 +01:00
EcAC(:) = 0d0
Ec(:,:) = 0d0
2020-09-24 11:56:06 +02:00
! Compute (singlet) RPA screening
isp_W = 1
EcRPA = 0d0
call phLR(isp_W,.true.,TDA_W,eta,nBas,nC,nO,nV,nR,nS,1d0,eW,ERI,EcRPA,OmRPA,XpY_RPA,XmY_RPA)
2023-07-04 10:51:03 +02:00
call GW_excitation_density(nBas,nC,nO,nR,nS,ERI,XpY_RPA,rho_RPA)
2020-09-24 11:56:06 +02:00
2023-07-17 15:17:02 +02:00
call BSE_static_kernel_KA(eta,nBas,nC,nO,nV,nR,nS,1d0,ERI,OmRPA,rho_RPA,KA)
call BSE_static_kernel_KB(eta,nBas,nC,nO,nV,nR,nS,1d0,ERI,OmRPA,rho_RPA,KB)
2020-01-17 15:31:38 +01:00
! Singlet manifold
2020-09-24 11:56:06 +02:00
if(singlet) then
2020-01-17 15:31:38 +01:00
ispin = 1
write(*,*) '--------------'
write(*,*) 'Singlet states'
write(*,*) '--------------'
write(*,*)
write(*,*) '-----------------------------------------------------------------------------------'
write(*,'(2X,A15,1X,A30,1X,A30)') 'lambda','Ec(lambda)','Tr(K x P_lambda)'
write(*,*) '-----------------------------------------------------------------------------------'
do iAC=1,nAC
lambda = rAC(iAC)
if(doXBS) then
call phLR(isp_W,.true.,TDA_W,eta,nBas,nC,nO,nV,nR,nS,lambda,eW,ERI,EcRPA,OmRPA,XpY_RPA,XmY_RPA)
2023-07-04 10:51:03 +02:00
call GW_excitation_density(nBas,nC,nO,nR,nS,ERI,XpY_RPA,rho_RPA)
2020-01-17 15:31:38 +01:00
2023-07-17 15:17:02 +02:00
call BSE_static_kernel_KA(eta,nBas,nC,nO,nV,nR,nS,lambda,ERI,OmRPA,rho_RPA,KA)
call BSE_static_kernel_KB(eta,nBas,nC,nO,nV,nR,nS,lambda,ERI,OmRPA,rho_RPA,KB)
2020-01-17 15:31:38 +01:00
end if
2023-07-17 15:17:02 +02:00
call linear_response_BSE(ispin,dRPA,TDA,BSE,eta,nBas,nC,nO,nV,nR,nS,lambda,e,ERI,KA,KB,EcAC(ispin),Om,XpY,XmY)
2020-01-17 15:31:38 +01:00
2023-07-17 15:17:02 +02:00
call phACFDT_correlation_energy(ispin,exchange_kernel,nBas,nC,nO,nV,nR,nS,ERI,XpY,XmY,Ec(iAC,ispin))
2020-01-17 15:31:38 +01:00
write(*,'(2X,F15.6,1X,F30.15,1X,F30.15)') lambda,EcAC(ispin),Ec(iAC,ispin)
end do
EcAC(ispin) = 0.5d0*dot_product(wAC,Ec(:,ispin))
2020-09-24 11:56:06 +02:00
if(exchange_kernel) EcAC(ispin) = 0.5d0*EcAC(ispin)
2020-01-17 15:31:38 +01:00
write(*,*) '-----------------------------------------------------------------------------------'
write(*,'(2X,A50,1X,F15.6)') ' Ec(AC) via Gauss-Legendre quadrature:',EcAC(ispin)
write(*,*) '-----------------------------------------------------------------------------------'
write(*,*)
end if
! Triplet manifold
2020-09-24 11:56:06 +02:00
if(triplet) then
2020-01-17 15:31:38 +01:00
ispin = 2
write(*,*) '--------------'
write(*,*) 'Triplet states'
write(*,*) '--------------'
write(*,*)
write(*,*) '-----------------------------------------------------------------------------------'
write(*,'(2X,A15,1X,A30,1X,A30)') 'lambda','Ec(lambda)','Tr(K x P_lambda)'
write(*,*) '-----------------------------------------------------------------------------------'
do iAC=1,nAC
lambda = rAC(iAC)
2020-01-28 09:00:11 +01:00
if(doXBS) then
2020-01-17 15:31:38 +01:00
call phLR(isp_W,.true.,TDA_W,eta,nBas,nC,nO,nV,nR,nS,lambda,eW,ERI,EcRPA,OmRPA,XpY_RPA,XmY_RPA)
2023-07-04 10:51:03 +02:00
call GW_excitation_density(nBas,nC,nO,nR,nS,ERI,XpY_RPA,rho_RPA)
2020-01-17 15:31:38 +01:00
2023-07-17 15:17:02 +02:00
call BSE_static_kernel_KA(eta,nBas,nC,nO,nV,nR,nS,lambda,ERI,OmRPA,rho_RPA,KA)
call BSE_static_kernel_KB(eta,nBas,nC,nO,nV,nR,nS,lambda,ERI,OmRPA,rho_RPA,KB)
2020-01-28 09:00:11 +01:00
end if
2020-01-17 15:31:38 +01:00
2023-07-17 15:17:02 +02:00
call linear_response_BSE(ispin,dRPA,TDA,BSE,eta,nBas,nC,nO,nV,nR,nS,lambda,e,ERI,KA,KB,EcAC(ispin),Om,XpY,XmY)
2020-01-17 15:31:38 +01:00
2023-07-17 15:17:02 +02:00
call phACFDT_correlation_energy(ispin,exchange_kernel,nBas,nC,nO,nV,nR,nS,ERI,XpY,XmY,Ec(iAC,ispin))
2020-01-17 15:31:38 +01:00
2020-01-28 09:00:11 +01:00
write(*,'(2X,F15.6,1X,F30.15,1X,F30.15)') lambda,EcAC(ispin),Ec(iAC,ispin)
2020-01-17 15:31:38 +01:00
end do
EcAC(ispin) = 0.5d0*dot_product(wAC,Ec(:,ispin))
2020-09-24 11:56:06 +02:00
if(exchange_kernel) EcAC(ispin) = 1.5d0*EcAC(ispin)
2020-01-17 15:31:38 +01:00
write(*,*) '-----------------------------------------------------------------------------------'
write(*,'(2X,A50,1X,F15.6)') ' Ec(AC) via Gauss-Legendre quadrature:',EcAC(ispin)
write(*,*) '-----------------------------------------------------------------------------------'
write(*,*)
end if
end subroutine