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

OMP for phLR

This commit is contained in:
AbdAmmar 2024-12-01 11:22:29 +01:00
parent da4f8df0f8
commit 6cfe4a5dec
2 changed files with 163 additions and 60 deletions

View File

@ -25,6 +25,9 @@ subroutine phLR_A(ispin,dRPA,nBas,nC,nO,nV,nR,nS,lambda,e,ERI,Aph)
double precision,external :: Kronecker_delta double precision,external :: Kronecker_delta
integer :: i,j,a,b,ia,jb integer :: i,j,a,b,ia,jb
integer :: nn,jb0
logical :: i_eq_j
double precision :: ct1,ct2
! Output variables ! Output variables
@ -39,22 +42,49 @@ subroutine phLR_A(ispin,dRPA,nBas,nC,nO,nV,nR,nS,lambda,e,ERI,Aph)
if(ispin == 1) then if(ispin == 1) then
ia = 0 nn = nBas - nR - nO
do i=nC+1,nO ct1 = 2d0 * lambda
do a=nO+1,nBas-nR ct2 = - (1d0 - delta_dRPA) * lambda
ia = ia + 1 !$OMP PARALLEL DEFAULT(NONE) &
jb = 0 !$OMP PRIVATE (i, a, j, b, i_eq_j, ia, jb0, jb) &
do j=nC+1,nO !$OMP SHARED (nC, nO, nR, nBas, nn, ct1, ct2, e, ERI, Aph)
do b=nO+1,nBas-nR !$OMP DO COLLAPSE(2)
jb = jb + 1 do i = nC+1, nO
do a = nO+1, nBas-nR
Aph(ia,jb) = (e(a) - e(i))*Kronecker_delta(i,j)*Kronecker_delta(a,b) & ia = a - nO + (i - nC - 1) * nn
+ 2d0*lambda*ERI(i,b,a,j) - (1d0 - delta_dRPA)*lambda*ERI(i,b,j,a)
end do do j = nC+1, nO
end do i_eq_j = i == j
end do jb0 = (j - nC - 1) * nn - nO
end do
do b = nO+1, nBas-nR
jb = b + jb0
Aph(ia,jb) = ct1 * ERI(i,b,a,j) + ct2 * ERI(i,b,j,a)
if(i_eq_j) then
if(a == b) Aph(ia,jb) = Aph(ia,jb) + e(a) - e(i)
endif
enddo
enddo
enddo
enddo
!$OMP END DO
!$OMP END PARALLEL
!ia = 0
!do i=nC+1,nO
! do a=nO+1,nBas-nR
! ia = ia + 1
! jb = 0
! do j=nC+1,nO
! do b=nO+1,nBas-nR
! jb = jb + 1
! Aph(ia,jb) = (e(a) - e(i))*Kronecker_delta(i,j)*Kronecker_delta(a,b) &
! + 2d0*lambda*ERI(i,b,a,j) - (1d0 - delta_dRPA)*lambda*ERI(i,b,j,a)
! end do
! end do
! end do
!end do
end if end if
@ -62,22 +92,48 @@ subroutine phLR_A(ispin,dRPA,nBas,nC,nO,nV,nR,nS,lambda,e,ERI,Aph)
if(ispin == 2) then if(ispin == 2) then
ia = 0 nn = nBas - nR - nO
do i=nC+1,nO ct2 = - (1d0 - delta_dRPA) * lambda
do a=nO+1,nBas-nR !$OMP PARALLEL DEFAULT(NONE) &
ia = ia + 1 !$OMP PRIVATE (i, a, j, b, i_eq_j, ia, jb0, jb) &
jb = 0 !$OMP SHARED (nC, nO, nR, nBas, nn, ct2, e, ERI, Aph)
do j=nC+1,nO !$OMP DO COLLAPSE(2)
do b=nO+1,nBas-nR do i = nC+1, nO
jb = jb + 1 do a = nO+1, nBas-nR
ia = a - nO + (i - nC - 1) * nn
Aph(ia,jb) = (e(a) - e(i))*Kronecker_delta(i,j)*Kronecker_delta(a,b) &
- (1d0 - delta_dRPA)*lambda*ERI(i,b,j,a)
end do do j = nC+1, nO
end do i_eq_j = i == j
end do jb0 = (j - nC - 1) * nn - nO
end do
do b = nO+1, nBas-nR
jb = b + jb0
Aph(ia,jb) = ct2 * ERI(i,b,j,a)
if(i_eq_j) then
if(a == b) Aph(ia,jb) = Aph(ia,jb) + e(a) - e(i)
endif
enddo
enddo
enddo
enddo
!$OMP END DO
!$OMP END PARALLEL
! ia = 0
! do i=nC+1,nO
! do a=nO+1,nBas-nR
! ia = ia + 1
! jb = 0
! do j=nC+1,nO
! do b=nO+1,nBas-nR
! jb = jb + 1
! Aph(ia,jb) = (e(a) - e(i))*Kronecker_delta(i,j)*Kronecker_delta(a,b) &
! - (1d0 - delta_dRPA)*lambda*ERI(i,b,j,a)
! end do
! end do
! end do
! end do
end if end if

View File

@ -17,6 +17,8 @@ subroutine phLR_B(ispin,dRPA,nBas,nC,nO,nV,nR,nS,lambda,ERI,Bph)
double precision :: delta_dRPA double precision :: delta_dRPA
integer :: i,j,a,b,ia,jb integer :: i,j,a,b,ia,jb
integer :: nn,jb0
double precision :: ct1,ct2
! Output variables ! Output variables
@ -31,21 +33,44 @@ subroutine phLR_B(ispin,dRPA,nBas,nC,nO,nV,nR,nS,lambda,ERI,Bph)
if(ispin == 1) then if(ispin == 1) then
ia = 0 nn = nBas - nR - nO
do i=nC+1,nO ct1 = 2d0 * lambda
do a=nO+1,nBas-nR ct2 = - (1d0 - delta_dRPA) * lambda
ia = ia + 1 !$OMP PARALLEL DEFAULT(NONE) &
jb = 0 !$OMP PRIVATE (i, a, j, b, ia, jb0, jb) &
do j=nC+1,nO !$OMP SHARED (nC, nO, nR, nBas, nn, ct1, ct2, ERI, Bph)
do b=nO+1,nBas-nR !$OMP DO COLLAPSE(2)
jb = jb + 1 do i = nC+1, nO
do a = nO+1, nBas-nR
Bph(ia,jb) = 2d0*lambda*ERI(i,j,a,b) - (1d0 - delta_dRPA)*lambda*ERI(i,j,b,a) ia = a - nO + (i - nC - 1) * nn
end do do j = nC+1, nO
end do jb0 = (j - nC - 1) * nn - nO
end do
end do do b = nO+1, nBas-nR
jb = b + jb0
Bph(ia,jb) = ct1 * ERI(i,j,a,b) + ct2 * ERI(i,j,b,a)
enddo
enddo
enddo
enddo
!$OMP END DO
!$OMP END PARALLEL
!ia = 0
!do i=nC+1,nO
! do a=nO+1,nBas-nR
! ia = ia + 1
! jb = 0
! do j=nC+1,nO
! do b=nO+1,nBas-nR
! jb = jb + 1
! Bph(ia,jb) = 2d0*lambda*ERI(i,j,a,b) - (1d0 - delta_dRPA)*lambda*ERI(i,j,b,a)
! end do
! end do
! end do
!end do
end if end if
@ -53,21 +78,43 @@ subroutine phLR_B(ispin,dRPA,nBas,nC,nO,nV,nR,nS,lambda,ERI,Bph)
if(ispin == 2) then if(ispin == 2) then
ia = 0 nn = nBas - nR - nO
do i=nC+1,nO ct2 = - (1d0 - delta_dRPA) * lambda
do a=nO+1,nBas-nR !$OMP PARALLEL DEFAULT(NONE) &
ia = ia + 1 !$OMP PRIVATE (i, a, j, b, ia, jb0, jb) &
jb = 0 !$OMP SHARED (nC, nO, nR, nBas, nn, ct2, ERI, Bph)
do j=nC+1,nO !$OMP DO COLLAPSE(2)
do b=nO+1,nBas-nR do i = nC+1, nO
jb = jb + 1 do a = nO+1, nBas-nR
ia = a - nO + (i - nC - 1) * nn
Bph(ia,jb) = - (1d0 - delta_dRPA)*lambda*ERI(i,j,b,a)
do j = nC+1, nO
end do jb0 = (j - nC - 1) * nn - nO
end do
end do do b = nO+1, nBas-nR
end do jb = b + jb0
Bph(ia,jb) = ct2 * ERI(i,j,b,a)
enddo
enddo
enddo
enddo
!$OMP END DO
!$OMP END PARALLEL
! ia = 0
! do i=nC+1,nO
! do a=nO+1,nBas-nR
! ia = ia + 1
! jb = 0
! do j=nC+1,nO
! do b=nO+1,nBas-nR
! jb = jb + 1
! Bph(ia,jb) = - (1d0 - delta_dRPA)*lambda*ERI(i,j,b,a)
! end do
! end do
! end do
! end do
end if end if