get_threshold_based_SIs Subroutine

private subroutine get_threshold_based_SIs(ref_buf, refs_found)

Arguments

Type IntentOptional Attributes Name
integer(kind=n_int), intent(out) :: ref_buf(0:NIfTot,maxNRefs)
integer, intent(out) :: refs_found

Contents


Source Code

    subroutine get_threshold_based_SIs(ref_buf, refs_found)
        implicit none
        integer(n_int), intent(out) :: ref_buf(0:NIfTot, maxNRefs)
        integer, intent(out) :: refs_found
        integer(n_int), allocatable :: tmp(:, :)

!       integer :: i, nBlocks
        integer(int64) :: i
        integer :: nBlocks

        integer, parameter :: blockSize = 5000
        real(dp) :: sgn(lenof_sign)
        real(dp) :: repAvSgn

        ref_buf = 0
        refs_found = 0
        nBlocks = 1
        repAvSgn = 0.0_dp

        allocate(tmp(0:NIfTot, blockSize))
        tmp = 0

        do i = 1, TotWalkers
            call extract_sign(CurrentDets(:, i), sgn)
            ! either compare the sum of the signed or unsigned walker numbers to the
            ! initiator threshold
            if (tSignedRepAv) then
                repAvSgn = abs(sum(sgn) / inum_runs)
            else
                repAvSgn = av_pop(sgn)
            end if
            if ((repAvSgn >= NoTypeN)) then
                refs_found = refs_found + 1

                ! If the temporary is full, resize it
                if (refs_found > nBlocks * blockSize) then
                    nBlocks = nBlocks + 1
                    call resize_ilut_list(tmp, (nBlocks - 1) * blockSize, nBlocks * blockSize)
                end if

                ! add the possible SI to the temporary
                tmp(:, refs_found) = CurrentDets(:, i)
            end if
        end do

        ! we only keep at most maxNRefs determinants
        if (refs_found > maxNRefs .and. NoTypeN > 1) then
            write(stdout, '(A,I5,A,I5,A,I5,A)') "On proc ", iProcIndex, " found ", refs_found, &
                " SIs, which is more than the maximum of ", maxNRefs, " - truncating"
            ! in case we found more, take the maxNRefs with the highest population
            call sort(tmp(0:NIfTot, 1:refs_found), sign_gt, sign_lt)
            ref_buf(:, 1:maxNRefs) = tmp(:, 1:maxNRefs)
            refs_found = maxNRefs
        else
            ! else, take all of them
            ref_buf(:, 1:refs_found) = tmp(:, 1:refs_found)
        end if

        deallocate(tmp)
    end subroutine get_threshold_based_SIs