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

168 lines
5.9 KiB
Fortran
Raw Normal View History

2020-10-07 22:51:30 +02:00
subroutine URPA(TDA,doACFDT,exchange_kernel,spin_conserved,spin_flip,eta,nBas,nC,nO,nV,nR,nS,ENuc,EUHF, &
2020-10-05 23:00:56 +02:00
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)
2020-09-24 16:39:15 +02:00
double precision,intent(in) :: eta
2020-09-23 09:46:44 +02:00
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
integer :: ispin
integer :: nS_aa,nS_bb,nS_sc
double precision,allocatable :: Omega_sc(:)
double precision,allocatable :: XpY_sc(:,:)
double precision,allocatable :: XmY_sc(:,:)
integer :: nS_ab,nS_ba,nS_sf
double precision,allocatable :: Omega_sf(:)
double precision,allocatable :: XpY_sf(:,:)
double precision,allocatable :: XmY_sf(:,:)
double precision :: rho_sc,rho_sf
double precision :: EcRPA(nspin)
double precision :: EcAC(nspin)
! Hello world
write(*,*)
write(*,*)'**************************************************************'
write(*,*)'| Unrestricted direct random phase approximation calculation |'
write(*,*)'**************************************************************'
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
EcRPA(:) = 0d0
2020-10-07 22:51:30 +02:00
EcAC(:) = 0d0
2020-09-23 09:46:44 +02:00
! Spin-conserved transitions
if(spin_conserved) then
ispin = 1
! Memory allocation
nS_aa = nS(1)
nS_bb = nS(2)
nS_sc = nS_aa + nS_bb
allocate(Omega_sc(nS_sc),XpY_sc(nS_sc,nS_sc),XmY_sc(nS_sc,nS_sc))
2020-09-24 16:39:15 +02:00
call unrestricted_linear_response(ispin,.true.,TDA,.false.,eta,nBas,nC,nO,nV,nR,nS_aa,nS_bb,nS_sc,nS_sc,1d0,e, &
2020-10-05 23:00:56 +02:00
ERI_aaaa,ERI_aabb,ERI_bbbb,Omega_sc,rho_sc,EcRPA(ispin),Omega_sc,XpY_sc,XmY_sc)
2020-09-23 09:46:44 +02:00
call print_excitation('URPA ',5,nS_sc,Omega_sc)
2020-10-05 16:58:19 +02:00
call print_unrestricted_transition_vectors(ispin,nBas,nC,nO,nV,nR,nS,nS_aa,nS_bb,nS_sc,dipole_int_aa,dipole_int_bb, &
c,S,Omega_sc,XpY_sc,XmY_sc)
2020-09-23 09:46:44 +02:00
2020-10-07 22:51:30 +02:00
deallocate(Omega_sc,XpY_sc,XmY_sc)
2020-09-23 09:46:44 +02:00
endif
! Spin-flip transitions
if(spin_flip) then
ispin = 2
! Memory allocation
nS_ab = (nO(1) - nC(1))*(nV(2) - nR(2))
nS_ba = (nO(2) - nC(2))*(nV(1) - nR(1))
nS_sf = nS_ab + nS_ba
allocate(Omega_sf(nS_sf),XpY_sf(nS_sf,nS_sf),XmY_sf(nS_sf,nS_sf))
2020-09-30 15:02:48 +02:00
call unrestricted_linear_response(ispin,.true.,TDA,.false.,eta,nBas,nC,nO,nV,nR,nS_ab,nS_ba,nS_sf,nS_sf,1d0,e, &
2020-10-05 23:00:56 +02:00
ERI_aaaa,ERI_aabb,ERI_bbbb,Omega_sf,rho_sf,EcRPA(ispin),Omega_sf,XpY_sf,XmY_sf)
2020-09-23 09:46:44 +02:00
call print_excitation('URPA ',6,nS_sf,Omega_sf)
2020-10-05 16:58:19 +02:00
call print_unrestricted_transition_vectors(ispin,nBas,nC,nO,nV,nR,nS,nS_ab,nS_ba,nS_sf,dipole_int_aa,dipole_int_bb, &
c,S,Omega_sf,XpY_sf,XmY_sf)
2020-09-23 09:46:44 +02:00
2020-10-07 22:51:30 +02:00
deallocate(Omega_sf,XpY_sf,XmY_sf)
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
2020-10-07 22:51:30 +02:00
call unrestricted_ACFDT(exchange_kernel,.false.,.true.,.false.,TDA,.false.,spin_conserved,spin_flip,eta, &
2020-10-07 18:26:24 +02:00
nBas,nC,nO,nV,nR,nS,ERI_aaaa,ERI_aabb,ERI_bbbb,e,e,EcAC)
2020-09-23 09:46:44 +02:00
2021-04-01 22:04:23 +02:00
if(exchange_kernel) then
EcAC(1) = 0.5d0*EcAC(1)
EcAC(2) = 1.5d0*EcAC(2)
end if
2020-10-07 18:26:24 +02:00
write(*,*)
write(*,*)'-------------------------------------------------------------------------------'
2020-10-07 22:51:30 +02:00
write(*,'(2X,A50,F20.10)') 'AC@URPA correlation energy (spin-conserved) =',EcAC(1)
write(*,'(2X,A50,F20.10)') 'AC@URPA correlation energy (spin-flip) =',EcAC(2)
write(*,'(2X,A50,F20.10)') 'AC@URPA correlation energy =',EcAC(1) + EcAC(2)
write(*,'(2X,A50,F20.10)') 'AC@URPA total energy =',ENuc + EUHF + EcAC(1) + EcAC(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
2020-10-07 22:51:30 +02:00
end subroutine URPA