communicate_rdm_spawn_t Subroutine

public subroutine communicate_rdm_spawn_t(spawn, rdm_recv)

Arguments

Type IntentOptional Attributes Name
type(rdm_spawn_t), intent(inout) :: spawn
type(rdm_list_t), intent(inout) :: rdm_recv

Contents


Source Code

    subroutine communicate_rdm_spawn_t(spawn, rdm_recv)

        ! Perform communication of RDM elements in the spawn object, to the
        ! rdm_recv object. The hash table and free_slots array for the spawn
        ! object are then reset at the end of this routine. The newly
        ! received spawnings will be added to rdm_recv *without* overwriting
        ! elements currently in the list (which is useful for situations where
        ! multiple communications are required, due to limited space in the
        ! spawning array).

        use hash, only: clear_hash_table
        use Parallel_neci, only: MPIAlltoAll, MPIAlltoAllv

        type(rdm_spawn_t), intent(inout) :: spawn
        type(rdm_list_t), intent(inout) :: rdm_recv

        integer :: iproc, nelements_old, new_nelements, ierr, i
        integer(MPIArg) :: send_sizes(0:nProcessors - 1), recv_sizes(0:nProcessors - 1)
        integer(MPIArg) :: send_displs(0:nProcessors - 1), recv_displs(0:nProcessors - 1)

        nelements_old = rdm_recv%nelements

        ! How many rows of data to send to each processor.
        do iproc = 0, nProcessors - 1
            send_sizes(iproc) = int(spawn%free_slots(iproc) - spawn%init_free_slots(iproc), MPIArg)
            ! The displacement of the beginning of each processor's section of the
            ! free_slots array, relative to the first element of this array.
            send_displs(iproc) = int(spawn%init_free_slots(iproc) - 1, MPIArg)
        end do

        ! this does not work with some compilers
        !send_displs = int(spawn%init_free_slots(0:nProcessors-1) - 1, MPIArg)

        call MPIAlltoAll(send_sizes, 1, recv_sizes, 1, ierr)

        recv_displs(0) = 0
        do iproc = 1, nProcessors - 1
            recv_displs(iproc) = recv_displs(iproc - 1) + recv_sizes(iproc - 1)
        end do

        ! The total number of RDM elements in the list after the receive.
        new_nelements = rdm_recv%nelements + int(sum(recv_sizes))

        ! If we don't have enough memory in the receiving list, try
        ! reallocating it to be big enough.
        if (new_nelements > rdm_recv%max_nelements) then
            call try_rdm_list_realloc(rdm_recv, new_nelements, .true.)
        end if

        ! Update the number of valid RDM elements in the received list.
        rdm_recv%nelements = new_nelements

        send_sizes = send_sizes * size(spawn%rdm_send%elements, 1)
        recv_sizes = recv_sizes * size(spawn%rdm_send%elements, 1)
        send_displs = send_displs * size(spawn%rdm_send%elements, 1)
        recv_displs = recv_displs * size(spawn%rdm_send%elements, 1)

        rdm_recv%elements(:, nelements_old + 1:) = 0
        ! Perform the communication.
        call MPIAlltoAllv(spawn%rdm_send%elements, send_sizes, send_displs, &
                          rdm_recv%elements(:, nelements_old + 1:), &
                          recv_sizes, recv_displs, ierr)

        ! Now we can reset the free_slots array and reset the hash table.
        do iproc = 0, nProcessors - 1
            spawn%free_slots(iproc) = spawn%init_free_slots(iproc)
        end do
        call clear_hash_table(spawn%rdm_send%hash_table)

    end subroutine communicate_rdm_spawn_t