1 /* splitFrontMtx.c */
2
3 #include "../spoolesMPI.h"
4
5 /*--------------------------------------------------------------------*/
6 typedef struct _msg Msg ;
7 struct _msg {
8 int rowid ;
9 int colid ;
10 int nbytes ;
11 SubMtx *mtx ;
12 MPI_Request req ;
13 Msg *next ;
14 } ;
15 static void split ( FrontMtx *frontmtx, SolveMap *solvemap, char cflag,
16 int stats[], int msglvl, FILE *msgFile, int firsttag, MPI_Comm comm ) ;
17 /*--------------------------------------------------------------------*/
18 /*
19 -------------------------------------------------------------
20 purpose -- after the factorization has been computed and the
21 front matrices have been split into submatrices, and after
22 a solve map object has been computed, the submatrices
23 are sent to the process that owns them.
24
25 frontmtx -- stores the factor matrix
26 solvemap -- stores the map from submatrices to processes
27 created -- 98may21, cca
28 -------------------------------------------------------------
29 */
30 void
FrontMtx_MPI_split(FrontMtx * frontmtx,SolveMap * solvemap,int stats[],int msglvl,FILE * msgFile,int firsttag,MPI_Comm comm)31 FrontMtx_MPI_split (
32 FrontMtx *frontmtx,
33 SolveMap *solvemap,
34 int stats[],
35 int msglvl,
36 FILE *msgFile,
37 int firsttag,
38 MPI_Comm comm
39 ) {
40 /*
41 --------------
42 check the data
43 --------------
44 */
45 if ( frontmtx == NULL || solvemap == NULL
46 || (msglvl > 0 && msgFile == NULL) ) {
47 fprintf(msgFile, "\n fatal error in FrontMtx_MPI_split()"
48 "\n frontmtx %p, solvemap %p, firsttag %d"
49 "\n stats %p, msglvl %d, msgFile %p"
50 "\n bad input\n",
51 frontmtx, solvemap, firsttag, stats, msglvl, msgFile) ;
52 exit(-1) ;
53 }
54 split(frontmtx, solvemap, 'U', stats,
55 msglvl, msgFile, firsttag, comm) ;
56 if ( FRONTMTX_IS_NONSYMMETRIC(frontmtx) ) {
57 split(frontmtx, solvemap, 'L', stats,
58 msglvl, msgFile, firsttag, comm) ;
59 }
60 return ; }
61
62 /*--------------------------------------------------------------------*/
63 /*
64 -------------------------------------------------------------
65 purpose -- after the factorization has been computed and the
66 front matrices have been split into submatrices, and after
67 a solve map object has been computed, the submatrices
68 are sent to the process that owns them.
69
70 frontmtx -- stores the factor matrix
71 solvemap -- stores the map from submatrices to processes
72 cflag -- 'U' --> split the upper matrix
73 'L' --> split the lower matrix
74
75 created -- 98may21, cca
76 -------------------------------------------------------------
77 */
78 static void
split(FrontMtx * frontmtx,SolveMap * solvemap,char cflag,int stats[],int msglvl,FILE * msgFile,int firsttag,MPI_Comm comm)79 split (
80 FrontMtx *frontmtx,
81 SolveMap *solvemap,
82 char cflag,
83 int stats[],
84 int msglvl,
85 FILE *msgFile,
86 int firsttag,
87 MPI_Comm comm
88 ) {
89 SubMtx *mtx ;
90 SubMtxManager *mtxmanager ;
91 int colid, count, destination, flag, incount, iproc,
92 J, JK, kk, K, left, maxcount, myid, nblock, nbytes,
93 nproc, offset, outcount, right, rowid, source, tag ;
94 int *colids, *inbuff, *incounts, *map, *offsets, *owners,
95 *outbuff, *outcounts, *rowids ;
96 int **p_inbuff, **p_outbuff ;
97 I2Ohash *hash ;
98 Msg *messages, *msg, *nextmsg, *recvhead, *sendhead ;
99 MPI_Status status ;
100
101 if ( msglvl > 1 ) {
102 fprintf(msgFile, "\n\n ##### inside FrontMtx_MPI_split()") ;
103 fflush(msgFile) ;
104 }
105 /*
106 --------------
107 check the data
108 --------------
109 */
110 if ( frontmtx == NULL || solvemap == NULL
111 || (msglvl > 0 && msgFile == NULL) ) {
112 fprintf(msgFile, "\n fatal error in FrontMtx_MPI_split()"
113 "\n frontmtx %p, solvemap %p, firsttag %d"
114 "\n stats %p, msglvl %d, msgFile %p"
115 "\n bad input\n",
116 frontmtx, solvemap, firsttag, stats, msglvl, msgFile) ;
117 exit(-1) ;
118 }
119 /*
120 ---------------------------------
121 get id of self and # of processes
122 ---------------------------------
123 */
124 MPI_Comm_rank(comm, &myid) ;
125 MPI_Comm_size(comm, &nproc) ;
126 mtxmanager = frontmtx->manager ;
127 if ( cflag == 'U' ) {
128 hash = frontmtx->upperhash ;
129 } else {
130 hash = frontmtx->lowerhash ;
131 }
132 if ( cflag == 'U' ) {
133 nblock = SolveMap_nblockUpper(solvemap) ;
134 owners = SolveMap_owners(solvemap) ;
135 rowids = SolveMap_rowidsUpper(solvemap) ;
136 colids = SolveMap_colidsUpper(solvemap) ;
137 map = SolveMap_mapUpper(solvemap) ;
138 } else if ( frontmtx->pivotingflag != 0 ) {
139 nblock = SolveMap_nblockLower(solvemap) ;
140 owners = SolveMap_owners(solvemap) ;
141 rowids = SolveMap_rowidsLower(solvemap) ;
142 colids = SolveMap_colidsLower(solvemap) ;
143 map = SolveMap_mapLower(solvemap) ;
144 } else {
145 nblock = SolveMap_nblockUpper(solvemap) ;
146 owners = SolveMap_owners(solvemap) ;
147 colids = SolveMap_rowidsUpper(solvemap) ;
148 rowids = SolveMap_colidsUpper(solvemap) ;
149 map = SolveMap_mapUpper(solvemap) ;
150 }
151 /*
152 ------------------------------------------------
153 step 1: compute the number of submatrices this
154 process will send to every other process
155 ------------------------------------------------
156 */
157 incounts = IVinit(nproc, 0) ;
158 outcounts = IVinit(nproc, 0) ;
159 if ( cflag == 'U' ) {
160 for ( JK = 0 ; JK < nblock ; JK++ ) {
161 J = rowids[JK] ;
162 K = colids[JK] ;
163 if ( owners[J] == myid
164 && (iproc = map[JK]) != myid
165 && (mtx = FrontMtx_upperMtx(frontmtx, J, K)) != NULL ) {
166 outcounts[iproc]++ ;
167 }
168 }
169 } else {
170 for ( JK = 0 ; JK < nblock ; JK++ ) {
171 K = rowids[JK] ;
172 J = colids[JK] ;
173 if ( owners[J] == myid
174 && (iproc = map[JK]) != myid
175 && (mtx = FrontMtx_lowerMtx(frontmtx, K, J)) != NULL ) {
176 outcounts[iproc]++ ;
177 }
178 }
179 }
180 /*
181 --------------------------------------------------------
182 step 2: do an all-to-all gather/scatter to get the
183 number of incoming submatrices from each process
184 --------------------------------------------------------
185 */
186 MPI_Alltoall((void *) outcounts, 1, MPI_INT,
187 (void *) incounts, 1, MPI_INT, comm) ;
188 if ( msglvl > 1 ) {
189 fprintf(msgFile, "\n\n incounts") ;
190 IVfprintf(msgFile, nproc, incounts) ;
191 fprintf(msgFile, "\n\n outcounts") ;
192 IVfprintf(msgFile, nproc, outcounts) ;
193 fflush(msgFile) ;
194 }
195 /*
196 -----------------------------------
197 step 3: allocate the buffer vectors
198 -----------------------------------
199 */
200 ALLOCATE(p_outbuff, int *, nproc) ;
201 ALLOCATE(p_inbuff, int *, nproc) ;
202 for ( iproc = 0 ; iproc < nproc ; iproc++ ) {
203 if ( outcounts[iproc] > 0 ) {
204 p_outbuff[iproc] = IVinit(3*outcounts[iproc], 0) ;
205 } else {
206 p_outbuff[iproc] = NULL ;
207 }
208 if ( incounts[iproc] > 0 ) {
209 p_inbuff[iproc] = IVinit(3*incounts[iproc], 0) ;
210 } else {
211 p_inbuff[iproc] = NULL ;
212 }
213 }
214 /*
215 -----------------------------------
216 step 4: fill the outbuffer vectors
217 ----------------------------------
218 */
219 offsets = IVinit(nproc, 0) ;
220 if ( cflag == 'U' ) {
221 for ( JK = kk = 0 ; JK < nblock ; JK++ ) {
222 J = rowids[JK] ;
223 K = colids[JK] ;
224 if ( owners[J] == myid
225 && (iproc = map[JK]) != myid
226 && (mtx = FrontMtx_upperMtx(frontmtx, J, K)) != NULL ) {
227 outbuff = p_outbuff[iproc] ;
228 kk = offsets[iproc] ;
229 outbuff[kk++] = J ;
230 outbuff[kk++] = K ;
231 outbuff[kk++] = SubMtx_nbytesInUse(mtx) ;
232 offsets[iproc] = kk ;
233 }
234 }
235 } else {
236 for ( JK = kk = 0 ; JK < nblock ; JK++ ) {
237 K = rowids[JK] ;
238 J = colids[JK] ;
239 if ( owners[J] == myid
240 && (iproc = map[JK]) != myid
241 && (mtx = FrontMtx_lowerMtx(frontmtx, K, J)) != NULL ) {
242 outbuff = p_outbuff[iproc] ;
243 kk = offsets[iproc] ;
244 outbuff[kk++] = K ;
245 outbuff[kk++] = J ;
246 outbuff[kk++] = SubMtx_nbytesInUse(mtx) ;
247 offsets[iproc] = kk ;
248 }
249 }
250 }
251 IVfree(offsets) ;
252 /*
253 -------------------------------------------------
254 step 5: send/receive the buffer vectors in stages
255 -------------------------------------------------
256 */
257 for ( offset = 1, tag = firsttag ; offset < nproc ; offset++, tag++ ) {
258 right = (myid + offset) % nproc ;
259 left = (nproc + myid - offset) % nproc ;
260 outcount = outcounts[right] ;
261 incount = incounts[left] ;
262 if ( msglvl > 1 ) {
263 fprintf(msgFile,
264 "\n ### process %d, send %d to right %d, recv %d from left %d",
265 myid, outcount, right, incount, left) ;
266 fflush(msgFile) ;
267 }
268 if ( outcount > 0 ) {
269 destination = right ;
270 stats[0]++ ;
271 stats[2] += 3*outcount*sizeof(int) ;
272 outbuff = p_outbuff[right] ;
273 } else {
274 destination = MPI_PROC_NULL ;
275 outbuff = NULL ;
276 }
277 if ( incount > 0 ) {
278 source = left ;
279 stats[1]++ ;
280 stats[3] += 3*incount*sizeof(int) ;
281 inbuff = p_inbuff[left] ;
282 } else {
283 source = MPI_PROC_NULL ;
284 inbuff = NULL ;
285 }
286 /*
287 -----------------
288 do a send/receive
289 -----------------
290 */
291 if ( msglvl > 1 ) {
292 fprintf(msgFile,
293 "\n ### Sendrecv: dest %d, tag %d, source %d, tag %d",
294 destination, tag, source, tag) ;
295 fflush(msgFile) ;
296 }
297 MPI_Sendrecv(outbuff, 3*outcount, MPI_INT, destination, tag,
298 inbuff, 3*incount, MPI_INT, source, tag,
299 comm, &status) ;
300 if ( source != MPI_PROC_NULL ) {
301 MPI_Get_count(&status, MPI_INT, &count) ;
302 if ( count != 3*incount ) {
303 fprintf(stderr,
304 "\n 1. fatal error in FrontMtx_MPI_split()"
305 "\n proc %d : source = %d, count = %d, incount = %d\n",
306 myid, source, count, 3*incount) ;
307 exit(-1) ;
308 }
309 }
310 if ( incount > 0 && msglvl > 1 ) {
311 fprintf(msgFile, "\n inbuffer from proc %d", left) ;
312 for ( kk = 0 ; kk < incount ; kk++ ) {
313 fprintf(msgFile, "\n %8d %8d %8d",
314 inbuff[3*kk], inbuff[3*kk+1], inbuff[3*kk+2]) ;
315 }
316 fflush(msgFile) ;
317 }
318 }
319 if ( msglvl > 1 ) {
320 fprintf(msgFile,
321 "\n ### outbuff[] and inbuff[] sent and received\n") ;
322 fflush(msgFile) ;
323 }
324 /*
325 ----------------------------------------------
326 step 6: compute the maximum number of messages
327 sent and received in any stage
328 ----------------------------------------------
329 */
330 maxcount = 0 ;
331 for ( offset = 1 ; offset < nproc ; offset++ ) {
332 right = (myid + offset) % nproc ;
333 left = (nproc + myid - offset) % nproc ;
334 outcount = outcounts[right] ;
335 incount = incounts[left] ;
336 if ( maxcount < incount + outcount ) {
337 maxcount = incount + outcount ;
338 }
339 }
340 if ( msglvl > 1 ) {
341 fprintf(msgFile, "\n\n maxcount = %d", maxcount) ;
342 fflush(msgFile) ;
343 }
344 ALLOCATE(messages, struct _msg, maxcount) ;
345 /*
346 ----------------------------------------------
347 step 8: send/receive the submatrices in stages
348 ----------------------------------------------
349 */
350 if ( msglvl > 1 ) {
351 fprintf(msgFile, "\n\n starting to send/receive the submatrices") ;
352 fflush(msgFile) ;
353 }
354 tag = firsttag ;
355 for ( offset = 1 ; offset < nproc ; offset++ ) {
356 right = (myid + offset) % nproc ;
357 left = (nproc + myid - offset) % nproc ;
358 outcount = outcounts[right] ;
359 incount = incounts[left] ;
360 if ( msglvl > 1 ) {
361 fprintf(msgFile,
362 "\n ### process %d, send %d to right %d, recv %d from left %d",
363 myid, outcount, right, incount, left) ;
364 fflush(msgFile) ;
365 }
366 msg = messages ;
367 /*
368 ---------------------------------------------------
369 post the receives for submatrices from process left
370 ---------------------------------------------------
371 */
372 inbuff = p_inbuff[left] ;
373 for ( JK = 0, recvhead = NULL ; JK < incount ; JK++, msg++ ) {
374 rowid = inbuff[3*JK] ;
375 colid = inbuff[3*JK+1] ;
376 nbytes = inbuff[3*JK+2] ;
377 if ( msglvl > 1 ) {
378 fprintf(msgFile,
379 "\n ready to receive mtx(%d,%d) with %d bytes from %d",
380 rowid, colid, nbytes, left) ;
381 fflush(msgFile) ;
382 }
383 mtx = SubMtxManager_newObjectOfSizeNbytes(mtxmanager, nbytes);
384 I2Ohash_insert(hash, rowid, colid, mtx) ;
385 msg->rowid = rowid ;
386 msg->colid = colid ;
387 msg->mtx = mtx ;
388 msg->nbytes = nbytes ;
389 msg->next = recvhead ;
390 recvhead = msg ;
391 if ( msglvl > 1 ) {
392 fprintf(msgFile,
393 "\n ### posting Irecv, left %d, tag %d, nbytes %d",
394 left, JK + firsttag, nbytes) ;
395 fflush(msgFile) ;
396 }
397 MPI_Irecv(SubMtx_workspace(mtx), nbytes, MPI_BYTE,
398 left, JK + firsttag, comm, &msg->req) ;
399 if ( msglvl > 1 ) {
400 fprintf(msgFile, ", Irecv posted") ;
401 fflush(msgFile) ;
402 }
403 stats[1]++ ;
404 stats[3] += nbytes ;
405 }
406 /*
407 -----------------------------------------------------
408 post the sends for submatrices going to process right
409 -----------------------------------------------------
410 */
411 outbuff = p_outbuff[right] ;
412 for ( JK = kk = 0, sendhead = NULL ; JK < outcount ; JK++, msg++ ) {
413 rowid = outbuff[3*JK] ;
414 colid = outbuff[3*JK+1] ;
415 nbytes = outbuff[3*JK+2] ;
416 I2Ohash_remove(hash, rowid, colid, (void *) &mtx) ;
417 if ( msglvl > 1 ) {
418 fprintf(msgFile,
419 "\n ready to send mtx(%d,%d) with %d bytes to %d",
420 rowid, colid, nbytes, right) ;
421 fflush(msgFile) ;
422 }
423 msg->rowid = rowid ;
424 msg->colid = colid ;
425 msg->mtx = mtx ;
426 msg->nbytes = nbytes ;
427 msg->next = sendhead ;
428 sendhead = msg ;
429 if ( msglvl > 1 ) {
430 fprintf(msgFile,
431 "\n ### posting Isend, right %d, tag %d, nbytes %d",
432 right, JK + firsttag, nbytes) ;
433 SubMtx_writeForHumanEye(msg->mtx, msgFile) ;
434 fprintf(msgFile, "\n first seven entries of buffer") ;
435 IVfprintf(msgFile, 7, (int *) SubMtx_workspace(msg->mtx)) ;
436 fflush(msgFile) ;
437 }
438 MPI_Isend(SubMtx_workspace(mtx), nbytes, MPI_BYTE,
439 right, JK + firsttag, comm, &msg->req) ;
440 if ( msglvl > 1 ) {
441 fprintf(msgFile, "\n Isend posted") ;
442 fflush(msgFile) ;
443 }
444 stats[0]++ ;
445 stats[2] += nbytes ;
446 }
447 /*
448 -----------------------------------------
449 loop while send/receives are not complete
450 -----------------------------------------
451 */
452 while ( sendhead != NULL || recvhead != NULL ) {
453 for ( msg = sendhead, sendhead = NULL ;
454 msg != NULL ;
455 msg = nextmsg ) {
456 nextmsg = msg->next ;
457 if ( msglvl > 1 ) {
458 fprintf(msgFile,
459 "\n testing sent msg %p, (%d,%d) with %d bytes",
460 msg, msg->rowid, msg->colid, msg->nbytes) ;
461 fflush(msgFile) ;
462 }
463 MPI_Test(&msg->req, &flag, &status) ;
464 if ( flag == 1 ) {
465 /*
466 ---------------------------------------------
467 message has been received, release the matrix
468 ---------------------------------------------
469 */
470 if ( msglvl > 1 ) {
471 fprintf(msgFile, ", send complete, releasing mtx") ;
472 fflush(msgFile) ;
473 }
474 SubMtxManager_releaseObject(mtxmanager, msg->mtx) ;
475 } else {
476 /*
477 ---------------------------------------------------
478 message has not been received, keep message in list
479 ---------------------------------------------------
480 */
481 msg->next = sendhead ;
482 sendhead = msg ;
483 }
484 }
485 for ( msg = recvhead, recvhead = NULL ;
486 msg != NULL ;
487 msg = nextmsg ) {
488 nextmsg = msg->next ;
489 if ( msglvl > 1 ) {
490 fprintf(msgFile,
491 "\n testing recv msg %p, (%d,%d) with %d bytes",
492 msg, msg->rowid, msg->colid, msg->nbytes) ;
493 fflush(msgFile) ;
494 }
495 MPI_Test(&msg->req, &flag, &status) ;
496 if ( flag == 1 ) {
497 /*
498 ------------------------------------------------
499 message has been received, initialize the matrix
500 ------------------------------------------------
501 */
502 if ( msglvl > 1 ) {
503 fprintf(msgFile, ", recv complete, initializing mtx") ;
504 fprintf(msgFile, "\n first seven entries of buffer") ;
505 IVfprintf(msgFile, 7, (int *) SubMtx_workspace(msg->mtx)) ;
506 fflush(msgFile) ;
507 }
508 SubMtx_initFromBuffer(msg->mtx) ;
509 if ( msglvl > 1 ) {
510 SubMtx_writeForHumanEye(msg->mtx, msgFile) ;
511 fflush(msgFile) ;
512 }
513 } else {
514 /*
515 ---------------------------------------------------
516 message has not been received, keep message in list
517 ---------------------------------------------------
518 */
519 msg->next = recvhead ;
520 recvhead = msg ;
521 }
522 }
523 }
524 }
525 /*
526 ------------------------
527 free the working storage
528 ------------------------
529 */
530 IVfree(incounts) ;
531 IVfree(outcounts) ;
532 for ( iproc = 0 ; iproc < nproc ; iproc++ ) {
533 if ( (inbuff = p_inbuff[iproc]) != NULL ) {
534 IVfree(inbuff) ;
535 }
536 if ( (outbuff = p_outbuff[iproc]) != NULL ) {
537 IVfree(outbuff) ;
538 }
539 }
540 FREE(p_inbuff) ;
541 FREE(p_outbuff) ;
542 FREE(messages) ;
543
544 return ; }
545
546 /*--------------------------------------------------------------------*/
547