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