calc_all_excitations Subroutine

public subroutine calc_all_excitations(ilut, non_zero_list, n_non_zero)

Arguments

Type IntentOptional Attributes Name
integer(kind=n_int), intent(in) :: ilut(0:niftot)
integer(kind=n_int), intent(out), allocatable :: non_zero_list(:,:)
integer, intent(out) :: n_non_zero

Contents

Source Code


Source Code

    subroutine calc_all_excitations(ilut, non_zero_list, n_non_zero)
        ! routine which calculates all the excitations for a given
        ! determinant.. have to move that to a different file later on..
        ! does not fit in here!
        use bit_reps, only: decode_bit_det
        use GenRandSymExcitNUMod, only: init_excit_gen_store
        use SymExcit3, only: CountExcitations3, GenExcitations3
        use DetBitOps, only: EncodeBitDet
        use sort_mod, only: sort

        integer(n_int), intent(in) :: ilut(0:niftot)
        integer(n_int), intent(out), allocatable :: non_zero_list(:, :)
        integer, intent(out) :: n_non_zero
        character(*), parameter :: this_routine = "calc_all_excitations"

        integer :: src_det(nel), nsing, ndoub, ndet, ex(2, 2), flag, det(nel), &
                   nexcit, i, cnt
        type(excit_gen_store_type) :: store
        logical :: found_all, par
        HElement_t(dp), allocatable :: hel_list(:)
        integer(n_int), allocatable :: det_list(:, :)
        logical, allocatable :: t_non_zero(:)

        ! Decode the determiant
        call decode_bit_det(src_det, ilut)

        ! Initialise
        call init_excit_gen_store(store)

        ! How many connected determinants are we expecting?
        call CountExcitations3(src_det, 3, nsing, ndoub)
        nexcit = nsing + ndoub
        allocate(det_list(0:NIfTot, nexcit))
        allocate(hel_list(nexcit))
        hel_list = 0.0_dp

        ! Loop through all of the possible excitations
        ndet = 0
        found_all = .false.
        ex = 0
        flag = 3
        call GenExcitations3(src_det, ilut, det, flag, ex, par, found_all, &
                             .false.)

        if (tGUGA) then
            call stop_all(this_routine, "modify get_helement for GUGA")
        end if

        do while (.not. found_all)
            ndet = ndet + 1
            call EncodeBitDet(det, det_list(:, ndet))

            hel_list(ndet) = get_helement(src_det, det, ilut, det_list(:, ndet))

            call GenExcitations3(src_det, ilut, det, flag, ex, par, &
                                 found_all, .false.)

        end do

        if (ndet /= nexcit) &
            call stop_all(this_routine, "Incorrect number of excitations found")

        ! now i should loop over the hel_list and check the actual number
        ! of non-zero excitations
        allocate(t_non_zero(nexcit))
        t_non_zero = .false.

        n_non_zero = 0

        do i = 1, ndet
            if (abs(hel_list(i)) > 1.0e-6) then
                n_non_zero = n_non_zero + 1
                t_non_zero(i) = .true.
            end if
        end do

        allocate(non_zero_list(0:niftot, n_non_zero))

        cnt = 0

        do i = 1, ndet
            if (t_non_zero(i)) then
                cnt = cnt + 1
                non_zero_list(:, cnt) = det_list(:, i)
            end if
        end do

        deallocate(det_list)
        deallocate(t_non_zero)
        deallocate(hel_list)

        ! now i have a list of all non-zero excitations..

        ! for now this creates all excitation, independent if the matrix
        ! element is actually zero..
        ! so also check matrix element!

        ! Sort the dets, so they are easy to find by binary searching
        call sort(non_zero_list, ilut_lt, ilut_gt)

    end subroutine calc_all_excitations