4
1
mirror of https://github.com/pfloos/quack synced 2024-06-01 19:05:27 +02:00
quack/src/QuAcK/linear_response_pp.f90

159 lines
4.4 KiB
Fortran
Raw Normal View History

2020-03-21 16:31:39 +01:00
subroutine linear_response_pp(ispin,ortho_eigvec,BSE,nBas,nC,nO,nV,nR,nOO,nVV, &
2020-03-14 23:00:44 +01:00
e,ERI,Omega1,X1,Y1,Omega2,X2,Y2,EcppRPA)
2019-10-05 22:06:25 +02:00
2020-03-14 23:00:44 +01:00
! Compute the p-p 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
2020-03-21 16:31:39 +01:00
logical,intent(in) :: ortho_eigvec
2019-10-05 22:06:25 +02:00
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
2020-03-14 23:00:44 +01:00
logical :: dump_matrices = .false.
2019-10-30 10:31:05 +01:00
integer :: ab,cd,ij,kl
2020-03-11 16:41:20 +01:00
integer :: p,q,r,s
2019-10-05 22:06:25 +02:00
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-22 23:03:01 +02:00
double precision,intent(out) :: EcppRPA
2019-10-05 23:09:20 +02:00
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-24 08:07:31 +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
2020-03-14 23:00:44 +01:00
! Dump ppRPA matrices
if(dump_matrices) then
open(unit=42,file='B.dat')
open(unit=43,file='C.dat')
open(unit=44,file='D.dat')
open(unit=45,file='ERI.dat')
open(unit=46,file='eps.dat')
do ab=1,nVV
do ij=1,nOO
if(abs(B(ab,ij)) > 1d-15) write(42,*) ab,ij,B(ab,ij)
end do
2020-03-11 16:41:20 +01:00
end do
2020-03-14 23:00:44 +01:00
do ab=1,nVV
do cd=1,nVV
if(abs(C(ab,cd)) > 1d-15) write(43,*) ab,cd,C(ab,cd)
end do
2020-03-11 16:41:20 +01:00
end do
2020-03-14 23:00:44 +01:00
do ij=1,nOO
do kl=1,nOO
if(abs(D(ij,kl)) > 1d-15) write(44,*) ij,kl,D(ij,kl)
2020-03-11 16:41:20 +01:00
end do
end do
2020-03-14 23:00:44 +01:00
do p=1,nBas
write(46,*) p,e(p)
do q=1,nBas
do r=1,nBas
do s=1,nBas
if(abs(ERI(p,q,r,s)) > 1d-15) write(45,*) p,q,r,s,ERI(p,q,r,s)
end do
end do
end do
2020-03-11 16:41:20 +01:00
end do
2020-03-14 23:00:44 +01:00
close(42)
close(43)
close(44)
close(45)
close(46)
2019-10-30 10:31:05 +01:00
2020-03-14 23:00:44 +01:00
end if
2019-10-28 16:34:09 +01: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-28 16:34:09 +01: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
! Split the various quantities in p-p and h-h parts
2019-10-05 22:06:25 +02:00
2020-03-21 16:31:39 +01:00
call sort_ppRPA(ortho_eigvec,nOO,nVV,Omega(:),Z(:,:),Omega1(:),X1(:,:),Y1(:,:),Omega2(:),X2(:,:),Y2(:,:))
2019-10-05 22:06:25 +02:00
! Compute the RPA correlation energy
2019-10-22 23:03:01 +02:00
! EcppRPA = 0.5d0*( sum(Omega1(:)) - sum(Omega2(:)) - trace_matrix(nVV,C(:,:)) - trace_matrix(nOO,D(:,:)) )
2020-03-21 16:31:39 +01:00
! print*,+sum(Omega1(:)),- trace_matrix(nVV,C(:,:))
! print*,-sum(Omega2(:)),- trace_matrix(nOO,D(:,:))
2019-10-22 23:03:01 +02:00
EcppRPA = +sum(Omega1(:)) - trace_matrix(nVV,C(:,:))
! EcppRPA = -sum(Omega2(:)) - trace_matrix(nOO,D(:,:))
2019-10-05 22:06:25 +02:00
2019-10-18 23:16:37 +02:00
! write(*,*)'X1'
! call matout(nVV,nVV,X1)
! write(*,*)'Y1'
! call matout(nVV,nOO,Y1)
! write(*,*)'X2'
! call matout(nOO,nVV,X2)
! write(*,*)'Y2'
! call matout(nOO,nOO,Y2)
2019-10-22 23:03:01 +02:00
! print*,'Ec(pp-RPA) = ',EcppRPA
2019-10-05 22:06:25 +02:00
2019-10-28 16:34:09 +01:00
! print*,'Eigenvalues'
! call matout(nOO+nVV,1,Omega)
! print*,'Eigenvectors'
! call matout(nOO+nVV,nOO+nVV,matmul(transpose(Z),Z))
2019-10-05 22:06:25 +02:00
end subroutine linear_response_pp