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