10
0
mirror of https://github.com/QuantumPackage/qp2.git synced 2024-07-22 18:57:36 +02:00

binary_search_cfg accelerated

This commit is contained in:
Anthony Scemama 2021-03-06 19:23:21 +01:00
parent 3750f5446d
commit 52c4fc1ce4
2 changed files with 65 additions and 33 deletions

View File

@ -748,34 +748,66 @@ subroutine binary_search_cfg(cfgInp,addcfg)
implicit none
BEGIN_DOC
! Documentation for binary_search
!
! Does a binary search to find
!
! Does a binary search to find
! the address of a configuration in a list of
! configurations.
END_DOC
integer(bit_kind), intent(in) :: cfgInp(N_int,2)
integer , intent(out) :: addcfg
integer :: i,j,k,r,l
integer*8 :: key, key2
logical :: found
!integer*8, allocatable :: bit_tmp(:)
!integer*8, external :: configuration_search_key
!allocate(bit_tmp(0:N_configuration))
!bit_tmp(0) = 0
do i=1,N_configuration
!bit_tmp(i) = configuration_search_key(psi_configuration(1,1,i),N_int)
found = .True.
do k=1,N_int
found = found .and. (psi_configuration(k,1,i) == cfgInp(k,1)) &
.and. (psi_configuration(k,2,i) == cfgInp(k,2))
enddo
if (found) then
addcfg = i
exit
logical :: found
integer :: l, r, j, k
integer*8 :: bit_tmp, key
integer*8, external :: configuration_search_key
key = configuration_search_key(cfgInp,N_int)
! Binary search
l = 0
r = N_configuration+1
j = shiftr(r-l,1)
do while (j>=1)
j = j+l
bit_tmp = configuration_search_key(psi_configuration(1,1,j),N_int)
if (bit_tmp == key) then
! Find 1st element which matches the key
if (j > 1) then
bit_tmp = configuration_search_key(psi_configuration(1,1,j-1),N_int)
do while (j>1 .and. bit_tmp == key)
j = j-1
bit_tmp = configuration_search_key(psi_configuration(1,1,j-1),N_int)
enddo
bit_tmp = key
endif
! Find correct element matching the key
do while (bit_tmp == key)
found = .True.
do k=1,N_int
found = found .and. (psi_configuration(k,1,j) == cfgInp(k,1))&
.and. (psi_configuration(k,2,j) == cfgInp(k,2))
enddo
if (found) then
addcfg = j
return
endif
j = j+1
bit_tmp = configuration_search_key(psi_configuration(1,1,j),N_int)
enddo
addcfg = -1
return
else if (bit_tmp > key) then
r = j
else
l = j
endif
j = shiftr(r-l,1)
enddo
addcfg = -1
return
end subroutine
BEGIN_PROVIDER [ integer, psi_configuration_to_psi_det, (2,N_configuration) ]
@ -812,7 +844,7 @@ end subroutine
psi_configuration_to_psi_det(2,k) = N_det
! Reorder determinants according to generation
! Reorder determinants according to generation
! --------------------------------------------
integer(bit_kind), allocatable :: dets(:,:,:)

View File

@ -306,23 +306,23 @@ subroutine davidson_diag_csf_hjj(dets_in,u_in,H_jj,energies,dim_in,sze,sze_csf,N
PROVIDE mo_two_e_integrals_in_map mo_integrals_map big_array_exchange_integrals
if ((sze > 100000).and.distributed_davidson) then
!call convertWFfromCSFtoDET(N_st_diag,U_csf(1,shift+1),U)
!call convertWFfromCSFtoDET(N_st_diag,W_csf(1,shift+1),W)
!call H_u_0_nstates_zmq (W,U,N_st_diag,sze)
!call convertWFfromDETtoCSF(N_st_diag,U,U_csf(1,shift+1))
!call convertWFfromDETtoCSF(N_st_diag,W,W_csf(1,shift+1))
!call calculate_sigma_vector_cfg_nst(W_csf(1,shift+1),U_csf(1,shift+1),N_st_diag,sze_csf,1,sze_csf,0,1)
! TODO : psi_det_size ? for psi_det
! call convertWFfromCSFtoDET(N_st_diag,U_csf(1,shift+1),U)
! call convertWFfromCSFtoDET(N_st_diag,W_csf(1,shift+1),W)
! call H_u_0_nstates_zmq (W,U,N_st_diag,sze)
! call convertWFfromDETtoCSF(N_st_diag,U,U_csf(1,shift+1))
! call convertWFfromDETtoCSF(N_st_diag,W,W_csf(1,shift+1))
! call calculate_sigma_vector_cfg_nst(W_csf(1,shift+1),U_csf(1,shift+1),N_st_diag,sze_csf,1,sze_csf,0,1)
! ! TODO : psi_det_size ? for psi_det
do kk=1,N_st_diag
call calculate_sigma_vector_cfg_nst(W_csf(1,shift+kk),U_csf(1,shift+kk),1,sze_csf,1,sze_csf,0,1)
enddo
else
!call convertWFfromCSFtoDET(N_st_diag,U_csf(1,shift+1),U)
!call convertWFfromCSFtoDET(N_st_diag,W_csf(1,shift+1),W)
!call H_u_0_nstates_openmp(W,U,N_st_diag,sze)
!call convertWFfromDETtoCSF(N_st_diag,U,U_csf(1,shift+1))
!call convertWFfromDETtoCSF(N_st_diag,W,W_csf(1,shift+1))
!call calculate_sigma_vector_cfg_nst(W_csf(1,shift+1),U_csf(1,shift+1),N_st_diag,sze_csf,1,sze_csf,0,1)
! call convertWFfromCSFtoDET(N_st_diag,U_csf(1,shift+1),U)
! call convertWFfromCSFtoDET(N_st_diag,W_csf(1,shift+1),W)
! call H_u_0_nstates_openmp(W,U,N_st_diag,sze)
! call convertWFfromDETtoCSF(N_st_diag,U,U_csf(1,shift+1))
! call convertWFfromDETtoCSF(N_st_diag,W,W_csf(1,shift+1))
! call calculate_sigma_vector_cfg_nst(W_csf(1,shift+1),U_csf(1,shift+1),N_st_diag,sze_csf,1,sze_csf,0,1)
do kk=1,N_st_diag
call calculate_sigma_vector_cfg_nst(W_csf(1,shift+kk),U_csf(1,shift+kk),1,sze_csf,1,sze_csf,0,1)
enddo