4
1
mirror of https://github.com/pfloos/quack synced 2024-12-22 20:35:36 +01:00

OpenMP -> DGEMM for ispin=2,4 in GTpp_excitation_density

This commit is contained in:
Abdallah Ammar 2024-08-21 09:03:40 +02:00
parent c9fa0470aa
commit 5a738f00b9

View File

@ -2,6 +2,11 @@ subroutine GTpp_excitation_density(ispin,nBas,nC,nO,nV,nR,nOO,nVV,ERI,X1,Y1,rho1
! Compute excitation densities for T-matrix self-energy ! Compute excitation densities for T-matrix self-energy
! TODO
! debug DGEMM for nC != 0
! and nR != 0
implicit none implicit none
@ -127,81 +132,132 @@ subroutine GTpp_excitation_density(ispin,nBas,nC,nO,nV,nR,nOO,nVV,ERI,X1,Y1,rho1
if(ispin == 2 .or. ispin == 4) then if(ispin == 2 .or. ispin == 4) then
!$OMP PARALLEL DEFAULT(NONE) & dim_1 = (nBas - nO) * (nBas - nO - 1) / 2
!$OMP PRIVATE(p, q, a, b, ab, c, d, cd, i, j, ij, k, l, kl) & dim_2 = nO * (nO - 1) / 2
!$OMP SHARED(nC, nBas, nR, nO, rho1, rho2, ERI, X1, Y1, X2, Y2)
allocate(ERI_1(nBas,nBas,dim_1), ERI_2(nBas,nBas,dim_2))
ERI_1 = 0.d0
ERI_2 = 0.d0
!$OMP PARALLEL DEFAULT(NONE) &
!$OMP PRIVATE(p, q, c, d, cd, k, l, kl) &
!$OMP SHARED(nC, nBas, nR, nO, ERI_1, ERI_2, ERI)
!$OMP DO COLLAPSE(2) !$OMP DO COLLAPSE(2)
do q = nC+1, nBas-nR do q = nC+1, nBas-nR
do p = nC+1, nBas-nR do p = nC+1, nBas-nR
cd = 0
do c = nO+1, nBas-nR
do d = c+1, nBas-nR
cd = cd + 1
ERI_1(p,q,cd) = ERI(p,q,c,d) - ERI(p,q,d,c)
enddo
enddo
kl = 0
do k = nC+1, nO
do l = k+1, nO
kl = kl + 1
ERI_2(p,q,kl) = ERI(p,q,k,l) - ERI(p,q,l,k)
end do
end do
enddo
enddo
!$OMP END DO
!$OMP END PARALLEL
ab = 0 call dgemm("N", "N", nBas*nBas, dim_1, dim_1, 1.d0, &
ERI_1(1,1,1), nBas*nBas, X1(1,1), dim_1, &
0.d0, rho1(1,1,1), nBas*nBas)
do a = nO+1, nBas-nR call dgemm("N", "N", nBas*nBas, dim_1, dim_2, 1.d0, &
do b = a+1, nBas-nR ERI_2(1,1,1), nBas*nBas, Y1(1,1), dim_2, &
1.d0, rho1(1,1,1), nBas*nBas)
ab = ab + 1 call dgemm("N", "N", nBas*nBas, dim_2, dim_1, 1.d0, &
ERI_1(1,1,1), nBas*nBas, X2(1,1), dim_1, &
0.d0, rho2(1,1,1), nBas*nBas)
cd = 0 call dgemm("N", "N", nBas*nBas, dim_2, dim_2, 1.d0, &
do c = nO+1, nBas-nR ERI_2(1,1,1), nBas*nBas, Y2(1,1), dim_2, &
do d = c+1, nBas-nR 1.d0, rho2(1,1,1), nBas*nBas)
cd = cd + 1 deallocate(ERI_1, ERI_2)
rho1(p,q,ab) = rho1(p,q,ab) &
+ (ERI(p,q,c,d) - ERI(p,q,d,c))*X1(cd,ab)
end do ! d
end do ! c
kl = 0 ! !$OMP PARALLEL DEFAULT(NONE) &
do k = nC+1, nO ! !$OMP PRIVATE(p, q, a, b, ab, c, d, cd, i, j, ij, k, l, kl) &
do l = k+1, nO ! !$OMP SHARED(nC, nBas, nR, nO, rho1, rho2, ERI, X1, Y1, X2, Y2)
! !$OMP DO COLLAPSE(2)
kl = kl + 1 ! do q = nC+1, nBas-nR
! do p = nC+1, nBas-nR
rho1(p,q,ab) = rho1(p,q,ab) & !
+ (ERI(p,q,k,l) - ERI(p,q,l,k))*Y1(kl,ab) ! ab = 0
end do ! l !
end do ! k ! do a = nO+1, nBas-nR
! do b = a+1, nBas-nR
end do ! b !
end do ! a ! ab = ab + 1
!
ij = 0 ! cd = 0
do i = nC+1, nO ! do c = nO+1, nBas-nR
do j = i+1, nO ! do d = c+1, nBas-nR
!
ij = ij + 1 ! cd = cd + 1
!
cd = 0 ! rho1(p,q,ab) = rho1(p,q,ab) &
! + (ERI(p,q,c,d) - ERI(p,q,d,c))*X1(cd,ab)
do c = nO+1, nBas-nR ! end do ! d
do d = c+1, nBas-nR ! end do ! c
!
cd = cd + 1 ! kl = 0
! do k = nC+1, nO
rho2(p,q,ij) = rho2(p,q,ij) & ! do l = k+1, nO
+ (ERI(p,q,c,d) - ERI(p,q,d,c))*X2(cd,ij) !
end do ! d ! kl = kl + 1
end do ! c !
! rho1(p,q,ab) = rho1(p,q,ab) &
kl = 0 ! + (ERI(p,q,k,l) - ERI(p,q,l,k))*Y1(kl,ab)
do k = nC+1, nO ! end do ! l
do l = k+1, nO ! end do ! k
!
kl = kl + 1 ! end do ! b
! end do ! a
rho2(p,q,ij) = rho2(p,q,ij) & !
+ (ERI(p,q,k,l) - ERI(p,q,l,k))*Y2(kl,ij) ! ij = 0
end do ! l ! do i = nC+1, nO
end do ! k ! do j = i+1, nO
!
end do ! j ! ij = ij + 1
end do ! i !
! cd = 0
end do ! p !
end do ! q ! do c = nO+1, nBas-nR
!$OMP END DO ! do d = c+1, nBas-nR
!$OMP END PARALLEL !
! cd = cd + 1
!
! rho2(p,q,ij) = rho2(p,q,ij) &
! + (ERI(p,q,c,d) - ERI(p,q,d,c))*X2(cd,ij)
! end do ! d
! end do ! c
!
! kl = 0
! do k = nC+1, nO
! do l = k+1, nO
!
! kl = kl + 1
!
! rho2(p,q,ij) = rho2(p,q,ij) &
! + (ERI(p,q,k,l) - ERI(p,q,l,k))*Y2(kl,ij)
! end do ! l
! end do ! k
!
! end do ! j
! end do ! i
!
! end do ! p
! end do ! q
! !$OMP END DO
! !$OMP END PARALLEL
end if end if
@ -209,9 +265,6 @@ subroutine GTpp_excitation_density(ispin,nBas,nC,nO,nV,nR,nOO,nVV,ERI,X1,Y1,rho1
! alpha-beta block ! alpha-beta block
!---------------------------------------------- !----------------------------------------------
! TODO
! debug for nC & nR
if(ispin == 3) then if(ispin == 3) then
dim_1 = (nBas - nO) * (nBas - nO) dim_1 = (nBas - nO) * (nBas - nO)