Actual source code: mmdense.c

  1: #define PETSCMAT_DLL

  3: /*
  4:    Support for the parallel dense matrix vector multiply
  5: */
 6:  #include ../src/mat/impls/dense/mpi/mpidense.h

 10: PetscErrorCode MatSetUpMultiply_MPIDense(Mat mat)
 11: {
 12:   Mat_MPIDense *mdn = (Mat_MPIDense*)mat->data;
 14:   IS           from,to;
 15:   Vec          gvec;

 18:   /* Create local vector that is used to scatter into */
 19:   VecCreateSeq(PETSC_COMM_SELF,mat->cmap->N,&mdn->lvec);

 21:   /* Create temporary index set for building scatter gather */
 22:   ISCreateStride(((PetscObject)mat)->comm,mat->cmap->N,0,1,&from);
 23:   ISCreateStride(PETSC_COMM_SELF,mat->cmap->N,0,1,&to);

 25:   /* Create temporary global vector to generate scatter context */
 26:   /* n    = mdn->cowners[mdn->rank+1] - mdn->cowners[mdn->rank]; */

 28:   VecCreateMPIWithArray(((PetscObject)mat)->comm,mdn->nvec,mat->cmap->N,PETSC_NULL,&gvec);

 30:   /* Generate the scatter context */
 31:   VecScatterCreate(gvec,from,mdn->lvec,to,&mdn->Mvctx);
 32:   PetscLogObjectParent(mat,mdn->Mvctx);
 33:   PetscLogObjectParent(mat,mdn->lvec);
 34:   PetscLogObjectParent(mat,from);
 35:   PetscLogObjectParent(mat,to);
 36:   PetscLogObjectParent(mat,gvec);

 38:   ISDestroy(to);
 39:   ISDestroy(from);
 40:   VecDestroy(gvec);
 41:   return(0);
 42: }

 44: EXTERN PetscErrorCode MatGetSubMatrices_MPIDense_Local(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat*);
 47: PetscErrorCode MatGetSubMatrices_MPIDense(Mat C,PetscInt ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submat[])
 48: {
 50:   PetscInt           nmax,nstages_local,nstages,i,pos,max_no;

 53:   /* Allocate memory to hold all the submatrices */
 54:   if (scall != MAT_REUSE_MATRIX) {
 55:     PetscMalloc((ismax+1)*sizeof(Mat),submat);
 56:   }
 57:   /* Determine the number of stages through which submatrices are done */
 58:   nmax          = 20*1000000 / (C->cmap->N * sizeof(PetscInt));
 59:   if (!nmax) nmax = 1;
 60:   nstages_local = ismax/nmax + ((ismax % nmax)?1:0);

 62:   /* Make sure every processor loops through the nstages */
 63:   MPI_Allreduce(&nstages_local,&nstages,1,MPIU_INT,MPI_MAX,((PetscObject)C)->comm);


 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:     MatGetSubMatrices_MPIDense_Local(C,max_no,isrow+pos,iscol+pos,scall,*submat+pos);
 71:     pos += max_no;
 72:   }
 73:   return(0);
 74: }
 75: /* -------------------------------------------------------------------------*/
 78: PetscErrorCode MatGetSubMatrices_MPIDense_Local(Mat C,PetscInt ismax,const IS isrow[],const 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;
 84:   PetscMPIInt    rank,size,tag0,tag1,idex,end,i;
 85:   PetscInt       N = C->cmap->N,rstart = C->rmap->rstart,count;
 86:   const PetscInt **irow,**icol,*irow_i;
 87:   PetscInt       *nrow,*ncol,*w1,*w3,*w4,*rtable,start;
 88:   PetscInt       **sbuf1,m,j,k,l,ct1,**rbuf1,row,proc;
 89:   PetscInt       nrqs,msz,**ptr,*ctr,*pa,*tmp,bsz,nrqr;
 90:   PetscInt       is_no,jmax,**rmap,*rmap_i;
 91:   PetscInt       ctr_j,*sbuf1_j,*rbuf1_i;
 92:   MPI_Request    *s_waits1,*r_waits1,*s_waits2,*r_waits2;
 93:   MPI_Status     *r_status1,*r_status2,*s_status1,*s_status2;
 94:   MPI_Comm       comm;
 95:   PetscScalar    **rbuf2,**sbuf2;
 96:   PetscTruth     sorted;

 99:   comm   = ((PetscObject)C)->comm;
100:   tag0   = ((PetscObject)C)->tag;
101:   size   = c->size;
102:   rank   = c->rank;
103:   m      = C->rmap->N;
104: 
105:   /* Get some new tags to keep the communication clean */
106:   PetscObjectGetNewTag((PetscObject)C,&tag1);

108:     /* Check if the col indices are sorted */
109:   for (i=0; i<ismax; i++) {
110:     ISSorted(isrow[i],&sorted);
111:     if (!sorted) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"ISrow is not sorted");
112:     ISSorted(iscol[i],&sorted);
113:     if (!sorted) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IScol is not sorted");
114:   }

116:   PetscMalloc(2*ismax*(sizeof(PetscInt*)+sizeof(PetscInt)) + (m+1)*sizeof(PetscInt),&irow);
117:   icol   = irow + ismax;
118:   nrow   = (PetscInt*)(icol + ismax);
119:   ncol   = nrow + ismax;
120:   rtable = ncol + ismax;

122:   for (i=0; i<ismax; i++) {
123:     ISGetIndices(isrow[i],&irow[i]);
124:     ISGetIndices(iscol[i],&icol[i]);
125:     ISGetLocalSize(isrow[i],&nrow[i]);
126:     ISGetLocalSize(iscol[i],&ncol[i]);
127:   }

129:   /* Create hash table for the mapping :row -> proc*/
130:   for (i=0,j=0; i<size; i++) {
131:     jmax = C->rmap->range[i+1];
132:     for (; j<jmax; j++) {
133:       rtable[j] = i;
134:     }
135:   }

137:   /* evaluate communication - mesg to who,length of mesg, and buffer space
138:      required. Based on this, buffers are allocated, and data copied into them*/
139:   PetscMalloc(size*4*sizeof(PetscInt),&w1); /* mesg size */
140:   w3     = w1 + 2*size;      /* no of IS that needs to be sent to proc i */
141:   w4     = w3 + size;      /* temp work space used in determining w1,  w3 */
142:   PetscMemzero(w1,size*3*sizeof(PetscInt)); /* initialize work vector*/
143:   for (i=0; i<ismax; i++) {
144:     PetscMemzero(w4,size*sizeof(PetscInt)); /* initialize work vector*/
145:     jmax   = nrow[i];
146:     irow_i = irow[i];
147:     for (j=0; j<jmax; j++) {
148:       row  = irow_i[j];
149:       proc = rtable[row];
150:       w4[proc]++;
151:     }
152:     for (j=0; j<size; j++) {
153:       if (w4[j]) { w1[2*j] += w4[j];  w3[j]++;}
154:     }
155:   }
156: 
157:   nrqs       = 0;              /* no of outgoing messages */
158:   msz        = 0;              /* total mesg length (for all procs) */
159:   w1[2*rank] = 0;              /* no mesg sent to self */
160:   w3[rank]   = 0;
161:   for (i=0; i<size; i++) {
162:     if (w1[2*i])  { w1[2*i+1] = 1; nrqs++;} /* there exists a message to proc i */
163:   }
164:   PetscMalloc((nrqs+1)*sizeof(PetscInt),&pa); /*(proc -array)*/
165:   for (i=0,j=0; i<size; i++) {
166:     if (w1[2*i]) { pa[j] = i; j++; }
167:   }

169:   /* Each message would have a header = 1 + 2*(no of IS) + data */
170:   for (i=0; i<nrqs; i++) {
171:     j       = pa[i];
172:     w1[2*j] += w1[2*j+1] + 2* w3[j];
173:     msz     += w1[2*j];
174:   }
175:   /* Do a global reduction to determine how many messages to expect*/
176:   PetscMaxSum(comm,w1,&bsz,&nrqr);

178:   /* Allocate memory for recv buffers . Prob none if nrqr = 0 ???? */
179:   PetscMalloc((nrqr+1)*sizeof(PetscInt*) + nrqr*bsz*sizeof(PetscInt),&rbuf1);
180:   rbuf1[0] = (PetscInt*)(rbuf1 + nrqr);
181:   for (i=1; i<nrqr; ++i) rbuf1[i] = rbuf1[i-1] + bsz;
182: 
183:   /* Post the receives */
184:   PetscMalloc((nrqr+1)*sizeof(MPI_Request),&r_waits1);
185:   for (i=0; i<nrqr; ++i) {
186:     MPI_Irecv(rbuf1[i],bsz,MPIU_INT,MPI_ANY_SOURCE,tag0,comm,r_waits1+i);
187:   }

189:   /* Allocate Memory for outgoing messages */
190:   PetscMalloc(2*size*sizeof(PetscInt*) + 2*msz*sizeof(PetscInt)+ size*sizeof(PetscInt),&sbuf1);
191:   ptr   = sbuf1 + size;   /* Pointers to the data in outgoing buffers */
192:   PetscMemzero(sbuf1,2*size*sizeof(PetscInt*));
193:   /* allocate memory for outgoing data + buf to receive the first reply */
194:   tmp   = (PetscInt*)(ptr + size);
195:   ctr   = tmp + 2*msz;

197:   {
198:     PetscInt *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:   }

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:     PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(PetscInt));
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:     PetscMemzero(ctr,size*sizeof(PetscInt));
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:   }

241:   /*  Now  post the sends */
242:   PetscMalloc((nrqs+1)*sizeof(MPI_Request),&s_waits1);
243:   for (i=0; i<nrqs; ++i) {
244:     j = pa[i];
245:     MPI_Isend(sbuf1[j],w1[2*j],MPIU_INT,j,tag0,comm,s_waits1+i);
246:   }

248:   /* Post recieves to capture the row_data from other procs */
249:   PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits2);
250:   PetscMalloc((nrqs+1)*sizeof(PetscScalar*),&rbuf2);
251:   for (i=0; i<nrqs; i++) {
252:     j        = pa[i];
253:     count    = (w1[2*j] - (2*sbuf1[j][0] + 1))*N;
254:     PetscMalloc((count+1)*sizeof(PetscScalar),&rbuf2[i]);
255:     MPI_Irecv(rbuf2[i],count,MPIU_SCALAR,j,tag1,comm,r_waits2+i);
256:   }

258:   /* Receive messages(row_nos) and then, pack and send off the rowvalues
259:      to the correct processors */

261:   PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits2);
262:   PetscMalloc((nrqr+1)*sizeof(MPI_Status),&r_status1);
263:   PetscMalloc((nrqr+1)*sizeof(PetscScalar*),&sbuf2);
264: 
265:   {
266:     PetscScalar *sbuf2_i,*v_start;
267:     PetscInt         s_proc;
268:     for (i=0; i<nrqr; ++i) {
269:       MPI_Waitany(nrqr,r_waits1,&idex,r_status1+i);
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:       MPI_Get_count(r_status1+i,MPIU_INT,&end);
276:       /* allocate memory sufficinet to hold all the row values */
277:       PetscMalloc((end-start)*N*sizeof(PetscScalar),&sbuf2[idex]);
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->rmap->n;
286:         }
287:       }
288:       /* Now send off the data */
289:       MPI_Isend(sbuf2[idex],(end-start)*N,MPIU_SCALAR,s_proc,tag1,comm,s_waits2+i);
290:     }
291:   }
292:   /* End Send-Recv of IS + row_numbers */
293:   PetscFree(r_status1);
294:   PetscFree(r_waits1);
295:   PetscMalloc((nrqs+1)*sizeof(MPI_Status),&s_status1);
296:   if (nrqs) {MPI_Waitall(nrqs,s_waits1,s_status1);}
297:   PetscFree(s_status1);
298:   PetscFree(s_waits1);

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]->rmap->n != nrow[i]) || (submats[i]->cmap->n != ncol[i])) {
305:         SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size");
306:       }
307:       PetscMemzero(mat->v,submats[i]->rmap->n*submats[i]->cmap->n*sizeof(PetscScalar));
308:       submats[i]->factor = C->factor;
309:     }
310:   } else {
311:     for (i=0; i<ismax; i++) {
312:       MatCreate(PETSC_COMM_SELF,submats+i);
313:       MatSetSizes(submats[i],nrow[i],ncol[i],nrow[i],ncol[i]);
314:       MatSetType(submats[i],((PetscObject)A)->type_name);
315:       MatSeqDenseSetPreallocation(submats[i],PETSC_NULL);
316:     }
317:   }
318: 
319:   /* Assemble the matrices */
320:   {
321:     PetscInt         col;
322:     PetscScalar *imat_v,*mat_v,*imat_vi,*mat_vi;
323: 
324:     for (i=0; i<ismax; i++) {
325:       mat       = (Mat_SeqDense*)submats[i]->data;
326:       mat_v     = a->v;
327:       imat_v    = mat->v;
328:       irow_i    = irow[i];
329:       m         = nrow[i];
330:       for (j=0; j<m; j++) {
331:         row      = irow_i[j] ;
332:         proc     = rtable[row];
333:         if (proc == rank) {
334:           row      = row - rstart;
335:           mat_vi   = mat_v + row;
336:           imat_vi  = imat_v + j;
337:           for (k=0; k<ncol[i]; k++) {
338:             col = icol[i][k];
339:             imat_vi[k*m] = mat_vi[col*C->rmap->n];
340:           }
341:         }
342:       }
343:     }
344:   }

346:   /* Create row map-> This maps c->row to submat->row for each submat*/
347:   /* this is a very expensive operation wrt memory usage */
348:   PetscMalloc((1+ismax)*sizeof(PetscInt*)+ ismax*C->rmap->N*sizeof(PetscInt),&rmap);
349:   rmap[0] = (PetscInt*)(rmap + ismax);
350:   PetscMemzero(rmap[0],ismax*C->rmap->N*sizeof(PetscInt));
351:   for (i=1; i<ismax; i++) { rmap[i] = rmap[i-1] + C->rmap->N;}
352:   for (i=0; i<ismax; i++) {
353:     rmap_i = rmap[i];
354:     irow_i = irow[i];
355:     jmax   = nrow[i];
356:     for (j=0; j<jmax; j++) {
357:       rmap_i[irow_i[j]] = j;
358:     }
359:   }
360: 
361:   /* Now Receive the row_values and assemble the rest of the matrix */
362:   PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status2);

364:   {
365:     PetscInt    is_max,tmp1,col,*sbuf1_i,is_sz;
366:     PetscScalar *rbuf2_i,*imat_v,*imat_vi;
367: 
368:     for (tmp1=0; tmp1<nrqs; tmp1++) { /* For each message */
369:       MPI_Waitany(nrqs,r_waits2,&i,r_status2+tmp1);
370:       /* Now dig out the corresponding sbuf1, which contains the IS data_structure */
371:       sbuf1_i = sbuf1[pa[i]];
372:       is_max  = sbuf1_i[0];
373:       ct1     = 2*is_max+1;
374:       rbuf2_i = rbuf2[i];
375:       for (j=1; j<=is_max; j++) { /* For each IS belonging to the message */
376:         is_no     = sbuf1_i[2*j-1];
377:         is_sz     = sbuf1_i[2*j];
378:         mat       = (Mat_SeqDense*)submats[is_no]->data;
379:         imat_v    = mat->v;
380:         rmap_i    = rmap[is_no];
381:         m         = nrow[is_no];
382:         for (k=0; k<is_sz; k++,rbuf2_i+=N) {  /* For each row */
383:           row      = sbuf1_i[ct1]; ct1++;
384:           row      = rmap_i[row];
385:           imat_vi  = imat_v + row;
386:           for (l=0; l<ncol[is_no]; l++) { /* For each col */
387:             col = icol[is_no][l];
388:             imat_vi[l*m] = rbuf2_i[col];
389:           }
390:         }
391:       }
392:     }
393:   }
394:   /* End Send-Recv of row_values */
395:   PetscFree(r_status2);
396:   PetscFree(r_waits2);
397:   PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status2);
398:   if (nrqr) {MPI_Waitall(nrqr,s_waits2,s_status2);}
399:   PetscFree(s_status2);
400:   PetscFree(s_waits2);

402:   /* Restore the indices */
403:   for (i=0; i<ismax; i++) {
404:     ISRestoreIndices(isrow[i],irow+i);
405:     ISRestoreIndices(iscol[i],icol+i);
406:   }

408:   /* Destroy allocated memory */
409:   PetscFree(irow);
410:   PetscFree(w1);
411:   PetscFree(pa);


414:   for (i=0; i<nrqs; ++i) {
415:     PetscFree(rbuf2[i]);
416:   }
417:   PetscFree(rbuf2);
418:   PetscFree(sbuf1);
419:   PetscFree(rbuf1);

421:   for (i=0; i<nrqr; ++i) {
422:     PetscFree(sbuf2[i]);
423:   }

425:   PetscFree(sbuf2);
426:   PetscFree(rmap);

428:   for (i=0; i<ismax; i++) {
429:     MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);
430:     MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);
431:   }

433:   return(0);
434: }