gen_a_orb_cum_list_guga_mol Subroutine

private subroutine gen_a_orb_cum_list_guga_mol(csf_i, occ_orbs, cum_arr, tgt_orb, pgen)

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)
integer, intent(in), optional :: tgt_orb
real(kind=dp), intent(out), optional :: pgen

Contents


Source Code

    subroutine gen_a_orb_cum_list_guga_mol(csf_i, occ_orbs, cum_arr, tgt_orb, pgen)
        ! subroutine to generate the molecular cumullative probability
        ! distribution. there are no (atleast until now) addiditonal restrictions
        ! or some spin alignement restrictions to generate this list..
        type(CSF_Info_t), intent(in) :: csf_i
        integer, intent(in) :: occ_orbs(2)
        real(dp), intent(out) :: cum_arr(nSpatOrbs)
        integer, intent(in), optional :: tgt_orb
        real(dp), intent(out), optional :: pgen
        character(*), parameter :: this_routine = "gen_a_orb_cum_list_guga_mol"

        integer :: orb
        real(dp) :: cum_sum
        ! have to think about the ordering of the two electronic indices
        ! in the FCIDUMP integral choosing..
        ! since there is a relative sign, but i have not yet figured out if
        ! i can predetermine this relative sign..

        cum_sum = 0.0_dp

        if (present(tgt_orb)) then
            ASSERT(present(pgen))
            do orb = 1, nSpatOrbs
                if (csf_i%stepvector(orb) /= 3) then
                    cum_sum = cum_sum + get_guga_integral_contrib(occ_orbs, orb, -1)
                end if
                cum_arr(orb) = cum_sum

                if (orb == tgt_orb) then

                end if
            end do

        else
            do orb = 1, nSpatOrbs
                ! only check non-double occupied orbitals
                if (csf_i%stepvector(orb) /= 3) then
                    cum_sum = cum_sum + get_guga_integral_contrib(occ_orbs, orb, -1)
                end if
                cum_arr(orb) = cum_sum
            end do
        end if

    end subroutine gen_a_orb_cum_list_guga_mol