xref: /petsc/src/sys/utils/mpimesg.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
1e5c89e4eSSatish Balay 
2c6db04a5SJed Brown #include <petscsys.h>        /*I  "petscsys.h"  I*/
376eed172SJunchao Zhang #include <petsc/private/mpiutils.h>
4e5c89e4eSSatish Balay 
5e5c89e4eSSatish Balay /*@C
6e5c89e4eSSatish Balay   PetscGatherNumberOfMessages -  Computes the number of messages a node expects to receive
7e5c89e4eSSatish Balay 
8d083f849SBarry Smith   Collective
9e5c89e4eSSatish Balay 
10e5c89e4eSSatish Balay   Input Parameters:
11e5c89e4eSSatish Balay + comm     - Communicator
12e5c89e4eSSatish Balay . iflags   - an array of integers of length sizeof(comm). A '1' in ilengths[i] represent a
130298fd71SBarry Smith              message from current node to ith node. Optionally NULL
14e5c89e4eSSatish Balay - ilengths - Non zero ilengths[i] represent a message to i of length ilengths[i].
150298fd71SBarry Smith              Optionally NULL.
16e5c89e4eSSatish Balay 
17e5c89e4eSSatish Balay   Output Parameters:
18e5c89e4eSSatish Balay . nrecvs    - number of messages received
19e5c89e4eSSatish Balay 
20e5c89e4eSSatish Balay   Level: developer
21e5c89e4eSSatish Balay 
22e5c89e4eSSatish Balay   Notes:
23e5c89e4eSSatish Balay   With this info, the correct message lengths can be determined using
24e5c89e4eSSatish Balay   PetscGatherMessageLengths()
25e5c89e4eSSatish Balay 
26e5c89e4eSSatish Balay   Either iflags or ilengths should be provided.  If iflags is not
270298fd71SBarry Smith   provided (NULL) it can be computed from ilengths. If iflags is
28e5c89e4eSSatish Balay   provided, ilengths is not required.
29e5c89e4eSSatish Balay 
30e5c89e4eSSatish Balay .seealso: PetscGatherMessageLengths()
31e5c89e4eSSatish Balay @*/
327087cfbeSBarry Smith PetscErrorCode  PetscGatherNumberOfMessages(MPI_Comm comm,const PetscMPIInt iflags[],const PetscMPIInt ilengths[],PetscMPIInt *nrecvs)
33e5c89e4eSSatish Balay {
340298fd71SBarry Smith   PetscMPIInt    size,rank,*recv_buf,i,*iflags_local = NULL,*iflags_localm = NULL;
35e5c89e4eSSatish Balay 
36e5c89e4eSSatish Balay   PetscFunctionBegin;
375f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(comm,&size));
385f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(comm,&rank));
39e5c89e4eSSatish Balay 
405f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc2(size,&recv_buf,size,&iflags_localm));
41e5c89e4eSSatish Balay 
42e5c89e4eSSatish Balay   /* If iflags not provided, compute iflags from ilengths */
43e5c89e4eSSatish Balay   if (!iflags) {
44*28b400f6SJacob Faibussowitsch     PetscCheck(ilengths,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Either iflags or ilengths should be provided");
45e5c89e4eSSatish Balay     iflags_local = iflags_localm;
46e5c89e4eSSatish Balay     for (i=0; i<size; i++) {
47e5c89e4eSSatish Balay       if (ilengths[i]) iflags_local[i] = 1;
48e5c89e4eSSatish Balay       else iflags_local[i] = 0;
49e5c89e4eSSatish Balay     }
50a297a907SKarl Rupp   } else iflags_local = (PetscMPIInt*) iflags;
51e5c89e4eSSatish Balay 
52e5c89e4eSSatish Balay   /* Post an allreduce to determine the numer of messages the current node will receive */
535f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPIU_Allreduce(iflags_local,recv_buf,size,MPI_INT,MPI_SUM,comm));
54e5c89e4eSSatish Balay   *nrecvs = recv_buf[rank];
55e5c89e4eSSatish Balay 
565f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2(recv_buf,iflags_localm));
57e5c89e4eSSatish Balay   PetscFunctionReturn(0);
58e5c89e4eSSatish Balay }
59e5c89e4eSSatish Balay 
60e5c89e4eSSatish Balay /*@C
61e5c89e4eSSatish Balay   PetscGatherMessageLengths - Computes info about messages that a MPI-node will receive,
62e5c89e4eSSatish Balay   including (from-id,length) pairs for each message.
63e5c89e4eSSatish Balay 
64d083f849SBarry Smith   Collective
65e5c89e4eSSatish Balay 
66e5c89e4eSSatish Balay   Input Parameters:
67e5c89e4eSSatish Balay + comm      - Communicator
68e5c89e4eSSatish Balay . nsends    - number of messages that are to be sent.
69e5c89e4eSSatish Balay . nrecvs    - number of messages being received
70e5c89e4eSSatish Balay - ilengths  - an array of integers of length sizeof(comm)
71e5c89e4eSSatish Balay               a non zero ilengths[i] represent a message to i of length ilengths[i]
72e5c89e4eSSatish Balay 
73e5c89e4eSSatish Balay   Output Parameters:
74e5c89e4eSSatish Balay + onodes    - list of node-ids from which messages are expected
75e5c89e4eSSatish Balay - olengths  - corresponding message lengths
76e5c89e4eSSatish Balay 
77e5c89e4eSSatish Balay   Level: developer
78e5c89e4eSSatish Balay 
79e5c89e4eSSatish Balay   Notes:
80e5c89e4eSSatish Balay   With this info, the correct MPI_Irecv() can be posted with the correct
81e5c89e4eSSatish Balay   from-id, with a buffer with the right amount of memory required.
82e5c89e4eSSatish Balay 
83e5c89e4eSSatish Balay   The calling function deallocates the memory in onodes and olengths
84e5c89e4eSSatish Balay 
85c2916339SPierre Jolivet   To determine nrecvs, one can use PetscGatherNumberOfMessages()
86e5c89e4eSSatish Balay 
87e5c89e4eSSatish Balay .seealso: PetscGatherNumberOfMessages()
88e5c89e4eSSatish Balay @*/
897087cfbeSBarry Smith PetscErrorCode  PetscGatherMessageLengths(MPI_Comm comm,PetscMPIInt nsends,PetscMPIInt nrecvs,const PetscMPIInt ilengths[],PetscMPIInt **onodes,PetscMPIInt **olengths)
90e5c89e4eSSatish Balay {
916bfd7d4fSJunchao Zhang   PetscMPIInt    size,rank,tag,i,j;
920298fd71SBarry Smith   MPI_Request    *s_waits  = NULL,*r_waits = NULL;
930298fd71SBarry Smith   MPI_Status     *w_status = NULL;
94e5c89e4eSSatish Balay 
95e5c89e4eSSatish Balay   PetscFunctionBegin;
965f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(comm,&size));
975f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(comm,&rank));
985f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCommGetNewTag(comm,&tag));
99e5c89e4eSSatish Balay 
100e5c89e4eSSatish Balay   /* cannot use PetscMalloc3() here because in the call to MPI_Waitall() they MUST be contiguous */
1015f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc2(nrecvs+nsends,&r_waits,nrecvs+nsends,&w_status));
102e5c89e4eSSatish Balay   s_waits = r_waits+nrecvs;
103e5c89e4eSSatish Balay 
104e5c89e4eSSatish Balay   /* Post the Irecv to get the message length-info */
1055f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nrecvs,olengths));
106e5c89e4eSSatish Balay   for (i=0; i<nrecvs; i++) {
1075f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Irecv((*olengths)+i,1,MPI_INT,MPI_ANY_SOURCE,tag,comm,r_waits+i));
108e5c89e4eSSatish Balay   }
109e5c89e4eSSatish Balay 
110e5c89e4eSSatish Balay   /* Post the Isends with the message length-info */
111e5c89e4eSSatish Balay   for (i=0,j=0; i<size; ++i) {
112e5c89e4eSSatish Balay     if (ilengths[i]) {
1135f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Isend((void*)(ilengths+i),1,MPI_INT,i,tag,comm,s_waits+j));
114e5c89e4eSSatish Balay       j++;
115e5c89e4eSSatish Balay     }
116e5c89e4eSSatish Balay   }
117e5c89e4eSSatish Balay 
118e5c89e4eSSatish Balay   /* Post waits on sends and receivs */
1195f80ce2aSJacob Faibussowitsch   if (nrecvs+nsends) CHKERRMPI(MPI_Waitall(nrecvs+nsends,r_waits,w_status));
120e5c89e4eSSatish Balay 
121e5c89e4eSSatish Balay   /* Pack up the received data */
1225f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nrecvs,onodes));
1236bfd7d4fSJunchao Zhang   for (i=0; i<nrecvs; ++i) {
1246bfd7d4fSJunchao Zhang     (*onodes)[i] = w_status[i].MPI_SOURCE;
1256bfd7d4fSJunchao Zhang #if defined(PETSC_HAVE_OMPI_MAJOR_VERSION)
1266bfd7d4fSJunchao Zhang     /* This line is a workaround for a bug in OpenMPI-2.1.1 distributed by Ubuntu-18.04.2 LTS.
1276bfd7d4fSJunchao Zhang        It happens in self-to-self MPI_Send/Recv using MPI_ANY_SOURCE for message matching. OpenMPI
1286bfd7d4fSJunchao Zhang        does not put correct value in recv buffer. See also
1296bfd7d4fSJunchao Zhang        https://lists.mcs.anl.gov/pipermail/petsc-dev/2019-July/024803.html
1306bfd7d4fSJunchao Zhang        https://www.mail-archive.com/users@lists.open-mpi.org//msg33383.html
1316bfd7d4fSJunchao Zhang      */
1326bfd7d4fSJunchao Zhang     if (w_status[i].MPI_SOURCE == rank) (*olengths)[i] = ilengths[rank];
1336bfd7d4fSJunchao Zhang #endif
1346bfd7d4fSJunchao Zhang   }
1355f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2(r_waits,w_status));
136e5c89e4eSSatish Balay   PetscFunctionReturn(0);
137e5c89e4eSSatish Balay }
138dd6ea824SBarry Smith 
13976eed172SJunchao Zhang /* Same as PetscGatherNumberOfMessages(), except using PetscInt for ilengths[] */
14076eed172SJunchao Zhang PetscErrorCode  PetscGatherNumberOfMessages_Private(MPI_Comm comm,const PetscMPIInt iflags[],const PetscInt ilengths[],PetscMPIInt *nrecvs)
14176eed172SJunchao Zhang {
14276eed172SJunchao Zhang   PetscMPIInt    size,rank,*recv_buf,i,*iflags_local = NULL,*iflags_localm = NULL;
14376eed172SJunchao Zhang 
14476eed172SJunchao Zhang   PetscFunctionBegin;
1455f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(comm,&size));
1465f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(comm,&rank));
14776eed172SJunchao Zhang 
1485f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc2(size,&recv_buf,size,&iflags_localm));
14976eed172SJunchao Zhang 
15076eed172SJunchao Zhang   /* If iflags not provided, compute iflags from ilengths */
15176eed172SJunchao Zhang   if (!iflags) {
152*28b400f6SJacob Faibussowitsch     PetscCheck(ilengths,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Either iflags or ilengths should be provided");
15376eed172SJunchao Zhang     iflags_local = iflags_localm;
15476eed172SJunchao Zhang     for (i=0; i<size; i++) {
15576eed172SJunchao Zhang       if (ilengths[i]) iflags_local[i] = 1;
15676eed172SJunchao Zhang       else iflags_local[i] = 0;
15776eed172SJunchao Zhang     }
15876eed172SJunchao Zhang   } else iflags_local = (PetscMPIInt*) iflags;
15976eed172SJunchao Zhang 
16076eed172SJunchao Zhang   /* Post an allreduce to determine the numer of messages the current node will receive */
1615f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPIU_Allreduce(iflags_local,recv_buf,size,MPI_INT,MPI_SUM,comm));
16276eed172SJunchao Zhang   *nrecvs = recv_buf[rank];
16376eed172SJunchao Zhang 
1645f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2(recv_buf,iflags_localm));
16576eed172SJunchao Zhang   PetscFunctionReturn(0);
16676eed172SJunchao Zhang }
16776eed172SJunchao Zhang 
16876eed172SJunchao Zhang /* Same as PetscGatherMessageLengths(), except using PetscInt for message lengths */
16976eed172SJunchao Zhang PetscErrorCode  PetscGatherMessageLengths_Private(MPI_Comm comm,PetscMPIInt nsends,PetscMPIInt nrecvs,const PetscInt ilengths[],PetscMPIInt **onodes,PetscInt **olengths)
17076eed172SJunchao Zhang {
17176eed172SJunchao Zhang   PetscMPIInt    size,rank,tag,i,j;
17276eed172SJunchao Zhang   MPI_Request    *s_waits  = NULL,*r_waits = NULL;
17376eed172SJunchao Zhang   MPI_Status     *w_status = NULL;
17476eed172SJunchao Zhang 
17576eed172SJunchao Zhang   PetscFunctionBegin;
1765f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(comm,&size));
1775f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(comm,&rank));
1785f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCommGetNewTag(comm,&tag));
17976eed172SJunchao Zhang 
18076eed172SJunchao Zhang   /* cannot use PetscMalloc3() here because in the call to MPI_Waitall() they MUST be contiguous */
1815f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc2(nrecvs+nsends,&r_waits,nrecvs+nsends,&w_status));
18276eed172SJunchao Zhang   s_waits = r_waits+nrecvs;
18376eed172SJunchao Zhang 
18476eed172SJunchao Zhang   /* Post the Irecv to get the message length-info */
1855f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nrecvs,olengths));
18676eed172SJunchao Zhang   for (i=0; i<nrecvs; i++) {
1875f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Irecv((*olengths)+i,1,MPIU_INT,MPI_ANY_SOURCE,tag,comm,r_waits+i));
18876eed172SJunchao Zhang   }
18976eed172SJunchao Zhang 
19076eed172SJunchao Zhang   /* Post the Isends with the message length-info */
19176eed172SJunchao Zhang   for (i=0,j=0; i<size; ++i) {
19276eed172SJunchao Zhang     if (ilengths[i]) {
1935f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Isend((void*)(ilengths+i),1,MPIU_INT,i,tag,comm,s_waits+j));
19476eed172SJunchao Zhang       j++;
19576eed172SJunchao Zhang     }
19676eed172SJunchao Zhang   }
19776eed172SJunchao Zhang 
19876eed172SJunchao Zhang   /* Post waits on sends and receivs */
1995f80ce2aSJacob Faibussowitsch   if (nrecvs+nsends) CHKERRMPI(MPI_Waitall(nrecvs+nsends,r_waits,w_status));
20076eed172SJunchao Zhang 
20176eed172SJunchao Zhang   /* Pack up the received data */
2025f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nrecvs,onodes));
20376eed172SJunchao Zhang   for (i=0; i<nrecvs; ++i) {
20476eed172SJunchao Zhang     (*onodes)[i] = w_status[i].MPI_SOURCE;
20576eed172SJunchao Zhang     if (w_status[i].MPI_SOURCE == rank) (*olengths)[i] = ilengths[rank]; /* See comments in PetscGatherMessageLengths */
20676eed172SJunchao Zhang   }
2075f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2(r_waits,w_status));
20876eed172SJunchao Zhang   PetscFunctionReturn(0);
20976eed172SJunchao Zhang }
21076eed172SJunchao Zhang 
211e5c89e4eSSatish Balay /*@C
212e5c89e4eSSatish Balay   PetscGatherMessageLengths2 - Computes info about messages that a MPI-node will receive,
213e5c89e4eSSatish Balay   including (from-id,length) pairs for each message. Same functionality as PetscGatherMessageLengths()
214e5c89e4eSSatish Balay   except it takes TWO ilenths and output TWO olengths.
215e5c89e4eSSatish Balay 
216d083f849SBarry Smith   Collective
217e5c89e4eSSatish Balay 
218e5c89e4eSSatish Balay   Input Parameters:
219e5c89e4eSSatish Balay + comm      - Communicator
220e5c89e4eSSatish Balay . nsends    - number of messages that are to be sent.
221e5c89e4eSSatish Balay . nrecvs    - number of messages being received
2226b867d5aSJose E. Roman . ilengths1 - first array of integers of length sizeof(comm)
2236b867d5aSJose E. Roman - ilengths2 - second array of integers of length sizeof(comm)
224e5c89e4eSSatish Balay 
225e5c89e4eSSatish Balay   Output Parameters:
226e5c89e4eSSatish Balay + onodes    - list of node-ids from which messages are expected
2276b867d5aSJose E. Roman . olengths1 - first corresponding message lengths
2286b867d5aSJose E. Roman - olengths2 - second  message lengths
229e5c89e4eSSatish Balay 
230e5c89e4eSSatish Balay   Level: developer
231e5c89e4eSSatish Balay 
232e5c89e4eSSatish Balay   Notes:
233e5c89e4eSSatish Balay   With this info, the correct MPI_Irecv() can be posted with the correct
234e5c89e4eSSatish Balay   from-id, with a buffer with the right amount of memory required.
235e5c89e4eSSatish Balay 
236e5c89e4eSSatish Balay   The calling function deallocates the memory in onodes and olengths
237e5c89e4eSSatish Balay 
238c2916339SPierre Jolivet   To determine nrecvs, one can use PetscGatherNumberOfMessages()
239e5c89e4eSSatish Balay 
240e5c89e4eSSatish Balay .seealso: PetscGatherMessageLengths() and PetscGatherNumberOfMessages()
241e5c89e4eSSatish Balay @*/
2427087cfbeSBarry Smith PetscErrorCode  PetscGatherMessageLengths2(MPI_Comm comm,PetscMPIInt nsends,PetscMPIInt nrecvs,const PetscMPIInt ilengths1[],const PetscMPIInt ilengths2[],PetscMPIInt **onodes,PetscMPIInt **olengths1,PetscMPIInt **olengths2)
243e5c89e4eSSatish Balay {
2440298fd71SBarry Smith   PetscMPIInt    size,tag,i,j,*buf_s = NULL,*buf_r = NULL,*buf_j = NULL;
2450298fd71SBarry Smith   MPI_Request    *s_waits  = NULL,*r_waits = NULL;
2460298fd71SBarry Smith   MPI_Status     *w_status = NULL;
247e5c89e4eSSatish Balay 
248e5c89e4eSSatish Balay   PetscFunctionBegin;
2495f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(comm,&size));
2505f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCommGetNewTag(comm,&tag));
251e5c89e4eSSatish Balay 
2523bf92927SBarry Smith   /* cannot use PetscMalloc5() because r_waits and s_waits must be contiguous for the call to MPI_Waitall() */
2535f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc4(nrecvs+nsends,&r_waits,2*nrecvs,&buf_r,2*nsends,&buf_s,nrecvs+nsends,&w_status));
254e5c89e4eSSatish Balay   s_waits = r_waits + nrecvs;
255e5c89e4eSSatish Balay 
256e5c89e4eSSatish Balay   /* Post the Irecv to get the message length-info */
2575f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nrecvs+1,olengths1));
2585f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nrecvs+1,olengths2));
259e5c89e4eSSatish Balay   for (i=0; i<nrecvs; i++) {
260e5c89e4eSSatish Balay     buf_j = buf_r + (2*i);
2615f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Irecv(buf_j,2,MPI_INT,MPI_ANY_SOURCE,tag,comm,r_waits+i));
262e5c89e4eSSatish Balay   }
263e5c89e4eSSatish Balay 
264e5c89e4eSSatish Balay   /* Post the Isends with the message length-info */
265e5c89e4eSSatish Balay   for (i=0,j=0; i<size; ++i) {
266e5c89e4eSSatish Balay     if (ilengths1[i]) {
267e5c89e4eSSatish Balay       buf_j    = buf_s + (2*j);
268e5c89e4eSSatish Balay       buf_j[0] = *(ilengths1+i);
269e5c89e4eSSatish Balay       buf_j[1] = *(ilengths2+i);
2705f80ce2aSJacob Faibussowitsch       CHKERRMPI(MPI_Isend(buf_j,2,MPI_INT,i,tag,comm,s_waits+j));
271e5c89e4eSSatish Balay       j++;
272e5c89e4eSSatish Balay     }
273e5c89e4eSSatish Balay   }
2742c71b3e2SJacob Faibussowitsch   PetscCheckFalse(j != nsends,PETSC_COMM_SELF,PETSC_ERR_PLIB,"j %d not equal to expected number of sends %d",j,nsends);
275e5c89e4eSSatish Balay 
276e5c89e4eSSatish Balay   /* Post waits on sends and receivs */
2775f80ce2aSJacob Faibussowitsch   if (nrecvs+nsends) CHKERRMPI(MPI_Waitall(nrecvs+nsends,r_waits,w_status));
278e5c89e4eSSatish Balay 
279e5c89e4eSSatish Balay   /* Pack up the received data */
2805f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nrecvs+1,onodes));
281e5c89e4eSSatish Balay   for (i=0; i<nrecvs; ++i) {
282e5c89e4eSSatish Balay     (*onodes)[i]    = w_status[i].MPI_SOURCE;
283e5c89e4eSSatish Balay     buf_j           = buf_r + (2*i);
284e5c89e4eSSatish Balay     (*olengths1)[i] = buf_j[0];
285e5c89e4eSSatish Balay     (*olengths2)[i] = buf_j[1];
286e5c89e4eSSatish Balay   }
287e5c89e4eSSatish Balay 
2885f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree4(r_waits,buf_r,buf_s,w_status));
289e5c89e4eSSatish Balay   PetscFunctionReturn(0);
290e5c89e4eSSatish Balay }
291e5c89e4eSSatish Balay 
292e5c89e4eSSatish Balay /*
293e5c89e4eSSatish Balay 
294a5b23f4aSJose E. Roman   Allocate a buffer sufficient to hold messages of size specified in olengths.
295e5c89e4eSSatish Balay   And post Irecvs on these buffers using node info from onodes
296e5c89e4eSSatish Balay 
297e5c89e4eSSatish Balay  */
2987087cfbeSBarry Smith PetscErrorCode  PetscPostIrecvInt(MPI_Comm comm,PetscMPIInt tag,PetscMPIInt nrecvs,const PetscMPIInt onodes[],const PetscMPIInt olengths[],PetscInt ***rbuf,MPI_Request **r_waits)
299e5c89e4eSSatish Balay {
300c05d87d6SBarry Smith   PetscInt       **rbuf_t,i,len = 0;
301e5c89e4eSSatish Balay   MPI_Request    *r_waits_t;
302e5c89e4eSSatish Balay 
303e5c89e4eSSatish Balay   PetscFunctionBegin;
304e5c89e4eSSatish Balay   /* compute memory required for recv buffers */
305e5c89e4eSSatish Balay   for (i=0; i<nrecvs; i++) len += olengths[i];  /* each message length */
306e5c89e4eSSatish Balay 
307e5c89e4eSSatish Balay   /* allocate memory for recv buffers */
3085f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nrecvs+1,&rbuf_t));
3095f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(len,&rbuf_t[0]));
310e5c89e4eSSatish Balay   for (i=1; i<nrecvs; ++i) rbuf_t[i] = rbuf_t[i-1] + olengths[i-1];
311e5c89e4eSSatish Balay 
312e5c89e4eSSatish Balay   /* Post the receives */
3135f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nrecvs,&r_waits_t));
314e5c89e4eSSatish Balay   for (i=0; i<nrecvs; ++i) {
3155f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Irecv(rbuf_t[i],olengths[i],MPIU_INT,onodes[i],tag,comm,r_waits_t+i));
316e5c89e4eSSatish Balay   }
317e5c89e4eSSatish Balay 
318e5c89e4eSSatish Balay   *rbuf    = rbuf_t;
319e5c89e4eSSatish Balay   *r_waits = r_waits_t;
320e5c89e4eSSatish Balay   PetscFunctionReturn(0);
321e5c89e4eSSatish Balay }
322e5c89e4eSSatish Balay 
3237087cfbeSBarry Smith PetscErrorCode  PetscPostIrecvScalar(MPI_Comm comm,PetscMPIInt tag,PetscMPIInt nrecvs,const PetscMPIInt onodes[],const PetscMPIInt olengths[],PetscScalar ***rbuf,MPI_Request **r_waits)
324e5c89e4eSSatish Balay {
325052f0c41SBarry Smith   PetscMPIInt    i;
326e5c89e4eSSatish Balay   PetscScalar    **rbuf_t;
327e5c89e4eSSatish Balay   MPI_Request    *r_waits_t;
328c05d87d6SBarry Smith   PetscInt       len = 0;
329e5c89e4eSSatish Balay 
330fe28d99cSBarry Smith   PetscFunctionBegin;
331e5c89e4eSSatish Balay   /* compute memory required for recv buffers */
332e5c89e4eSSatish Balay   for (i=0; i<nrecvs; i++) len += olengths[i];  /* each message length */
333e5c89e4eSSatish Balay 
334e5c89e4eSSatish Balay   /* allocate memory for recv buffers */
3355f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nrecvs+1,&rbuf_t));
3365f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(len,&rbuf_t[0]));
337e5c89e4eSSatish Balay   for (i=1; i<nrecvs; ++i) rbuf_t[i] = rbuf_t[i-1] + olengths[i-1];
338e5c89e4eSSatish Balay 
339e5c89e4eSSatish Balay   /* Post the receives */
3405f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nrecvs,&r_waits_t));
341e5c89e4eSSatish Balay   for (i=0; i<nrecvs; ++i) {
3425f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Irecv(rbuf_t[i],olengths[i],MPIU_SCALAR,onodes[i],tag,comm,r_waits_t+i));
343e5c89e4eSSatish Balay   }
344e5c89e4eSSatish Balay 
345e5c89e4eSSatish Balay   *rbuf    = rbuf_t;
346e5c89e4eSSatish Balay   *r_waits = r_waits_t;
347e5c89e4eSSatish Balay   PetscFunctionReturn(0);
348e5c89e4eSSatish Balay }
349