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

132 lines
4.1 KiB
Fortran
Raw Normal View History

subroutine phRPAx(TDA,doACFDT,exchange_kernel,singlet,triplet,eta,nBas,nC,nO,nV,nR,nS,ENuc,ERHF,ERI,dipole_int,eHF)
2020-01-28 09:00:11 +01:00
! Perform random phase approximation calculation with exchange (aka TDHF)
implicit none
include 'parameters.h'
include 'quadrature.h'
! Input variables
2020-09-24 16:39:15 +02:00
logical,intent(in) :: TDA
2020-01-28 09:00:11 +01:00
logical,intent(in) :: doACFDT
logical,intent(in) :: exchange_kernel
2020-09-23 09:46:44 +02:00
logical,intent(in) :: singlet
logical,intent(in) :: triplet
2021-11-10 22:41:40 +01:00
double precision,intent(in) :: eta
2020-01-28 09:00:11 +01: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
double precision,intent(in) :: ENuc
double precision,intent(in) :: ERHF
2020-09-24 11:56:06 +02:00
double precision,intent(in) :: eHF(nBas)
2020-01-28 09:00:11 +01:00
double precision,intent(in) :: ERI(nBas,nBas,nBas,nBas)
double precision,intent(in) :: dipole_int(nBas,nBas,ncart)
2020-01-28 09:00:11 +01:00
! Local variables
integer :: ispin
double precision,allocatable :: Om(:)
double precision,allocatable :: XpY(:,:)
double precision,allocatable :: XmY(:,:)
2020-01-28 09:00:11 +01:00
double precision :: EcRPAx(nspin)
double precision :: EcAC(nspin)
! Hello world
write(*,*)
write(*,*)'***********************************************************'
write(*,*)'| Random phase approximation calculation with exchange |'
write(*,*)'***********************************************************'
write(*,*)
2020-09-30 09:59:18 +02:00
! TDA
if(TDA) then
write(*,*) 'Tamm-Dancoff approximation activated!'
write(*,*) ' => RPAx + TDA = CIS '
write(*,*)
end if
2020-01-28 09:00:11 +01:00
! Initialization
EcRPAx(:) = 0d0
EcAC(:) = 0d0
! Memory allocation
allocate(Om(nS),XpY(nS,nS),XmY(nS,nS))
2020-01-28 09:00:11 +01:00
! Singlet manifold
2020-09-23 09:46:44 +02:00
if(singlet) then
2020-01-28 09:00:11 +01:00
ispin = 1
call phLR(ispin,.false.,TDA,eta,nBas,nC,nO,nV,nR,nS,1d0,eHF,ERI,EcRPAx(ispin),Om,XpY,XmY)
call print_excitation('RPAx@HF ',ispin,nS,Om)
call print_transition_vectors(.true.,nBas,nC,nO,nV,nR,nS,dipole_int,Om,XpY,XmY)
2020-06-05 22:34:32 +02:00
2020-01-28 09:00:11 +01:00
endif
! Triplet manifold
2020-09-23 09:46:44 +02:00
if(triplet) then
2020-01-28 09:00:11 +01:00
ispin = 2
call phLR(ispin,.false.,TDA,eta,nBas,nC,nO,nV,nR,nS,1d0,eHF,ERI,EcRPAx(ispin),Om,XpY,XmY)
call print_excitation('RPAx@HF ',ispin,nS,Om)
call print_transition_vectors(.false.,nBas,nC,nO,nV,nR,nS,dipole_int,Om,XpY,XmY)
2020-01-28 09:00:11 +01:00
endif
2021-11-10 14:47:26 +01:00
! if(exchange_kernel) then
2020-01-28 09:00:11 +01:00
EcRPAx(1) = 0.5d0*EcRPAx(1)
EcRPAx(2) = 1.5d0*EcRPAx(2)
2021-11-10 14:47:26 +01:00
! end if
2020-01-28 09:00:11 +01:00
write(*,*)
write(*,*)'-------------------------------------------------------------------------------'
write(*,'(2X,A50,F20.10)') 'Tr@RPAx correlation energy (singlet) =',EcRPAx(1)
write(*,'(2X,A50,F20.10)') 'Tr@RPAx correlation energy (triplet) =',EcRPAx(2)
write(*,'(2X,A50,F20.10)') 'Tr@RPAx correlation energy =',EcRPAx(1) + EcRPAx(2)
write(*,'(2X,A50,F20.10)') 'Tr@RPAx total energy =',ENuc + ERHF + EcRPAx(1) + EcRPAx(2)
write(*,*)'-------------------------------------------------------------------------------'
write(*,*)
! deallocate memory
deallocate(Om,XpY,XmY)
2020-01-28 09:00:11 +01:00
! Compute the correlation energy via the adiabatic connection
if(doACFDT) then
write(*,*) '-------------------------------------------------------'
write(*,*) 'Adiabatic connection version of RPAx correlation energy'
write(*,*) '-------------------------------------------------------'
write(*,*)
2020-09-24 16:39:15 +02:00
call ACFDT(exchange_kernel,.false.,.false.,.false.,TDA,.false.,singlet,triplet,eta, &
2020-09-24 11:56:06 +02:00
nBas,nC,nO,nV,nR,nS,ERI,eHF,eHF,EcAC)
2020-01-28 09:00:11 +01:00
write(*,*)
write(*,*)'-------------------------------------------------------------------------------'
write(*,'(2X,A50,F20.10)') 'AC@RPAx correlation energy (singlet) =',EcAC(1)
write(*,'(2X,A50,F20.10)') 'AC@RPAx correlation energy (triplet) =',EcAC(2)
write(*,'(2X,A50,F20.10)') 'AC@RPAx correlation energy =',EcAC(1) + EcAC(2)
write(*,'(2X,A50,F20.10)') 'AC@RPAx total energy =',ENuc + ERHF + EcAC(1) + EcAC(2)
write(*,*)'-------------------------------------------------------------------------------'
write(*,*)
end if
end subroutine