comms_gatherv_real_3 Subroutine

private subroutine comms_gatherv_real_3(array, localcount, rootglobalarray, counts, displs)

Gather real data to root node (for arrays of rank 3)

Arguments

Type IntentOptional AttributesName
real(kind=dp), intent(inout), dimension(:, :, :):: array

local array for sending data

integer, intent(in) :: localcount

localcount elements will be sent to the root node

real(kind=dp), intent(inout), dimension(:, :, :):: rootglobalarray

array on the root node to which data will be sent

integer, intent(in), dimension(num_nodes):: counts

how data should be partitioned, see MPI documentation or function comms_array_split

integer, intent(in), dimension(num_nodes):: displs

Calls

proc~~comms_gatherv_real_3~~CallsGraph proc~comms_gatherv_real_3 comms_gatherv_real_3 dcopy dcopy proc~comms_gatherv_real_3->dcopy

Called by

proc~~comms_gatherv_real_3~~CalledByGraph proc~comms_gatherv_real_3 comms_gatherv_real_3 interface~comms_gatherv comms_gatherv interface~comms_gatherv->proc~comms_gatherv_real_3 proc~k_path k_path proc~k_path->interface~comms_gatherv proc~k_slice k_slice proc~k_slice->interface~comms_gatherv

Contents

Source Code


Source Code

  subroutine comms_gatherv_real_3(array, localcount, rootglobalarray, counts, displs)
    !! Gather real data to root node (for arrays of rank 3)
    implicit none

    real(kind=dp), dimension(:, :, :), intent(inout) :: array
    !! local array for sending data
    integer, intent(in)                            :: localcount
    !! localcount elements will be sent to the root node
    real(kind=dp), dimension(:, :, :), intent(inout) :: rootglobalarray
    !! array on the root node to which data will be sent
    integer, dimension(num_nodes), intent(in)      :: counts
    !! how data should be partitioned, see MPI documentation or
    !! function comms_array_split
    integer, dimension(num_nodes), intent(in)      :: displs

#ifdef MPI
    integer :: error

    call MPI_gatherv(array, localcount, MPI_double_precision, rootglobalarray, counts, &
                     displs, MPI_double_precision, root_id, mpi_comm_world, error)

    if (error .ne. MPI_success) then
      call io_error('Error in comms_gatherv_real_3')
    end if

#else
    call dcopy(localcount, array, 1, rootglobalarray, 1)
#endif

    return

  end subroutine comms_gatherv_real_3