10
1
mirror of https://github.com/pfloos/quack synced 2025-01-03 01:55:57 +01:00

working on X-Fock with s8-ERI

This commit is contained in:
Abdallah Ammar 2024-12-10 15:35:13 +01:00
parent 995bcfd236
commit 995378de74
2 changed files with 116 additions and 80 deletions

View File

@ -45,126 +45,164 @@ subroutine exchange_matrix_AO_basis_hpc(nBas, ERI_size, P, ERI_chem, K)
double precision, intent(out) :: K(nBas,nBas) double precision, intent(out) :: K(nBas,nBas)
integer :: mu, nu, la, si integer :: mu, nu, la, si
integer :: nunu, lala, nula, lasi, numu integer :: nunu, nula, lanu, lasi, nusi, sinu
integer*8 :: nunununu, nunulala, nununula, nunulasi integer :: numu, mumu, mula, lamu, musi, simu
integer*8 :: lalanunu, lasinunu, numulala, lalanumu integer*8 :: nunununu, nulanula, lanulanu, nulanusi
integer*8 :: numunula, numulasi, lasinumu, nununumu integer*8 :: munulasi, lanunusi, lanusinu, numumumu
integer*8 :: munu0, munu integer*8 :: nulamula, nulalamu, lanulamu, nulamusi
integer*8 :: sila0, sila integer*8 :: nulasimu, lanumusi, lanusimu
integer*8 :: munulasi0, munulasi
integer*8, external :: Yoshimine_4ind integer*8, external :: Yoshimine_4ind
do nu = 1, nBas do nu = 1, nBas
munulasi = Yoshimine_4ind(nu, nu, nu, nu) nunu = (nu * (nu - 1)) / 2 + nu
K(nu,nu) = -P(nu,nu) * ERI_chem(munulasi) nunununu = (nunu * (nunu - 1)) / 2 + nunu
K(nu,nu) = -P(nu,nu) * ERI_chem(nunununu)
do la = 1, nu - 1 do la = 1, nu - 1
munulasi = Yoshimine_4ind(nu, la, nu, la) nula = (nu * (nu - 1)) / 2 + la
K(nu,nu) = K(nu,nu) - P(la,la) * ERI_chem(munulasi) nulanula = (nula * (nula - 1)) / 2 + nula
K(nu,nu) = K(nu,nu) - P(la,la) * ERI_chem(nulanula)
enddo enddo
do la = nu + 1, nBas do la = nu + 1, nBas
munulasi = Yoshimine_4ind(nu, la, nu, la) lanu = (la * (la - 1)) / 2 + nu
K(nu,nu) = K(nu,nu) - P(la,la) * ERI_chem(munulasi) lanulanu = (lanu * (lanu - 1)) / 2 + lanu
K(nu,nu) = K(nu,nu) - P(la,la) * ERI_chem(lanulanu)
enddo enddo
do la = 1, nu do la = 1, nu
nula = (nu * (nu - 1)) / 2 + la
do si = 1, la - 1 do si = 1, la - 1
munulasi = Yoshimine_4ind(nu, la, nu, si) nusi = (nu * (nu - 1)) / 2 + si
K(nu,nu) = K(nu,nu) - 2.d0 * P(si,la) * ERI_chem(munulasi) nulanusi = (nula * (nula - 1)) / 2 + nusi
K(nu,nu) = K(nu,nu) - 2.d0 * P(si,la) * ERI_chem(nulanusi)
enddo enddo
enddo enddo
do la = nu + 1, nBas do la = nu + 1, nBas
lanu = (la * (la - 1)) / 2 + nu
do si = 1, nu do si = 1, nu
munulasi = Yoshimine_4ind(nu, la, nu, si) nusi = (nu * (nu - 1)) / 2 + si
K(nu,nu) = K(nu,nu) - 2.d0 * P(si,la) * ERI_chem(munulasi) lanunusi = (lanu * (lanu - 1)) / 2 + nusi
K(nu,nu) = K(nu,nu) - 2.d0 * P(si,la) * ERI_chem(lanunusi)
enddo enddo
enddo
do la = nu + 1, nBas
do si = nu + 1, la - 1 do si = nu + 1, la - 1
munulasi = Yoshimine_4ind(nu, la, nu, si) sinu = (si * (si - 1)) / 2 + nu
K(nu,nu) = K(nu,nu) - 2.d0 * P(si,la) * ERI_chem(munulasi) lanusinu = (lanu * (lanu - 1)) / 2 + sinu
K(nu,nu) = K(nu,nu) - 2.d0 * P(si,la) * ERI_chem(lanusinu)
enddo enddo
enddo enddo
do mu = 1, nu - 1 do mu = 1, nu - 1
K(mu,nu) = 0.d0 numu = (nu * (nu - 1)) / 2 + mu
do la = 1, mu mumu = (mu * (mu - 1)) / 2 + mu
munulasi = Yoshimine_4ind(nu, la, mu, la) numumumu = (numu * (numu - 1)) / 2 + mumu
K(mu,nu) = K(mu,nu) - P(la,la) * ERI_chem(munulasi) K(mu,nu) = - P(mu,mu) * ERI_chem(numumumu)
do la = 1, mu - 1
mula = (mu * (mu - 1)) / 2 + la
nula = (nu * (nu - 1)) / 2 + la
nulamula = (nula * (nula - 1)) / 2 + mula
K(mu,nu) = K(mu,nu) - P(la,la) * ERI_chem(nulamula)
enddo enddo
do la = mu + 1, nu do la = mu + 1, nu
munulasi = Yoshimine_4ind(nu, la, mu, la) lamu = (la * (la - 1)) / 2 + mu
K(mu,nu) = K(mu,nu) - P(la,la) * ERI_chem(munulasi) nula = (nu * (nu - 1)) / 2 + la
nulalamu = (nula * (nula - 1)) / 2 + lamu
K(mu,nu) = K(mu,nu) - P(la,la) * ERI_chem(nulalamu)
enddo enddo
do la = nu + 1, nBas do la = nu + 1, nBas
munulasi = Yoshimine_4ind(nu, la, mu, la) lamu = (la * (la - 1)) / 2 + mu
K(mu,nu) = K(mu,nu) - P(la,la) * ERI_chem(munulasi) lanu = (la * (la - 1)) / 2 + nu
lanulamu = (lanu * (lanu - 1)) / 2 + lamu
K(mu,nu) = K(mu,nu) - P(la,la) * ERI_chem(lanulamu)
enddo enddo
do la = 1, mu do la = 1, mu
nula = (nu * (nu - 1)) / 2 + la
do si = 1, la - 1 do si = 1, la - 1
munulasi = Yoshimine_4ind(nu, la, si, mu) musi = (mu * (mu - 1)) / 2 + si
K(mu,nu) = K(mu,nu) - P(si,la) * ERI_chem(munulasi) nulamusi = (nula * (nula - 1)) / 2 + musi
K(mu,nu) = K(mu,nu) - P(si,la) * ERI_chem(nulamusi)
enddo enddo
enddo enddo
do la = mu + 1, nu do la = mu + 1, nu
nula = (nu * (nu - 1)) / 2 + la
do si = 1, mu do si = 1, mu
munulasi = Yoshimine_4ind(nu, la, si, mu) musi = (mu * (mu - 1)) / 2 + si
K(mu,nu) = K(mu,nu) - P(si,la) * ERI_chem(munulasi) nulamusi = (nula * (nula - 1)) / 2 + musi
K(mu,nu) = K(mu,nu) - P(si,la) * ERI_chem(nulamusi)
enddo enddo
do si = mu + 1, la - 1 do si = mu + 1, la - 1
munulasi = Yoshimine_4ind(nu, la, si, mu) simu = (si * (si - 1)) / 2 + mu
K(mu,nu) = K(mu,nu) - P(si,la) * ERI_chem(munulasi) nulasimu = (nula * (nula - 1)) / 2 + simu
K(mu,nu) = K(mu,nu) - P(si,la) * ERI_chem(nulasimu)
enddo enddo
enddo enddo
do la = nu + 1, nBas do la = nu + 1, nBas
lanu = (la * (la - 1)) / 2 + nu
do si = 1, mu do si = 1, mu
munulasi = Yoshimine_4ind(nu, la, si, mu) musi = (mu * (mu - 1)) / 2 + si
K(mu,nu) = K(mu,nu) - P(si,la) * ERI_chem(munulasi) lanumusi = (lanu * (lanu - 1)) / 2 + musi
K(mu,nu) = K(mu,nu) - P(si,la) * ERI_chem(lanumusi)
enddo enddo
do si = mu + 1, la - 1 do si = mu + 1, la - 1
munulasi = Yoshimine_4ind(nu, la, si, mu) simu = (si * (si - 1)) / 2 + mu
K(mu,nu) = K(mu,nu) - P(si,la) * ERI_chem(munulasi) lanusimu = (lanu * (lanu - 1)) / 2 + simu
K(mu,nu) = K(mu,nu) - P(si,la) * ERI_chem(lanusimu)
enddo enddo
enddo enddo
do la = 1, mu !TODO
do si = la + 1, mu ! do la = 1, mu
munulasi = Yoshimine_4ind(nu, la, si, mu) ! nula = (nu * (nu - 1)) / 2 + la
K(mu,nu) = K(mu,nu) - P(si,la) * ERI_chem(munulasi) ! do si = la + 1, mu
enddo ! musi = (mu * (mu - 1)) / 2 + si
do si = mu + 1, nBas ! nulamusi = (nula * (nula - 1)) / 2 + musi
munulasi = Yoshimine_4ind(nu, la, si, mu) ! !nulamusi = Yoshimine_4ind(nu, la, si, mu)
K(mu,nu) = K(mu,nu) - P(si,la) * ERI_chem(munulasi) ! K(mu,nu) = K(mu,nu) - P(si,la) * ERI_chem(nulamusi)
enddo ! enddo
enddo ! do si = mu + 1, nBas
do la = mu + 1, nu ! simu = (si * (si - 1)) / 2 + mu
do si = la + 1, nBas ! nulasimu = (nula * (nula - 1)) / 2 + simu
munulasi = Yoshimine_4ind(nu, la, si, mu) ! !nulasimu = Yoshimine_4ind(nu, la, si, mu)
K(mu,nu) = K(mu,nu) - P(si,la) * ERI_chem(munulasi) ! K(mu,nu) = K(mu,nu) - P(si,la) * ERI_chem(nulasimu)
enddo ! enddo
enddo ! enddo
do la = nu + 1, nBas ! do la = mu + 1, nu
do si = la + 1, nBas ! nula = (nu * (nu - 1)) / 2 + la
munulasi = Yoshimine_4ind(nu, la, si, mu) ! do si = la + 1, nu
K(mu,nu) = K(mu,nu) - P(si,la) * ERI_chem(munulasi) ! simu = (si * (si - 1)) / 2 + mu
enddo ! nulasimu = (nula * (nula - 1)) / 2 + simu
enddo ! !nulasimu = Yoshimine_4ind(nu, la, si, mu)
! K(mu,nu) = K(mu,nu) - P(si,la) * ERI_chem(nulasimu)
!do la = 1, nBas ! enddo
! do si = la + 1, nBas ! do si = nu + 1, nBas
! munulasi = Yoshimine_4ind(nu, la, mu, si) ! simu = (si * (si - 1)) / 2 + mu
! munulasi = Yoshimine_4ind(nu, la, si, mu)
! K(mu,nu) = K(mu,nu) - P(si,la) * ERI_chem(munulasi) ! K(mu,nu) = K(mu,nu) - P(si,la) * ERI_chem(munulasi)
! enddo ! enddo
! enddo ! enddo
! do la = nu + 1, nBas
! lanu = (la * (la - 1)) / 2 + nu
! do si = la + 1, mu
! simu = (si * (si - 1)) / 2 + mu
! munulasi = Yoshimine_4ind(nu, la, si, mu)
! K(mu,nu) = K(mu,nu) - P(si,la) * ERI_chem(munulasi)
! enddo
! do si = mu + 1, nBas
! musi = (mu * (mu - 1)) / 2 + si
! munulasi = Yoshimine_4ind(nu, la, si, mu)
! K(mu,nu) = K(mu,nu) - P(si,la) * ERI_chem(munulasi)
! enddo
! enddo
enddo ! mu enddo ! mu
enddo ! nu enddo ! nu

View File

@ -914,11 +914,9 @@ integer*8 function Yoshimine_4ind(a, b, c, d)
implicit none implicit none
integer, intent(in) :: a, b, c, d integer, intent(in) :: a, b, c, d
integer*8, external :: Yoshimine_2ind integer*8, external :: Yoshimine_2ind
integer*8 :: ab, cd
ab = Yoshimine_2ind(a, b) Yoshimine_4ind = Yoshimine_2ind(Yoshimine_2ind(a, b), &
cd = Yoshimine_2ind(c, d) Yoshimine_2ind(c, d))
Yoshimine_4ind = Yoshimine_2ind(ab, cd)
return return
end end
@ -928,7 +926,7 @@ end
integer*8 function Yoshimine_2ind(a, b) integer*8 function Yoshimine_2ind(a, b)
implicit none implicit none
integer, intent(in) :: a, b integer*8, intent(in) :: a, b
if(a > b) then if(a > b) then
Yoshimine_2ind = (a * (a - 1)) / 2 + b Yoshimine_2ind = (a * (a - 1)) / 2 + b