4
1
mirror of https://github.com/pfloos/quack synced 2024-11-06 22:24:03 +01:00
quack/src/QuAcK/exchange_matrix_AO_basis.f90

34 lines
659 B
Fortran
Raw Normal View History

2019-03-19 10:13:33 +01:00
subroutine exchange_matrix_AO_basis(nBas,P,G,K)
! Compute exchange matrix in the AO basis
implicit none
include 'parameters.h'
! Input variables
integer,intent(in) :: nBas
double precision,intent(in) :: P(nBas,nBas)
double precision,intent(in) :: G(nBas,nBas,nBas,nBas)
! Local variables
integer :: mu,nu,la,si
! Output variables
double precision,intent(out) :: K(nBas,nBas)
K = 0d0
do nu=1,nBas
do si=1,nBas
do la=1,nBas
do mu=1,nBas
2019-03-19 11:21:34 +01:00
K(mu,nu) = K(mu,nu) - P(la,si)*G(mu,la,si,nu)
2019-03-19 10:13:33 +01:00
enddo
enddo
enddo
enddo
end subroutine exchange_matrix_AO_basis