xref: /petsc/src/sys/utils/mpimesg.c (revision b2566f29c2b6470df769aa9f7deb9e2726b0959e)
1e5c89e4eSSatish Balay 
2c6db04a5SJed Brown #include <petscsys.h>        /*I  "petscsys.h"  I*/
3e5c89e4eSSatish Balay 
4e5c89e4eSSatish Balay 
5e5c89e4eSSatish Balay #undef __FUNCT__
6e5c89e4eSSatish Balay #define __FUNCT__ "PetscGatherNumberOfMessages"
7e5c89e4eSSatish Balay /*@C
8e5c89e4eSSatish Balay   PetscGatherNumberOfMessages -  Computes the number of messages a node expects to receive
9e5c89e4eSSatish Balay 
10e5c89e4eSSatish Balay   Collective on MPI_Comm
11e5c89e4eSSatish Balay 
12e5c89e4eSSatish Balay   Input Parameters:
13e5c89e4eSSatish Balay + comm     - Communicator
14e5c89e4eSSatish Balay . iflags   - an array of integers of length sizeof(comm). A '1' in ilengths[i] represent a
150298fd71SBarry Smith              message from current node to ith node. Optionally NULL
16e5c89e4eSSatish Balay - ilengths - Non zero ilengths[i] represent a message to i of length ilengths[i].
170298fd71SBarry Smith              Optionally NULL.
18e5c89e4eSSatish Balay 
19e5c89e4eSSatish Balay   Output Parameters:
20e5c89e4eSSatish Balay . nrecvs    - number of messages received
21e5c89e4eSSatish Balay 
22e5c89e4eSSatish Balay   Level: developer
23e5c89e4eSSatish Balay 
24e5c89e4eSSatish Balay   Concepts: mpi utility
25e5c89e4eSSatish Balay 
26e5c89e4eSSatish Balay   Notes:
27e5c89e4eSSatish Balay   With this info, the correct message lengths can be determined using
28e5c89e4eSSatish Balay   PetscGatherMessageLengths()
29e5c89e4eSSatish Balay 
30e5c89e4eSSatish Balay   Either iflags or ilengths should be provided.  If iflags is not
310298fd71SBarry Smith   provided (NULL) it can be computed from ilengths. If iflags is
32e5c89e4eSSatish Balay   provided, ilengths is not required.
33e5c89e4eSSatish Balay 
34e5c89e4eSSatish Balay .seealso: PetscGatherMessageLengths()
35e5c89e4eSSatish Balay @*/
367087cfbeSBarry Smith PetscErrorCode  PetscGatherNumberOfMessages(MPI_Comm comm,const PetscMPIInt iflags[],const PetscMPIInt ilengths[],PetscMPIInt *nrecvs)
37e5c89e4eSSatish Balay {
380298fd71SBarry Smith   PetscMPIInt    size,rank,*recv_buf,i,*iflags_local = NULL,*iflags_localm = NULL;
39e5c89e4eSSatish Balay   PetscErrorCode ierr;
40e5c89e4eSSatish Balay 
41e5c89e4eSSatish Balay   PetscFunctionBegin;
42e5c89e4eSSatish Balay   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
43e5c89e4eSSatish Balay   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
44e5c89e4eSSatish Balay 
45dcca6d9dSJed Brown   ierr = PetscMalloc2(size,&recv_buf,size,&iflags_localm);CHKERRQ(ierr);
46e5c89e4eSSatish Balay 
47e5c89e4eSSatish Balay   /* If iflags not provided, compute iflags from ilengths */
48e5c89e4eSSatish Balay   if (!iflags) {
49e32f2f54SBarry Smith     if (!ilengths) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Either iflags or ilengths should be provided");
50e5c89e4eSSatish Balay     iflags_local = iflags_localm;
51e5c89e4eSSatish Balay     for (i=0; i<size; i++) {
52e5c89e4eSSatish Balay       if (ilengths[i]) iflags_local[i] = 1;
53e5c89e4eSSatish Balay       else iflags_local[i] = 0;
54e5c89e4eSSatish Balay     }
55a297a907SKarl Rupp   } else iflags_local = (PetscMPIInt*) iflags;
56e5c89e4eSSatish Balay 
57e5c89e4eSSatish Balay   /* Post an allreduce to determine the numer of messages the current node will receive */
58*b2566f29SBarry Smith   ierr    = MPIU_Allreduce(iflags_local,recv_buf,size,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr);
59e5c89e4eSSatish Balay   *nrecvs = recv_buf[rank];
60e5c89e4eSSatish Balay 
61e5c89e4eSSatish Balay   ierr = PetscFree2(recv_buf,iflags_localm);CHKERRQ(ierr);
62e5c89e4eSSatish Balay   PetscFunctionReturn(0);
63e5c89e4eSSatish Balay }
64e5c89e4eSSatish Balay 
65e5c89e4eSSatish Balay 
66e5c89e4eSSatish Balay #undef __FUNCT__
67e5c89e4eSSatish Balay #define __FUNCT__ "PetscGatherMessageLengths"
68e5c89e4eSSatish Balay /*@C
69e5c89e4eSSatish Balay   PetscGatherMessageLengths - Computes info about messages that a MPI-node will receive,
70e5c89e4eSSatish Balay   including (from-id,length) pairs for each message.
71e5c89e4eSSatish Balay 
72e5c89e4eSSatish Balay   Collective on MPI_Comm
73e5c89e4eSSatish Balay 
74e5c89e4eSSatish Balay   Input Parameters:
75e5c89e4eSSatish Balay + comm      - Communicator
76e5c89e4eSSatish Balay . nsends    - number of messages that are to be sent.
77e5c89e4eSSatish Balay . nrecvs    - number of messages being received
78e5c89e4eSSatish Balay - ilengths  - an array of integers of length sizeof(comm)
79e5c89e4eSSatish Balay               a non zero ilengths[i] represent a message to i of length ilengths[i]
80e5c89e4eSSatish Balay 
81e5c89e4eSSatish Balay 
82e5c89e4eSSatish Balay   Output Parameters:
83e5c89e4eSSatish Balay + onodes    - list of node-ids from which messages are expected
84e5c89e4eSSatish Balay - olengths  - corresponding message lengths
85e5c89e4eSSatish Balay 
86e5c89e4eSSatish Balay   Level: developer
87e5c89e4eSSatish Balay 
88e5c89e4eSSatish Balay   Concepts: mpi utility
89e5c89e4eSSatish Balay 
90e5c89e4eSSatish Balay   Notes:
91e5c89e4eSSatish Balay   With this info, the correct MPI_Irecv() can be posted with the correct
92e5c89e4eSSatish Balay   from-id, with a buffer with the right amount of memory required.
93e5c89e4eSSatish Balay 
94e5c89e4eSSatish Balay   The calling function deallocates the memory in onodes and olengths
95e5c89e4eSSatish Balay 
96e5c89e4eSSatish Balay   To determine nrecevs, one can use PetscGatherNumberOfMessages()
97e5c89e4eSSatish Balay 
98e5c89e4eSSatish Balay .seealso: PetscGatherNumberOfMessages()
99e5c89e4eSSatish Balay @*/
1007087cfbeSBarry Smith PetscErrorCode  PetscGatherMessageLengths(MPI_Comm comm,PetscMPIInt nsends,PetscMPIInt nrecvs,const PetscMPIInt ilengths[],PetscMPIInt **onodes,PetscMPIInt **olengths)
101e5c89e4eSSatish Balay {
102e5c89e4eSSatish Balay   PetscErrorCode ierr;
103e5c89e4eSSatish Balay   PetscMPIInt    size,tag,i,j;
1040298fd71SBarry Smith   MPI_Request    *s_waits  = NULL,*r_waits = NULL;
1050298fd71SBarry Smith   MPI_Status     *w_status = NULL;
106e5c89e4eSSatish Balay 
107e5c89e4eSSatish Balay   PetscFunctionBegin;
108e5c89e4eSSatish Balay   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
109e5c89e4eSSatish Balay   ierr = PetscCommGetNewTag(comm,&tag);CHKERRQ(ierr);
110e5c89e4eSSatish Balay 
111e5c89e4eSSatish Balay   /* cannot use PetscMalloc3() here because in the call to MPI_Waitall() they MUST be contiguous */
112dcca6d9dSJed Brown   ierr    = PetscMalloc2(nrecvs+nsends,&r_waits,nrecvs+nsends,&w_status);CHKERRQ(ierr);
113e5c89e4eSSatish Balay   s_waits = r_waits+nrecvs;
114e5c89e4eSSatish Balay 
115e5c89e4eSSatish Balay   /* Post the Irecv to get the message length-info */
116785e854fSJed Brown   ierr = PetscMalloc1(nrecvs,olengths);CHKERRQ(ierr);
117e5c89e4eSSatish Balay   for (i=0; i<nrecvs; i++) {
118e5c89e4eSSatish Balay     ierr = MPI_Irecv((*olengths)+i,1,MPI_INT,MPI_ANY_SOURCE,tag,comm,r_waits+i);CHKERRQ(ierr);
119e5c89e4eSSatish Balay   }
120e5c89e4eSSatish Balay 
121e5c89e4eSSatish Balay   /* Post the Isends with the message length-info */
122e5c89e4eSSatish Balay   for (i=0,j=0; i<size; ++i) {
123e5c89e4eSSatish Balay     if (ilengths[i]) {
124300a7f5bSBarry Smith       ierr = MPI_Isend((void*)(ilengths+i),1,MPI_INT,i,tag,comm,s_waits+j);CHKERRQ(ierr);
125e5c89e4eSSatish Balay       j++;
126e5c89e4eSSatish Balay     }
127e5c89e4eSSatish Balay   }
128e5c89e4eSSatish Balay 
129e5c89e4eSSatish Balay   /* Post waits on sends and receivs */
130e5c89e4eSSatish Balay   if (nrecvs+nsends) {ierr = MPI_Waitall(nrecvs+nsends,r_waits,w_status);CHKERRQ(ierr);}
131e5c89e4eSSatish Balay 
132e5c89e4eSSatish Balay   /* Pack up the received data */
133785e854fSJed Brown   ierr = PetscMalloc1(nrecvs,onodes);CHKERRQ(ierr);
134a297a907SKarl Rupp   for (i=0; i<nrecvs; ++i) (*onodes)[i] = w_status[i].MPI_SOURCE;
135e5c89e4eSSatish Balay   ierr = PetscFree2(r_waits,w_status);CHKERRQ(ierr);
136e5c89e4eSSatish Balay   PetscFunctionReturn(0);
137e5c89e4eSSatish Balay }
138dd6ea824SBarry Smith 
139e5c89e4eSSatish Balay #undef __FUNCT__
140e5c89e4eSSatish Balay #define __FUNCT__ "PetscGatherMessageLengths2"
141e5c89e4eSSatish Balay /*@C
142e5c89e4eSSatish Balay   PetscGatherMessageLengths2 - Computes info about messages that a MPI-node will receive,
143e5c89e4eSSatish Balay   including (from-id,length) pairs for each message. Same functionality as PetscGatherMessageLengths()
144e5c89e4eSSatish Balay   except it takes TWO ilenths and output TWO olengths.
145e5c89e4eSSatish Balay 
146e5c89e4eSSatish Balay   Collective on MPI_Comm
147e5c89e4eSSatish Balay 
148e5c89e4eSSatish Balay   Input Parameters:
149e5c89e4eSSatish Balay + comm      - Communicator
150e5c89e4eSSatish Balay . nsends    - number of messages that are to be sent.
151e5c89e4eSSatish Balay . nrecvs    - number of messages being received
152e5c89e4eSSatish Balay - ilengths1, ilengths2 - array of integers of length sizeof(comm)
153e5c89e4eSSatish Balay               a non zero ilengths[i] represent a message to i of length ilengths[i]
154e5c89e4eSSatish Balay 
155e5c89e4eSSatish Balay   Output Parameters:
156e5c89e4eSSatish Balay + onodes    - list of node-ids from which messages are expected
157e5c89e4eSSatish Balay - olengths1, olengths2 - corresponding message lengths
158e5c89e4eSSatish Balay 
159e5c89e4eSSatish Balay   Level: developer
160e5c89e4eSSatish Balay 
161e5c89e4eSSatish Balay   Concepts: mpi utility
162e5c89e4eSSatish Balay 
163e5c89e4eSSatish Balay   Notes:
164e5c89e4eSSatish Balay   With this info, the correct MPI_Irecv() can be posted with the correct
165e5c89e4eSSatish Balay   from-id, with a buffer with the right amount of memory required.
166e5c89e4eSSatish Balay 
167e5c89e4eSSatish Balay   The calling function deallocates the memory in onodes and olengths
168e5c89e4eSSatish Balay 
169e5c89e4eSSatish Balay   To determine nrecevs, one can use PetscGatherNumberOfMessages()
170e5c89e4eSSatish Balay 
171e5c89e4eSSatish Balay .seealso: PetscGatherMessageLengths() and PetscGatherNumberOfMessages()
172e5c89e4eSSatish Balay @*/
1737087cfbeSBarry Smith PetscErrorCode  PetscGatherMessageLengths2(MPI_Comm comm,PetscMPIInt nsends,PetscMPIInt nrecvs,const PetscMPIInt ilengths1[],const PetscMPIInt ilengths2[],PetscMPIInt **onodes,PetscMPIInt **olengths1,PetscMPIInt **olengths2)
174e5c89e4eSSatish Balay {
175e5c89e4eSSatish Balay   PetscErrorCode ierr;
1760298fd71SBarry Smith   PetscMPIInt    size,tag,i,j,*buf_s = NULL,*buf_r = NULL,*buf_j = NULL;
1770298fd71SBarry Smith   MPI_Request    *s_waits  = NULL,*r_waits = NULL;
1780298fd71SBarry Smith   MPI_Status     *w_status = NULL;
179e5c89e4eSSatish Balay 
180e5c89e4eSSatish Balay   PetscFunctionBegin;
181e5c89e4eSSatish Balay   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
182e5c89e4eSSatish Balay   ierr = PetscCommGetNewTag(comm,&tag);CHKERRQ(ierr);
183e5c89e4eSSatish Balay 
184e5c89e4eSSatish Balay   /* cannot use PetscMalloc5() because r_waits and s_waits must be contiquous for the call to MPI_Waitall() */
185dcca6d9dSJed Brown   ierr = PetscMalloc4(nrecvs+nsends,&r_waits,2*nrecvs,&buf_r,2*nsends,&buf_s,nrecvs+nsends,&w_status);CHKERRQ(ierr);
186e5c89e4eSSatish Balay   s_waits = r_waits + nrecvs;
187e5c89e4eSSatish Balay 
188e5c89e4eSSatish Balay   /* Post the Irecv to get the message length-info */
189854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs+1,olengths1);CHKERRQ(ierr);
190854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs+1,olengths2);CHKERRQ(ierr);
191e5c89e4eSSatish Balay   for (i=0; i<nrecvs; i++) {
192e5c89e4eSSatish Balay     buf_j = buf_r + (2*i);
193e5c89e4eSSatish Balay     ierr  = MPI_Irecv(buf_j,2,MPI_INT,MPI_ANY_SOURCE,tag,comm,r_waits+i);CHKERRQ(ierr);
194e5c89e4eSSatish Balay   }
195e5c89e4eSSatish Balay 
196e5c89e4eSSatish Balay   /* Post the Isends with the message length-info */
197e5c89e4eSSatish Balay   for (i=0,j=0; i<size; ++i) {
198e5c89e4eSSatish Balay     if (ilengths1[i]) {
199e5c89e4eSSatish Balay       buf_j    = buf_s + (2*j);
200e5c89e4eSSatish Balay       buf_j[0] = *(ilengths1+i);
201e5c89e4eSSatish Balay       buf_j[1] = *(ilengths2+i);
202e5c89e4eSSatish Balay       ierr = MPI_Isend(buf_j,2,MPI_INT,i,tag,comm,s_waits+j);CHKERRQ(ierr);
203e5c89e4eSSatish Balay       j++;
204e5c89e4eSSatish Balay     }
205e5c89e4eSSatish Balay   }
206f327f304SBarry Smith   if (j != nsends) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"j %d not equal to expected number of sends %d\n",j,nsends);
207e5c89e4eSSatish Balay 
208e5c89e4eSSatish Balay   /* Post waits on sends and receivs */
209e5c89e4eSSatish Balay   if (nrecvs+nsends) {ierr = MPI_Waitall(nrecvs+nsends,r_waits,w_status);CHKERRQ(ierr);}
210e5c89e4eSSatish Balay 
211e5c89e4eSSatish Balay 
212e5c89e4eSSatish Balay   /* Pack up the received data */
213854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs+1,onodes);CHKERRQ(ierr);
214e5c89e4eSSatish Balay   for (i=0; i<nrecvs; ++i) {
215e5c89e4eSSatish Balay     (*onodes)[i]    = w_status[i].MPI_SOURCE;
216e5c89e4eSSatish Balay     buf_j           = buf_r + (2*i);
217e5c89e4eSSatish Balay     (*olengths1)[i] = buf_j[0];
218e5c89e4eSSatish Balay     (*olengths2)[i] = buf_j[1];
219e5c89e4eSSatish Balay   }
220e5c89e4eSSatish Balay 
221e5c89e4eSSatish Balay   ierr = PetscFree4(r_waits,buf_r,buf_s,w_status);CHKERRQ(ierr);
222e5c89e4eSSatish Balay   PetscFunctionReturn(0);
223e5c89e4eSSatish Balay }
224e5c89e4eSSatish Balay 
225e5c89e4eSSatish Balay /*
226e5c89e4eSSatish Balay 
227e5c89e4eSSatish Balay   Allocate a bufffer sufficient to hold messages of size specified in olengths.
228e5c89e4eSSatish Balay   And post Irecvs on these buffers using node info from onodes
229e5c89e4eSSatish Balay 
230e5c89e4eSSatish Balay  */
231e5c89e4eSSatish Balay #undef __FUNCT__
232e5c89e4eSSatish Balay #define __FUNCT__ "PetscPostIrecvInt"
2337087cfbeSBarry Smith PetscErrorCode  PetscPostIrecvInt(MPI_Comm comm,PetscMPIInt tag,PetscMPIInt nrecvs,const PetscMPIInt onodes[],const PetscMPIInt olengths[],PetscInt ***rbuf,MPI_Request **r_waits)
234e5c89e4eSSatish Balay {
235e5c89e4eSSatish Balay   PetscErrorCode ierr;
236c05d87d6SBarry Smith   PetscInt       **rbuf_t,i,len = 0;
237e5c89e4eSSatish Balay   MPI_Request    *r_waits_t;
238e5c89e4eSSatish Balay 
239e5c89e4eSSatish Balay   PetscFunctionBegin;
240e5c89e4eSSatish Balay   /* compute memory required for recv buffers */
241e5c89e4eSSatish Balay   for (i=0; i<nrecvs; i++) len += olengths[i];  /* each message length */
242e5c89e4eSSatish Balay 
243e5c89e4eSSatish Balay   /* allocate memory for recv buffers */
244854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs+1,&rbuf_t);CHKERRQ(ierr);
245785e854fSJed Brown   ierr = PetscMalloc1(len,&rbuf_t[0]);CHKERRQ(ierr);
246e5c89e4eSSatish Balay   for (i=1; i<nrecvs; ++i) rbuf_t[i] = rbuf_t[i-1] + olengths[i-1];
247e5c89e4eSSatish Balay 
248e5c89e4eSSatish Balay   /* Post the receives */
249785e854fSJed Brown   ierr = PetscMalloc1(nrecvs,&r_waits_t);CHKERRQ(ierr);
250e5c89e4eSSatish Balay   for (i=0; i<nrecvs; ++i) {
251e5c89e4eSSatish Balay     ierr = MPI_Irecv(rbuf_t[i],olengths[i],MPIU_INT,onodes[i],tag,comm,r_waits_t+i);CHKERRQ(ierr);
252e5c89e4eSSatish Balay   }
253e5c89e4eSSatish Balay 
254e5c89e4eSSatish Balay   *rbuf    = rbuf_t;
255e5c89e4eSSatish Balay   *r_waits = r_waits_t;
256e5c89e4eSSatish Balay   PetscFunctionReturn(0);
257e5c89e4eSSatish Balay }
258e5c89e4eSSatish Balay 
259e5c89e4eSSatish Balay #undef __FUNCT__
260e5c89e4eSSatish Balay #define __FUNCT__ "PetscPostIrecvScalar"
2617087cfbeSBarry Smith PetscErrorCode  PetscPostIrecvScalar(MPI_Comm comm,PetscMPIInt tag,PetscMPIInt nrecvs,const PetscMPIInt onodes[],const PetscMPIInt olengths[],PetscScalar ***rbuf,MPI_Request **r_waits)
262e5c89e4eSSatish Balay {
263e5c89e4eSSatish Balay   PetscErrorCode ierr;
264052f0c41SBarry Smith   PetscMPIInt    i;
265e5c89e4eSSatish Balay   PetscScalar    **rbuf_t;
266e5c89e4eSSatish Balay   MPI_Request    *r_waits_t;
267c05d87d6SBarry Smith   PetscInt       len = 0;
268e5c89e4eSSatish Balay 
269fe28d99cSBarry Smith   PetscFunctionBegin;
270e5c89e4eSSatish Balay   /* compute memory required for recv buffers */
271e5c89e4eSSatish Balay   for (i=0; i<nrecvs; i++) len += olengths[i];  /* each message length */
272e5c89e4eSSatish Balay 
273e5c89e4eSSatish Balay   /* allocate memory for recv buffers */
274854ce69bSBarry Smith   ierr = PetscMalloc1(nrecvs+1,&rbuf_t);CHKERRQ(ierr);
275785e854fSJed Brown   ierr = PetscMalloc1(len,&rbuf_t[0]);CHKERRQ(ierr);
276e5c89e4eSSatish Balay   for (i=1; i<nrecvs; ++i) rbuf_t[i] = rbuf_t[i-1] + olengths[i-1];
277e5c89e4eSSatish Balay 
278e5c89e4eSSatish Balay   /* Post the receives */
279785e854fSJed Brown   ierr = PetscMalloc1(nrecvs,&r_waits_t);CHKERRQ(ierr);
280e5c89e4eSSatish Balay   for (i=0; i<nrecvs; ++i) {
281e5c89e4eSSatish Balay     ierr = MPI_Irecv(rbuf_t[i],olengths[i],MPIU_SCALAR,onodes[i],tag,comm,r_waits_t+i);CHKERRQ(ierr);
282e5c89e4eSSatish Balay   }
283e5c89e4eSSatish Balay 
284e5c89e4eSSatish Balay   *rbuf    = rbuf_t;
285e5c89e4eSSatish Balay   *r_waits = r_waits_t;
286e5c89e4eSSatish Balay   PetscFunctionReturn(0);
287e5c89e4eSSatish Balay }
288