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

340 lines
11 KiB
Fortran
Raw Normal View History

2023-11-13 18:36:09 +01:00
subroutine RG0T0pp(dotest,doACFDT,exchange_kernel,doXBS,dophBSE,TDA_T,TDA,dBSE,dTDA,doppBSE,singlet,triplet, &
linearize,eta,regularize,nBas,nC,nO,nV,nR,nS,ENuc,ERHF,ERI,dipole_int,eHF)
2019-10-16 18:14:47 +02:00
2020-03-14 23:00:44 +01:00
! Perform one-shot calculation with a T-matrix self-energy (G0T0)
2019-10-16 18:14:47 +02:00
implicit none
include 'parameters.h'
! Input variables
2023-11-13 18:36:09 +01:00
logical,intent(in) :: dotest
2020-03-21 21:09:48 +01:00
logical,intent(in) :: doACFDT
logical,intent(in) :: exchange_kernel
logical,intent(in) :: doXBS
2023-07-21 12:30:29 +02:00
logical,intent(in) :: dophBSE
logical,intent(in) :: doppBSE
2021-10-15 22:32:22 +02:00
logical,intent(in) :: TDA_T
2020-03-21 21:09:48 +01:00
logical,intent(in) :: TDA
2020-06-14 21:20:01 +02:00
logical,intent(in) :: dBSE
2020-06-14 13:04:16 +02:00
logical,intent(in) :: dTDA
logical,intent(in) :: singlet
logical,intent(in) :: triplet
2020-03-21 21:09:48 +01:00
logical,intent(in) :: linearize
2019-10-16 18:14:47 +02:00
double precision,intent(in) :: eta
2021-12-17 11:41:40 +01:00
logical,intent(in) :: regularize
2019-10-16 18:14:47 +02:00
2020-03-21 21:09:48 +01:00
integer,intent(in) :: nBas
integer,intent(in) :: nC
integer,intent(in) :: nO
integer,intent(in) :: nV
integer,intent(in) :: nR
integer,intent(in) :: nS
2019-10-16 18:14:47 +02:00
double precision,intent(in) :: ENuc
double precision,intent(in) :: ERHF
double precision,intent(in) :: eHF(nBas)
2023-07-30 22:01:44 +02:00
double precision,intent(in) :: ERI(nBas,nBas,nBas,nBas)
double precision,intent(in) :: dipole_int(nBas,nBas,ncart)
2019-10-16 18:14:47 +02:00
! Local variables
2023-11-29 16:20:53 +01:00
logical :: print_T = .false.
2019-10-18 23:16:37 +02:00
integer :: ispin
2020-04-13 11:33:48 +02:00
integer :: iblock
2023-07-27 19:17:20 +02:00
integer :: nOOs,nOOt
integer :: nVVs,nVVt
2019-10-16 18:14:47 +02:00
double precision :: EcRPA(nspin)
double precision :: EcBSE(nspin)
2022-01-06 13:48:15 +01:00
double precision :: EcGM
double precision,allocatable :: Bpp(:,:)
double precision,allocatable :: Cpp(:,:)
double precision,allocatable :: Dpp(:,:)
2023-07-27 19:17:20 +02:00
double precision,allocatable :: Om1s(:),Om1t(:)
double precision,allocatable :: X1s(:,:),X1t(:,:)
double precision,allocatable :: Y1s(:,:),Y1t(:,:)
double precision,allocatable :: rho1s(:,:,:),rho1t(:,:,:)
double precision,allocatable :: Om2s(:),Om2t(:)
double precision,allocatable :: X2s(:,:),X2t(:,:)
double precision,allocatable :: Y2s(:,:),Y2t(:,:)
double precision,allocatable :: rho2s(:,:,:),rho2t(:,:,:)
2023-07-12 14:13:45 +02:00
double precision,allocatable :: Sig(:)
2019-10-18 23:16:37 +02:00
double precision,allocatable :: Z(:)
2023-07-17 13:35:24 +02:00
double precision,allocatable :: eGT(:)
2023-06-30 16:47:26 +02:00
double precision,allocatable :: eGTlin(:)
2019-10-16 18:14:47 +02:00
2020-03-11 16:41:20 +01:00
! Output variables
2019-10-16 18:14:47 +02:00
! Hello world
write(*,*)
2023-11-13 18:36:09 +01:00
write(*,*)'*********************************'
write(*,*)'* Restricted G0T0pp Calculation *'
write(*,*)'*********************************'
2019-10-16 18:14:47 +02:00
write(*,*)
2023-08-23 16:04:00 +02:00
! TDA for T
if(TDA_T) then
2023-11-29 16:20:53 +01:00
write(*,*) 'Tamm-Dancoff approximation activated for pp T-matrix!'
2023-08-23 16:04:00 +02:00
write(*,*)
end if
! TDA
if(TDA) then
write(*,*) 'Tamm-Dancoff approximation activated!'
write(*,*)
end if
2020-04-09 14:31:50 +02:00
! Dimensions of the pp-RPA linear reponse matrices
2019-10-16 18:14:47 +02:00
2023-07-31 14:47:09 +02:00
! nOOs = nO*(nO + 1)/2
! nVVs = nV*(nV + 1)/2
2023-07-27 19:17:20 +02:00
nOOs = nO*nO
nVVs = nV*nV
2020-04-13 11:33:48 +02:00
2023-07-27 19:17:20 +02:00
nOOt = nO*(nO - 1)/2
nVVt = nV*(nV - 1)/2
2019-10-18 23:16:37 +02:00
2019-10-16 18:14:47 +02:00
! Memory allocation
2023-07-27 19:17:20 +02:00
allocate(Om1s(nVVs),X1s(nVVs,nVVs),Y1s(nOOs,nVVs), &
Om2s(nOOs),X2s(nVVs,nOOs),Y2s(nOOs,nOOs), &
rho1s(nBas,nBas,nVVs),rho2s(nBas,nBas,nOOs), &
Om1t(nVVt),X1t(nVVt,nVVt),Y1t(nOOt,nVVt), &
Om2t(nOOt),X2t(nVVt,nOOt),Y2t(nOOt,nOOt), &
rho1t(nBas,nBas,nVVt),rho2t(nBas,nBas,nOOt), &
Sig(nBas),Z(nBas),eGT(nBas),eGTlin(nBas))
2019-10-16 18:14:47 +02:00
2019-10-20 08:18:58 +02:00
!----------------------------------------------
2020-04-13 11:33:48 +02:00
! alpha-beta block
2019-10-20 08:18:58 +02:00
!----------------------------------------------
2020-04-13 11:33:48 +02:00
ispin = 1
2023-07-31 14:47:09 +02:00
! iblock = 1
2022-01-09 23:12:21 +01:00
iblock = 3
2019-10-20 08:18:58 +02:00
2019-10-16 18:14:47 +02:00
! Compute linear response
2023-07-27 19:17:20 +02:00
allocate(Bpp(nVVs,nOOs),Cpp(nVVs,nVVs),Dpp(nOOs,nOOs))
2023-07-30 22:01:44 +02:00
if(.not.TDA_T) call ppLR_B(iblock,nBas,nC,nO,nV,nR,nOOs,nVVs,1d0,ERI,Bpp)
call ppLR_C(iblock,nBas,nC,nO,nV,nR,nVVs,1d0,eHF,ERI,Cpp)
call ppLR_D(iblock,nBas,nC,nO,nV,nR,nOOs,1d0,eHF,ERI,Dpp)
2023-07-27 19:17:20 +02:00
call ppLR(TDA_T,nOOs,nVVs,Bpp,Cpp,Dpp,Om1s,X1s,Y1s,Om2s,X2s,Y2s,EcRPA(ispin))
deallocate(Bpp,Cpp,Dpp)
2019-10-16 18:14:47 +02:00
2023-11-29 16:20:53 +01:00
if(print_T) call print_excitation_energies('ppRPA@RHF','2p (alpha-beta)',nVVs,Om1s)
if(print_T) call print_excitation_energies('ppRPA@RHF','2h (alpha-beta)',nOOs,Om2s)
2019-10-20 08:18:58 +02:00
!----------------------------------------------
2020-04-13 11:33:48 +02:00
! alpha-alpha block
2019-10-20 08:18:58 +02:00
!----------------------------------------------
2020-04-13 11:33:48 +02:00
ispin = 2
2023-07-31 14:47:09 +02:00
! iblock = 2
2020-04-13 11:33:48 +02:00
iblock = 4
2019-10-20 08:18:58 +02:00
2019-10-28 16:34:09 +01:00
! Compute linear response
2019-10-16 18:14:47 +02:00
2023-07-27 19:17:20 +02:00
allocate(Bpp(nVVt,nOOt),Cpp(nVVt,nVVt),Dpp(nOOt,nOOt))
2023-07-30 22:01:44 +02:00
if(.not.TDA_T) call ppLR_B(iblock,nBas,nC,nO,nV,nR,nOOt,nVVt,1d0,ERI,Bpp)
call ppLR_C(iblock,nBas,nC,nO,nV,nR,nVVt,1d0,eHF,ERI,Cpp)
call ppLR_D(iblock,nBas,nC,nO,nV,nR,nOOt,1d0,eHF,ERI,Dpp)
2023-07-27 19:17:20 +02:00
call ppLR(TDA_T,nOOt,nVVt,Bpp,Cpp,Dpp,Om1t,X1t,Y1t,Om2t,X2t,Y2t,EcRPA(ispin))
deallocate(Bpp,Cpp,Dpp)
2019-10-20 08:18:58 +02:00
2023-11-29 16:20:53 +01:00
if(print_T) call print_excitation_energies('ppRPA@RHF','2p (alpha-alpha)',nVVt,Om1t)
if(print_T) call print_excitation_energies('ppRPA@RHF','2h (alpha-alpha)',nOOt,Om2t)
2019-10-20 08:18:58 +02:00
!----------------------------------------------
2023-07-31 14:47:09 +02:00
! Compute excitation densities
2019-10-20 08:18:58 +02:00
!----------------------------------------------
2019-10-16 18:14:47 +02:00
2023-07-31 14:47:09 +02:00
! iblock = 1
2022-01-09 23:12:21 +01:00
iblock = 3
2023-07-30 22:01:44 +02:00
call GTpp_excitation_density(iblock,nBas,nC,nO,nV,nR,nOOs,nVVs,ERI,X1s,Y1s,rho1s,X2s,Y2s,rho2s)
2020-04-09 14:31:50 +02:00
2023-07-31 14:47:09 +02:00
! iblock = 2
2022-01-07 15:14:00 +01:00
iblock = 4
2023-07-30 22:01:44 +02:00
call GTpp_excitation_density(iblock,nBas,nC,nO,nV,nR,nOOt,nVVt,ERI,X1t,Y1t,rho1t,X2t,Y2t,rho2t)
2020-04-09 14:31:50 +02:00
2023-07-31 14:47:09 +02:00
!----------------------------------------------
! Compute T-matrix version of the self-energy
!----------------------------------------------
2023-08-01 17:15:44 +02:00
if(regularize) then
call GTpp_regularization(nBas,nC,nO,nV,nR,nOOs,nVVs,eHF,Om1s,rho1s,Om2s,rho2s)
call GTpp_regularization(nBas,nC,nO,nV,nR,nOOt,nVVt,eHF,Om1t,rho1t,Om2t,rho2t)
2023-08-01 17:15:44 +02:00
end if
2023-08-01 17:00:03 +02:00
2023-07-27 19:17:20 +02:00
call GTpp_self_energy_diag(eta,nBas,nC,nO,nV,nR,nOOs,nVVs,nOOt,nVVt,eHF,Om1s,rho1s,Om2s,rho2s, &
Om1t,rho1t,Om2t,rho2t,EcGM,Sig,Z)
2021-02-15 22:00:19 +01:00
2019-10-20 08:18:58 +02:00
!----------------------------------------------
2019-10-16 18:14:47 +02:00
! Solve the quasi-particle equation
2019-10-20 08:18:58 +02:00
!----------------------------------------------
2023-06-30 16:47:26 +02:00
2023-09-07 14:01:58 +02:00
eGTlin(:) = eHF(:) + Z(:)*Sig(:)
2020-03-21 21:09:48 +01:00
if(linearize) then
2023-06-30 16:47:26 +02:00
write(*,*) ' *** Quasiparticle energies obtained by linearization *** '
write(*,*)
2023-09-07 14:01:58 +02:00
eGT(:) = eGTlin(:)
2020-03-21 21:09:48 +01:00
else
2023-06-30 16:47:26 +02:00
2023-11-27 10:17:14 +01:00
write(*,*) ' *** Quasiparticle energies obtained by root search *** '
2023-06-30 16:47:26 +02:00
write(*,*)
2023-07-27 19:17:20 +02:00
call GTpp_QP_graph(eta,nBas,nC,nO,nV,nR,nOOs,nVVs,nOOt,nVVt,eHF,Om1s,rho1s,Om2s,rho2s, &
2023-10-02 21:41:37 +02:00
Om1t,rho1t,Om2t,rho2t,eGTlin,eHF,eGT,Z)
2020-03-21 21:09:48 +01:00
end if
2019-10-16 18:14:47 +02:00
2023-10-04 08:41:54 +02:00
! call GTpp_plot_self_energy(nBas,nC,nO,nV,nR,nOOs,nVVs,nOOt,nVVt,eHF,eGT,Om1s,rho1s,Om2s,rho2s, &
2023-09-07 22:28:47 +02:00
! Om1t,rho1t,Om2t,rho2t)
2019-10-20 08:18:58 +02:00
!----------------------------------------------
2019-10-16 18:14:47 +02:00
! Dump results
2019-10-20 08:18:58 +02:00
!----------------------------------------------
2019-10-16 18:14:47 +02:00
2020-04-13 14:19:14 +02:00
! Compute the ppRPA correlation energy
2020-04-13 23:05:19 +02:00
ispin = 1
2023-07-31 14:47:09 +02:00
! iblock = 1
2022-01-09 23:12:21 +01:00
iblock = 3
2022-01-08 13:59:45 +01:00
2023-11-29 16:20:53 +01:00
! allocate(Bpp(nVVs,nOOs),Cpp(nVVs,nVVs),Dpp(nOOs,nOOs))
2023-11-29 16:20:53 +01:00
! if(.not.TDA_T) call ppLR_B(iblock,nBas,nC,nO,nV,nR,nOOs,nVVs,1d0,ERI,Bpp)
! call ppLR_C(iblock,nBas,nC,nO,nV,nR,nVVs,1d0,eGT,ERI,Cpp)
! call ppLR_D(iblock,nBas,nC,nO,nV,nR,nOOs,1d0,eGT,ERI,Dpp)
2023-11-29 16:20:53 +01:00
! call ppLR(TDA_T,nOOs,nVVs,Bpp,Cpp,Dpp,Om1s,X1s,Y1s,Om2s,X2s,Y2s,EcRPA(ispin))
2023-11-29 16:20:53 +01:00
! deallocate(Bpp,Cpp,Dpp)
2020-04-13 23:05:19 +02:00
ispin = 2
2023-07-31 14:47:09 +02:00
! iblock = 2
2020-04-13 23:05:19 +02:00
iblock = 4
2022-09-09 21:48:50 +02:00
2023-11-29 16:20:53 +01:00
! allocate(Bpp(nVVt,nOOt),Cpp(nVVt,nVVt),Dpp(nOOt,nOOt))
2023-11-29 16:20:53 +01:00
! if(.not.TDA_T) call ppLR_B(iblock,nBas,nC,nO,nV,nR,nOOt,nVVt,1d0,ERI,Bpp)
! call ppLR_C(iblock,nBas,nC,nO,nV,nR,nVVt,1d0,eGT,ERI,Cpp)
! call ppLR_D(iblock,nBas,nC,nO,nV,nR,nOOt,1d0,eGT,ERI,Dpp)
2023-11-29 16:20:53 +01:00
! call ppLR(TDA_T,nOOt,nVVt,Bpp,Cpp,Dpp,Om1t,X1t,Y1t,Om2t,X2t,Y2t,EcRPA(ispin))
2023-11-29 16:20:53 +01:00
! deallocate(Bpp,Cpp,Dpp)
2022-01-02 10:24:30 +01:00
2023-11-29 16:20:53 +01:00
! EcRPA(1) = EcRPA(1) - EcRPA(2)
! EcRPA(2) = 3d0*EcRPA(2)
2020-04-13 14:19:14 +02:00
2023-11-20 17:42:50 +01:00
call print_RG0T0pp(nBas,nO,eHF,ENuc,ERHF,Sig,Z,eGT,EcGM,EcRPA)
2020-04-13 14:19:14 +02:00
2020-03-21 21:09:48 +01:00
! Perform BSE calculation
2023-07-21 12:30:29 +02:00
if(dophBSE) then
2020-03-21 21:09:48 +01:00
2023-07-27 19:17:20 +02:00
call GTpp_phBSE(TDA_T,TDA,dBSE,dTDA,singlet,triplet,eta,nBas,nC,nO,nV,nR,nS,nOOs,nVVs,nOOt,nVVt, &
Om1s,X1s,Y1s,Om2s,X2s,Y2s,rho1s,rho2s,Om1t,X1t,Y1t,Om2t,X2t,Y2t,rho1t,rho2t, &
2023-07-30 22:01:44 +02:00
ERI,dipole_int,eHF,eGT,EcBSE)
2020-03-21 21:09:48 +01:00
if(exchange_kernel) then
2023-07-21 12:30:29 +02:00
EcBSE(1) = 0.5d0*EcBSE(1)
EcBSE(2) = 1.5d0*EcBSE(1)
2020-03-21 21:09:48 +01:00
end if
write(*,*)
write(*,*)'-------------------------------------------------------------------------------'
2023-07-21 12:30:29 +02:00
write(*,'(2X,A50,F20.10,A3)') 'Tr@phBSE@G0T0pp correlation energy (singlet) =',EcBSE(1),' au'
write(*,'(2X,A50,F20.10,A3)') 'Tr@phBSE@G0T0pp correlation energy (triplet) =',EcBSE(2),' au'
2023-11-13 18:36:09 +01:00
write(*,'(2X,A50,F20.10,A3)') 'Tr@phBSE@G0T0pp correlation energy =',sum(EcBSE),' au'
write(*,'(2X,A50,F20.10,A3)') 'Tr@phBSE@G0T0pp total energy =',ENuc + ERHF + sum(EcBSE),' au'
2020-03-21 21:09:48 +01:00
write(*,*)'-------------------------------------------------------------------------------'
write(*,*)
! Compute the BSE correlation energy via the adiabatic connection
if(doACFDT) then
2023-07-21 12:30:29 +02:00
write(*,*) '--------------------------------------------------------'
write(*,*) 'Adiabatic connection version of phBSE correlation energy'
write(*,*) '--------------------------------------------------------'
2020-03-21 21:09:48 +01:00
write(*,*)
if(doXBS) then
write(*,*) '*** scaled screening version (XBS) ***'
write(*,*)
end if
2023-07-21 12:30:29 +02:00
call GTpp_phACFDT(exchange_kernel,doXBS,.false.,TDA_T,TDA,dophBSE,singlet,triplet,eta,nBas,nC,nO,nV,nR,nS, &
2023-07-27 19:17:20 +02:00
nOOs,nVVs,nOOt,nVVt,Om1s,X1s,Y1s,Om2s,X2s,Y2s,rho1s,rho2s,Om1t,X1t,Y1t, &
2023-07-30 22:01:44 +02:00
Om2t,X2t,Y2t,rho1t,rho2t,ERI,eHF,eGT,EcBSE)
2020-03-21 21:09:48 +01:00
if(exchange_kernel) then
2023-07-21 12:30:29 +02:00
EcBSE(1) = 0.5d0*EcBSE(1)
EcBSE(2) = 1.5d0*EcBSE(2)
2020-03-21 21:09:48 +01:00
end if
write(*,*)
write(*,*)'-------------------------------------------------------------------------------'
2023-11-29 16:20:53 +01:00
write(*,'(2X,A50,F20.10,A3)') 'AC@phBSE@G0T0pp@RHF correlation energy (singlet) = ',EcBSE(1),' au'
write(*,'(2X,A50,F20.10,A3)') 'AC@phBSE@G0T0pp@RHF correlation energy (triplet) = ',EcBSE(2),' au'
write(*,'(2X,A50,F20.10,A3)') 'AC@phBSE@G0T0pp@RHF correlation energy = ',sum(EcBSE),' au'
write(*,'(2X,A50,F20.10,A3)') 'AC@phBSE@G0T0pp@RHF total energy = ',ENuc + ERHF + sum(EcBSE),' au'
2020-03-21 21:09:48 +01:00
write(*,*)'-------------------------------------------------------------------------------'
write(*,*)
end if
end if
2023-07-21 12:30:29 +02:00
if(doppBSE) then
2022-09-09 21:48:50 +02:00
2023-07-27 19:17:20 +02:00
call GTpp_ppBSE(TDA_T,TDA,dBSE,dTDA,singlet,triplet,eta,nBas,nC,nO,nV,nR,nOOs,nVVs,nOOt,nVVt, &
Om1s,X1s,Y1s,Om2s,X2s,Y2s,rho1s,rho2s,Om1t,X1t,Y1t,Om2t,X2t,Y2t,rho1t,rho2t, &
2023-07-30 22:01:44 +02:00
ERI,dipole_int,eHF,eGT,EcBSE)
2022-09-09 21:48:50 +02:00
write(*,*)
write(*,*)'-------------------------------------------------------------------------------'
2023-11-29 16:20:53 +01:00
write(*,'(2X,A50,F20.10,A3)') 'Tr@ppBSE@G0T0pp@RHF correlation energy (singlet) = ',EcBSE(1),' au'
write(*,'(2X,A50,F20.10,A3)') 'Tr@ppBSE@G0T0pp@RHF correlation energy (triplet) = ',EcBSE(2),' au'
write(*,'(2X,A50,F20.10,A3)') 'Tr@ppBSE@G0T0pp@RHF correlation energy = ',sum(EcBSE),' au'
write(*,'(2X,A50,F20.10,A3)') 'Tr@ppBSE@G0T0pp@RHF total energy = ',ENuc + ERHF + sum(EcBSE),' au'
2022-09-09 21:48:50 +02:00
write(*,*)'-------------------------------------------------------------------------------'
write(*,*)
end if
2023-11-13 18:36:09 +01:00
! Testing zone
if(dotest) then
2023-11-14 14:31:27 +01:00
call dump_test_value('R','G0T0pp correlation energy',sum(EcRPA))
call dump_test_value('R','G0T0pp HOMO energy',eGT(nO))
call dump_test_value('R','G0T0pp LUMO energy',eGT(nO+1))
2023-11-13 18:36:09 +01:00
end if
2023-07-04 10:32:47 +02:00
end subroutine