proc_most_populated_states Subroutine

public subroutine proc_most_populated_states(n_keep, run, largest_walkers, opt_source, opt_source_size, norm)

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: n_keep
integer, intent(in) :: run
integer(kind=n_int), intent(out) :: largest_walkers(0:NIfTot,n_keep)
integer(kind=n_int), intent(in), optional :: opt_source(0:,1:)
integer(kind=int64), intent(in), optional :: opt_source_size
real(kind=dp), intent(out), optional :: norm

Contents


Source Code

    subroutine proc_most_populated_states(n_keep, run, &
                                          largest_walkers, opt_source, opt_source_size, norm)
        ! Return the most populated states in CurrentDets on *this* processor only.
        ! Also return the norm of these states, if requested.

        use DetBitOps, only: sign_lt, sign_gt
        use sort_mod, only: sort

        integer(int64), intent(in), optional :: opt_source_size
        integer(n_int), intent(in), optional :: opt_source(0:, 1:)
        integer, intent(in) :: n_keep, run
        integer(n_int), intent(out) :: largest_walkers(0:NIfTot, n_keep)
        real(dp), intent(out), optional :: norm
        integer :: i, j, smallest_pos, part_type
        real(dp) :: smallest_sign, sign_curr_real
        real(dp), dimension(lenof_sign) :: sign_curr, low_sign
        character(*), parameter :: this_routine = "return_most_populated_states"

        integer(n_int), allocatable :: loc_source(:, :)
        integer(int64) :: source_size

        if (present(opt_source)) then
            ASSERT(present(opt_source_size))

            source_size = opt_source_size
            ! ask Kai if I have to allocate
            allocate(loc_source(0:niftot, 1:source_size), &
                      source=opt_source(0:NIfTot, 1:source_size))
!             loc_source => opt_source

        else
            source_size = TotWalkers
            allocate(loc_source(0:niftot, 1:source_size), &
                      source=CurrentDets(0:NIfTot, 1:source_size))
!             loc_source => CurrentDets
        end if

        largest_walkers = 0_n_int
        smallest_sign = 0.0_dp
        smallest_pos = 1
        if (present(norm)) norm = 0.0_dp

        ! Run through all walkers on process.
        do i = 1, int(source_size)
            call extract_sign(loc_source(:, i), sign_curr)

            sign_curr_real = core_space_weight(sign_curr, run)
            if (present(norm)) norm = norm + (sign_curr_real**2.0)

            ! Is this determinant more populated than the smallest? First in
            ! the list is always the smallest.
            if (sign_curr_real > smallest_sign) then
                largest_walkers(:, smallest_pos) = loc_source(:, i)

                ! Instead of resorting, just find new smallest sign and position.
                call extract_sign(largest_walkers(:, 1), low_sign)

                smallest_sign = core_space_weight(low_sign, run)

                smallest_pos = 1
                do j = 2, n_keep
                    call extract_sign(largest_walkers(:, j), low_sign)
                    sign_curr_real = core_space_weight(low_sign, run)
                    if (sign_curr_real < smallest_sign .or. all(largest_walkers(:, j) == 0_n_int)) then
                        smallest_pos = j
                        smallest_sign = sign_curr_real
                    end if
                end do

            end if

        end do

        call sort(largest_walkers(:, 1:n_keep), sign_lt_run, sign_gt_run)

    contains

        pure function sign_lt_run(ilutI, ilutJ) result(bLt)

            ! This is a comparison function between two bit strings of length
            ! 0:NIfTot, and will return true if absolute value of the sign of
            ! ilutI is less than ilutJ

            integer(n_int), intent(in) :: iLutI(0:), iLutJ(0:)
            logical :: bLt
            real(dp) :: SignI(lenof_sign), SignJ(lenof_sign)

            call extract_sign(ilutI, SignI)
            call extract_sign(ilutJ, SignJ)

            bLt = core_space_weight(signI, run) < core_space_weight(signJ, run)

        end function sign_lt_run

        pure function sign_gt_run(ilutI, ilutJ) result(bGt)

            ! This is a comparison function between two bit strings of length
            ! 0:NIfTot, and will return true if the abs sign of ilutI is greater
            ! than ilutJ

            integer(n_int), intent(in) :: iLutI(0:), iLutJ(0:)
            logical :: bGt
            real(dp) :: SignI(lenof_sign), SignJ(lenof_sign)

            call extract_sign(ilutI, SignI)
            call extract_sign(ilutJ, SignJ)

            bGt = core_space_weight(signI, run) > core_space_weight(signJ, run)

        end function sign_gt_run

    end subroutine proc_most_populated_states