2019-05-23 09:53:23 +02:00
|
|
|
subroutine TDHF(singlet_manifold,triplet_manifold,nBas,nC,nO,nV,nR,nS,ENuc,ERHF,ERI,e)
|
2019-03-19 10:13:33 +01:00
|
|
|
|
|
|
|
! Perform random phase approximation calculation
|
|
|
|
|
|
|
|
implicit none
|
|
|
|
include 'parameters.h'
|
2020-01-13 23:08:03 +01:00
|
|
|
include 'quadrature.h'
|
2019-03-19 10:13:33 +01:00
|
|
|
|
|
|
|
! Input variables
|
|
|
|
|
2019-05-23 09:53:23 +02:00
|
|
|
logical,intent(in) :: singlet_manifold
|
|
|
|
logical,intent(in) :: triplet_manifold
|
|
|
|
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) :: ENuc
|
|
|
|
double precision,intent(in) :: ERHF
|
|
|
|
double precision,intent(in) :: e(nBas)
|
|
|
|
double precision,intent(in) :: ERI(nBas,nBas,nBas,nBas)
|
2019-03-19 10:13:33 +01:00
|
|
|
|
|
|
|
! Local variables
|
|
|
|
|
2019-05-23 09:53:23 +02:00
|
|
|
logical :: dRPA
|
|
|
|
logical :: TDA
|
|
|
|
logical :: BSE
|
2019-03-19 10:13:33 +01:00
|
|
|
integer :: ispin
|
2019-05-23 09:53:23 +02:00
|
|
|
double precision,allocatable :: Omega(:,:)
|
|
|
|
double precision,allocatable :: XpY(:,:,:)
|
2019-03-19 10:13:33 +01:00
|
|
|
|
|
|
|
double precision :: rho
|
2019-05-23 09:53:23 +02:00
|
|
|
double precision :: EcRPA(nspin)
|
2019-03-19 10:13:33 +01:00
|
|
|
|
2020-01-13 23:08:03 +01:00
|
|
|
logical :: AC
|
|
|
|
integer :: iAC
|
|
|
|
double precision :: lambda
|
|
|
|
double precision,allocatable :: EcACRPA(:,:)
|
|
|
|
double precision,allocatable :: EcAC(:,:)
|
|
|
|
|
2019-03-19 10:13:33 +01:00
|
|
|
! Hello world
|
|
|
|
|
|
|
|
write(*,*)
|
|
|
|
write(*,*)'************************************************'
|
|
|
|
write(*,*)'| Time-dependent Hartree-Fock calculation |'
|
|
|
|
write(*,*)'************************************************'
|
|
|
|
write(*,*)
|
|
|
|
|
2019-05-23 09:53:23 +02:00
|
|
|
! Initialization
|
|
|
|
|
|
|
|
EcRPA(:) = 0d0
|
|
|
|
|
2019-03-19 10:13:33 +01:00
|
|
|
! Switch on exchange for TDHF
|
|
|
|
|
|
|
|
dRPA = .false.
|
|
|
|
|
|
|
|
! Switch off Tamm-Dancoff approximation for TDHF
|
|
|
|
|
|
|
|
TDA = .false.
|
|
|
|
|
|
|
|
! Switch off Bethe-Salpeter equation for TDHF
|
|
|
|
|
|
|
|
BSE = .false.
|
|
|
|
|
|
|
|
! Memory allocation
|
|
|
|
|
|
|
|
allocate(Omega(nS,nspin),XpY(nS,nS,nspin))
|
|
|
|
|
2020-01-13 23:08:03 +01:00
|
|
|
AC = .true.
|
|
|
|
allocate(EcACRPA(nAC,nspin),EcAC(nAC,nspin))
|
|
|
|
|
2019-03-19 10:13:33 +01:00
|
|
|
! Singlet manifold
|
|
|
|
|
|
|
|
if(singlet_manifold) then
|
|
|
|
|
|
|
|
ispin = 1
|
|
|
|
|
2020-01-13 23:08:03 +01:00
|
|
|
call linear_response(ispin,dRPA,TDA,BSE,nBas,nC,nO,nV,nR,nS,1d0,e,ERI,rho, &
|
2019-10-05 22:06:25 +02:00
|
|
|
EcRPA(ispin),Omega(:,ispin),XpY(:,:,ispin))
|
2019-03-19 10:13:33 +01:00
|
|
|
call print_excitation('TDHF ',ispin,nS,Omega(:,ispin))
|
|
|
|
|
|
|
|
endif
|
|
|
|
|
|
|
|
! Triplet manifold
|
|
|
|
|
|
|
|
if(triplet_manifold) then
|
|
|
|
|
|
|
|
ispin = 2
|
|
|
|
|
2020-01-13 23:08:03 +01:00
|
|
|
call linear_response(ispin,dRPA,TDA,BSE,nBas,nC,nO,nV,nR,nS,1d0,e,ERI,rho, &
|
2019-10-05 22:06:25 +02:00
|
|
|
EcRPA(ispin),Omega(:,ispin),XpY(:,:,ispin))
|
2019-03-19 10:13:33 +01:00
|
|
|
call print_excitation('TDHF ',ispin,nS,Omega(:,ispin))
|
|
|
|
|
|
|
|
endif
|
|
|
|
|
2019-05-23 09:53:23 +02:00
|
|
|
write(*,*)
|
|
|
|
write(*,*)'-------------------------------------------------------------------------------'
|
|
|
|
write(*,'(2X,A40,F15.6)') 'RPA@TDHF correlation energy (singlet) =',EcRPA(1)
|
|
|
|
write(*,'(2X,A40,F15.6)') 'RPA@TDHF correlation energy (triplet) =',EcRPA(2)
|
|
|
|
write(*,'(2X,A40,F15.6)') 'RPA@TDHF correlation energy =',EcRPA(1) + EcRPA(2)
|
|
|
|
write(*,'(2X,A40,F15.6)') 'RPA@TDHF total energy =',ENuc + ERHF + EcRPA(1) + EcRPA(2)
|
|
|
|
write(*,*)'-------------------------------------------------------------------------------'
|
|
|
|
write(*,*)
|
|
|
|
|
2020-01-13 23:08:03 +01:00
|
|
|
! Compute the correlation energy via the adiabatic connection
|
|
|
|
|
|
|
|
if(AC) then
|
|
|
|
|
|
|
|
write(*,*) '------------------------------------------------------'
|
|
|
|
write(*,*) 'Adiabatic connection version of RPA correlation energy'
|
|
|
|
write(*,*) '------------------------------------------------------'
|
|
|
|
write(*,*)
|
|
|
|
|
|
|
|
if(singlet_manifold) then
|
|
|
|
|
|
|
|
ispin = 1
|
|
|
|
EcACRPA(:,ispin) = 0d0
|
|
|
|
|
|
|
|
write(*,*) '--------------'
|
|
|
|
write(*,*) 'Singlet states'
|
|
|
|
write(*,*) '--------------'
|
|
|
|
write(*,*)
|
|
|
|
|
|
|
|
write(*,*) '-----------------------------------------------------------------------------------'
|
|
|
|
write(*,'(2X,A15,1X,A30,1X,A30)') 'lambda','EcRPA(lambda)','Tr(V x P_lambda)'
|
|
|
|
write(*,*) '-----------------------------------------------------------------------------------'
|
|
|
|
|
|
|
|
do iAC=1,nAC
|
|
|
|
|
|
|
|
lambda = rAC(iAC)
|
|
|
|
|
|
|
|
! call linear_response(ispin,dRPA,TDA,.false.,nBas,nC,nO,nV,nR,nS,lambda,eHF,ERI, &
|
|
|
|
! rho(:,:,:,ispin),EcRPA(ispin),Omega(:,ispin),XpY(:,:,ispin))
|
|
|
|
! call excitation_density(nBas,nC,nO,nR,nS,ERI,XpY(:,:,ispin),rho(:,:,:,ispin))
|
|
|
|
|
|
|
|
call linear_response(ispin,dRPA,TDA,BSE,nBas,nC,nO,nV,nR,nS,lambda,e,ERI, &
|
|
|
|
rho,EcACRPA(iAC,ispin),Omega(:,ispin),XpY(:,:,ispin))
|
|
|
|
|
|
|
|
call Ec_AC(ispin,dRPA,nBas,nC,nO,nV,nR,nS,ERI,XpY(:,:,ispin),EcAC(iAC,ispin))
|
|
|
|
|
|
|
|
write(*,'(2X,F15.6,1X,F30.15,1X,F30.15)') lambda,EcACRPA(iAC,ispin),EcAC(iAC,ispin)
|
|
|
|
|
|
|
|
end do
|
|
|
|
|
|
|
|
write(*,*) '-----------------------------------------------------------------------------------'
|
|
|
|
write(*,'(2X,A50,1X,F15.6)') ' Ec(RPA) via Gauss-Legendre quadrature:',0.5d0*dot_product(wAC,EcAC(:,ispin))
|
|
|
|
write(*,*) '-----------------------------------------------------------------------------------'
|
|
|
|
write(*,*)
|
|
|
|
|
|
|
|
end if
|
|
|
|
|
|
|
|
if(triplet_manifold) then
|
|
|
|
|
|
|
|
ispin = 2
|
|
|
|
EcACRPA(:,ispin) = 0d0
|
|
|
|
|
|
|
|
write(*,*) '--------------'
|
|
|
|
write(*,*) 'Triplet states'
|
|
|
|
write(*,*) '--------------'
|
|
|
|
write(*,*)
|
|
|
|
|
|
|
|
write(*,*) '-----------------------------------------------------------------------------------'
|
|
|
|
write(*,'(2X,A15,1X,A30,1X,A30)') 'lambda','EcRPA(lambda)','Tr(V x P_lambda)'
|
|
|
|
write(*,*) '-----------------------------------------------------------------------------------'
|
|
|
|
|
|
|
|
do iAC=1,nAC
|
|
|
|
|
|
|
|
lambda = rAC(iAC)
|
|
|
|
|
|
|
|
! call linear_response(ispin,dRPA,TDA,.false.,nBas,nC,nO,nV,nR,nS,lambda,eHF,ERI, &
|
|
|
|
! rho(:,:,:,ispin),EcRPA(ispin),Omega(:,ispin),XpY(:,:,ispin))
|
|
|
|
! call excitation_density(nBas,nC,nO,nR,nS,ERI,XpY(:,:,ispin),rho(:,:,:,ispin))
|
|
|
|
|
|
|
|
call linear_response(ispin,dRPA,TDA,BSE,nBas,nC,nO,nV,nR,nS,lambda,e,ERI, &
|
|
|
|
rho,EcACRPA(iAC,ispin),Omega(:,ispin),XpY(:,:,ispin))
|
|
|
|
|
|
|
|
call Ec_AC(ispin,dRPA,nBas,nC,nO,nV,nR,nS,ERI,XpY(:,:,ispin),EcAC(iAC,ispin))
|
|
|
|
|
|
|
|
write(*,'(2X,F15.6,1X,F30.15,1X,F30.15)') lambda,EcACRPA(iAC,ispin),EcAC(iAC,ispin)
|
|
|
|
|
|
|
|
end do
|
|
|
|
|
|
|
|
write(*,*) '-----------------------------------------------------------------------------------'
|
|
|
|
write(*,'(2X,A50,1X,F15.6)') ' Ec(RPA) via Gauss-Legendre quadrature:',0.5d0*dot_product(wAC,EcAC(:,ispin))
|
|
|
|
write(*,*) '-----------------------------------------------------------------------------------'
|
|
|
|
write(*,*)
|
|
|
|
|
|
|
|
end if
|
|
|
|
|
|
|
|
end if
|
|
|
|
|
2019-03-19 10:13:33 +01:00
|
|
|
end subroutine TDHF
|