That shouldn&#39;t be the case; the request at index should get set to MPI_REQUEST_NULL and not hurt anything:<div><br><div>**** Test program; Only uses two nodes as written!****</div><div><div>#include &lt;iostream&gt;</div>
<div><br></div><div>#include &lt;mpi.h&gt;</div><div><br></div><div>int main(int argc, char * argv[])</div><div>{</div><div>  MPI_Init(&amp;argc, &amp;argv);</div><div>  int rank;</div><div>  MPI_Comm_rank(MPI_COMM_WORLD, &amp;rank);</div>
<div>  </div><div>#define BUFSIZE (0x1 &lt;&lt; 20) // 1M-element</div><div>#define COUNT 4</div><div>  if (rank == 1)</div><div>  {</div><div>    int * buffer = (int *) malloc(BUFSIZE * sizeof(int));</div><div>    for (int i=0; i &lt; COUNT; i++)</div>
<div>      MPI_Send(buffer, BUFSIZE, MPI_INT, 0, (i*3)%COUNT, MPI_COMM_WORLD);</div><div>  }</div><div>  else if (!rank)</div><div>  {</div><div>    int * buffers[COUNT];</div><div>    MPI_Request requests[COUNT];</div><div>
    for (int i=0; i &lt; COUNT; i++)</div><div>    {</div><div>      buffers[i] = (int *) malloc(BUFSIZE * sizeof(int));</div><div>      MPI_Irecv(buffers[i], BUFSIZE, MPI_INT, 1, i, MPI_COMM_WORLD, requests + i);</div><div>
    }</div><div>    int rec = 0;</div><div>    while (rec != MPI_UNDEFINED)</div><div>    {</div><div>      std::cout &lt;&lt; &quot;Waiting... &quot;;</div><div>      if(MPI_Waitany(COUNT, requests, &amp;rec, MPI_STATUS_IGNORE))</div>
<div>      {</div><div>        std::cout &lt;&lt; &quot;Error in wait?&quot; &lt;&lt; std::endl;</div><div>        break;</div><div>      }</div><div>      else if (rec != MPI_UNDEFINED)</div><div>        std::cout &lt;&lt; &quot;Received message &quot; &lt;&lt; rec;</div>
<div>      else </div><div>        std::cout &lt;&lt; &quot;None left!        &quot;;</div><div><br></div><div>      std::cout &lt;&lt; &quot; Requests: &quot;;</div><div>      for (int i=0; i &lt; COUNT; i++)</div><div>        if (requests[i] == MPI_REQUEST_NULL)</div>
<div>          std::cout &lt;&lt; &quot;N&quot;; // NULL</div><div>        else</div><div>          std::cout &lt;&lt; &quot;V&quot;; // VALID</div><div>      std::cout &lt;&lt; std::endl;</div><div>    }</div><div>  }</div>
<div><br></div><div>  MPI_Finalize();</div><div>  return 0;</div><div>}</div></div><div>*** Execution output</div><div><div><br></div><div><div>Waiting... Received message 0 Requests: NVVV</div><div>Waiting... Received message 3 Requests: NVVN</div>
<div>Waiting... Received message 2 Requests: NVNN</div><div>Waiting... Received message 1 Requests: NNNN</div><div>Waiting... None left!         Requests: NNNN</div></div><div><br></div><div> -Eric</div><br><div class="gmail_quote">
On Mon, Jan 3, 2011 at 10:19 PM, Xiao Li <span dir="ltr">&lt;<a href="mailto:shinelee.thewise@gmail.com">shinelee.thewise@gmail.com</a>&gt;</span> wrote:<br><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex;">
Hi Eric,<div><br></div><div>The if statement may make some error, I think it should be altered as .</div><div><br></div><div><div> if ++trunk_sent[index] != M do         </div><div>         MPI_Irecv(buffer[index], index, requests[index])        </div>

<div> else</div></div><div>        //remove the MPIRequest object of finished process, or else it might be halt forever</div><div>        remove requests[index]</div><div> end</div><div> </div><div>cheers</div><div>Xiao<br>

<br><div class="gmail_quote">On Mon, Jan 3, 2011 at 4:24 PM, Xiao Li <span dir="ltr">&lt;<a href="mailto:shinelee.thewise@gmail.com" target="_blank">shinelee.thewise@gmail.com</a>&gt;</span> wrote:<br><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex">

<div>Hi Eric,</div><div><br></div><div>You are right. An extra MPI_Irecv will be executed at end. Thanks for your comment.</div><div><br></div><div>cheers</div><div>Xiao<div><div></div><div><br><br><div class="gmail_quote">

On Mon, Jan 3, 2011 at 4:16 PM, Eric A. Borisch <span dir="ltr">&lt;<a href="mailto:eborisch@ieee.org" target="_blank">eborisch@ieee.org</a>&gt;</span> wrote:<br>
<blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex">Looks about right... I&#39;m assuming there is a &lt;do actual work here&gt; to be inserted between the MPI_Waitany and MPI_Irecv within the N*M-sized loop.... and I count from 0 rather than 1 by force of habit... :)<div>


<br>
</div><div>I think the logic below will end with one extra attempted MPI_Irecv than desired; perhaps change</div><div><div><br></div><div><div>if trunk_sent[index] != M do         </div><div>   MPI_Irecv(buffer[index], index, requests[index])</div>



<div>   trunk_sent[index]++</div><div>end</div></div></div><div><div><br></div><div>to</div><div><br></div><div><div><div><div>if ++trunk_sent[index] != M do         </div><div>   MPI_Irecv(buffer[index], index, requests[index])</div>



</div><div>end</div></div><div></div></div><div><br></div><div><font color="#888888"> -Eric</font><div><div></div><div><br><br><div class="gmail_quote">On Mon, Jan 3, 2011 at 3:00 PM, Xiao Li <span dir="ltr">&lt;<a href="mailto:shinelee.thewise@gmail.com" target="_blank">shinelee.thewise@gmail.com</a>&gt;</span> wrote:<br>



<blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex">Hi Eric,<div><br></div><div>Thanks for your detailed suggestion. After read MPI documents, I propose the following algorithm,</div>



<div><br></div><div>//begin collecting data for the first trunk</div><div>for i=1 to N do</div>
<div>    MPI_Irecv(buffer[i], i, requests[i])</div><div>end</div><div>//set data sending counter</div><div>for i=1 to N do</div><div>    trunk_sent[i] = 0</div><div>end</div><div>//begin collecting data</div><div>for i=1 to N*M do</div>




<div>    MPI_Waitany(N, requests, &amp;index, &amp;status)</div><div>    if trunk_sent[index] != M do         </div><div>         MPI_Irecv(buffer[index], index, requests[index])</div><div>         trunk_sent[index]++</div>




<div>    end</div><div>end</div><div><br></div><div>May I know what is your opinion of this algorithm?</div><div><br></div><div>cheers</div><div>Xiao<div><div></div><div><br><br><div class="gmail_quote">On Mon, Jan 3, 2011 at 3:31 PM, Eric A. Borisch <span dir="ltr">&lt;<a href="mailto:eborisch@ieee.org" target="_blank">eborisch@ieee.org</a>&gt;</span> wrote:<br>




<blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex">Xiao,<div><br></div><div>You should be able to get by with just N buffers, one for each client. After you have processed the i-th iteration for client n, re-issue an MPI_Irecv with the same buffer. This will match up with the next MPI_Send from client n. You don&#39;t have to worry about synchronizing -- the MPI_Irecv does not need to be posted before the MPI_Send. (But the MPI_Send won&#39;t complete until it has been, of course...)</div>





<div><br></div><div>You could always roll your own sockets, but MPI does a nice job of managing connections and messages for you. In addition, MPI can be used fairly efficiently on a wide range of interconnects, from shared memory to Infiniband with little to no change on the user&#39;s part.</div>





<div><div><div><br></div><div>In addition, you could likely improve performance in MPI by having two sets (call them A and B) of buffers to send from on each worker; one is in the &quot;send&quot; state (let&#39;s call this one A, started with an MPI_Isend after it was initially filled) while you&#39;re filling B. After B is filled, initiate a new MPI_Isend (very quick) on B and then wait for A&#39;s first send (MPI_Wait) to complete. Once the first send on A is completed, you can start populating A with the next iteration&#39;s output, initiate A&#39;s send, wait for B&#39;s send to complete, and the cycle begins again.</div>





<div><br></div><div>This approach allows you to overlay communication and computation times, and still works with the MPI_Waitany() approach to harvesting completed jobs in first-completed order on the master. This is an almost trivial thing to implement in MPI, but achieving it with sockets requires (IMHO) much more programmer overhead...</div>





<div><div><br></div></div></div><div>Just my 2c.</div><div><br></div><div><font color="#888888"> Eric</font><div><div></div><div><br><br><div class="gmail_quote">On Mon, Jan 3, 2011 at 1:24 PM, Xiao Li <span dir="ltr">&lt;<a href="mailto:shinelee.thewise@gmail.com" target="_blank">shinelee.thewise@gmail.com</a>&gt;</span> wrote:<br>





<blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex">Hi Eric,<div><br></div><div>Assume I have N workers and M trunks of sending data for each worker respectively, then I have to create N*M data buffer for MPI_Irecv usage. Is this method too costly?</div>





<div><br></div><div>
Or If I write raw socket programming, is that better? Just like traditional client/server socket programming? Master listens on port and spawn a new thread to accept worker&#39;s data storage request?</div><div><br></div>






<div>cheers</div><div>Xiao<div><div></div><div><br><br><div class="gmail_quote">On Mon, Jan 3, 2011 at 2:13 PM, Eric A. Borisch <span dir="ltr">&lt;<a href="mailto:eborisch@ieee.org" target="_blank">eborisch@ieee.org</a>&gt;</span> wrote:<br>





<blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex">
Look at the documentation for MPI_Irecv and MPI_Testany ... these should help you do what you want.<div><br></div><div> Eric<br><div class="gmail_quote"><br></div><div class="gmail_quote"><div><div></div><div>On Mon, Jan 3, 2011 at 12:45 PM, Xiao Li <span dir="ltr">&lt;<a href="mailto:shinelee.thewise@gmail.com" target="_blank">shinelee.thewise@gmail.com</a>&gt;</span> wrote:<br>







</div></div><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><div><div></div><div>Hi MPICH2 people,<div><br></div><div>Now, I have a application that composed of single master and many workers. The application requirement is very simple: workers finish some jobs and send data to master and master store these data into files separately. I can simply use MPI_Send on worker side to send data to master. But master does not know the data sending sequence. Some worker go fast while some are slow. More specifically, suppose there are 5 workers, then the data send sequence may be 1,3,4,5,2 or 2,5,4,1,3. If I just write a for loop for(i=1 to 5) on master side with MPI_Recv to get data, the master and some faster worker have to wait for a long time. I know MPI_Gather can implement this. But I am not sure is MPI_Gather works parallelly or just a sequential MPI_Recv? Another issue is my data is extremely large, more than 1GB data needed to be sent to master. If I divide the data into pieces, I do not think MPI_Gather can work. I also tried to think about raw socket programming, but I do not think it is a good practice. Would you give me some suggestion please?</div>








<div><br></div><div>cheers</div><div>Xiao</div>
<br></div></div>_______________________________________________<br>
mpich-discuss mailing list<br>
<a href="mailto:mpich-discuss@mcs.anl.gov" target="_blank">mpich-discuss@mcs.anl.gov</a><br>
<a href="https://lists.mcs.anl.gov/mailman/listinfo/mpich-discuss" target="_blank">https://lists.mcs.anl.gov/mailman/listinfo/mpich-discuss</a><br><br></blockquote></div>
</div>
<br>_______________________________________________<br>
mpich-discuss mailing list<br>
<a href="mailto:mpich-discuss@mcs.anl.gov" target="_blank">mpich-discuss@mcs.anl.gov</a><br>
<a href="https://lists.mcs.anl.gov/mailman/listinfo/mpich-discuss" target="_blank">https://lists.mcs.anl.gov/mailman/listinfo/mpich-discuss</a><br>
<br></blockquote></div><br></div></div></div>
<br>_______________________________________________<br>
mpich-discuss mailing list<br>
<a href="mailto:mpich-discuss@mcs.anl.gov" target="_blank">mpich-discuss@mcs.anl.gov</a><br>
<a href="https://lists.mcs.anl.gov/mailman/listinfo/mpich-discuss" target="_blank">https://lists.mcs.anl.gov/mailman/listinfo/mpich-discuss</a><br><br></blockquote></div>
</div></div></div></div>
<br>_______________________________________________<br>
mpich-discuss mailing list<br>
<a href="mailto:mpich-discuss@mcs.anl.gov" target="_blank">mpich-discuss@mcs.anl.gov</a><br>
<a href="https://lists.mcs.anl.gov/mailman/listinfo/mpich-discuss" target="_blank">https://lists.mcs.anl.gov/mailman/listinfo/mpich-discuss</a><br>
<br></blockquote></div><br></div></div></div>
<br>_______________________________________________<br>
mpich-discuss mailing list<br>
<a href="mailto:mpich-discuss@mcs.anl.gov" target="_blank">mpich-discuss@mcs.anl.gov</a><br>
<a href="https://lists.mcs.anl.gov/mailman/listinfo/mpich-discuss" target="_blank">https://lists.mcs.anl.gov/mailman/listinfo/mpich-discuss</a><br>
<br></blockquote></div><br>
</div></div></div></div>
<br>_______________________________________________<br>
mpich-discuss mailing list<br>
<a href="mailto:mpich-discuss@mcs.anl.gov" target="_blank">mpich-discuss@mcs.anl.gov</a><br>
<a href="https://lists.mcs.anl.gov/mailman/listinfo/mpich-discuss" target="_blank">https://lists.mcs.anl.gov/mailman/listinfo/mpich-discuss</a><br>
<br></blockquote></div><br></div></div></div>
</blockquote></div><br></div>
<br>_______________________________________________<br>
mpich-discuss mailing list<br>
<a href="mailto:mpich-discuss@mcs.anl.gov">mpich-discuss@mcs.anl.gov</a><br>
<a href="https://lists.mcs.anl.gov/mailman/listinfo/mpich-discuss" target="_blank">https://lists.mcs.anl.gov/mailman/listinfo/mpich-discuss</a><br>
<br></blockquote></div><br>
</div></div>