4
1
mirror of https://github.com/pfloos/quack synced 2024-06-02 11:25:28 +02:00
quack/src/RPA/crRRPA.f90

141 lines
4.4 KiB
Fortran
Raw Permalink Normal View History

2023-11-14 13:50:53 +01:00
subroutine crRRPA(dotest,TDA,doACFDT,exchange_kernel,singlet,triplet,nBas,nC,nO,nV,nR,nS,ENuc,ERHF,ERI,dipole_int,eHF)
2021-11-10 09:42:30 +01:00
2021-11-10 14:47:26 +01:00
! Crossed-ring channel of the random phase approximation
2021-11-10 09:42:30 +01:00
implicit none
include 'parameters.h'
include 'quadrature.h'
! Input variables
2023-11-11 23:00:00 +01:00
logical,intent(in) :: dotest
2021-11-10 09:42:30 +01:00
logical,intent(in) :: TDA
logical,intent(in) :: doACFDT
logical,intent(in) :: exchange_kernel
logical,intent(in) :: singlet
logical,intent(in) :: triplet
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
2023-11-14 13:50:53 +01:00
double precision,intent(in) :: ERHF
double precision,intent(in) :: eHF(nBas)
2021-11-10 09:42:30 +01:00
double precision,intent(in) :: ERI(nBas,nBas,nBas,nBas)
double precision,intent(in) :: dipole_int(nBas,nBas,ncart)
! Local variables
integer :: ispin
2023-07-18 11:53:38 +02:00
logical :: dRPA
double precision,allocatable :: Aph(:,:)
double precision,allocatable :: Bph(:,:)
double precision,allocatable :: Om(:)
double precision,allocatable :: XpY(:,:)
double precision,allocatable :: XmY(:,:)
2023-07-31 16:01:12 +02:00
double precision :: EcRPA(nspin)
2021-11-10 09:42:30 +01:00
! Hello world
write(*,*)
2023-11-14 13:50:53 +01:00
write(*,*)'*********************************'
write(*,*)'* Restricted cr-RPA Calculation *'
write(*,*)'*********************************'
2021-11-10 09:42:30 +01:00
write(*,*)
! TDA
if(TDA) then
write(*,*) 'Tamm-Dancoff approximation activated!'
write(*,*)
end if
! Initialization
2023-11-14 13:50:53 +01:00
dRPA = .false.
2023-07-31 16:01:12 +02:00
EcRPA(:) = 0d0
2021-11-10 09:42:30 +01:00
! Memory allocation
2023-11-20 21:34:08 +01:00
allocate(Om(nS),XpY(nS,nS),XmY(nS,nS),Aph(nS,nS),Bph(nS,nS))
2021-11-10 09:42:30 +01:00
! Singlet manifold
if(singlet) then
ispin = 1
2023-11-14 13:50:53 +01:00
call phLR_A(ispin,dRPA,nBas,nC,nO,nV,nR,nS,-1d0,eHF,ERI,Aph)
2023-07-18 11:53:38 +02:00
if(.not.TDA) call phLR_B(ispin,dRPA,nBas,nC,nO,nV,nR,nS,-1d0,ERI,Bph)
2023-07-31 16:01:12 +02:00
call phLR(TDA,nS,Aph,Bph,EcRPA(ispin),Om,XpY,XmY)
2023-11-22 10:07:23 +01:00
call print_excitation_energies('crRPA@RHF','singlet',nS,Om)
2023-07-28 14:14:35 +02:00
call phLR_transition_vectors(.true.,nBas,nC,nO,nV,nR,nS,dipole_int,Om,XpY,XmY)
2021-11-10 09:42:30 +01:00
2023-12-03 18:47:30 +01:00
end if
2021-11-10 09:42:30 +01:00
! Triplet manifold
if(triplet) then
ispin = 2
2023-11-14 13:50:53 +01:00
call phLR_A(ispin,dRPA,nBas,nC,nO,nV,nR,nS,-1d0,eHF,ERI,Aph)
2023-07-18 11:53:38 +02:00
if(.not.TDA) call phLR_B(ispin,dRPA,nBas,nC,nO,nV,nR,nS,-1d0,ERI,Bph)
2023-07-31 16:01:12 +02:00
call phLR(TDA,nS,Aph,Bph,EcRPA(ispin),Om,XpY,XmY)
2023-11-22 14:54:05 +01:00
call print_excitation_energies('crRPA@RHF','triplet',nS,Om)
2023-07-28 14:14:35 +02:00
call phLR_transition_vectors(.false.,nBas,nC,nO,nV,nR,nS,dipole_int,Om,XpY,XmY)
2021-11-10 09:42:30 +01:00
2023-12-03 18:47:30 +01:00
end if
2021-11-10 09:42:30 +01:00
2023-07-18 11:53:38 +02:00
if(exchange_kernel) then
2021-11-10 09:42:30 +01:00
2023-07-31 16:01:12 +02:00
EcRPA(1) = 0.5d0*EcRPA(1)
EcRPA(2) = 1.5d0*EcRPA(2)
2021-11-10 09:42:30 +01:00
2023-07-18 11:53:38 +02:00
end if
2021-11-10 09:42:30 +01:00
write(*,*)
write(*,*)'-------------------------------------------------------------------------------'
write(*,'(2X,A50,F20.10,A3)') 'Tr@crRRPA correlation energy (singlet) =',EcRPA(1),' au'
write(*,'(2X,A50,F20.10,A3)') 'Tr@crRRPA correlation energy (triplet) =',EcRPA(2),' au'
write(*,'(2X,A50,F20.10,A3)') 'Tr@crRRPA correlation energy =',sum(EcRPA),' au'
2023-11-14 13:50:53 +01:00
write(*,'(2X,A50,F20.10,A3)') 'Tr@crRRPA total energy =',ENuc + ERHF + sum(EcRPA),' au'
2021-11-10 09:42:30 +01:00
write(*,*)'-------------------------------------------------------------------------------'
write(*,*)
! Compute the correlation energy via the adiabatic connection
2021-11-10 22:41:40 +01:00
if(doACFDT) then
2021-11-10 09:42:30 +01:00
2021-11-10 22:41:40 +01:00
write(*,*) '-------------------------------------------------------'
write(*,*) 'Adiabatic connection version of crRPA correlation energy'
write(*,*) '-------------------------------------------------------'
write(*,*)
2021-11-10 09:42:30 +01:00
2023-11-14 13:50:53 +01:00
call crACFDT(exchange_kernel,dRPA,TDA,singlet,triplet,nBas,nC,nO,nV,nR,nS,ERI,eHF,EcRPA)
2021-11-10 09:42:30 +01:00
2021-11-10 22:41:40 +01:00
write(*,*)
write(*,*)'-------------------------------------------------------------------------------'
write(*,'(2X,A50,F20.10,A3)') 'AC@crRRPA correlation energy (singlet) =',EcRPA(1),' au'
write(*,'(2X,A50,F20.10,A3)') 'AC@crRRPA correlation energy (triplet) =',EcRPA(2),' au'
write(*,'(2X,A50,F20.10,A3)') 'AC@crRRPA correlation energy =',sum(EcRPA),' au'
2023-11-14 13:50:53 +01:00
write(*,'(2X,A50,F20.10,A3)') 'AC@crRRPA total energy =',ENuc + ERHF + sum(EcRPA),' au'
2021-11-10 22:41:40 +01:00
write(*,*)'-------------------------------------------------------------------------------'
write(*,*)
2021-11-10 09:42:30 +01:00
2021-11-10 22:41:40 +01:00
end if
2021-11-10 09:42:30 +01:00
if(dotest) then
2023-11-14 13:50:53 +01:00
call dump_test_value('R','crRPA correlation energy',sum(EcRPA))
end if
2023-07-17 15:17:02 +02:00
end subroutine