send.c 4.7 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250
  1. #include "mpiP.h"
  2. /*
  3. * SENDING
  4. *
  5. */
  6. static int mpi_match_recv(void *r, void *tag)
  7. {
  8. return( ((Req *)r)->tag == MPI_ANY_TAG ||
  9. ((Req *)r)->tag == *((int *)tag) );
  10. }
  11. /*
  12. *
  13. */
  14. FC_FUNC( mpi_isend , MPI_ISEND )(void *buf, int *count, int *datatype,
  15. int *dest, int *tag, int *comm, int *req, int *ierror)
  16. {
  17. *ierror=MPI_Isend(buf,*count,*datatype,*dest,*tag,
  18. *comm, (void *)req);
  19. }
  20. int MPI_Isend(void *buf, int count, MPI_Datatype datatype,
  21. int dest, int tag, MPI_Comm comm, MPI_Request *request)
  22. {
  23. pListitem match;
  24. Comm *mycomm;
  25. Req *rreq, *sreq;
  26. mycomm=mpi_handle_to_ptr(comm); /* (Comm *)comm; */
  27. #ifdef INFO
  28. fflush(stdout);
  29. fprintf(stderr,"MPI_Isend: Comm=%d tag=%d count=%d type=%d\n",
  30. mycomm->num,tag,count,datatype);
  31. #endif
  32. if (dest!=0 && dest!=MPI_PROC_NULL)
  33. {
  34. fprintf(stderr,"MPI_Isend: send to %d\n",dest);
  35. abort();
  36. }
  37. mpi_alloc_handle(request,(void **) &sreq);
  38. if (dest==MPI_PROC_NULL)
  39. {
  40. sreq->complete=1;
  41. return(MPI_SUCCESS);
  42. }
  43. if ( match=AP_list_search_func(mycomm->recvlist,mpi_match_recv,&tag) )
  44. {
  45. rreq=(Req *)AP_listitem_data(match);
  46. AP_list_delete_item(mycomm->recvlist,match);
  47. memcpy(rreq->buf,buf,count * datatype);
  48. rreq->complete=1;
  49. rreq->source=0;
  50. rreq->tag=tag; /* in case rreq->tag was MPI_ANY_TAG */
  51. sreq->complete=1;
  52. #ifdef DEBUG
  53. printf("Completion(send) value=%d tag=%d\n",
  54. *((int *)buf),rreq->tag);
  55. #endif
  56. return(MPI_SUCCESS);
  57. }
  58. sreq->buf=buf;
  59. sreq->tag=tag;
  60. sreq->complete=0;
  61. sreq->listitem=AP_list_append(mycomm->sendlist,sreq);
  62. #ifdef INFO
  63. print_list(mycomm->sendlist,"sendlist for comm ",mycomm->num);
  64. #endif
  65. return(MPI_SUCCESS);
  66. }
  67. /*********/
  68. FC_FUNC(mpi_send, MPI_SEND) ( void *buf, int *count, int *datatype,
  69. int *dest, int *tag, int *comm, int *ierror)
  70. {
  71. *ierror=MPI_Send(buf, *count, *datatype, *dest, *tag, *comm);
  72. }
  73. int MPI_Send(void* buf, int count, MPI_Datatype datatype,
  74. int dest, int tag, MPI_Comm comm)
  75. {
  76. MPI_Request request;
  77. MPI_Status status;
  78. #ifdef INFO
  79. fflush(stdout);
  80. fprintf(stderr,"MPI_Send: ");
  81. #endif
  82. MPI_Isend(buf,count,datatype,dest,tag,comm,&request);
  83. MPI_Wait(&request,&status);
  84. return(MPI_SUCCESS);
  85. }
  86. /*********/
  87. FC_FUNC(mpi_ssend, MPI_SSEND) ( void *buf, int *count, int *datatype,
  88. int *dest, int *tag, int *comm, int *ierror)
  89. {
  90. *ierror=MPI_Send(buf, *count, *datatype, *dest, *tag, *comm);
  91. }
  92. int MPI_Ssend(void* buf, int count, MPI_Datatype datatype,
  93. int dest, int tag, MPI_Comm comm)
  94. {
  95. return(MPI_Send(buf,count,datatype,dest,tag,comm));
  96. }
  97. /*********/
  98. FC_FUNC(mpi_rsend, MPI_RSEND) ( void *buf, int *count, int *datatype,
  99. int *dest, int *tag, int *comm, int *ierror)
  100. {
  101. *ierror=MPI_Send(buf, *count, *datatype, *dest, *tag, *comm);
  102. }
  103. int MPI_Rsend(void* buf, int count, MPI_Datatype datatype,
  104. int dest, int tag, MPI_Comm comm)
  105. {
  106. return(MPI_Send(buf,count,datatype,dest,tag,comm));
  107. }
  108. /*********/
  109. FC_FUNC( mpi_irsend , MPI_IRSEND )(void *buf, int *count, int *datatype,
  110. int *dest, int *tag, int *comm, int *req, int *ierror)
  111. {
  112. *ierror=MPI_Irsend(buf,*count,*datatype,*dest,*tag,
  113. *comm, (void *)req);
  114. }
  115. int MPI_Irsend(void *buf, int count, MPI_Datatype datatype,
  116. int dest, int tag, MPI_Comm comm, MPI_Request *request)
  117. {
  118. MPI_Status status;
  119. Req *req;
  120. MPI_Isend(buf,count,datatype,dest,tag,comm,request);
  121. /* Ready mode implied a receive must already be posted,
  122. * so the Isend should have completed already.
  123. * Can't use MPI_Test here for the error check because
  124. * it would clear the request prematurely.
  125. */
  126. req=mpi_handle_to_ptr(*request);
  127. if ( !req->complete )
  128. {
  129. fprintf(stderr,"MPI_Irsend: no matching receive found\n");
  130. abort();
  131. }
  132. return(MPI_SUCCESS);
  133. }
  134. /*********/
  135. FC_FUNC(mpi_sendrecv, MPI_SENDRECV) (
  136. void *sendbuf, int *sendcount, int *sendtype, int *dest, int *sendtag,
  137. void *recvbuf, int *recvcount, int *recvtype, int *source, int *recvtag,
  138. int *comm, int *status,
  139. int *ierror)
  140. {
  141. *ierror=MPI_Sendrecv(sendbuf, *sendcount, *sendtype, *dest, *sendtag,
  142. recvbuf, *recvcount, *recvtype, *source, *recvtag,
  143. *comm, (MPI_Status *)status);
  144. }
  145. int MPI_Sendrecv(void* sendbuf, int sendcount, MPI_Datatype sendtype,
  146. int dest, int sendtag,
  147. void *recvbuf, int recvcount, MPI_Datatype recvtype,
  148. int source, int recvtag,
  149. MPI_Comm comm, MPI_Status *status)
  150. {
  151. MPI_Request request;
  152. MPI_Irecv(recvbuf, recvcount, recvtype, source, recvtag, comm, &request);
  153. MPI_Send(sendbuf, sendcount, sendtype, dest, sendtag, comm);
  154. MPI_Wait(&request,status);
  155. return(MPI_SUCCESS);
  156. }