Skip to content

Commit

Permalink
OMP for phLR
Browse files Browse the repository at this point in the history
  • Loading branch information
AbdAmmar committed Dec 1, 2024
1 parent da4f8df commit 6cfe4a5
Show file tree
Hide file tree
Showing 2 changed files with 165 additions and 62 deletions.
120 changes: 88 additions & 32 deletions src/LR/phLR_A.f90
Original file line number Diff line number Diff line change
Expand Up @@ -25,6 +25,9 @@ subroutine phLR_A(ispin,dRPA,nBas,nC,nO,nV,nR,nS,lambda,e,ERI,Aph)
double precision,external :: Kronecker_delta

integer :: i,j,a,b,ia,jb
integer :: nn,jb0
logical :: i_eq_j
double precision :: ct1,ct2

! Output variables

Expand All @@ -39,45 +42,98 @@ subroutine phLR_A(ispin,dRPA,nBas,nC,nO,nV,nR,nS,lambda,e,ERI,Aph)

if(ispin == 1) then

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
nn = nBas - nR - nO
ct1 = 2d0 * lambda
ct2 = - (1d0 - delta_dRPA) * lambda
!$OMP PARALLEL DEFAULT(NONE) &
!$OMP PRIVATE (i, a, j, b, i_eq_j, ia, jb0, jb) &
!$OMP SHARED (nC, nO, nR, nBas, nn, ct1, ct2, e, ERI, Aph)
!$OMP DO COLLAPSE(2)
do i = nC+1, nO
do a = nO+1, nBas-nR
ia = a - nO + (i - nC - 1) * nn

do j = nC+1, nO
i_eq_j = i == j
jb0 = (j - nC - 1) * nn - nO

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

! Build A matrix for triplet manifold

if(ispin == 2) then

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
nn = nBas - nR - nO
ct2 = - (1d0 - delta_dRPA) * lambda
!$OMP PARALLEL DEFAULT(NONE) &
!$OMP PRIVATE (i, a, j, b, i_eq_j, ia, jb0, jb) &
!$OMP SHARED (nC, nO, nR, nBas, nn, ct2, e, ERI, Aph)
!$OMP DO COLLAPSE(2)
do i = nC+1, nO
do a = nO+1, nBas-nR
ia = a - nO + (i - nC - 1) * nn

do j = nC+1, nO
i_eq_j = i == j
jb0 = (j - nC - 1) * nn - nO

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

Expand Down
107 changes: 77 additions & 30 deletions src/LR/phLR_B.f90
Original file line number Diff line number Diff line change
Expand Up @@ -17,6 +17,8 @@ subroutine phLR_B(ispin,dRPA,nBas,nC,nO,nV,nR,nS,lambda,ERI,Bph)
double precision :: delta_dRPA

integer :: i,j,a,b,ia,jb
integer :: nn,jb0
double precision :: ct1,ct2

! Output variables

Expand All @@ -31,43 +33,88 @@ subroutine phLR_B(ispin,dRPA,nBas,nC,nO,nV,nR,nS,lambda,ERI,Bph)

if(ispin == 1) then

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
nn = nBas - nR - nO
ct1 = 2d0 * lambda
ct2 = - (1d0 - delta_dRPA) * lambda
!$OMP PARALLEL DEFAULT(NONE) &
!$OMP PRIVATE (i, a, j, b, ia, jb0, jb) &
!$OMP SHARED (nC, nO, nR, nBas, nn, ct1, ct2, ERI, Bph)
!$OMP DO COLLAPSE(2)
do i = nC+1, nO
do a = nO+1, nBas-nR
ia = a - nO + (i - nC - 1) * nn

do j = nC+1, nO
jb0 = (j - nC - 1) * nn - nO

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

! Build B matrix for triplet manifold

if(ispin == 2) then

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
nn = nBas - nR - nO
ct2 = - (1d0 - delta_dRPA) * lambda
!$OMP PARALLEL DEFAULT(NONE) &
!$OMP PRIVATE (i, a, j, b, ia, jb0, jb) &
!$OMP SHARED (nC, nO, nR, nBas, nn, ct2, ERI, Bph)
!$OMP DO COLLAPSE(2)
do i = nC+1, nO
do a = nO+1, nBas-nR
ia = a - nO + (i - nC - 1) * nn

do j = nC+1, nO
jb0 = (j - nC - 1) * nn - nO

do b = nO+1, nBas-nR
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

Expand Down

0 comments on commit 6cfe4a5

Please sign in to comment.