create_cum_list_rs_hubbard Subroutine

public subroutine create_cum_list_rs_hubbard(ilutI, src, neighbors, cum_arr, cum_sum, tgt, cpt)

Arguments

Type IntentOptional Attributes Name
integer(kind=n_int), intent(in) :: ilutI(0:NIfTot)
integer, intent(in) :: src
integer, intent(in) :: neighbors(:)
real(kind=dp), intent(out), allocatable :: cum_arr(:)
real(kind=dp), intent(out) :: cum_sum
integer, intent(in), optional :: tgt
real(kind=dp), intent(out), optional :: cpt

Contents


Source Code

    subroutine create_cum_list_rs_hubbard(ilutI, src, neighbors, cum_arr, cum_sum, &
                                          tgt, cpt)
        integer(n_int), intent(in) :: ilutI(0:NIfTot)
        integer, intent(in) :: src, neighbors(:)
        real(dp), intent(out) :: cum_sum
        real(dp), intent(out), allocatable :: cum_arr(:)
        integer, intent(in), optional :: tgt
        real(dp), intent(out), optional :: cpt
#ifdef DEBUG_
        character(*), parameter :: this_routine = "create_cum_list_rs_hubbard"
#endif

        real(dp) :: elem
        integer :: i, nI(nel), ex(2), ex2(2), nJ(nel)

        ASSERT(IsOcc(ilutI, src))

        call decode_bit_det(nI, ilutI)

        ex(1) = src

        allocate(cum_arr(size(neighbors)))
        cum_arr = 0.0_dp
        cum_sum = 0.0_dp
        if (present(tgt)) then
            ASSERT(present(cpt))
            do i = 1, ubound(neighbors, 1)
                elem = 0.0_dp
                ASSERT(is_beta(src) .eqv. is_beta(neighbors(i)))
                if (IsNotOcc(ilutI, neighbors(i))) then
                    ! change the order of determinants to reflect
                    ! non-hermiticity correctly
                    ! old implo:
                    ex(2) = neighbors(i)
                    call swap_excitations(nI, ex, nJ, ex2)
                    elem = abs(get_offdiag_helement_rs_hub(nJ, ex2, .false.))

                end if
                cum_sum = cum_sum + elem
                if (neighbors(i) == tgt) then
                    cpt = elem
                end if
            end do
            if (cum_sum < EPS) then
                cpt = 0.0
            else
                cpt = cpt / cum_sum
            end if
        else
            do i = 1, ubound(neighbors, 1)
                elem = 0.0_dp
                ASSERT(is_beta(src) .eqv. is_beta(neighbors(i)))
                if (IsNotOcc(ilutI, neighbors(i))) then
                    ex(2) = neighbors(i)
                    call swap_excitations(nI, ex, nJ, ex2)
                    elem = abs(get_offdiag_helement_rs_hub(nJ, ex2, .false.))
                end if
                cum_sum = cum_sum + elem
                cum_arr(i) = cum_sum
            end do
        end if

    end subroutine create_cum_list_rs_hubbard