4
1
mirror of https://github.com/pfloos/quack synced 2024-06-26 07:02:17 +02:00
quack/src/QuAcK/linear_response_pp.f90

101 lines
3.2 KiB
Fortran
Raw Normal View History

2019-10-07 22:31:45 +02:00
subroutine linear_response_pp(ispin,BSE,nBas,nC,nO,nV,nR,nOO,nVV,e,ERI,Omega1,X1,Y1,Omega2,X2,Y2,Ec_ppRPA)
2019-10-05 22:06:25 +02:00
2019-10-06 20:08:38 +02:00
! Compute the p-p channel of the linear response: see Scueria et al. JCP 139, 104113 (2013)
2019-10-05 22:06:25 +02:00
implicit none
include 'parameters.h'
! Input variables
logical,intent(in) :: BSE
2019-10-06 22:35:36 +02:00
integer,intent(in) :: ispin,nBas,nC,nO,nV,nR
integer,intent(in) :: nOO
integer,intent(in) :: nVV
2019-10-05 22:06:25 +02:00
double precision,intent(in) :: e(nBas)
double precision,intent(in) :: ERI(nBas,nBas,nBas,nBas)
! Local variables
double precision :: trace_matrix
double precision,allocatable :: B(:,:)
double precision,allocatable :: C(:,:)
double precision,allocatable :: D(:,:)
double precision,allocatable :: M(:,:)
2019-10-06 20:08:38 +02:00
double precision,allocatable :: Z(:,:)
2019-10-06 22:35:36 +02:00
double precision,allocatable :: Omega(:)
2019-10-05 22:06:25 +02:00
! Output variables
2019-10-06 22:35:36 +02:00
double precision,intent(out) :: Omega1(nVV)
double precision,intent(out) :: X1(nVV,nVV)
double precision,intent(out) :: Y1(nOO,nVV)
double precision,intent(out) :: Omega2(nOO)
double precision,intent(out) :: X2(nVV,nOO)
double precision,intent(out) :: Y2(nOO,nOO)
2019-10-05 23:09:20 +02:00
double precision,intent(out) :: Ec_ppRPA
2019-10-05 22:06:25 +02:00
! Memory allocation
2019-10-06 22:35:36 +02:00
allocate(B(nVV,nOO),C(nVV,nVV),D(nOO,nOO),M(nOO+nVV,nOO+nVV),Z(nOO+nVV,nOO+nVV),Omega(nOO+nVV))
2019-10-05 22:06:25 +02:00
! Build B, C and D matrices for the pp channel
2019-10-07 22:31:45 +02:00
call linear_response_B_pp(ispin,nBas,nC,nO,nV,nR,nOO,nVV,e,ERI,B)
call linear_response_C_pp(ispin,nBas,nC,nO,nV,nR,nOO,nVV,e,ERI,C)
call linear_response_D_pp(ispin,nBas,nC,nO,nV,nR,nOO,nVV,e,ERI,D)
2019-10-05 22:06:25 +02:00
!------------------------------------------------------------------------
! Solve the p-p eigenproblem
!------------------------------------------------------------------------
!
! | C -B | | X1 X2 | | w1 0 | | X1 X2 |
! | | | | = | | | |
! | Bt -D | | Y1 Y2 | | 0 w2 | | Y1 Y2 |
!
! Diagonal blocks
2019-10-06 20:08:38 +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)
2019-10-05 22:06:25 +02:00
! Off-diagonal blocks
2019-10-05 23:09:20 +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))
2019-10-05 22:06:25 +02:00
2019-10-05 23:09:20 +02:00
! print*, 'pp-RPA matrix'
! call matout(nOO+nVV,nOO+nVV,M(:,:))
2019-10-05 22:06:25 +02:00
! Diagonalize the p-h matrix
2019-10-06 20:08:38 +02:00
Z(:,:) = M(:,:)
2019-10-15 17:56:27 +02:00
call diagonalize_general_matrix(nOO+nVV,M(:,:),Omega(:),Z(:,:))
2019-10-05 22:06:25 +02:00
2019-10-15 23:13:00 +02:00
! write(*,*) 'pp-RPA excitation energies'
! call matout(nOO+nVV,1,Omega(:))
! write(*,*)
2019-10-05 22:06:25 +02:00
2019-10-15 23:13:00 +02:00
! call matout(nOO+nVV,nOO+nVV,matmul(transpose(Z(:,:)),matmul(W(:,:),Z(:,:))))
! write(*,*)
2019-10-05 22:06:25 +02:00
2019-10-15 23:13:00 +02:00
! Split the various quantities in p-p and h-h parts
2019-10-05 22:06:25 +02:00
2019-10-15 23:13:00 +02:00
call sort_ppRPA(nOO,nVV,Omega(:),Z(:,:),Omega1(:),X1(:,:),Y1(:,:),Omega2(:),X2(:,:),Y2(:,:))
! Omega1(:) = Omega(nOO+1:nOO+nVV)
! Omega2(:) = Omega(1:nOO)
2019-10-06 20:08:38 +02:00
2019-10-15 23:13:00 +02:00
! X1(:,:) = Z(nOO+1:nOO+nVV,nOO+1:nOO+nVV)
! Y1(:,:) = Z( 1:nOO ,nOO+1:nOO+nVV)
! X2(:,:) = Z(nOO+1:nOO+nVV, 1:nOO )
! Y2(:,:) = Z( 1:nOO ,nOO+1:nOO+nVV)
2019-10-05 22:06:25 +02:00
! Compute the RPA correlation energy
2019-10-15 23:20:43 +02:00
! Ec_ppRPA = 0.5d0*( sum(Omega1(:)) - sum(Omega2(:)) - trace_matrix(nVV,C(:,:)) - trace_matrix(nOO,D(:,:)) )
Ec_ppRPA = +sum(Omega1(:)) - trace_matrix(nVV,C(:,:))
! Ec_ppRPA = -sum(Omega2(:)) - trace_matrix(nOO,D(:,:))
2019-10-05 22:06:25 +02:00
print*,'Ec(pp-RPA) = ',Ec_ppRPA
end subroutine linear_response_pp