4
1
mirror of https://github.com/pfloos/quack synced 2024-11-08 07:03:59 +01:00
quack/src/RPA/phURPA.f90

167 lines
5.5 KiB
Fortran
Raw Normal View History

2023-07-23 11:16:42 +02:00
subroutine phURPA(TDA,doACFDT,exchange_kernel,spin_conserved,spin_flip,nBas,nC,nO,nV,nR,nS,ENuc,EUHF, &
ERI_aaaa,ERI_aabb,ERI_bbbb,dipole_int_aa,dipole_int_bb,e,c,S)
2020-09-23 09:46:44 +02:00
! Perform random phase approximation calculation with exchange (aka TDHF) in the unrestricted formalism
implicit none
include 'parameters.h'
include 'quadrature.h'
! Input variables
2020-09-24 16:39:15 +02:00
logical,intent(in) :: TDA
2020-09-23 09:46:44 +02:00
logical,intent(in) :: doACFDT
logical,intent(in) :: exchange_kernel
logical,intent(in) :: spin_conserved
logical,intent(in) :: spin_flip
integer,intent(in) :: nBas
integer,intent(in) :: nC(nspin)
integer,intent(in) :: nO(nspin)
integer,intent(in) :: nV(nspin)
integer,intent(in) :: nR(nspin)
integer,intent(in) :: nS(nspin)
double precision,intent(in) :: ENuc
double precision,intent(in) :: EUHF
double precision,intent(in) :: e(nBas,nspin)
2020-10-05 16:58:19 +02:00
double precision,intent(in) :: c(nBas,nBas,nspin)
double precision,intent(in) :: S(nBas,nBas)
2020-09-23 09:46:44 +02:00
double precision,intent(in) :: ERI_aaaa(nBas,nBas,nBas,nBas)
double precision,intent(in) :: ERI_aabb(nBas,nBas,nBas,nBas)
double precision,intent(in) :: ERI_bbbb(nBas,nBas,nBas,nBas)
2020-09-28 22:58:58 +02:00
double precision,intent(in) :: dipole_int_aa(nBas,nBas,ncart)
double precision,intent(in) :: dipole_int_bb(nBas,nBas,ncart)
2020-09-23 09:46:44 +02:00
! Local variables
logical :: dRPA
2020-09-23 09:46:44 +02:00
integer :: ispin
integer :: nSa,nSb,nSt
double precision,allocatable :: Aph(:,:)
double precision,allocatable :: Bph(:,:)
double precision,allocatable :: Om(:)
double precision,allocatable :: XpY(:,:)
double precision,allocatable :: XmY(:,:)
2020-09-23 09:46:44 +02:00
double precision :: rho
2020-09-23 09:46:44 +02:00
double precision :: EcRPA(nspin)
! Hello world
write(*,*)
write(*,*)'***********************************'
write(*,*)'* Unrestricted ph-RPA Calculation *'
write(*,*)'***********************************'
2020-09-23 09:46:44 +02:00
write(*,*)
2020-09-30 09:59:18 +02:00
! TDA
if(TDA) then
write(*,*) 'Tamm-Dancoff approximation activated!'
write(*,*)
end if
2020-09-23 09:46:44 +02:00
! Initialization
dRPA = .true.
2020-09-23 09:46:44 +02:00
EcRPA(:) = 0d0
! Spin-conserved transitions
if(spin_conserved) then
ispin = 1
! Memory allocation
nSa = nS(1)
nSb = nS(2)
nSt = nSa + nSb
2020-09-23 09:46:44 +02:00
allocate(Aph(nSt,nSt),Bph(nSt,nSt),Om(nSt),XpY(nSt,nSt),XmY(nSt,nSt))
2020-09-23 09:46:44 +02:00
call phULR_A(ispin,dRPA,nBas,nC,nO,nV,nR,nSa,nSb,nSt,1d0,e,ERI_aaaa,ERI_aabb,ERI_bbbb,Aph)
if(.not.TDA) call phULR_B(ispin,dRPA,nBas,nC,nO,nV,nR,nSa,nSb,nSt,1d0,ERI_aaaa,ERI_aabb,ERI_bbbb,Bph)
2020-09-23 09:46:44 +02:00
call phULR(TDA,nSa,nSb,nSt,Aph,Bph,EcRPA(ispin),Om,XpY,XmY)
call print_excitation_energies('phRPA@UHF',5,nSt,Om)
call phULR_transition_vectors(ispin,nBas,nC,nO,nV,nR,nS,nSa,nSb,nSt,dipole_int_aa,dipole_int_bb,c,S,Om,XpY,XmY)
deallocate(Aph,Bph,Om,XpY,XmY)
2020-09-23 09:46:44 +02:00
endif
! Spin-flip transitions
if(spin_flip) then
ispin = 2
! Memory allocation
nSa = (nO(1) - nC(1))*(nV(2) - nR(2))
nSb = (nO(2) - nC(2))*(nV(1) - nR(1))
nSt = nSa + nSb
allocate(Om(nSt),XpY(nSt,nSt),XmY(nSt,nSt))
2020-09-23 09:46:44 +02:00
call phULR_A(ispin,dRPA,nBas,nC,nO,nV,nR,nSa,nSb,nSt,1d0,e,ERI_aaaa,ERI_aabb,ERI_bbbb,Aph)
if(.not.TDA) call phULR_B(ispin,dRPA,nBas,nC,nO,nV,nR,nSa,nSb,nSt,1d0,ERI_aaaa,ERI_aabb,ERI_bbbb,Bph)
2020-09-23 09:46:44 +02:00
call phULR(TDA,nSa,nSa,nSt,Aph,Bph,EcRPA(ispin),Om,XpY,XmY)
call print_excitation_energies('phRPA@UHF',6,nSt,Om)
call phULR_transition_vectors(ispin,nBas,nC,nO,nV,nR,nS,nSa,nSb,nSt,dipole_int_aa,dipole_int_bb,c,S,Om,XpY,XmY)
2020-09-23 09:46:44 +02:00
deallocate(Aph,Bph,Om,XpY,XmY)
2020-10-07 22:51:30 +02:00
2020-09-23 09:46:44 +02:00
endif
2020-10-07 18:26:24 +02:00
if(exchange_kernel) then
2020-09-23 09:46:44 +02:00
2020-10-07 18:26:24 +02:00
EcRPA(1) = 0.5d0*EcRPA(1)
EcRPA(2) = 1.5d0*EcRPA(2)
2020-09-23 09:46:44 +02:00
2020-10-07 18:26:24 +02:00
end if
2020-09-23 09:46:44 +02:00
write(*,*)
write(*,*)'-------------------------------------------------------------------------------'
write(*,'(2X,A50,F20.10)') 'Tr@URPA correlation energy (spin-conserved) =',EcRPA(1)
write(*,'(2X,A50,F20.10)') 'Tr@URPA correlation energy (spin-flip) =',EcRPA(2)
write(*,'(2X,A50,F20.10)') 'Tr@URPA correlation energy =',EcRPA(1) + EcRPA(2)
write(*,'(2X,A50,F20.10)') 'Tr@URPA total energy =',ENuc + EUHF + EcRPA(1) + EcRPA(2)
write(*,*)'-------------------------------------------------------------------------------'
write(*,*)
! Compute the correlation energy via the adiabatic connection
2020-10-07 18:26:24 +02:00
if(doACFDT) then
2020-09-23 09:46:44 +02:00
2020-10-07 22:51:30 +02:00
write(*,*) '---------------------------------------------------------'
write(*,*) ' Adiabatic connection version of URPA correlation energy '
write(*,*) '---------------------------------------------------------'
2020-10-07 18:26:24 +02:00
write(*,*)
2020-09-23 09:46:44 +02:00
2023-07-23 11:16:42 +02:00
call phUACFDT(exchange_kernel,.false.,.true.,.false.,TDA,.false.,spin_conserved,spin_flip, &
nBas,nC,nO,nV,nR,nS,ERI_aaaa,ERI_aabb,ERI_bbbb,e,e,EcRPA)
2020-09-23 09:46:44 +02:00
2021-04-01 22:04:23 +02:00
if(exchange_kernel) then
2023-07-22 19:41:45 +02:00
EcRPA(1) = 0.5d0*EcRPA(1)
EcRPA(2) = 1.5d0*EcRPA(2)
2021-04-01 22:04:23 +02:00
end if
2020-10-07 18:26:24 +02:00
write(*,*)
write(*,*)'-------------------------------------------------------------------------------'
2023-07-22 19:41:45 +02:00
write(*,'(2X,A50,F20.10)') 'AC@URPA correlation energy (spin-conserved) =',EcRPA(1)
write(*,'(2X,A50,F20.10)') 'AC@URPA correlation energy (spin-flip) =',EcRPA(2)
write(*,'(2X,A50,F20.10)') 'AC@URPA correlation energy =',EcRPA(1) + EcRPA(2)
write(*,'(2X,A50,F20.10)') 'AC@URPA total energy =',ENuc + EUHF + EcRPA(1) + EcRPA(2)
2020-10-07 18:26:24 +02:00
write(*,*)'-------------------------------------------------------------------------------'
write(*,*)
2020-09-23 09:46:44 +02:00
2020-10-07 18:26:24 +02:00
end if
2020-09-23 09:46:44 +02:00
2023-07-22 19:41:45 +02:00
end subroutine