10
1
mirror of https://github.com/pfloos/quack synced 2024-10-20 23:08:23 +02:00
QuAcK/src/LR/ppLR.f90

92 lines
2.7 KiB
Fortran
Raw Normal View History

subroutine ppLR(TDA,nOO,nVV,B,C,D,Om1,X1,Y1,Om2,X2,Y2,EcRPA)
2019-10-05 22:06:25 +02:00
2023-07-18 09:59:25 +02:00
! Compute the pp channel of the linear response: see Scuseria et al. JCP 139, 104113 (2013)
2019-10-05 22:06:25 +02:00
implicit none
include 'parameters.h'
! Input variables
2021-10-18 22:05:26 +02:00
logical,intent(in) :: TDA
2019-10-06 22:35:36 +02:00
integer,intent(in) :: nOO
integer,intent(in) :: nVV
2023-07-17 22:31:28 +02:00
double precision,intent(in) :: B(nVV,nOO)
double precision,intent(in) :: C(nVV,nVV)
double precision,intent(in) :: D(nOO,nOO)
2019-10-05 22:06:25 +02:00
! Local variables
double precision :: trace_matrix
2020-03-24 12:28:56 +01:00
double precision :: EcRPA1
double precision :: EcRPA2
2019-10-05 22:06:25 +02:00
double precision,allocatable :: M(:,:)
2019-10-06 20:08:38 +02:00
double precision,allocatable :: Z(:,:)
2023-07-17 22:31:28 +02:00
double precision,allocatable :: Om(:)
2019-10-05 22:06:25 +02:00
! Output variables
2023-07-17 22:31:28 +02:00
double precision,intent(out) :: Om1(nVV)
2019-10-06 22:35:36 +02:00
double precision,intent(out) :: X1(nVV,nVV)
double precision,intent(out) :: Y1(nOO,nVV)
2023-07-17 22:31:28 +02:00
double precision,intent(out) :: Om2(nOO)
2019-10-06 22:35:36 +02:00
double precision,intent(out) :: X2(nVV,nOO)
double precision,intent(out) :: Y2(nOO,nOO)
2020-03-24 12:28:56 +01:00
double precision,intent(out) :: EcRPA
2019-10-05 23:09:20 +02:00
2019-10-05 22:06:25 +02:00
! Memory allocation
2023-07-17 22:31:28 +02:00
allocate(M(nOO+nVV,nOO+nVV),Z(nOO+nVV,nOO+nVV),Om(nOO+nVV))
2021-10-18 22:05:26 +02:00
!-------------------------------------------------!
! Solve the p-p eigenproblem !
!-------------------------------------------------!
! !
2021-11-06 14:36:54 +01:00
! | C B | | X1 X2 | | w1 0 | | X1 X2 | !
2021-10-18 22:05:26 +02:00
! | | | | = | | | | !
2021-11-06 14:36:54 +01:00
! | -Bt -D | | Y1 Y2 | | 0 w2 | | Y1 Y2 | !
2021-10-18 22:05:26 +02:00
! !
!-------------------------------------------------!
if(TDA) then
2019-10-05 22:06:25 +02:00
2021-10-18 22:05:26 +02:00
X1(:,:) = +C(:,:)
Y1(:,:) = 0d0
2023-07-17 22:31:28 +02:00
if(nVV > 0) call diagonalize_matrix(nVV,X1,Om1)
2022-09-09 21:48:50 +02:00
2021-10-18 22:05:26 +02:00
X2(:,:) = 0d0
Y2(:,:) = -D(:,:)
2023-07-17 22:31:28 +02:00
if(nOO > 0) call diagonalize_matrix(nOO,Y2,Om2)
2019-10-05 22:06:25 +02:00
2021-10-18 22:05:26 +02:00
else
2019-10-05 22:06:25 +02:00
2021-10-18 22:05:26 +02:00
! Diagonal blocks
2019-10-05 22:06:25 +02:00
2021-10-18 22:05:26 +02:00
M( 1:nVV , 1:nVV) = + C(1:nVV,1:nVV)
M(nVV+1:nVV+nOO,nVV+1:nVV+nOO) = - D(1:nOO,1:nOO)
2020-03-14 23:00:44 +01:00
2021-10-18 22:05:26 +02:00
! Off-diagonal blocks
2020-03-14 23:00:44 +01:00
2021-10-18 22:05:26 +02:00
M( 1:nVV ,nVV+1:nOO+nVV) = - B(1:nVV,1:nOO)
M(nVV+1:nOO+nVV, 1:nVV) = + transpose(B(1:nVV,1:nOO))
2022-08-17 14:32:14 +02:00
! Diagonalize the p-p matrix
2019-10-05 22:06:25 +02:00
2023-07-17 22:31:28 +02:00
if(nOO+nVV > 0) call diagonalize_general_matrix(nOO+nVV,M,Om,Z)
2019-10-05 22:06:25 +02:00
2021-10-18 22:05:26 +02:00
! Split the various quantities in p-p and h-h parts
2020-04-16 17:02:01 +02:00
2023-07-18 09:59:25 +02:00
call sort_ppRPA(nOO,nVV,Om,Z,Om1,X1,Y1,Om2,X2,Y2)
2019-10-05 22:06:25 +02:00
2021-10-18 22:05:26 +02:00
end if
2020-04-16 17:02:01 +02:00
2019-10-05 22:06:25 +02:00
! Compute the RPA correlation energy
2023-07-18 09:59:25 +02:00
EcRPA = 0.5d0*( sum(Om1) - sum(Om2) - trace_matrix(nVV,C) - trace_matrix(nOO,D) )
EcRPA1 = +sum(Om1) - trace_matrix(nVV,C)
EcRPA2 = -sum(Om2) - trace_matrix(nOO,D)
2022-08-17 14:32:14 +02:00
2020-03-25 21:10:21 +01:00
if(abs(EcRPA - EcRPA1) > 1d-6 .or. abs(EcRPA - EcRPA2) > 1d-6) &
2020-03-24 12:28:56 +01:00
print*,'!!! Issue in pp-RPA linear reponse calculation RPA1 != RPA2 !!!'
2019-10-05 22:06:25 +02:00
end subroutine