xref: /petsc/src/mat/impls/dense/mpi/mmdense.c (revision 99dcc5326efaeeb96f7fc5515ee10e547e6b62bd)
1 #ifdef PETSC_RCS_HEADER
2 static char vcid[] = "$Id: mmdense.c,v 1.15 1998/09/09 21:10:00 balay Exp bsmith $";
3 #endif
4 
5 /*
6    Support for the parallel dense matrix vector multiply
7 */
8 #include "src/mat/impls/dense/mpi/mpidense.h"
9 #include "src/vec/vecimpl.h"
10 
11 #undef __FUNC__
12 #define __FUNC__ "MatSetUpMultiply_MPIDense"
13 int MatSetUpMultiply_MPIDense(Mat mat)
14 {
15   Mat_MPIDense *mdn = (Mat_MPIDense *) mat->data;
16   int          ierr,n;
17   IS           tofrom;
18   Vec          gvec;
19 
20   PetscFunctionBegin;
21   /* Create local vector that is used to scatter into */
22   ierr = VecCreateSeq(PETSC_COMM_SELF,mdn->N,&mdn->lvec); CHKERRQ(ierr);
23 
24   /* Create temporary index set for building scatter gather */
25   ierr = ISCreateStride(PETSC_COMM_SELF,mdn->N,0,1,&tofrom); CHKERRQ(ierr);
26 
27   /* Create temporary global vector to generate scatter context */
28   n    = mdn->cowners[mdn->rank+1] - mdn->cowners[mdn->rank];
29   ierr = VecCreateMPI(mat->comm,n,mdn->N,&gvec); CHKERRQ(ierr);
30 
31   /* Generate the scatter context */
32   ierr = VecScatterCreate(gvec,tofrom,mdn->lvec,tofrom,&mdn->Mvctx); CHKERRQ(ierr);
33   PLogObjectParent(mat,mdn->Mvctx);
34   PLogObjectParent(mat,mdn->lvec);
35   PLogObjectParent(mat,tofrom);
36   PLogObjectParent(mat,gvec);
37 
38   ierr = ISDestroy(tofrom); CHKERRQ(ierr);
39   ierr = VecDestroy(gvec); CHKERRQ(ierr);
40   PetscFunctionReturn(0);
41 }
42 
43 extern int MatGetSubMatrices_MPIDense_Local(Mat,int,IS*,IS*,MatReuse,Mat*);
44 #undef __FUNC__
45 #define __FUNC__ "MatGetSubMatrices_MPIDense"
46 int MatGetSubMatrices_MPIDense(Mat C,int ismax,IS *isrow,IS *iscol,
47                              MatReuse scall,Mat **submat)
48 {
49   Mat_MPIDense  *c = (Mat_MPIDense *) C->data;
50   int           nmax,nstages_local,nstages,i,pos,max_no,ierr;
51 
52   PetscFunctionBegin;
53   /* Allocate memory to hold all the submatrices */
54   if (scall != MAT_REUSE_MATRIX) {
55     *submat = (Mat *)PetscMalloc((ismax+1)*sizeof(Mat));CHKPTRQ(*submat);
56   }
57   /* Determine the number of stages through which submatrices are done */
58   nmax          = 20*1000000 / (c->N * sizeof(int));
59   if (nmax == 0) nmax = 1;
60   nstages_local = ismax/nmax + ((ismax % nmax)?1:0);
61 
62   /* Make sure every processor loops through the nstages */
63   ierr = MPI_Allreduce(&nstages_local,&nstages,1,MPI_INT,MPI_MAX,C->comm);CHKERRQ(ierr);
64 
65 
66   for ( i=0,pos=0; i<nstages; i++ ) {
67     if (pos+nmax <= ismax) max_no = nmax;
68     else if (pos == ismax) max_no = 0;
69     else                   max_no = ismax-pos;
70     ierr = MatGetSubMatrices_MPIDense_Local(C,max_no,isrow+pos,iscol+pos,scall,*submat+pos);CHKERRQ(ierr);
71     pos += max_no;
72   }
73   PetscFunctionReturn(0);
74 }
75 /* -------------------------------------------------------------------------*/
76 #undef __FUNC__
77 #define __FUNC__ "MatGetSubMatrices_MPIDense_Local"
78 int MatGetSubMatrices_MPIDense_Local(Mat C,int ismax,IS *isrow,IS *iscol,MatReuse scall,Mat *submats)
79 {
80   Mat_MPIDense  *c = (Mat_MPIDense *) C->data;
81   Mat            A = c->A;
82   Mat_SeqDense  *a = (Mat_SeqDense*)A->data, *mat;
83   int           N = c->N, rstart = c->rstart,count;
84   int           **irow,**icol,*nrow,*ncol,*w1,*w2,*w3,*w4,*rtable,start,end,size;
85   int           **sbuf1, rank, m,i,j,k,l,ct1,ierr, **rbuf1,row,proc;
86   int           nrqs, msz, **ptr,index,*ctr,*pa,*tmp,bsz,nrqr;
87   int           is_no,jmax,*irow_i,**rmap,*rmap_i;
88   int           len,ctr_j,*sbuf1_j,*rbuf1_i;
89   int           tag0,tag1;
90   MPI_Request   *s_waits1,*r_waits1,*s_waits2,*r_waits2;
91   MPI_Status    *r_status1,*r_status2,*s_status1,*s_status2;
92   MPI_Comm      comm;
93   Scalar        **rbuf2,**sbuf2;
94 
95   PetscFunctionBegin;
96   comm   = C->comm;
97   tag0   = C->tag;
98   size   = c->size;
99   rank   = c->rank;
100   m      = c->M;
101 
102   /* Get some new tags to keep the communication clean */
103   ierr = PetscObjectGetNewTag((PetscObject)C,&tag1); CHKERRQ(ierr);
104 
105     /* Check if the col indices are sorted */
106   for ( i=0; i<ismax; i++ ) {
107     ierr = ISSorted(isrow[i],(PetscTruth*)&j);
108     if (!j) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"ISrow is not sorted");
109     ierr = ISSorted(iscol[i],(PetscTruth*)&j);
110     if (!j) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"IScol is not sorted");
111   }
112 
113   len    =  2*ismax*(sizeof(int *)+sizeof(int)) + (m+1)*sizeof(int);
114   irow   = (int **)PetscMalloc(len); CHKPTRQ(irow);
115   icol   = irow + ismax;
116   nrow   = (int *) (icol + ismax);
117   ncol   = nrow + ismax;
118   rtable = ncol + ismax;
119 
120   for ( i=0; i<ismax; i++ ) {
121     ierr = ISGetIndices(isrow[i],&irow[i]);  CHKERRQ(ierr);
122     ierr = ISGetIndices(iscol[i],&icol[i]);  CHKERRQ(ierr);
123     ierr = ISGetSize(isrow[i],&nrow[i]);  CHKERRQ(ierr);
124     ierr = ISGetSize(iscol[i],&ncol[i]);  CHKERRQ(ierr);
125   }
126 
127   /* Create hash table for the mapping :row -> proc*/
128   for ( i=0,j=0; i<size; i++ ) {
129     jmax = c->rowners[i+1];
130     for ( ; j<jmax; j++ ) {
131       rtable[j] = i;
132     }
133   }
134 
135   /* evaluate communication - mesg to who, length of mesg, and buffer space
136      required. Based on this, buffers are allocated, and data copied into them*/
137   w1     = (int *)PetscMalloc(size*4*sizeof(int));CHKPTRQ(w1); /* mesg size */
138   w2     = w1 + size;      /* if w2[i] marked, then a message to proc i*/
139   w3     = w2 + size;      /* no of IS that needs to be sent to proc i */
140   w4     = w3 + size;      /* temp work space used in determining w1, w2, w3 */
141   PetscMemzero(w1,size*3*sizeof(int)); /* initialize work vector*/
142   for ( i=0; i<ismax; i++ ) {
143     PetscMemzero(w4,size*sizeof(int)); /* initialize work vector*/
144     jmax   = nrow[i];
145     irow_i = irow[i];
146     for ( j=0; j<jmax; j++ ) {
147       row  = irow_i[j];
148       proc = rtable[row];
149       w4[proc]++;
150     }
151     for ( j=0; j<size; j++ ) {
152       if (w4[j]) { w1[j] += w4[j];  w3[j]++;}
153     }
154   }
155 
156   nrqs     = 0;              /* no of outgoing messages */
157   msz      = 0;              /* total mesg length (for all procs) */
158   w1[rank] = 0;              /* no mesg sent to self */
159   w3[rank] = 0;
160   for ( i=0; i<size; i++ ) {
161     if (w1[i])  { w2[i] = 1; nrqs++;} /* there exists a message to proc i */
162   }
163   pa = (int *)PetscMalloc((nrqs+1)*sizeof(int));CHKPTRQ(pa); /*(proc -array)*/
164   for ( i=0, j=0; i<size; i++ ) {
165     if (w1[i]) { pa[j] = i; j++; }
166   }
167 
168   /* Each message would have a header = 1 + 2*(no of IS) + data */
169   for ( i=0; i<nrqs; i++ ) {
170     j     = pa[i];
171     w1[j] += w2[j] + 2* w3[j];
172     msz   += w1[j];
173   }
174   /* Do a global reduction to determine how many messages to expect*/
175   {
176     int *rw1, *rw2;
177     rw1   = (int *)PetscMalloc(2*size*sizeof(int)); CHKPTRQ(rw1);
178     rw2   = rw1+size;
179     ierr  = MPI_Allreduce(w1, rw1, size, MPI_INT, MPI_MAX, comm);CHKERRQ(ierr);
180     bsz   = rw1[rank];
181     ierr  = MPI_Allreduce(w2, rw2, size, MPI_INT, MPI_SUM, comm);CHKERRQ(ierr);
182     nrqr  = rw2[rank];
183     PetscFree(rw1);
184   }
185 
186   /* Allocate memory for recv buffers . Prob none if nrqr = 0 ???? */
187   len      = (nrqr+1)*sizeof(int*) + nrqr*bsz*sizeof(int);
188   rbuf1    = (int**) PetscMalloc(len);  CHKPTRQ(rbuf1);
189   rbuf1[0] = (int *) (rbuf1 + nrqr);
190   for ( i=1; i<nrqr; ++i ) rbuf1[i] = rbuf1[i-1] + bsz;
191 
192   /* Post the receives */
193   r_waits1 = (MPI_Request *) PetscMalloc((nrqr+1)*sizeof(MPI_Request));CHKPTRQ(r_waits1);
194   for ( i=0; i<nrqr; ++i ) {
195     ierr = MPI_Irecv(rbuf1[i],bsz,MPI_INT,MPI_ANY_SOURCE,tag0,comm,r_waits1+i);CHKERRQ(ierr);
196   }
197 
198   /* Allocate Memory for outgoing messages */
199   len      = 2*size*sizeof(int*) + 2*msz*sizeof(int) + size*sizeof(int);
200   sbuf1    = (int **)PetscMalloc(len); CHKPTRQ(sbuf1);
201   ptr      = sbuf1 + size;   /* Pointers to the data in outgoing buffers */
202   PetscMemzero(sbuf1,2*size*sizeof(int*));
203   /* allocate memory for outgoing data + buf to receive the first reply */
204   tmp      = (int *) (ptr + size);
205   ctr      = tmp + 2*msz;
206 
207   {
208     int *iptr = tmp,ict = 0;
209     for ( i=0; i<nrqs; i++ ) {
210       j         = pa[i];
211       iptr     += ict;
212       sbuf1[j]  = iptr;
213       ict       = w1[j];
214     }
215   }
216 
217   /* Form the outgoing messages */
218   /* Initialize the header space */
219   for ( i=0; i<nrqs; i++ ) {
220     j           = pa[i];
221     sbuf1[j][0] = 0;
222     PetscMemzero(sbuf1[j]+1, 2*w3[j]*sizeof(int));
223     ptr[j]      = sbuf1[j] + 2*w3[j] + 1;
224   }
225 
226   /* Parse the isrow and copy data into outbuf */
227   for ( i=0; i<ismax; i++ ) {
228     PetscMemzero(ctr,size*sizeof(int));
229     irow_i = irow[i];
230     jmax   = nrow[i];
231     for ( j=0; j<jmax; j++ ) {  /* parse the indices of each IS */
232       row  = irow_i[j];
233       proc = rtable[row];
234       if (proc != rank) { /* copy to the outgoing buf*/
235         ctr[proc]++;
236         *ptr[proc] = row;
237         ptr[proc]++;
238       }
239     }
240     /* Update the headers for the current IS */
241     for ( j=0; j<size; j++ ) { /* Can Optimise this loop too */
242       if ((ctr_j = ctr[j])) {
243         sbuf1_j        = sbuf1[j];
244         k              = ++sbuf1_j[0];
245         sbuf1_j[2*k]   = ctr_j;
246         sbuf1_j[2*k-1] = i;
247       }
248     }
249   }
250 
251   /*  Now  post the sends */
252   s_waits1 = (MPI_Request *) PetscMalloc((nrqs+1)*sizeof(MPI_Request));CHKPTRQ(s_waits1);
253   for ( i=0; i<nrqs; ++i ) {
254     j = pa[i];
255     ierr = MPI_Isend( sbuf1[j], w1[j], MPI_INT, j, tag0, comm, s_waits1+i);CHKERRQ(ierr);
256   }
257 
258   /* Post recieves to capture the row_data from other procs */
259   r_waits2 = (MPI_Request *) PetscMalloc((nrqs+1)*sizeof(MPI_Request));CHKPTRQ(r_waits2);
260   rbuf2    = (Scalar**)PetscMalloc((nrqs+1)*sizeof(Scalar *));CHKPTRQ(rbuf2);
261   for ( i=0; i<nrqs; i++ ) {
262     j        = pa[i];
263     count    = (w1[j] - (2*sbuf1[j][0] + 1))*N;
264     rbuf2[i] = (Scalar *)PetscMalloc((count+1)*sizeof(Scalar)); CHKPTRQ(rbuf2[i]);
265     ierr     = MPI_Irecv( rbuf2[i], count, MPIU_SCALAR, j, tag1, comm, r_waits2+i);CHKERRQ(ierr);
266   }
267 
268   /* Receive messages(row_nos) and then, pack and send off the rowvalues
269      to the correct processors */
270 
271   s_waits2  = (MPI_Request *) PetscMalloc((nrqr+1)*sizeof(MPI_Request));CHKPTRQ(s_waits2);
272   r_status1 = (MPI_Status *) PetscMalloc((nrqr+1)*sizeof(MPI_Status));CHKPTRQ(r_status1);
273   sbuf2       = (Scalar **) PetscMalloc((nrqr+1)*sizeof(Scalar*));CHKPTRQ(sbuf2);
274 
275   {
276     Scalar *sbuf2_i,*v_start;
277     int    s_proc;
278     for ( i=0; i<nrqr; ++i ) {
279       ierr = MPI_Waitany(nrqr, r_waits1, &index, r_status1+i);CHKERRQ(ierr);
280       s_proc          = r_status1[i].MPI_SOURCE; /* send processor */
281       rbuf1_i         = rbuf1[index]; /* Actual message from s_proc */
282       /* no of rows = end - start; since start is array index[], 0index, whel end
283          is length of the buffer - which is 1index */
284       start           = 2*rbuf1_i[0] + 1;
285       MPI_Get_count(r_status1+i,MPI_INT, &end);
286       /* allocate memory sufficinet to hold all the row values */
287       sbuf2[index] = (Scalar *)PetscMalloc((end-start)*N*sizeof(Scalar));CHKPTRQ(sbuf2[index]);
288       sbuf2_i      = sbuf2[index];
289       /* Now pack the data */
290       for ( j=start; j<end; j++ ) {
291         row = rbuf1_i[j] - rstart;
292         v_start = a->v + row;
293         for ( k=0; k<N; k++ ) {
294           sbuf2_i[0] = v_start[0];
295           sbuf2_i++; v_start+=a->m;
296         }
297       }
298       /* Now send off the data */
299       ierr = MPI_Isend(sbuf2[index],(end-start)*N,MPIU_SCALAR,s_proc,tag1,comm,s_waits2+i); CHKERRQ(ierr);
300     }
301   }
302   /* End Send-Recv of IS + row_numbers */
303   PetscFree(r_status1); PetscFree(r_waits1);
304   s_status1 = (MPI_Status *) PetscMalloc((nrqs+1)*sizeof(MPI_Status));CHKPTRQ(s_status1);
305   ierr      = MPI_Waitall(nrqs,s_waits1,s_status1);CHKERRQ(ierr);
306   PetscFree(s_status1); PetscFree(s_waits1);
307 
308   /* Create the submatrices */
309   if (scall == MAT_REUSE_MATRIX) {
310     for ( i=0; i<ismax; i++ ) {
311       mat = (Mat_SeqDense *)(submats[i]->data);
312       if ((mat->m != nrow[i]) || (mat->n != ncol[i])) {
313         SETERRQ(PETSC_ERR_ARG_SIZ,0,"Cannot reuse matrix. wrong size");
314       }
315       PetscMemzero(mat->v,mat->m*mat->n*sizeof(Scalar));
316       submats[i]->factor = C->factor;
317     }
318   } else {
319     for ( i=0; i<ismax; i++ ) {
320       ierr = MatCreateSeqDense(PETSC_COMM_SELF,nrow[i],ncol[i],PETSC_NULL,submats+i);CHKERRQ(ierr);
321     }
322   }
323 
324   /* Assemble the matrices */
325   {
326     int    col;
327     Scalar *imat_v,*mat_v,*imat_vi,*mat_vi;
328 
329     for ( i=0; i<ismax; i++ ) {
330       mat       = (Mat_SeqDense *) submats[i]->data;
331       mat_v     = a->v;
332       imat_v    = mat->v;
333       irow_i    = irow[i];
334       m         = nrow[i];
335       for ( j=0; j<m; j++ ) {
336         row      = irow_i[j] ;
337         proc     = rtable[row];
338         if (proc == rank) {
339           row      = row - rstart;
340           mat_vi   = mat_v + row;
341           imat_vi  = imat_v + j;
342           for ( k=0; k<ncol[i]; k++ ) {
343             col = icol[i][k];
344             imat_vi[k*m] = mat_vi[col*a->m];
345           }
346         }
347       }
348     }
349   }
350 
351   /* Create row map. This maps c->row to submat->row for each submat*/
352   /* this is a very expensive operation wrt memory usage */
353   len     = (1+ismax)*sizeof(int *) + ismax*c->M*sizeof(int);
354   rmap    = (int **)PetscMalloc(len); CHKPTRQ(rmap);
355   rmap[0] = (int *)(rmap + ismax);
356   PetscMemzero(rmap[0],ismax*c->M*sizeof(int));
357   for ( i=1; i<ismax; i++ ) { rmap[i] = rmap[i-1] + c->M;}
358   for ( i=0; i<ismax; i++ ) {
359     rmap_i = rmap[i];
360     irow_i = irow[i];
361     jmax   = nrow[i];
362     for ( j=0; j<jmax; j++ ) {
363       rmap_i[irow_i[j]] = j;
364     }
365   }
366 
367   /* Now Receive the row_values and assemble the rest of the matrix */
368 
369   r_status2 = (MPI_Status *) PetscMalloc((nrqs+1)*sizeof(MPI_Status));CHKPTRQ(r_status2);
370 
371   {
372     int    is_max,tmp1,col,*sbuf1_i,is_sz;
373     Scalar *rbuf2_i,*imat_v,*imat_vi;
374 
375     for ( tmp1=0; tmp1<nrqs; tmp1++ ) { /* For each message */
376       ierr = MPI_Waitany(nrqs,r_waits2,&i,r_status2+tmp1);CHKERRQ(ierr);
377       /* Now dig out the corresponding sbuf1, which contains the IS data_structure */
378       sbuf1_i = sbuf1[pa[i]];
379       is_max  = sbuf1_i[0];
380       ct1     = 2*is_max+1;
381       rbuf2_i = rbuf2[i];
382       for ( j=1; j<=is_max; j++ ) { /* For each IS belonging to the message */
383         is_no     = sbuf1_i[2*j-1];
384         is_sz     = sbuf1_i[2*j];
385         mat       = (Mat_SeqDense *) submats[is_no]->data;
386         imat_v    = mat->v;
387         rmap_i    = rmap[is_no];
388         m         = nrow[is_no];
389         for ( k=0; k<is_sz; k++,rbuf2_i+=N) {  /* For each row */
390           row      = sbuf1_i[ct1]; ct1++;
391           row      = rmap_i[row];
392           imat_vi  = imat_v + row;
393           for ( l=0; l<ncol[is_no]; l++ ) { /* For each col */
394             col = icol[is_no][l];
395             imat_vi[l*m] = rbuf2_i[col];
396           }
397         }
398       }
399     }
400   }
401   /* End Send-Recv of row_values */
402   PetscFree(r_status2); PetscFree(r_waits2);
403   s_status2 = (MPI_Status *) PetscMalloc((nrqr+1)*sizeof(MPI_Status));CHKPTRQ(s_status1);
404   ierr      = MPI_Waitall(nrqr,s_waits2,s_status2); CHKERRQ(ierr);
405   PetscFree(s_status2); PetscFree(s_waits2);
406 
407   /* Restore the indices */
408   for ( i=0; i<ismax; i++ ) {
409     ierr = ISRestoreIndices(isrow[i], irow+i); CHKERRQ(ierr);
410     ierr = ISRestoreIndices(iscol[i], icol+i); CHKERRQ(ierr);
411   }
412 
413   /* Destroy allocated memory */
414   PetscFree(irow);
415   PetscFree(w1);
416   PetscFree(pa);
417 
418 
419   for ( i=0; i<nrqs; ++i ) {
420     PetscFree(rbuf2[i]);
421   }
422   PetscFree(rbuf2);
423   PetscFree(sbuf1);
424   PetscFree(rbuf1);
425 
426   for ( i=0; i<nrqr; ++i ) {
427     PetscFree(sbuf2[i]);
428   }
429 
430   PetscFree(sbuf2);
431   PetscFree(rmap);
432 
433   ierr = PetscObjectRestoreNewTag((PetscObject)C,&tag1); CHKERRQ(ierr);
434 
435   for ( i=0; i<ismax; i++ ) {
436     ierr = MatAssemblyBegin(submats[i], MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
437     ierr = MatAssemblyEnd(submats[i], MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
438   }
439 
440   PetscFunctionReturn(0);
441 }
442 
443