4
1
mirror of https://github.com/pfloos/quack synced 2024-07-12 14:14:40 +02:00
quack/src/RPA/ppRPA.f90

150 lines
5.0 KiB
Fortran
Raw Normal View History

2023-07-18 09:59:25 +02:00
subroutine ppRPA(TDA,doACFDT,singlet,triplet,nBas,nC,nO,nV,nR,ENuc,EHF,ERI,dipole_int,e)
2019-10-05 23:09:20 +02:00
2023-07-18 09:59:25 +02:00
! Perform ppRPA calculation
2019-10-05 23:09:20 +02:00
implicit none
include 'parameters.h'
! Input variables
2021-10-18 22:05:26 +02:00
logical,intent(in) :: TDA
2021-11-10 22:41:40 +01:00
logical,intent(in) :: doACFDT
2021-10-18 22:05:26 +02:00
logical,intent(in) :: singlet
logical,intent(in) :: triplet
2019-10-05 23:09:20 +02:00
integer,intent(in) :: nBas
integer,intent(in) :: nC
integer,intent(in) :: nO
integer,intent(in) :: nV
integer,intent(in) :: nR
double precision,intent(in) :: ENuc
2023-07-18 09:59:25 +02:00
double precision,intent(in) :: EHF
2019-10-05 23:09:20 +02:00
double precision,intent(in) :: e(nBas)
double precision,intent(in) :: ERI(nBas,nBas,nBas,nBas)
2022-08-17 15:52:13 +02:00
double precision,intent(in) :: dipole_int(nBas,nBas,ncart)
2019-10-05 23:09:20 +02:00
! Local variables
integer :: ispin
2022-08-17 15:52:13 +02:00
integer :: nOO
integer :: nVV
2023-07-18 09:59:25 +02:00
double precision,allocatable :: Bpp(:,:)
double precision,allocatable :: Cpp(:,:)
double precision,allocatable :: Dpp(:,:)
double precision,allocatable :: Om1(:)
2022-08-17 15:52:13 +02:00
double precision,allocatable :: X1(:,:)
double precision,allocatable :: Y1(:,:)
double precision,allocatable :: Om2(:)
2022-08-17 15:52:13 +02:00
double precision,allocatable :: X2(:,:)
double precision,allocatable :: Y2(:,:)
2019-10-05 23:09:20 +02:00
2023-07-18 09:59:25 +02:00
double precision :: EcRPA(nspin)
2019-10-05 23:09:20 +02:00
! Hello world
write(*,*)
2019-10-05 23:32:17 +02:00
write(*,*)'****************************************'
write(*,*)'| particle-particle RPA calculation |'
write(*,*)'****************************************'
2019-10-05 23:09:20 +02:00
write(*,*)
! Initialization
2023-07-18 09:59:25 +02:00
EcRPA(:) = 0d0
2019-10-05 23:09:20 +02:00
2022-08-17 15:52:13 +02:00
! Singlet manifold
2019-10-05 23:09:20 +02:00
2022-08-17 15:52:13 +02:00
if(singlet) then
2019-10-07 22:31:45 +02:00
2022-09-09 21:48:50 +02:00
write(*,*) '****************'
write(*,*) '*** Singlets ***'
write(*,*) '****************'
write(*,*)
2022-08-17 15:52:13 +02:00
ispin = 1
2019-10-07 22:31:45 +02:00
2022-08-17 15:52:13 +02:00
nOO = nO*(nO+1)/2
nVV = nV*(nV+1)/2
2021-11-10 22:41:40 +01:00
2023-07-18 09:59:25 +02:00
allocate(Om1(nVV),X1(nVV,nVV),Y1(nOO,nVV),Om2(nOO),X2(nVV,nOO),Y2(nOO,nOO), &
Bpp(nVV,nOO),Cpp(nVV,nVV),Dpp(nOO,nOO))
2021-11-10 22:41:40 +01:00
2023-07-18 09:59:25 +02:00
if(.not.TDA) call ppLR_B(ispin,nBas,nC,nO,nV,nR,nOO,nVV,1d0,ERI,Bpp)
call ppLR_C(ispin,nBas,nC,nO,nV,nR,nVV,1d0,e,ERI,Cpp)
call ppLR_D(ispin,nBas,nC,nO,nV,nR,nOO,1d0,e,ERI,Dpp)
call ppLR(TDA,nOO,nVV,Bpp,Cpp,Dpp,Om1,X1,Y1,Om2,X2,Y2,EcRPA(ispin))
2019-10-07 22:31:45 +02:00
! call print_transition_vectors_pp(.true.,nBas,nC,nO,nV,nR,nOO,nVV,dipole_int,Om1,X1,Y1,Om2,X2,Y2)
2019-10-07 22:31:45 +02:00
2023-07-28 14:35:14 +02:00
call print_excitation_energies('ppRPA@HF (N+2)',ispin,nVV,Om1)
call print_excitation_energies('ppRPA@HF (N-2)',ispin,nOO,Om2)
2022-09-09 21:48:50 +02:00
2023-07-18 09:59:25 +02:00
deallocate(Om1,X1,Y1,Om2,X2,Y2,Bpp,Cpp,Dpp)
2019-10-07 22:31:45 +02:00
2019-10-05 23:09:20 +02:00
endif
! Triplet manifold
2021-10-18 22:05:26 +02:00
if(triplet) then
2019-10-05 23:09:20 +02:00
2022-09-09 21:48:50 +02:00
write(*,*) '****************'
write(*,*) '*** Triplets ***'
write(*,*) '****************'
write(*,*)
2019-10-05 23:09:20 +02:00
ispin = 2
2022-08-17 15:52:13 +02:00
nOO = nO*(nO-1)/2
nVV = nV*(nV-1)/2
2023-07-18 09:59:25 +02:00
allocate(Om1(nVV),X1(nVV,nVV),Y1(nOO,nVV),Om2(nOO),X2(nVV,nOO),Y2(nOO,nOO), &
Bpp(nVV,nOO),Cpp(nVV,nVV),Dpp(nOO,nOO))
if(.not.TDA) call ppLR_B(ispin,nBas,nC,nO,nV,nR,nOO,nVV,1d0,ERI,Bpp)
call ppLR_C(ispin,nBas,nC,nO,nV,nR,nVV,1d0,e,ERI,Cpp)
call ppLR_D(ispin,nBas,nC,nO,nV,nR,nOO,1d0,e,ERI,Dpp)
2022-08-17 15:52:13 +02:00
2023-07-18 09:59:25 +02:00
call ppLR(TDA,nOO,nVV,Bpp,Cpp,Dpp,Om1,X1,Y1,Om2,X2,Y2,EcRPA(ispin))
2022-08-17 15:52:13 +02:00
! call print_transition_vectors_pp(.false.,nBas,nC,nO,nV,nR,nOO,nVV,dipole_int,Om1,X1,Y1,Om2,X2,Y2)
2019-10-07 22:31:45 +02:00
2023-07-28 14:35:14 +02:00
call print_excitation_energies('ppRPA@HF (N+2)',ispin,nVV,Om1)
call print_excitation_energies('ppRPA@HF (N-2)',ispin,nOO,Om2)
2022-09-09 21:48:50 +02:00
2023-07-18 09:59:25 +02:00
deallocate(Om1,X1,Y1,Om2,X2,Y2,Bpp,Cpp,Dpp)
2019-10-07 22:31:45 +02:00
2019-10-05 23:09:20 +02:00
endif
write(*,*)
write(*,*)'-------------------------------------------------------------------------------'
2023-07-18 09:59:25 +02:00
write(*,'(2X,A50,F20.10)') 'Tr@ppRPA correlation energy (singlet) =',EcRPA(1)
write(*,'(2X,A50,F20.10)') 'Tr@ppRPA correlation energy (triplet) =',3d0*EcRPA(2)
write(*,'(2X,A50,F20.10)') 'Tr@ppRPA correlation energy =',EcRPA(1) + 3d0*EcRPA(2)
write(*,'(2X,A50,F20.10)') 'Tr@ppRPA total energy =',ENuc + EHF + EcRPA(1) + 3d0*EcRPA(2)
2019-10-05 23:09:20 +02:00
write(*,*)'-------------------------------------------------------------------------------'
write(*,*)
2021-11-10 22:41:40 +01:00
! Compute the correlation energy via the adiabatic connection
if(doACFDT) then
2023-07-18 09:59:25 +02:00
write(*,*) '--------------------------------------------------------'
write(*,*) 'Adiabatic connection version of ppRPA correlation energy'
write(*,*) '--------------------------------------------------------'
2021-11-10 22:41:40 +01:00
write(*,*)
2023-07-31 16:01:12 +02:00
call ppACFDT(TDA,singlet,triplet,nBas,nC,nO,nV,nR,ERI,e,EcRPA)
2021-11-10 22:41:40 +01:00
write(*,*)
write(*,*)'-------------------------------------------------------------------------------'
2023-07-31 16:01:12 +02:00
write(*,'(2X,A50,F20.10,A3)') 'AC@ppRPA correlation energy (singlet) =',EcRPA(1),' au'
write(*,'(2X,A50,F20.10,A3)') 'AC@ppRPA correlation energy (triplet) =',EcRPA(2),' au'
write(*,'(2X,A50,F20.10,A3)') 'AC@ppRPA correlation energy =',EcRPA(1) + EcRPA(2),' au'
write(*,'(2X,A50,F20.10,A3)') 'AC@ppRPA total energy =',ENuc + EHF + EcRPA(1) + EcRPA(2),' au'
2021-11-10 22:41:40 +01:00
write(*,*)'-------------------------------------------------------------------------------'
write(*,*)
end if
2023-07-18 09:59:25 +02:00
end subroutine