4
1
mirror of https://github.com/pfloos/quack synced 2024-06-14 17:25:33 +02:00
quack/src/QuAcK/linear_response.f90

99 lines
2.2 KiB
Fortran
Raw Normal View History

2020-01-08 10:17:19 +01:00
subroutine linear_response(ispin,dRPA,TDA,BSE,nBas,nC,nO,nV,nR,nS,lambda,e,ERI,rho,EcRPA,Omega,XpY)
2019-03-19 10:13:33 +01:00
! Compute linear response
implicit none
include 'parameters.h'
! Input variables
logical,intent(in) :: dRPA,TDA,BSE
integer,intent(in) :: ispin,nBas,nC,nO,nV,nR,nS
2020-01-08 10:17:19 +01:00
double precision,intent(in) :: lambda
double precision,intent(in) :: e(nBas)
double precision,intent(in) :: rho(nBas,nBas,nS)
double precision,intent(in) :: ERI(nBas,nBas,nBas,nBas)
2019-03-19 10:13:33 +01:00
! Local variables
double precision :: trace_matrix
double precision,allocatable :: A(:,:),B(:,:),ApB(:,:),AmB(:,:),AmBSq(:,:),Z(:,:)
! Output variables
double precision,intent(out) :: EcRPA
double precision,intent(out) :: Omega(nS),XpY(nS,nS)
! Memory allocation
allocate(A(nS,nS),B(nS,nS),ApB(nS,nS),AmB(nS,nS),AmBSq(nS,nS),Z(nS,nS))
! Build A and B matrices
2020-01-08 10:17:19 +01:00
call linear_response_A_matrix(ispin,dRPA,nBas,nC,nO,nV,nR,nS,lambda,e,ERI,A)
if(BSE) call Bethe_Salpeter_A_matrix(nBas,nC,nO,nV,nR,nS,lambda,ERI,Omega,rho,A)
2019-03-19 10:13:33 +01:00
! Tamm-Dancoff approximation
B = 0d0
if(.not. TDA) then
2020-01-08 10:17:19 +01:00
call linear_response_B_matrix(ispin,dRPA,nBas,nC,nO,nV,nR,nS,lambda,ERI,B)
if(BSE) call Bethe_Salpeter_B_matrix(nBas,nC,nO,nV,nR,nS,lambda,ERI,Omega,rho,B)
2019-03-19 10:13:33 +01:00
endif
! Build A + B and A - B matrices
ApB = A + B
2019-04-29 09:43:33 +02:00
AmB = A - B
2019-03-19 10:13:33 +01:00
2019-05-23 09:53:23 +02:00
! print*,'A'
! call matout(nS,nS,A)
! print*,'B'
! call matout(nS,nS,B)
2019-03-19 10:13:33 +01:00
! print*,'A+B'
! call matout(nS,nS,ApB)
! print*,'A-B'
! call matout(nS,nS,AmB)
! Diagonalize TD-HF matrix
call diagonalize_matrix(nS,AmB,Omega)
if(minval(Omega) < 0d0) &
call print_warning('You may have instabilities in linear response!!')
call ADAt(nS,AmB,sqrt(Omega),AmBSq)
Z = matmul(AmBSq,matmul(ApB,AmBSq))
2019-04-29 09:43:33 +02:00
! print*,'Z'
! call matout(nS,nS,Z)
2019-03-19 10:13:33 +01:00
call diagonalize_matrix(nS,Z,Omega)
if(minval(Omega) < 0d0) &
call print_warning('You may have instabilities in linear response!!')
Omega = sqrt(Omega)
XpY = matmul(transpose(Z),AmBSq)
2019-10-28 16:34:09 +01:00
call DA(nS,1d0/sqrt(abs(Omega)),XpY)
2019-03-19 10:13:33 +01:00
2019-04-29 09:43:33 +02:00
! print*,'X+Y'
! call matout(nS,nS,XpY)
2019-03-19 10:13:33 +01:00
! print*,'RPA excitations'
! call matout(nS,1,Omega)
! Compute the RPA correlation energy
EcRPA = 0.5d0*(sum(Omega) - trace_matrix(nS,A))
2019-10-24 08:07:31 +02:00
! print*,'EcRPA = ',EcRPA
2019-05-23 09:53:23 +02:00
2019-03-19 10:13:33 +01:00
end subroutine linear_response