exchange the communication buffers with a non-blocking mpi communication using preposted irecv and isend with a waitall
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(tem_communication_type), | intent(inout) | :: | send | |||
type(tem_communication_type), | intent(inout) | :: | recv | |||
real(kind=rk), | intent(inout) | :: | state(*) | |||
integer, | intent(in) | :: | message_flag | |||
real(kind=rk), | intent(in), | optional | :: | send_state(*) | ||
integer, | intent(in) | :: | comm |
mpi communicator |
subroutine comm_isend_irecv_real( send, recv, state, message_flag, & & send_state, comm ) ! -------------------------------------------------------------------- ! type( tem_communication_type ), intent(inout) :: send type( tem_communication_type ), intent(inout) :: recv real(kind=rk), intent(inout) :: state(*) !< state vector to update integer, intent(in) :: message_flag real(kind=rk), intent(in), optional :: send_state(*) !< data to send !> mpi communicator integer, intent(in) :: comm ! -------------------------------------------------------------------- ! ! request handle for messages ! @todo request handle array could exist during complete code runtime ! integer :: rq_handle( recv%nprocs + send%nprocs ) integer :: status( mpi_status_size, max(recv%nprocs, send%nprocs) ) integer :: ierr !< error flag integer :: iproc, ival integer :: nsendvals, nrecvvals ! -------------------------------------------------------------------- ! if (present(send_state)) then do iproc = 1, send%nprocs ! fill communication message nsendvals = send%buf_real( iproc )%nvals !$nec ivdep do ival = 1, nsendvals send%buf_real( iproc )%val( ival ) & & = send_state( send%buf_real( iproc )%pos( ival ) ) end do end do else do iproc = 1, send%nprocs ! fill communication message nsendvals = send%buf_real( iproc )%nvals !$nec ivdep do ival = 1, nsendvals send%buf_real( iproc )%val( ival ) & & = state( send%buf_real( iproc )%pos( ival ) ) end do end do end if do iproc = 1, recv%nprocs ! start receive communications call mpi_irecv( & & recv%buf_real( iproc )%val, & ! me & recv%buf_real( iproc )%nvals, & ! me size & rk_mpi, & ! data type & recv%proc(iproc), & ! target me & message_flag, & ! flag & comm, & ! communicator & recv%rqhandle(iproc), & ! handle & ierr ) ! error status enddo ! start the sending communications do iproc = 1, send%nprocs call mpi_isend( & & send%buf_real( iproc )%val, & ! buffer & send%buf_real( iproc )%nvals, & ! count & rk_mpi, & ! data type & send%proc(iproc), & ! target & message_flag, & ! tag & comm, & ! communicator & send%rqhandle( iproc ), & ! handle & ierr ) ! error status enddo ! iproc ! wait for receive buffer to be ready if ( recv%nprocs /= 0 ) then call mpi_waitall(recv%nprocs, & ! count & recv%rqhandle, & ! request handles & status, & ! mpi status & ierr ) ! error status end if ! now values from recv me can be copied to the actual state array do iproc = 1, recv%nprocs nrecvvals = recv%buf_real( iproc )%nvals !$nec ivdep do ival = 1, nrecvvals ! write the values from the recv me to the state array ! to positions given in recvme%pos state( recv%buf_real( iproc )%pos( ival ) ) & & = recv%buf_real( iproc )%val( ival ) end do end do ! wait for send buffer to be ready if ( send%nprocs /= 0 ) then call mpi_waitall(send%nprocs, & ! count & send%rqhandle, & ! request handles & status, & ! mpi status & ierr ) ! error status end if end subroutine comm_isend_irecv_real