gen_ab_cum_list_3 Subroutine

private subroutine gen_ab_cum_list_3(csf_i, occ_orbs, cum_arr, excit_arr, orb_arr)

Arguments

Type IntentOptional Attributes Name
type(CSF_Info_t), intent(in) :: csf_i
integer, intent(in) :: occ_orbs(2)
real(kind=dp), intent(out) :: cum_arr(nSpatOrbs)
type(ExcitationInformation_t), intent(out) :: excit_arr(nSpatOrbs)
integer, intent(out) :: orb_arr(nSpatOrbs)

Contents

Source Code


Source Code

    subroutine gen_ab_cum_list_3(csf_i, occ_orbs, cum_arr, excit_arr, orb_arr)
        ! specific routine, when 2 already picked orbtitals are from same
        ! spatial orbital
        type(CSF_Info_t), intent(in) :: csf_i
        integer, intent(in) :: occ_orbs(2)
        real(dp), intent(out) :: cum_arr(nSpatOrbs)
        integer, intent(out) :: orb_arr(nSpatOrbs)
        type(ExcitationInformation_t), intent(out) :: excit_arr(nSpatOrbs)

        integer :: ki(3), kj(3), ka(3), kb(3), a, b, en, st, i, z_ind
        real(dp) :: cum_sum, contrib
        type(ExcitationInformation_t) :: excitInfo
        logical :: tSwitch

        z_ind = ubound(kPointToBasisFn, 3)

        ki = G1(occ_orbs(1))%k
        kj = G1(occ_orbs(2))%k

        i = gtID(occ_orbs(1))
        ! redo this!
        cum_sum = 0.0_dp
        orb_arr = 0
        b = 0

        ! GUGA restrictions change depending on n(i), n(j)
        ! I == J -> n(i) = 3
        ! no restrictions on a
        ! todo!! stupid! if i = j => ki = kj => ka = kb! only have to loop
        ! once over orbitals!
        ! NO! 2ki = ka + kb -> kann trotzdem noch mehrere möglichkeiten geben!
        ! change that! (but maybe reuse it for 3 3 case
        ! aber vl. doch spatial orbitals...
        do a = 1, i - 1

            contrib = 0.0_dp
            excitInfo%valid = .false.
            tSwitch = .true.

            if (csf_i%stepvector(a) /= 3) then
                ! determine fiting b by k-point restrictions ->
                ! and only allow b > a
                ka = G1(2 * a)%k
                kb = ki + kj - ka

                call MomPbcSym(kb, nBasisMax)

                ! is b allowed by the size of space? -> TODO: what does that mean?
                ! this is specific for the UEG.. so i do not need it for the
                ! Hubbard model i am implementing right now..

                ! there is no "spin" restriction in the guga case
                ! so both possible b orbs have to be checked
                ! its actually NOT possible that ka = ki !! so do not have
                ! to check that case!
                b = gtID(KPointToBasisFn(kb(1), kb(2), z_ind, 1))

                ! this should work as we pick beta orbital first

                ! orb_a should not be equal orb_b due to k point restrictions
                ! hm... it could happen i guess.. but i have to check it
                ! since i did not include this in the possibilities
                ! below..
                ! ok so i just realised this case can happen...
                ! so i have to fix that and take all the possible
                ! combinations into account..
                if (a == b) then
                    ! then orb a has to be empty!
                    ! actually this type of excitation is not even possible
                    ! with momentum conservation.. todo!
                    if (csf_i%stepvector(a) == 0) then
                        !todo: remove this possibility then!
                        ! _RR_(ab) > ^RR^(ij)
                        excitInfo = assign_excitInfo_values_double( &
                                    excit_type%fullstart_stop_alike, &
                                    gen_type%R, gen_type%R, gen_type%R, gen_type%R, gen_type%R, &
                                    a, i, a, i, a, a, i, i, 0, 2, 1.0_dp, 1.0_dp)

                        ! use uniform, since all the integrals are equal
                        ! anyway
                        ! if (a==b) there is only one entry in the cum-list
                        ! which leads to that excitation, compared to the
                        ! usual 2 for other types of excitations..
                        ! so i could consider multiplying it by 2 just
                        ! to make the pgens more homogenous.. since matrix
                        ! elements usually will be the same i guess..
                        ! no, do it by making it twice as big later..
                        ! NO! i have to do it here or otherwise the actual
                        ! probability of picking this orbital is not twice as
                        ! high..
                        contrib = 2.0_dp

                    end if
                else
                    if (csf_i%stepvector(b) /= 3) then
                        ! only then its a possible excitation
                        ! hm... i guess i should check if there is a possible
                        ! switch in this case too.. if both a and b have the
                        ! same stepvector there needs to be a switch possible
                        ! at some place to lead to a valid excitation..
                        ! i have not yet done that in the molecular excitation
                        ! generator and i do not know why.. was there a reason?

                        st = min(a, b)
                        en = max(a, b)

                        if (csf_i%stepvector(a) == csf_i%stepvector(b)) then
                            if (csf_i%stepvector(a) == 1 .and. &
                                count_alpha_orbs_ij(csf_i, st, en) == 0) tSwitch = .false.
                            if (csf_i%stepvector(a) == 2 .and. &
                                count_beta_orbs_ij(csf_i, st, en) == 0) tSwitch = .false.
                        end if

                        if (tSwitch) then
                            contrib = 1.0_dp

                            ! then have to determine the order of orbitals
                            if (b > i) then
                                ! _R(a) > ^RL_(ij) > ^L(b)
                                excitInfo = assign_excitInfo_values_double( &
                                            excit_type%single_overlap_R_to_L, &
                                            gen_type%L, gen_type%R, gen_type%R, gen_type%R, gen_type%L, &
                                            b, i, a, i, a, i, i, b, &
                                            0, 2, 1.0_dp, 1.0_dp, 1)

                            else
                                ! only the extrema count
                                ! if both have the same stepvalue i need to check
                                ! if there is a switch possible i guess..
                                ! _R(min) > _RR(max) > ^RR^(ij)
                                excitInfo = assign_excitInfo_values_double( &
                                            excit_type%fullstop_raising, &
                                            gen_type%R, gen_type%R, gen_type%R, gen_type%R, gen_type%R, &
                                            st, i, en, i, st, en, i, i, 0, 2, 1.0_dp, 1.0_dp)
                            end if
                        end if
                    end if
                end if
            end if
            cum_sum = cum_sum + contrib
            cum_arr(a) = cum_sum
            excit_arr(a) = excitInfo
            orb_arr(a) = b
        end do
        ! zero out orbital (i)
        if (i == 1) then
            cum_arr(i) = 0.0_dp
        else
            cum_arr(i) = cum_arr(i - 1)
        end if

        do a = i + 1, nSpatOrbs

            contrib = 0.0_dp
            excitInfo%valid = .false.
            tSwitch = .true.

            if (csf_i%stepvector(a) /= 3) then
                ka = G1(2 * a)%k
                kb = ki + kj - ka

                call MomPbcSym(kb, nBasisMax)

                ! there is no "spin" restriction in the guga case
                ! so both possible b orbs have to be checked
                ! its actually NOT possible that ka = ki !! so do not have
                ! to check that case!
                b = gtID(KPointToBasisFn(kb(1), kb(2), z_ind, 1))

                ! this should work as we pick beta orbital first
                ! as i realised this possibility below is possible!
                if (a == b) then
                    !todo: i think with momentum conservation this below is not
                    ! even possible, if really not -> remove it
                    if (csf_i%stepvector(a) == 0) then
                        ! _LL_(ij) > ^LL^(ab)
                        excitInfo = assign_excitInfo_values_double( &
                                    excit_type%fullstart_stop_alike, &
                                    gen_type%L, gen_type%L, gen_type%L, gen_type%L, gen_type%L, &
                                    a, i, a, i, i, i, a, a, 0, 2, 1.0_dp, 1.0_dp)

                        ! same convention as above since only one entry
                        ! for this kind of ab combination in list
                        ! no! just generally add up the 2 combinations later
                        ! on at the end of the orbital picking!
                        contrib = 2.0_dp

                    end if
                else
                    if (csf_i%stepvector(b) /= 3) then
                        ! only then its a possible excitation

                        st = min(a, b)
                        en = max(a, b)

                        if (csf_i%stepvector(a) == csf_i%stepvector(b)) then
                            if (csf_i%stepvector(a) == 1 .and. &
                                count_alpha_orbs_ij(csf_i, st, en) == 0) tSwitch = .false.
                            if (csf_i%stepvector(a) == 2 .and. &
                                count_beta_orbs_ij(csf_i, st, en) == 0) tSwitch = .false.
                        end if

                        if (tSwitch) then
                            contrib = 1.0_dp

                            ! now we know a is bigger than (i)
                            if (b < i) then
                                ! _R(b) > ^RL_(ij) > ^L(a)
                                excitInfo = assign_excitInfo_values_double( &
                                            excit_type%single_overlap_R_to_L, &
                                            gen_type%L, gen_type%R, gen_type%R, gen_type%R, gen_type%L, &
                                            a, i, b, i, b, i, i, a, &
                                            0, 2, 1.0_dp, 1.0_dp, 1)

                            else
                                ! only extrema count
                                ! _LL_(ij) > ^LL(min) > ^L(max)
                                excitInfo = assign_excitInfo_values_double( &
                                            excit_type%fullstart_lowering, &
                                            gen_type%L, gen_type%L, gen_type%L, gen_type%L, gen_type%L, &
                                            st, i, en, i, i, i, st, en, 0, 2, 1.0_dp, 1.0_dp)
                            end if
                        end if
                    end if
                end if
            end if
            ! update the cumsum
            cum_sum = cum_sum + contrib
            cum_arr(a) = cum_sum
            excit_arr(a) = excitInfo
            orb_arr(a) = b
        end do

        excitInfo%valid = .false.
        excit_arr(i) = excitInfo

    end subroutine gen_ab_cum_list_3