singles_search_guga Subroutine

private subroutine singles_search_guga(recvcounts, recvdisps)

Arguments

Type IntentOptional Attributes Name
integer(kind=MPIArg), intent(in) :: recvcounts(nProcessors)
integer(kind=MPIArg), intent(in) :: recvdisps(nProcessors)

Contents

Source Code


Source Code

    subroutine singles_search_guga(recvcounts, recvdisps)
        ! make a tailored GUGA search for explicit single excitations

        ! this routine is very similar to Sing_SearchOccDets in the rdm_explicit
        ! module. see there for comments
        integer(MPIArg), intent(in) :: recvcounts(nProcessors), recvdisps(nProcessors)

        integer :: i, NoDets, StartDets, nJ(nel), PartInd, FlagsDj
        integer :: j
        integer(int_rdm) :: rdm_ind
        integer(n_int) :: ilutJ(0:GugaBits%len_tot), ilutI(0:GugaBits%len_tot)
        real(dp) :: mat_ele, sign_i(lenof_sign), sign_j(lenof_sign)
        logical :: tDetFound

        do i = 1, nProcessors

            NoDets = recvcounts(i) / (GugaBits%len_tot + 1)
            StartDets = (recvdisps(i) / (GugaBits%len_tot + 1)) + 1

            if (NoDets > 1) then

                ! extract the determinant, the matrix element and the
                ! RDM index from the single excitation array

                ! the convention is: i encode the coupling coefficient
                ! in the x0 matrix element and the sign of Di in the
                ! x1 element. and the combined RDM index in the
                ! Delta B value position!
                ilutI = Sing_ExcDjs2(:, StartDets)

                sign_i = extract_matrix_element(ilutI, 2)

                do j = StartDets + 1, (NoDets + StartDets - 1)

                    ! apparently D_i is in the first spot and all
                    ! excitations come here afterwards..

                    ilutJ = Sing_ExcDjs2(:, j)

                    call BinSearchParts_rdm(ilutJ, 1, int(TotWalkers), PartInd, tDetFound)

                    if (tDetFound) then

                        call extract_bit_rep(CurrentDets(:, PartInd), nJ, sign_j, FlagsDj)

                        mat_ele = extract_matrix_element(Sing_ExcDjs2(:, j), 1)
                        rdm_ind = extract_rdm_ind(Sing_ExcDjs2(:, j))

                        if (RDMExcitLevel == 1) then
                            call fill_sings_1rdm_guga(one_rdms, sign_i, sign_j, &
                                                      mat_ele, rdm_ind)

                        else if (RDMExcitLevel == 3) then
                            call fill_sings_2rdm_guga(two_rdm_spawn, ilutI, &
                                                      ilutJ, sign_i, sign_j, mat_ele, rdm_ind)

                        end if
                    end if
                end do
            end if
        end do

    end subroutine singles_search_guga