4
1
mirror of https://github.com/pfloos/quack synced 2024-11-08 15:13:53 +01:00
quack/src/QuAcK/TDHF.f90

97 lines
2.7 KiB
Fortran
Raw Normal View History

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'
! 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
! 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))
! Singlet manifold
if(singlet_manifold) then
ispin = 1
2019-10-05 22:06:25 +02:00
call linear_response(ispin,dRPA,TDA,BSE,nBas,nC,nO,nV,nR,nS,e,ERI,rho, &
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
2019-10-05 22:06:25 +02:00
call linear_response(ispin,dRPA,TDA,BSE,nBas,nC,nO,nV,nR,nS,e,ERI,rho, &
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(*,*)
2019-03-19 10:13:33 +01:00
end subroutine TDHF