return_haloCounts Subroutine

private subroutine return_haloCounts(sendbuffer, recvbuffer, comm)

Report the actually existing elements, which were requested as halos from remote

first report the number of actual halos

Arguments

Type IntentOptional Attributes Name
type(tem_communication_type), intent(in) :: sendbuffer

send buffer

type(tem_communication_type), intent(inout) :: recvbuffer

recv buffer

integer, intent(in) :: comm

MPI communicator


Calls

proc~~return_halocounts~~CallsGraph proc~return_halocounts return_haloCounts mpi_irecv mpi_irecv proc~return_halocounts->mpi_irecv mpi_isend mpi_isend proc~return_halocounts->mpi_isend mpi_waitall mpi_waitall proc~return_halocounts->mpi_waitall

Called by

proc~~return_halocounts~~CalledByGraph proc~return_halocounts return_haloCounts proc~communicate_elements communicate_elements proc~communicate_elements->proc~return_halocounts proc~tem_find_allelements tem_find_allElements proc~tem_find_allelements->proc~communicate_elements proc~tem_create_leveldesc tem_create_levelDesc proc~tem_create_leveldesc->proc~tem_find_allelements proc~tem_dimbydim_construction tem_dimByDim_construction proc~tem_dimbydim_construction->proc~tem_create_leveldesc

Source Code

  subroutine return_haloCounts( sendbuffer, recvbuffer, comm )
    ! ---------------------------------------------------------------------------
    !> send buffer
    type(tem_communication_type), intent(in) :: sendbuffer
    !> recv buffer
    type(tem_communication_type), intent(inout) :: recvbuffer
    !> MPI communicator
    integer, intent(in) :: comm
    ! ---------------------------------------------------------------------------
    integer :: iProc ! process iterator
    integer :: iErr  ! error flag
    integer, parameter :: message_tag = 23 ! arbitrary message flag
    integer, allocatable :: rq_handle(:)
    integer, allocatable :: status(:,:)
    integer :: nCommunications ! amount of communciations
    ! ---------------------------------------------------------------------------

    nCommunications = sendbuffer%nProcs + recvbuffer%nProcs
    allocate( rq_handle( nCommunications ) )
    allocate( status( mpi_status_size, nCommunications ) )
    ! ---------------------------------------------------------------------
    !! first report the number of actual halos
    rq_handle(:) = MPI_REQUEST_NULL

    do iProc = 1, recvbuffer%nProcs
      ! Receive the number of actual elements
      call mpi_irecv( recvbuffer%nElemsProc( iProc ),                          &
       &              1,                                                   &
       &              mpi_integer,                                             &
       &              recvbuffer%proc(iProc),                                  &
       &              message_tag,                                             &
       &              comm,                                               &
       &              rq_handle(iProc),                                        &
       &              iErr  )
    end do ! iProc

    ! I send the number of actual existing halo elements to the requester.
    do iProc = 1, sendbuffer%nProcs
      call mpi_isend( sendbuffer%nElemsProc( iProc ),                          &
       &              1,                                                   &
       &              mpi_integer,                                             &
       &              sendbuffer%proc( iProc ),                                &
       &              message_tag,                                             &
       &              comm,                                               &
       &              rq_handle(iProc + recvbuffer%nProcs),                    &
       &              iErr )
    end do ! iProc

    call mpi_waitall( nCommunications, rq_handle, status, iErr)
    ! Now we know the number of how actually existent cells, to receive
    ! allocate the recv buffers and communicate
    ! ---------------------------------------------------------------------

  end subroutine return_haloCounts