gen_all_singles_rs_hub_hop_transcorr Subroutine

private subroutine gen_all_singles_rs_hub_hop_transcorr(nI, n_excits, det_list)

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: nI(nel)
integer, intent(out) :: n_excits
integer(kind=n_int), intent(out), allocatable :: det_list(:,:)

Contents


Source Code

    subroutine gen_all_singles_rs_hub_hop_transcorr(nI, n_excits, det_list)
        integer, intent(in) :: nI(nel)
        integer, intent(out) :: n_excits
        integer(n_int), intent(out), allocatable :: det_list(:, :)
        integer(n_int) :: ilut(0:NIfTot), ilutJ(0:NIfTot)
        integer :: n_bound, i, src, ex(2, 2), pos, a
        integer(n_int), allocatable :: temp_list(:, :)
        real(dp) :: elem

        call EncodeBitDet(nI, ilut)

        n_excits = 1

        n_bound = nel * (nbasis - nel)
        allocate(temp_list(0:NIfTot, n_bound))
        temp_list = 0_n_int
        ex = 0

        ! loop over electrons
        do i = 1, nel
            ! but now loop over all orbitals
            src = nI(i)
            ex(1, 1) = src

            do a = 1, nbasis

                ! only same-spin excitations possible
                if (same_spin(src, a) .and. IsNotOcc(ilut, a)) then
                    ex(2, 1) = a
                    elem = abs(get_helement_lattice(nI, 1, ex, .false.))

                    if (elem > EPS) then

                        ilutJ = make_ilutJ(ilut, ex, 1)

                        ! actually a search is not really necessary.. since
                        ! all the single excitations are unique.. but
                        ! just to be sure
                        pos = binary_search_ilut(temp_list(:, 1:n_excits), ilutJ, nifd + 1)

                        if (pos < 0) then

                            temp_list(:, n_excits) = ilutJ
                            call sort(temp_list(:, 1:n_excits), ilut_lt, ilut_gt)
                            n_excits = n_excits + 1
                            ! damn.. i have to sort everytime i guess..
                        end if
                    end if
                end if
            end do
        end do

        n_excits = n_excits - 1
        allocate(det_list(0:NIfTot, n_excits), source=temp_list(:, 1:n_excits))

        call sort(det_list, ilut_lt, ilut_gt)

    end subroutine gen_all_singles_rs_hub_hop_transcorr