10
1
mirror of https://github.com/pfloos/quack synced 2025-01-10 13:08:19 +01:00
QuAcK/src/GW/GW_ppBSE_static_kernel_C.f90

371 lines
8.7 KiB
Fortran
Raw Normal View History

2023-07-27 10:11:35 +02:00
subroutine GW_ppBSE_static_kernel_C(ispin,eta,nBas,nC,nO,nV,nR,nS,nVV,lambda,ERI,Om,rho,KC)
2022-08-17 14:32:14 +02:00
! Compute the VVVV block of the static screening W for the pp-BSE
implicit none
include 'parameters.h'
! Input variables
integer,intent(in) :: ispin
integer,intent(in) :: nBas
integer,intent(in) :: nC
integer,intent(in) :: nO
integer,intent(in) :: nV
integer,intent(in) :: nR
integer,intent(in) :: nS
integer,intent(in) :: nVV
double precision,intent(in) :: eta
double precision,intent(in) :: lambda
double precision,intent(in) :: ERI(nBas,nBas,nBas,nBas)
2023-07-21 13:04:29 +02:00
double precision,intent(in) :: Om(nS)
2022-08-17 14:32:14 +02:00
double precision,intent(in) :: rho(nBas,nBas,nS)
! Local variables
2022-08-17 16:53:59 +02:00
double precision,external :: Kronecker_delta
2022-08-17 14:32:14 +02:00
double precision :: chi
double precision :: eps
2024-08-22 15:23:08 +02:00
double precision :: tmp_ab, lambda4, eta2
2022-08-17 14:32:14 +02:00
integer :: a,b,c,d,ab,cd,m
2024-08-22 15:23:08 +02:00
integer :: a0, aa
double precision, allocatable :: Om_tmp(:)
double precision, allocatable :: tmp_m(:,:,:)
double precision, allocatable :: tmp(:,:,:,:)
2022-08-17 14:32:14 +02:00
! Output variables
2023-07-24 15:01:54 +02:00
double precision,intent(out) :: KC(nVV,nVV)
2022-08-17 14:32:14 +02:00
!---------------!
! Singlet block !
!---------------!
if(ispin == 1) then
2024-08-22 15:23:08 +02:00
a0 = nBas - nR - nO
lambda4 = 4.d0 * lambda
eta2 = eta * eta
allocate(tmp_m(nBas,nBas,nS))
allocate(tmp(nBas,nBas,nBas,nBas))
2024-08-22 15:23:08 +02:00
do m = 1, nS
eps = Om(m) / (Om(m)*Om(m) + eta2)
do c = 1, nBas
do a = 1, nBas
tmp_m(a,c,m) = eps * rho(a,c,m)
enddo
enddo
2024-08-22 15:23:08 +02:00
enddo
call dgemm("N", "T", nBas*nBas, nBas*nBas, nS, 1.d0, &
tmp_m(1,1,1), nBas*nBas, rho(1,1,1), nBas*nBas, &
0.d0, tmp(1,1,1,1), nBas*nBas)
deallocate(tmp_m)
!$OMP PARALLEL DEFAULT(NONE) &
!$OMP PRIVATE(a, b, aa, ab, c, d, cd, tmp_ab) &
!$OMP SHARED(nO, nBas, nR, nS, a0, lambda4, tmp, KC)
2024-08-22 15:23:08 +02:00
!$OMP DO
do a = nO+1, nBas-nR
aa = a0 * (a - nO - 1) - (a - nO - 1) * (a - nO) / 2 - nO
do b = a, nBas-nR
ab = aa + b
tmp_ab = lambda4
if(a .eq. b) then
tmp_ab = 0.7071067811865475d0 * lambda4
endif
2022-08-17 14:32:14 +02:00
cd = 0
2024-08-22 16:44:16 +02:00
do c = nO+1, nBas-nR
2024-08-22 15:23:08 +02:00
do d = c, nBas-nR
2022-09-09 21:48:50 +02:00
cd = cd + 1
KC(ab,cd) = -tmp_ab * (tmp(a,c,b,d) + tmp(a,d,b,c))
2024-08-22 15:23:08 +02:00
if(c .eq. d) then
KC(ab,cd) = 0.7071067811865475d0 * KC(ab,cd)
endif
enddo
enddo
enddo
enddo
!$OMP END DO
!$OMP END PARALLEL
deallocate(tmp)
2024-08-22 16:44:16 +02:00
! do a=nO+1,nBas-nR
! do b=a,nBas-nR
! ab = ab + 1
! cd = 0
! do c=nO+1,nBas-nR
! do d=c,nBas-nR
! cd = cd + 1
!
! chi = 0d0
! do m=1,nS
! eps = Om(m)**2 + eta**2
! chi = chi - rho(a,c,m)*rho(b,d,m)*Om(m)/eps &
! - rho(a,d,m)*rho(b,c,m)*Om(m)/eps
! end do
! --- --- ---
! OpenMP implementation
! --- --- ---
!
! a0 = nBas - nR - nO
! lambda4 = 4.d0 * lambda
! eta2 = eta * eta
!
! allocate(Om_tmp(nS))
!
! !$OMP PARALLEL DEFAULT(NONE) PRIVATE(m) SHARED(nS, eta2, Om, Om_tmp)
! !$OMP DO
! do m = 1, nS
! Om_tmp(m) = Om(m) / (Om(m)*Om(m) + eta2)
! enddo
! !$OMP END DO
! !$OMP END PARALLEL
!
! !$OMP PARALLEL DEFAULT(NONE) &
! !$OMP PRIVATE(a, b, aa, ab, c, d, cd, m, tmp_ab) &
! !$OMP SHARED(nO, nBas, nR, nS, a0, lambda4, Om_tmp, rho, KC)
! !$OMP DO
! do a = nO+1, nBas-nR
! aa = a0 * (a - nO - 1) - (a - nO - 1) * (a - nO) / 2 - nO
! do b = a, nBas-nR
! ab = aa + b
!
! tmp_ab = lambda4
! if(a .eq. b) then
! tmp_ab = 0.7071067811865475d0 * lambda4
! endif
!
! cd = 0
! do c = nO+1, nBas-nR
! do d = c, nBas-nR
! cd = cd + 1
!
! KC(ab,cd) = 0d0
! do m = 1, nS
! KC(ab,cd) = KC(ab,cd) - rho(a,c,m) * rho(b,d,m) * Om_tmp(m) &
! - rho(a,d,m) * rho(b,c,m) * Om_tmp(m)
! end do
!
! KC(ab,cd) = tmp_ab * KC(ab,cd)
! if(c .eq. d) then
! KC(ab,cd) = 0.7071067811865475d0 * KC(ab,cd)
! endif
! enddo
! enddo
! enddo
! enddo
! !$OMP END DO
! !$OMP END PARALLEL
!
! deallocate(Om_tmp)
! --- --- ---
! --- --- ---
! Naive implementation
! --- --- ---
!
2024-08-22 15:23:08 +02:00
! ab = 0
! do a=nO+1,nBas-nR
! do b=a,nBas-nR
! ab = ab + 1
! cd = 0
! do c=nO+1,nBas-nR
! do d=c,nBas-nR
! cd = cd + 1
!
! chi = 0d0
! do m=1,nS
! eps = Om(m)**2 + eta**2
! chi = chi - rho(a,c,m)*rho(b,d,m)*Om(m)/eps &
! - rho(a,d,m)*rho(b,c,m)*Om(m)/eps
! end do
!
! KC(ab,cd) = 4d0*lambda*chi/sqrt((1d0 + Kronecker_delta(a,b))*(1d0 + Kronecker_delta(c,d)))
!
! end do
! end do
! end do
! end do
! --- --- ---
2022-09-09 21:48:50 +02:00
end if
!---------------!
! Triplet block !
!---------------!
if(ispin == 2) then
2024-08-22 16:44:16 +02:00
a0 = nBas - nR - nO - 1
lambda4 = 4.d0 * lambda
eta2 = eta * eta
allocate(tmp_m(nBas,nBas,nS))
allocate(tmp(nBas,nBas,nBas,nBas))
2024-08-22 16:44:16 +02:00
do m = 1, nS
eps = Om(m) / (Om(m)*Om(m) + eta2)
do c = 1, nBas
do a = 1, nBas
tmp_m(a,c,m) = eps * rho(a,c,m)
enddo
enddo
2024-08-22 16:44:16 +02:00
enddo
call dgemm("N", "T", nBas*nBas, nBas*nBas, nS, 1.d0, &
tmp_m(1,1,1), nBas*nBas, rho(1,1,1), nBas*nBas, &
0.d0, tmp(1,1,1,1), nBas*nBas)
deallocate(tmp_m)
!$OMP PARALLEL DEFAULT(NONE) &
!$OMP PRIVATE(a, b, aa, ab, c, d, cd) &
!$OMP SHARED(nO, nBas, nR, nS, a0, lambda4, tmp, KC)
2024-08-22 16:44:16 +02:00
!$OMP DO
do a = nO+1, nBas-nR
aa = a0 * (a - nO - 1) - (a - nO - 1) * (a - nO) / 2 - nO - 1
do b = a+1, nBas-nR
ab = aa + b
2022-09-09 21:48:50 +02:00
cd = 0
2024-08-22 16:44:16 +02:00
do c = nO+1, nBas-nR
do d = c+1, nBas-nR
2022-08-17 14:32:14 +02:00
cd = cd + 1
KC(ab,cd) = lambda4 * (-tmp(a,c,b,d) + tmp(a,d,b,c))
2024-08-22 16:44:16 +02:00
enddo
enddo
enddo
enddo
!$OMP END DO
!$OMP END PARALLEL
deallocate(tmp)
! --- --- ---
! OpenMP implementation
! --- --- ---
!
! a0 = nBas - nR - nO - 1
! lambda4 = 4.d0 * lambda
! eta2 = eta * eta
!
! allocate(Om_tmp(nS))
!
! !$OMP PARALLEL DEFAULT(NONE) PRIVATE(m) SHARED(nS, eta2, Om, Om_tmp)
! !$OMP DO
! do m = 1, nS
! Om_tmp(m) = Om(m) / (Om(m)*Om(m) + eta2)
! enddo
! !$OMP END DO
! !$OMP END PARALLEL
!
! !$OMP PARALLEL DEFAULT(NONE) &
! !$OMP PRIVATE(a, b, aa, ab, c, d, cd, m) &
! !$OMP SHARED(nO, nBas, nR, nS, a0, lambda4, Om_tmp, rho, KC)
! !$OMP DO
! do a = nO+1, nBas-nR
! aa = a0 * (a - nO - 1) - (a - nO - 1) * (a - nO) / 2 - nO - 1
! do b = a+1, nBas-nR
! ab = aa + b
!
! cd = 0
! do c = nO+1, nBas-nR
! do d = c+1, nBas-nR
! cd = cd + 1
!
! KC(ab,cd) = 0d0
! do m = 1, nS
! KC(ab,cd) = KC(ab,cd) - rho(a,c,m) * rho(b,d,m) * Om_tmp(m) &
! + rho(a,d,m) * rho(b,c,m) * Om_tmp(m)
! end do
!
! KC(ab,cd) = lambda4 * KC(ab,cd)
! enddo
! enddo
! enddo
! enddo
! !$OMP END DO
! !$OMP END PARALLEL
!
! deallocate(Om_tmp)
! --- --- ---
2024-08-22 16:44:16 +02:00
! --- --- ---
! Naive implementation
! --- --- ---
!
2024-08-22 16:44:16 +02:00
! ab = 0
! do a=nO+1,nBas-nR
! do b=a+1,nBas-nR
! ab = ab + 1
! cd = 0
! do c=nO+1,nBas-nR
! do d=c+1,nBas-nR
! cd = cd + 1
!
! chi = 0d0
! do m=1,nS
! eps = Om(m)**2 + eta**2
! chi = chi - rho(a,c,m)*rho(b,d,m)*Om(m)/eps &
! + rho(a,d,m)*rho(b,c,m)*Om(m)/eps
! end do
!
! KC(ab,cd) = 4d0*lambda*chi
!
! end do
! end do
! end do
! end do
! --- --- ---
2022-08-17 14:32:14 +02:00
end if
!---------------!
2022-09-09 21:48:50 +02:00
! SpinOrb block !
2022-08-17 14:32:14 +02:00
!---------------!
2022-09-09 21:48:50 +02:00
if(ispin == 4) then
2022-08-17 14:32:14 +02:00
ab = 0
do a=nO+1,nBas-nR
2022-08-17 15:52:13 +02:00
do b=a+1,nBas-nR
2022-08-17 14:32:14 +02:00
ab = ab + 1
cd = 0
do c=nO+1,nBas-nR
2022-08-17 15:52:13 +02:00
do d=c+1,nBas-nR
2022-08-17 14:32:14 +02:00
cd = cd + 1
chi = 0d0
do m=1,nS
2023-07-21 13:04:29 +02:00
eps = Om(m)**2 + eta**2
2023-07-24 15:01:54 +02:00
chi = chi - rho(a,c,m)*rho(b,d,m)*Om(m)/eps &
+ rho(a,d,m)*rho(b,c,m)*Om(m)/eps
2023-12-03 18:47:30 +01:00
end do
2022-08-17 14:32:14 +02:00
2024-07-02 14:41:12 +02:00
KC(ab,cd) = 2d0*lambda*chi
2022-08-17 14:32:14 +02:00
end do
end do
end do
end do
end if
2023-07-20 15:54:24 +02:00
end subroutine