mirror of
https://github.com/pfloos/quack
synced 2025-01-08 04:15:56 +01:00
79 lines
2.3 KiB
Fortran
79 lines
2.3 KiB
Fortran
|
subroutine phGRPA(TDA,nBas,nC,nO,nV,nR,nS,ENuc,EHF,ERI,dipole_int,e)
|
||
|
|
||
|
! Perform a direct random phase approximation calculation
|
||
|
|
||
|
implicit none
|
||
|
include 'parameters.h'
|
||
|
include 'quadrature.h'
|
||
|
|
||
|
! Input variables
|
||
|
|
||
|
logical,intent(in) :: TDA
|
||
|
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) :: EHF
|
||
|
double precision,intent(in) :: e(nBas)
|
||
|
double precision,intent(in) :: ERI(nBas,nBas,nBas,nBas)
|
||
|
double precision,intent(in) :: dipole_int(nBas,nBas,ncart)
|
||
|
|
||
|
! Local variables
|
||
|
|
||
|
integer :: ispin
|
||
|
logical :: dRPA
|
||
|
double precision,allocatable :: Aph(:,:)
|
||
|
double precision,allocatable :: Bph(:,:)
|
||
|
double precision,allocatable :: Om(:)
|
||
|
double precision,allocatable :: XpY(:,:)
|
||
|
double precision,allocatable :: XmY(:,:)
|
||
|
|
||
|
double precision :: EcRPA
|
||
|
|
||
|
! Hello world
|
||
|
|
||
|
write(*,*)
|
||
|
write(*,*)'***********************************************'
|
||
|
write(*,*)'| Random-phase approximation calculation |'
|
||
|
write(*,*)'***********************************************'
|
||
|
write(*,*)
|
||
|
|
||
|
! TDA
|
||
|
|
||
|
if(TDA) then
|
||
|
write(*,*) 'Tamm-Dancoff approximation activated!'
|
||
|
write(*,*)
|
||
|
end if
|
||
|
|
||
|
! Initialization
|
||
|
|
||
|
dRPA = .true.
|
||
|
|
||
|
EcRPA = 0d0
|
||
|
|
||
|
! Memory allocation
|
||
|
|
||
|
allocate(Om(nS),XpY(nS,nS),XmY(nS,nS),Aph(nS,nS))
|
||
|
if(.not.TDA) allocate(Bph(nS,nS))
|
||
|
|
||
|
ispin = 3
|
||
|
|
||
|
call phLR_A(ispin,dRPA,nBas,nC,nO,nV,nR,nS,1d0,e,ERI,Aph)
|
||
|
if(.not.TDA) call phLR_B(ispin,dRPA,nBas,nC,nO,nV,nR,nS,1d0,ERI,Bph)
|
||
|
|
||
|
call phLR(TDA,nS,Aph,Bph,EcRPA,Om,XpY,XmY)
|
||
|
call print_excitation_energies('phRPA@GHF',ispin,nS,Om)
|
||
|
call phLR_transition_vectors(.true.,nBas,nC,nO,nV,nR,nS,dipole_int,Om,XpY,XmY)
|
||
|
|
||
|
write(*,*)
|
||
|
write(*,*)'-------------------------------------------------------------------------------'
|
||
|
write(*,'(2X,A50,F20.10)') 'Tr@phRPA correlation energy =',EcRPA
|
||
|
write(*,'(2X,A50,F20.10)') 'Tr@phRPA total energy =',ENuc + EHF + EcRPA
|
||
|
write(*,*)'-------------------------------------------------------------------------------'
|
||
|
write(*,*)
|
||
|
|
||
|
end subroutine
|