add_to_rdm_spawn_t Subroutine

public subroutine add_to_rdm_spawn_t(spawn, i, j, k, l, contrib_sign, spinfree, nearly_full)

Arguments

Type IntentOptional Attributes Name
type(rdm_spawn_t), intent(inout) :: spawn
integer, intent(in) :: i
integer, intent(in) :: j
integer, intent(in) :: k
integer, intent(in) :: l
real(kind=dp), intent(in) :: contrib_sign(spawn%rdm_send%sign_length)
logical, intent(in) :: spinfree
logical, intent(inout), optional :: nearly_full

Contents

Source Code


Source Code

    subroutine add_to_rdm_spawn_t(spawn, i, j, k, l, contrib_sign, spinfree, nearly_full)

        ! In/Out: rdm_spawn - the rdm_spawn_t object to which contributions will be added.
        ! In: i, j, k, l - orbitals labels for the RDM contribution, with i<j, k<l.
        ! In: contrib_sign - the sign (amplitude) of the contribution to be added.
        ! In: spinfree - is the RDM being created to be output directly in spinfree form?
        ! In/Out: nearly_full - make this logical true if we come close to filling a
        !     processes section of the spawning list.

        use hash, only: hash_table_lookup, add_hash_table_entry
        use SystemData, only: nbasis

        type(rdm_spawn_t), intent(inout) :: spawn
        integer, intent(in) :: i, j, k, l ! spin or spatial orbital
        real(dp), intent(in) :: contrib_sign(spawn%rdm_send%sign_length)
        logical, intent(in) :: spinfree
        logical, intent(inout), optional :: nearly_full

        integer(int_rdm) :: ijkl
        integer :: ij_compressed, proc, ind, hash_val, slots_left, kl_compressed
        real(dp) :: real_sign_old(spawn%rdm_send%sign_length), real_sign_new(spawn%rdm_send%sign_length)
        logical :: tSuccess
        character(*), parameter :: t_r = 'add_to_rdm_spawn_t'

        associate(rdm => spawn%rdm_send)

            ! Calculate combined RDM labels. A different ordering is used for
            ! outputting RDMs with and without spin. The following definitions
            ! will aid ordering via a sort operation later.
            if (tGUGA) then
                ijkl = contract_2_rdm_ind(i, j, k, l)
            else
                if (spinfree) then
                    call calc_combined_rdm_label(k, l, j, i, ijkl)
                else
                    call calc_combined_rdm_label(i, j, k, l, ijkl)
                end if
            end if

            ! Search to see if this RDM element is already in the RDM array.
            ! If it, tSuccess will be true and ind will hold the position of the
            ! entry in rdm%elements.
            call hash_table_lookup((/i, j, k, l/), (/ijkl/), 0, rdm%hash_table, rdm%elements, ind, hash_val, tSuccess)

            if (tSuccess) then
                ! Extract the existing sign.
                call extract_sign_rdm(rdm%elements(:, ind), real_sign_old)
                ! Update the total sign.
                real_sign_new = real_sign_old + contrib_sign
                ! Encode the new sign.
                call encode_sign_rdm(rdm%elements(:, ind), real_sign_new)
            else
                ! Determine the process label.
                if (spinfree .or. tGUGA) then
                    ! For spin-free case, we halve the number of labels. Also,
                    ! the last two labels are dominant in the ordering, so use
                    ! these instead, to allow writing out in the correct order.
                    kl_compressed = int(contract_1_rdm_ind(k, l))
                    proc = (kl_compressed - 1) * nProcessors / (nbasis**2 / 4)
                else
                    ! The following maps (p,q), with p<q, to single integers
                    ! with no gaps. It is benefical to have no gaps here, for
                    ! good load balancing. The final integers are ordered so
                    ! that p is dominant over q.
                    ij_compressed = nbasis * (i - 1) - i * (i - 1) / 2 + j - i
                    ! Calculate the process for the element.
                    proc = (ij_compressed - 1) * nProcessors / spawn%nrows
                end if

                ! Check that there is enough memory for the new spawned RDM entry.
                slots_left = spawn%init_free_slots(proc + 1) - spawn%free_slots(proc)

                if (slots_left <= 0) call try_rdm_spawn_realloc(spawn, proc, spinfree)

                if (present(nearly_full)) then
                    ! 10 chosen somewhat arbitrarily, although there are times
                    ! when we call this routine 8 times in a row, so best not
                    ! to make any smaller.
                    if (slots_left <= 10) nearly_full = .true.
                end if

                rdm%elements(0, spawn%free_slots(proc)) = ijkl
                call encode_sign_rdm(rdm%elements(:, spawn%free_slots(proc)), contrib_sign)

                call add_hash_table_entry(rdm%hash_table, spawn%free_slots(proc), hash_val)

                spawn%free_slots(proc) = spawn%free_slots(proc) + 1
            end if

        end associate

    end subroutine add_to_rdm_spawn_t