generate_sing_doub_determinants Subroutine

public subroutine generate_sing_doub_determinants(ilut_list, space_size, only_keep_conn)

Arguments

Type IntentOptional Attributes Name
integer(kind=n_int), intent(inout) :: ilut_list(0:,:)
integer, intent(inout) :: space_size
logical, intent(in) :: only_keep_conn

Contents


Source Code

    subroutine generate_sing_doub_determinants(ilut_list, space_size, only_keep_conn)

        ! In/Out: ilut_list - List of determinants generated.
        ! In/Out: space_size - Number of determinants in the generated space.
        !             If ilut_list is not empty on input and you want to keep
        !             the states already in it, then on input space_size should
        !             be equal to the number of states to be kept in ilut_list,
        !             and new states will be added in from space_size+1.
        !             Otherwise, space_size must equal 0 on input.
        !             On output space_size will equal the total number of
        !             generated plus what space_size was on input.

        integer(n_int), intent(inout) :: ilut_list(0:, :)
        integer, intent(inout) :: space_size
        logical, intent(in) :: only_keep_conn

        integer(n_int) :: ilut(0:NIfTot)
        integer :: nI(nel)
        integer :: excit(2, 2)
        integer :: nsing, ndoub, ex_flag
        logical :: tAllExcitFound, tParity
        HElement_t(dp) :: HEl

        type(ExcitGenSessionType) :: session

        ! Always generate both the single and double excitations.
        ex_flag = 3

        ! Start by adding the HF state.
        call add_state_to_space(ilutHF, ilut_list, space_size)

        if (tGUGA) then
            call stop_all("generate_sing_doub_determinants", &
                          "modify get_helement for GUGA")
        end if
        if (tKPntSym) then
            call enumerate_sing_doub_kpnt(ex_flag, only_keep_conn, nsing, ndoub, .true., ilut_list, space_size)
        else

            tAllExcitFound = .false.
            excit = 0

            do while (.true.)
                ! Generate the next determinant.
                if (tReltvy) then
                    call GenExcitations4(session, HFDet, nI, ex_flag, excit, tParity, tAllExcitFound, .false.)
                else
                    call GenExcitations3(HFDet, ilutHF, nI, ex_flag, excit, tParity, tAllExcitFound, .false.)
                end if
                if (tAllExcitFound) exit

                call EncodeBitDet(nI, ilut)
                ! If using a deterministic space connected to the Hartree-Fock
                ! then check that this determinant is actually connected to it!
                if (only_keep_conn) then
                    HEl = get_helement(HFDet, nI, ilutHF, ilut)
                    ! [W.D. 15.5.2017:]
                    ! is this still enough, even for Hamiltonians containing
                    ! complex entries??
                    ! and why is the cast to real(dp) done here??
!                     if (abs(real(HEl,dp)) < 1.e-12_dp) cycle
                    if (abs(HEl) < 1.e-12_dp) cycle
                end if
                call add_state_to_space(ilut, ilut_list, space_size, nI)
            end do
        end if

    end subroutine generate_sing_doub_determinants