create_cum_list_rs_hubbard_transcorr_single Subroutine

private subroutine create_cum_list_rs_hubbard_transcorr_single(nI, ilutI, src, cum_arr, cum_sum, tgt, p_orb)

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: nI(nel)
integer(kind=n_int), intent(in) :: ilutI(0:NIfTot)
integer, intent(in) :: src
real(kind=dp), intent(out) :: cum_arr(nBasis/2)
real(kind=dp), intent(out) :: cum_sum
integer, intent(in), optional :: tgt
real(kind=dp), intent(out), optional :: p_orb

Contents


Source Code

    subroutine create_cum_list_rs_hubbard_transcorr_single(nI, ilutI, src, &
                                                           cum_arr, cum_sum, tgt, p_orb)
        ! with transcorrelation use a different cum-list creator, due to
        ! longer range single excitations possible now.
        integer, intent(in) :: nI(nel), src
        integer(n_int), intent(in) :: ilutI(0:NIfTot)
        real(dp), intent(out) :: cum_arr(nBasis / 2), cum_sum
        integer, intent(in), optional :: tgt
        real(dp), intent(out), optional :: p_orb
#ifdef DEBUG_
        character(*), parameter :: this_routine = "create_cum_list_rs_hubbard_transcorr_single"
#endif
        integer :: spin, ex(2, maxExcit), nJ(nel), i, orb
        integer, allocatable :: ex2(:, :)
        real(dp) :: elem

        ASSERT(IsOcc(ilutI, src))

        cum_arr = 0.0_dp
        cum_sum = 0.0_dp

        ! 0.. alpha
        ! 1.. beta
        spin = get_spin(src) - 1

        ex = 0
        ex(1, 1) = src

        if (present(tgt)) then
            ASSERT(present(p_orb))
            ASSERT(same_spin(src, tgt))

            p_orb = 0.0_dp

            do i = 1, nBasis / 2
                elem = 0.0_dp

                ! take the same spin
                orb = 2 * i - spin

                ASSERT(same_spin(src, orb))

                if (IsNotOcc(ilutI, orb)) then

                    ! i am still not sure about the ordering of these weights..
                    ex(2, 1) = orb
                    call swap_excitations(nI, ex, nJ, ex2)
                    elem = abs(get_single_helem_rs_hub_transcorr(nJ, ex2(:, 1), .false.))
                    ! elem = abs(get_single_helem_rs_hub_transcorr(nI, ex(:,1), .false.))

                end if
                cum_sum = cum_sum + elem

                if (tgt == orb) then
                    p_orb = elem
                end if
            end do
            if (cum_sum < EPS) then
                p_orb = 0.0_dp
            else
                p_orb = p_orb / cum_sum
            end if
        else
            do i = 1, nBasis / 2
                elem = 0.0_dp
                orb = 2 * i - spin

                ASSERT(same_spin(src, orb))

                if (IsNotOcc(ilutI, orb)) then
                    ex(2, 1) = orb
                    call swap_excitations(nI, ex, nJ, ex2)
                    ! elem = abs(get_single_helem_rs_hub_transcorr(nI, ex(:,1), .false.))
                    elem = abs(get_single_helem_rs_hub_transcorr(nJ, ex2(:, 1), .false.))
                end if

                cum_sum = cum_sum + elem
                cum_arr(i) = cum_sum
            end do
        end if

    end subroutine create_cum_list_rs_hubbard_transcorr_single