1 /*  DenseMtx_gather.c  */
2 
3 #include "../spoolesMPI.h"
4 
5 /*--------------------------------------------------------------------*/
6 typedef struct _Msg Msg ;
7 struct _Msg {
8    int           id    ;
9    int           size  ;
10    double        *base ;
11    MPI_Request   req   ;
12    Msg           *next ;
13 } ;
14 /*--------------------------------------------------------------------*/
15 /*
16    --------------------------------------------------------------------
17    purpose -- to gather entries from a distributed X into Y .
18       what rows of X to be sent to other processors are
19       found in sendIVL. what rows of Y to be received
20       from other processors are found in recvIVL.
21 
22    Y  -- on return, contains the rows specified by recvIVL.
23       row indices of Y will be in ascending order.
24    X  -- this processor's part of the distributed partitioned
25       DenseMtx object. row indices of X are assumed to be
26       in ascending order.
27    sendIVL -- list jproc contains the global ids of rows in X
28       that need to be sent to processor jproc. note, lists are
29       assumed to be in ascending order and are local with respect to X.
30    recvIVL -- list jproc contains the global ids of rows in jproc's
31       part of X that will to be sent to this processor. note,
32       lists are assumed to be in ascending order and are local
33       with respect to Y.
34 
35    created -- 98jul31, cca
36    --------------------------------------------------------------------
37 */
38 void
DenseMtx_MPI_gatherRows(DenseMtx * Y,DenseMtx * X,IVL * sendIVL,IVL * recvIVL,int stats[],int msglvl,FILE * msgFile,int firsttag,MPI_Comm comm)39 DenseMtx_MPI_gatherRows (
40    DenseMtx   *Y,
41    DenseMtx   *X,
42    IVL        *sendIVL,
43    IVL        *recvIVL,
44    int        stats[],
45    int        msglvl,
46    FILE       *msgFile,
47    int        firsttag,
48    MPI_Comm   comm
49 ) {
50 double       *recvvec, *sendvec ;
51 int          flag, iirow, irow, jproc, jrow, lasttag, myid, ncolX,
52              ncolY, nproc, nrecv, nrowX, nrowY, nsend, nword, tag,
53              tagbound ;
54 int          *colindX, *colindY, *recvrowids,
55              *rowindX, *rowindY, *sendrowids ;
56 Msg          *msg, *nextmsg, *recvhead, *sendhead ;
57 MPI_Status   status ;
58 /*
59    ---------------
60    check the input
61    ---------------
62 */
63 if ( Y == NULL || X == NULL || sendIVL == NULL || recvIVL == NULL
64      || (msglvl > 0 && msgFile == NULL) ) {
65    fprintf(stderr, "\n fatal error in DenseMtx_MPI_gatherRows()"
66            "\n bad input\n") ;
67    fprintf(stderr, "\n Y %p, X %p, sendIVL %p, recvIVL %p",
68            Y, X, sendIVL, recvIVL) ;
69    exit(-1) ;
70 }
71 MPI_Comm_rank(comm, &myid) ;
72 MPI_Comm_size(comm, &nproc) ;
73 lasttag = firsttag + nproc*nproc ;
74 if ( lasttag > (tagbound = maxTagMPI(comm)) ) {
75    fprintf(stderr, "\n fatal error in DenseMtx_MPI_gatherRows()"
76           "\n lasttag = %d, tag_bound = %d", lasttag, tagbound) ;
77    exit(-1) ;
78 }
79 if ( DENSEMTX_IS_REAL(X) ) {
80    nword = 1 ;
81 } else if ( DENSEMTX_IS_COMPLEX(X) ) {
82    nword = 2 ;
83 } else {
84    fprintf(stderr, "\n fatal error in DenseMtx_MPI_gatherRows()"
85            "\n X->type = %d\n", X->type) ;
86    exit(-1) ;
87 }
88 DenseMtx_columnIndices(Y, &ncolY,  &colindY) ;
89 DenseMtx_rowIndices(Y, &nrowY, &rowindY) ;
90 DenseMtx_columnIndices(X, &ncolX,  &colindX) ;
91 DenseMtx_rowIndices(X, &nrowX, &rowindX) ;
92 if ( msglvl > 2 ) {
93    fprintf(msgFile, "\n\n sendIVL ") ;
94    IVL_writeForHumanEye(sendIVL, msgFile) ;
95    fprintf(msgFile, "\n\n recvIVL ") ;
96    IVL_writeForHumanEye(recvIVL, msgFile) ;
97    fflush(msgFile) ;
98 }
99 /*
100    ----------------------
101    load the internal rows
102    ----------------------
103 */
104 if ( msglvl > 2 ) {
105    fprintf(msgFile, "\n\n loading internal rows") ;
106    fflush(msgFile) ;
107 }
108 IVL_listAndSize(sendIVL, myid, &nsend, &sendrowids) ;
109 IVL_listAndSize(recvIVL, myid, &nrecv, &recvrowids) ;
110 for ( iirow = 0 ; iirow < nsend ; iirow++ ) {
111    irow = sendrowids[iirow] ;
112    jrow = recvrowids[iirow] ;
113    if ( msglvl > 2 ) {
114       fprintf(msgFile, "\n irow %d, jrow %d", irow, jrow) ;
115       fflush(msgFile) ;
116    }
117 /*
118    DenseMtx_copyRow(Y, jrow, X, irow) ;
119 */
120    DenseMtx_copyRowAndIndex(Y, jrow, X, irow) ;
121 }
122 if ( msglvl > 2 ) {
123    fprintf(msgFile, "\n\n after loading internal rows") ;
124    DenseMtx_writeForHumanEye(Y, msgFile) ;
125    fflush(msgFile) ;
126 }
127 /*
128    -------------------------------
129    post the sends and the receives
130    -------------------------------
131 */
132 recvhead = sendhead = NULL ;
133 for ( jproc = 0 ; jproc < nproc ; jproc++ ) {
134    if ( jproc != myid ) {
135       IVL_listAndSize(sendIVL, jproc, &nsend, &sendrowids) ;
136       IVL_listAndSize(recvIVL, jproc, &nrecv, &recvrowids) ;
137       if ( msglvl > 2 ) {
138          fprintf(msgFile, "\n\n jproc %d, nsend %d, nrecv %d",
139                  jproc, nsend, nrecv) ;
140          fflush(msgFile) ;
141       }
142       if ( nsend > 0 ) {
143 /*
144          -------------------------
145          create the message object
146          -------------------------
147 */
148          ALLOCATE(msg, struct _Msg, 1) ;
149          msg->id   = jproc ;
150          msg->size = nword * nsend * ncolY ;
151          msg->base = sendvec = DVinit(msg->size, 0.0) ;
152          msg->next = sendhead, sendhead = msg ;
153          tag       = firsttag + myid*nproc + jproc ;
154 /*
155          -----------------------------------
156          fill the buffer with matrix entries
157          -----------------------------------
158 */
159          for ( iirow = 0 ; iirow < nsend ; iirow++ ) {
160             irow = sendrowids[iirow] ;
161             DenseMtx_copyRowIntoVector(X, irow, sendvec) ;
162             sendvec += nword*ncolY ;
163          }
164          if ( msglvl > 2 ) {
165             fprintf(msgFile, "\n sendvec") ;
166             DVfprintf(msgFile, msg->size, msg->base) ;
167             fflush(msgFile) ;
168          }
169 /*
170          -------------
171          post the send
172          -------------
173 */
174          stats[0]++ ;
175          stats[2] += msg->size * sizeof(double) ;
176          if ( msglvl > 2 ) {
177             fprintf(msgFile,
178                     "\n posting Isend to %d, size %d, tag %d",
179                     jproc, msg->size, tag) ;
180             fflush(msgFile) ;
181          }
182          MPI_Isend(msg->base, msg->size, MPI_DOUBLE,
183                    jproc, tag, comm, &msg->req) ;
184 
185       }
186       if ( nrecv > 0 ) {
187 /*
188          -------------------------
189          create the message object
190          -------------------------
191 */
192          ALLOCATE(msg, struct _Msg, 1) ;
193          msg->id   = jproc ;
194          msg->size = nword * nrecv * ncolY ;
195          msg->base = (void *) DVinit(msg->size, 0.0) ;
196          msg->next = recvhead, recvhead = msg ;
197          tag       = firsttag + jproc*nproc + myid ;
198 /*
199          ----------------
200          post the receive
201          ----------------
202 */
203          stats[1]++ ;
204          stats[3] += msg->size * sizeof(double) ;
205          if ( msglvl > 2 ) {
206             fprintf(msgFile,
207                     "\n posting Irecv from %d, size %d, tag %d",
208                     jproc, msg->size, tag) ;
209             fflush(msgFile) ;
210          }
211          MPI_Irecv(msg->base, msg->size, MPI_DOUBLE,
212                    jproc, tag, comm, &msg->req) ;
213 
214       }
215    }
216 }
217 /*
218    -------------------------------------------
219    loop while there are messages to receive or
220    sent messages that have not been received
221    -------------------------------------------
222 */
223 while ( sendhead != NULL || recvhead != NULL ) {
224    for ( msg = sendhead, sendhead = NULL ;
225          msg != NULL ;
226          msg = nextmsg ) {
227       nextmsg = msg->next ;
228       if ( msglvl > 2 ) {
229          fprintf(msgFile, "\n msg %p to %d", msg, msg->id) ;
230          fflush(msgFile) ;
231       }
232       MPI_Test(&msg->req, &flag, &status) ;
233       if ( flag == 1 ) {
234          if ( msglvl > 2 ) {
235             fprintf(msgFile, ", received") ;
236             fflush(msgFile) ;
237          }
238          DVfree((double *) msg->base) ;
239          FREE(msg) ;
240       } else {
241          msg->next = sendhead, sendhead = msg ;
242       }
243    }
244    for ( msg = recvhead, recvhead = NULL ;
245          msg != NULL ;
246          msg = nextmsg ) {
247       nextmsg = msg->next ;
248       if ( msglvl > 2 ) {
249          fprintf(msgFile, "\n msg %p from %d", msg, msg->id) ;
250          fflush(msgFile) ;
251       }
252       MPI_Test(&msg->req, &flag, &status) ;
253       if ( flag == 1 ) {
254          jproc = msg->id ;
255          if ( msglvl > 2 ) {
256             fprintf(msgFile, ", received") ;
257             fflush(msgFile) ;
258          }
259          IVL_listAndSize(recvIVL, jproc, &nrecv, &recvrowids) ;
260          if ( msglvl > 2 ) {
261             fprintf(msgFile, "\n recvrowids") ;
262             IVfprintf(msgFile, nrecv, recvrowids) ;
263             fflush(msgFile) ;
264          }
265          recvvec = msg->base ;
266          if ( msglvl > 2 ) {
267             fprintf(msgFile, "\n recvvec") ;
268             DVfprintf(msgFile, nword*nrecv*ncolY, recvvec) ;
269             fflush(msgFile) ;
270          }
271          for ( iirow = 0 ; iirow < nrecv ; iirow++ ) {
272             irow = recvrowids[iirow] ;
273             DenseMtx_copyVectorIntoRow(Y, irow, recvvec) ;
274             recvvec += nword*ncolY ;
275          }
276          DVfree((double *) msg->base) ;
277          FREE(msg) ;
278       } else {
279          msg->next = recvhead, recvhead = msg ;
280       }
281    }
282 }
283 return ; }
284 
285 /*--------------------------------------------------------------------*/
286