generate_connection_normal Subroutine

public subroutine generate_connection_normal(nI, ilutI, nJ, ilutJ, ex_flag, excit, tAllExcitFound, hel, ncon)

Arguments

Type IntentOptional Attributes Name
integer :: nI(nel)
integer(kind=n_int), intent(in) :: ilutI(0:NIfTot)
integer, intent(out) :: nJ(nel)
integer(kind=n_int), intent(out) :: ilutJ(0:NIfTot)
integer, intent(inout) :: ex_flag
integer, intent(inout) :: excit(2,2)
logical, intent(inout) :: tAllExcitFound
real(kind=dp), intent(out), optional :: hel
integer, intent(inout), optional :: ncon

Contents


Source Code

    subroutine generate_connection_normal(nI, ilutI, nJ, ilutJ, ex_flag, excit, &
                                          tAllExcitFound, hel, ncon)

        use procedure_pointers, only: get_conn_helement
        use SymExcit4, only: GenExcitations4

        integer :: nI(nel)
        integer(n_int), intent(in) :: ilutI(0:NIfTot)
        integer, intent(out) :: nJ(nel)
        integer(n_int), intent(out) :: ilutJ(0:NIfTot)
        integer, intent(inout) :: ex_flag
        integer, intent(inout) :: excit(2, 2)
        logical, intent(inout) :: tAllExcitFound
        HElement_t(dp), optional, intent(out) :: hel
        integer, optional, intent(inout) :: ncon

        integer(n_int) :: ilut_tmp(0:NIfTot)
        integer :: ic
        logical :: tParity
        HElement_t(dp) :: hel_unused

        ! Generate the next determinant.
        if (tReltvy) then
            call GenExcitations4(session, nI, nJ, ex_flag, excit, tParity, tAllExcitFound, .false.)
        else
            call GenExcitations3(nI, ilutI, nJ, ex_flag, excit, tParity, &
                                 tAllExcitFound, .false.)
        end if
        if (tAllExcitFound) return

        ! Encode nJ in ilutJ.
        call EncodeBitDet(nJ, ilutJ)

        ! If using HPHFs, and if this isn't the allowed state for this HPHF, find
        ! the the allowed state and continue with this as ilut.
        if (tHPHF .and. (.not. TestClosedShellDet(ilutJ))) then
            if (.not. IsAllowedHPHF(ilutJ(0:NIfD))) then
                ilut_tmp = ilutJ
                call FindExcitBitDetSym(ilut_tmp, ilutJ)
            end if
        end if

        if (present(HEl)) then
            ! Map ex_flag values (1, 2 or 3) to 1, 2 and 1, respectively.
            ic = 2 - mod(ex_flag, 2)
            HEl = get_conn_helement(nI, nJ, ilutI, ilutJ, ic, excit, tParity, hel_unused)
        end if

        if (present(ncon)) ncon = ncon + 1

    end subroutine generate_connection_normal