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

88 lines
2.4 KiB
Fortran
Raw Normal View History

2023-11-14 13:50:53 +01:00
subroutine phGRPA(dotest,TDA,nBas,nC,nO,nV,nR,nS,ENuc,EGHF,ERI,dipole_int,eHF)
2023-10-26 18:00:16 +02:00
! Perform a direct random phase approximation calculation
implicit none
include 'parameters.h'
include 'quadrature.h'
! Input variables
2023-11-11 23:00:00 +01:00
logical,intent(in) :: dotest
2023-10-26 18:00:16 +02:00
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
2023-11-14 13:50:53 +01:00
double precision,intent(in) :: EGHF
double precision,intent(in) :: eHF(nBas)
2023-10-26 18:00:16 +02:00
double precision,intent(in) :: ERI(nBas,nBas,nBas,nBas)
double precision,intent(in) :: dipole_int(nBas,nBas,ncart)
! Local variables
integer :: ispin
logical :: dRPA
2023-11-22 14:54:05 +01:00
double precision :: lambda
2023-10-26 18:00:16 +02:00
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(*,*)
2023-11-11 23:00:00 +01:00
write(*,*)'**********************************'
write(*,*)'* Generalized ph-RPA Calculation |'
write(*,*)'**********************************'
2023-10-26 18:00:16 +02:00
write(*,*)
! TDA
if(TDA) then
write(*,*) 'Tamm-Dancoff approximation activated!'
write(*,*)
end if
! Initialization
2023-11-14 13:50:53 +01:00
dRPA = .true.
2023-10-26 18:00:16 +02:00
EcRPA = 0d0
2023-11-22 14:54:05 +01:00
lambda = 0d0
2023-10-26 18:00:16 +02: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))
2023-10-26 18:00:16 +02:00
ispin = 3
2023-11-22 14:54:05 +01:00
call phLR_A(ispin,dRPA,nBas,nC,nO,nV,nR,nS,lambda,eHF,ERI,Aph)
if(.not.TDA) call phLR_B(ispin,dRPA,nBas,nC,nO,nV,nR,nS,lambda,ERI,Bph)
2023-10-26 18:00:16 +02:00
call phLR(TDA,nS,Aph,Bph,EcRPA,Om,XpY,XmY)
2023-11-22 10:07:23 +01:00
call print_excitation_energies('phRPA@GHF','spinorbital',nS,Om)
2023-10-26 18:00:16 +02:00
call phLR_transition_vectors(.true.,nBas,nC,nO,nV,nR,nS,dipole_int,Om,XpY,XmY)
write(*,*)
write(*,*)'-------------------------------------------------------------------------------'
2023-11-14 13:50:53 +01:00
write(*,'(2X,A50,F20.10,A3)') 'Tr@phGRPA correlation energy = ',EcRPA,' au'
write(*,'(2X,A50,F20.10,A3)') 'Tr@phGRPA total energy = ',ENuc + EGHF + EcRPA,' au'
2023-10-26 18:00:16 +02:00
write(*,*)'-------------------------------------------------------------------------------'
write(*,*)
2023-11-11 23:00:00 +01:00
if(dotest) then
2023-11-14 13:50:53 +01:00
call dump_test_value('G','phRPA corrlation energy',EcRPA)
2023-11-11 23:00:00 +01:00
end if
2023-10-26 18:00:16 +02:00
end subroutine