generate_connected_space_kpnt Subroutine

public subroutine generate_connected_space_kpnt(original_space_size, original_space, connected_space_size, connected_space, tSinglesOnlyOpt)

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: original_space_size
integer(kind=n_int), intent(in) :: original_space(0:,:)
integer, intent(inout) :: connected_space_size
integer(kind=n_int), intent(out), optional :: connected_space(0:,:)
logical, intent(in), optional :: tSinglesOnlyOpt

Contents


Source Code

    subroutine generate_connected_space_kpnt(original_space_size, original_space, &
                                             connected_space_size, connected_space, tSinglesOnlyOpt)

        use SymExcit2, only: gensymexcitit2par_worker

        ! This is the same as generate_connected_space, but using the old excitations
        ! generators because the new ones don't work with the tKPntSym option.

        ! This routine either counts or generates all the determinants connected to those in
        ! original_space. If connected_space is not present then they will only be counted,
        ! else they will be stored in connected_space. If tSinglesOnlyOpt is present and
        ! also .true. then actually this routine will only generate all single excitations
        ! (regardless of the system being studied), otherwise it will generate all connected
        ! determinants.

        use guga_bitRepOps, only: convert_ilut_toGUGA, convert_ilut_toNECI
        use guga_excitations, only: actHamiltonian
        use SystemData, only: tGUGA

        integer :: nexcit, j
        integer(n_int) :: ilutG(0:nifguga)
        integer(n_int), allocatable :: excitations(:, :)

        integer, intent(in) :: original_space_size
        integer(n_int), intent(in) :: original_space(0:, :)
        integer, intent(inout) :: connected_space_size
        integer(n_int), optional, intent(out) :: connected_space(0:, :)
        logical, intent(in), optional :: tSinglesOnlyOpt

        integer(n_int) :: ilutJ(0:NIfTot)
        integer :: nI(nel), nJ(nel)
        integer :: i, excit(2, 2), ex_flag
        integer :: nStore(6)
        integer, allocatable :: excit_gen(:)
        logical :: tStoreConnSpace, tSinglesOnly, tTempUseBrill, tAllExcitFound
        character(*), parameter :: this_routine = "generate_connected_space_kpnt"
        integer :: n_excits
        integer(n_int), allocatable :: temp_dets(:, :)

        if (present(connected_space)) then
            tStoreConnSpace = .true.
        else
            tStoreConnSpace = .false.
        end if

        tSinglesOnly = .false.
        if (present(tSinglesOnlyOpt)) then
            if (tSinglesOnlyOpt) tSinglesOnly = .true.
        end if

        connected_space_size = 0

        ! Over all the states in the original list:
        do i = 1, original_space_size

            call decode_bit_det(nI, original_space(0:NIfTot, i))

            if (t_k_space_hubbard) then

                ! for every loop we have to save the excitations per
                ! do we have to check if the list is unique?? i guess i do
#ifdef CMPLX_
                call stop_all(this_routine, "not implemented for complex")
#else
                call gen_all_excits_k_space_hubbard(nI, n_excits, temp_dets)
#endif

                if (tStoreConnSpace) then
                    connected_space(0:nifd, connected_space_size + 1:connected_space_size + n_excits) &
                        = temp_dets(0:nifd, :)
                end if
                connected_space_size = connected_space_size + n_excits

                ! GUGA changes:
                ! my actHamiltonian routine seems to satisfy the k-point symmetry
                ! so it should be straight forward to implemt it here too
            else if (tGUGA) then
                if (tSinglesOnly) call stop_all(this_routine, "don't use tSinglesOnly with GUGA")

                call convert_ilut_toGUGA(original_space(:, i), ilutG)

                call actHamiltonian(ilutG, CSF_Info_t(ilutG), excitations, nexcit)

                if (tStoreConnSpace) then
                    do j = 1, nexcit
                        call convert_ilut_toNECI(excitations(:, j), &
                                                 connected_space(:, connected_space_size + j))
                    end do
                end if

                connected_space_size = connected_space_size + nexcit

                deallocate(excitations)
                call LogMemDealloc(this_routine, tag_excitations)

            else

                call init_generate_connected_space(nI, ex_flag, tAllExcitFound, excit, excit_gen, nstore, tTempUseBrill)

                if (tSinglesOnly) ex_flag = 1

                do while (.true.)

                    call generate_connection_kpnt(nI, original_space(:, i), nJ, &
                                                  ilutJ, ex_flag, tAllExcitFound, nStore, excit_gen, &
                                                  ncon=connected_space_size)

                    if (tStoreConnSpace) then
                        connected_space(0:NIfD, connected_space_size) = ilutJ(0:NIfD)
                    end if

                end do

                call generate_connection_kpnt(nI, original_space(:, i), nJ, ilutJ, ex_flag, tAllExcitFound, &
                                              nStore, excit_gen, ncon=connected_space_size)

                if (tAllExcitFound) exit

                if (tStoreConnSpace) connected_space(0:NIfD, connected_space_size) = ilutJ(0:NIfD)

                deallocate(excit_gen)

            end if
        end do

    end subroutine generate_connected_space_kpnt