gen_all_excits Subroutine

public subroutine gen_all_excits(GAS_spec, nI, n_excits, det_list, ic)

@brief Get all excitated determinants from det_I that are allowed under GAS constraints.

@param[in] GAS_spec GAS specification @param[in] nI Starting determinant @param[out] n_excits Number of determinants @param[out] det_list Allocatable array of determinants in ilut format @param[in] ic Optional input for excitation level (ic=1 => singles, ic=2 => doubles) If ommited generate all.

Arguments

Type IntentOptional Attributes Name
class(GASSpec_t), intent(in) :: GAS_spec
integer, intent(in) :: nI(:)
integer, intent(out) :: n_excits
integer(kind=n_int), intent(out), allocatable :: det_list(:,:)
integer, intent(in), optional :: ic

Contents

Source Code


Source Code

    subroutine gen_all_excits(GAS_spec, nI, n_excits, det_list, ic)
        class(GASSpec_t), intent(in) :: GAS_spec
        integer, intent(in) :: nI(:)
        integer, intent(out) :: n_excits
        integer(n_int), intent(out), allocatable :: det_list(:,:)
        integer, intent(in), optional :: ic

        integer, allocatable :: singles(:, :), doubles(:, :)
        integer :: i, j

        if (present(ic)) then
            select case(ic)
            case(1)
                singles = get_available_singles(GAS_spec, nI)
                allocate(doubles(0, 0))
            case(2)
                allocate(singles(0, 0))
                doubles = get_available_doubles(GAS_spec, nI)
            end select
        else
            singles = get_available_singles(GAS_spec, nI)
            doubles = get_available_doubles(GAS_spec, nI)
        end if

        n_excits = size(singles, 2) + size(doubles, 2)
        allocate(det_list(0:niftot, n_excits))
        j = 1
        do i = 1, size(singles, 2)
            call EncodeBitDet(singles(:, i), det_list(:, j))
            j = j + 1
        end do

        do i = 1, size(doubles, 2)
            call EncodeBitDet(doubles(:, i), det_list(:, j))
            j = j + 1
        end do

        call sort(det_list, ilut_lt, ilut_gt)
    end subroutine gen_all_excits