gen_cum_list_guga_single_3 Subroutine

private subroutine gen_cum_list_guga_single_3(nI, csf_i, orb_i, cc_i, cum_arr)

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: nI(nel)
type(CSF_Info_t), intent(in) :: csf_i
integer, intent(in) :: orb_i
integer, intent(in) :: cc_i
real(kind=dp), intent(out) :: cum_arr(OrbClassCount(cc_i))

Contents


Source Code

    subroutine gen_cum_list_guga_single_3(nI, csf_i, orb_i, cc_i, cum_arr)
        ! specific single orbital picker if stepvector of electron (i) is 3
        integer, intent(in) :: nI(nel)
        type(CSF_Info_t), intent(in) :: csf_i
        integer, intent(in) :: orb_i, cc_i
        real(dp), intent(out) :: cum_arr(OrbClassCount(cc_i))
        character(*), parameter :: this_routine = "gen_cum_list_guga_single_3"

        integer :: nOrb, i, label_index, j, n_id(nEl), id_i, s_orb, spin
        real(dp) :: cum_sum, hel
        ! in the case of a 3 there are actually no additional, restrictions

        nOrb = OrbClassCount(cc_i)
        label_index = SymLabelCounts2(1, cc_i)

        n_id = gtID(nI)
        id_i = gtID(orb_i)

        cum_sum = 0.0_dp

        if (is_beta(orb_i)) then
            spin = 1
        else
            spin = 0
        end if

        do i = 1, nOrb
            ! ich glaub es ist besser ueber spatial orbitals zu loopen,
            ! sonst ist das so doof mit dem ignorieren der spin symmetrie..
            s_orb = sym_label_list_spat(label_index + i - 1)

            if (s_orb == id_i) then
                cum_arr(i) = cum_sum
                cycle
            end if

            hel = 0.0_dp
            ! here i can exclude it just with check if its not occupied,
            ! since if it a 1 or 2 i will at some point
            ! new spat. orb. impl. : check if not three
            ! or do a select case on stepvector?
            select case (csf_i%stepvector(s_orb))

                ! with the predetermination of the stepvalue at (a) and since
                ! this routine only is called in the case of d(i) = 3 it makes
                ! it much easier to determine the sign of the two-particle
                ! contribution to the single excitation matrix element!
            case (0)

                hel = hel + abs(GetTMatEl(orb_i, 2 * s_orb - spin))

                ! here a contribution from orbital id_i

                if (.not. t_mixed_hubbard) then
                    hel = hel + abs(get_umat_el(id_i, id_i, s_orb, id_i))

                    ! do the loop over all the other electrons
                    ! (is this always symmetrie allowed?..)
                    do j = 1, nEl

                        ! todo: finish all contributions later for now only do
                        ! those which are the same for all
                        ! have to exclude both electrons at spatial orb i
                        if (n_id(j) == id_i) cycle
                        hel = hel + abs(get_umat_el(id_i, n_id(j), s_orb, n_id(j)))
                        hel = hel + abs(get_umat_el(id_i, n_id(j), n_id(j), s_orb))

                    end do
                end if

            case (1)

                hel = hel + abs(GetTMatEl(orb_i, 2 * s_orb - spin))

                if (.not. t_mixed_hubbard) then
                    ! now contribution for both start and end
                    hel = hel + abs(get_umat_el(id_i, id_i, s_orb, id_i))
                    hel = hel + abs(get_umat_el(id_i, s_orb, s_orb, s_orb))

                    ! do the loop over all the other electrons
                    ! (is this always symmetrie allowed?..)
                    do j = 1, nEl

                        ! todo: finish all contributions later for now only do
                        ! those which are the same for all
                        if (n_id(j) == id_i .or. n_id(j) == s_orb) cycle
                        hel = hel + abs(get_umat_el(id_i, n_id(j), s_orb, n_id(j)))
                        hel = hel + abs(get_umat_el(id_i, n_id(j), n_id(j), s_orb))

                    end do
                end if

            case (2)

                hel = hel + abs(GetTMatEl(orb_i, 2 * s_orb - spin))
                ! do the loop over all the other electrons
                ! (is this always symmetrie allowed?..)

                if (.not. t_mixed_hubbard) then
                    hel = hel + abs(get_umat_el(id_i, id_i, s_orb, id_i))
                    hel = hel + abs(get_umat_el(id_i, s_orb, s_orb, s_orb))

                    do j = 1, nEl

                        ! todo: finish all contributions later for now only do
                        ! those which are the same for all
                        if (n_id(j) == id_i .or. n_id(j) == s_orb) cycle

                        hel = hel + abs(get_umat_el(id_i, n_id(j), s_orb, n_id(j)))
                        hel = hel + abs(get_umat_el(id_i, n_id(j), n_id(j), s_orb))

                    end do
                end if

            case (3)
                ! do nothing actually

            case default
                ! error happend
                call stop_all(this_routine, "stepvalue /= {0,1,2,3}! something's wrong!")

            end select

            cum_sum = cum_sum + abs_l1(hel)
            cum_arr(i) = cum_sum

        end do

    end subroutine gen_cum_list_guga_single_3