Actual source code: mmdense.c

  2: /*
  3:    Support for the parallel dense matrix vector multiply
  4: */
  5: #include <../src/mat/impls/dense/mpi/mpidense.h>
  6: #include <petscblaslapack.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: }

 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:   PetscBool      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_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"ISrow is not sorted");
112:     ISSorted(iscol[i],&sorted);
113:     if (!sorted) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"IScol is not sorted");
114:   }

116:   PetscMalloc5(ismax,const PetscInt*,&irow,ismax,const PetscInt*,&icol,ismax,PetscInt,&nrow,ismax,PetscInt,&ncol,m,PetscInt,&rtable);
117:   for (i=0; i<ismax; i++) {
118:     ISGetIndices(isrow[i],&irow[i]);
119:     ISGetIndices(iscol[i],&icol[i]);
120:     ISGetLocalSize(isrow[i],&nrow[i]);
121:     ISGetLocalSize(iscol[i],&ncol[i]);
122:   }

124:   /* Create hash table for the mapping :row -> proc*/
125:   for (i=0,j=0; i<size; i++) {
126:     jmax = C->rmap->range[i+1];
127:     for (; j<jmax; j++) {
128:       rtable[j] = i;
129:     }
130:   }

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

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

172:   /* Allocate memory for recv buffers . Make sure rbuf1[0] exists by adding 1 to the buffer length */
173:   PetscMalloc((nrqr+1)*sizeof(PetscInt*),&rbuf1);
174:   PetscMalloc(nrqr*bsz*sizeof(PetscInt),&rbuf1[0]);
175:   for (i=1; i<nrqr; ++i) rbuf1[i] = rbuf1[i-1] + bsz;
176: 
177:   /* Post the receives */
178:   PetscMalloc((nrqr+1)*sizeof(MPI_Request),&r_waits1);
179:   for (i=0; i<nrqr; ++i) {
180:     MPI_Irecv(rbuf1[i],bsz,MPIU_INT,MPI_ANY_SOURCE,tag0,comm,r_waits1+i);
181:   }

183:   /* Allocate Memory for outgoing messages */
184:   PetscMalloc4(size,PetscInt*,&sbuf1,size,PetscInt*,&ptr,2*msz,PetscInt,&tmp,size,PetscInt,&ctr);
185:   PetscMemzero(sbuf1,size*sizeof(PetscInt*));
186:   PetscMemzero(ptr,size*sizeof(PetscInt*));
187:   {
188:     PetscInt *iptr = tmp,ict = 0;
189:     for (i=0; i<nrqs; i++) {
190:       j         = pa[i];
191:       iptr     += ict;
192:       sbuf1[j]  = iptr;
193:       ict       = w1[2*j];
194:     }
195:   }

197:   /* Form the outgoing messages */
198:   /* Initialize the header space */
199:   for (i=0; i<nrqs; i++) {
200:     j           = pa[i];
201:     sbuf1[j][0] = 0;
202:     PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(PetscInt));
203:     ptr[j]      = sbuf1[j] + 2*w3[j] + 1;
204:   }
205: 
206:   /* Parse the isrow and copy data into outbuf */
207:   for (i=0; i<ismax; i++) {
208:     PetscMemzero(ctr,size*sizeof(PetscInt));
209:     irow_i = irow[i];
210:     jmax   = nrow[i];
211:     for (j=0; j<jmax; j++) {  /* parse the indices of each IS */
212:       row  = irow_i[j];
213:       proc = rtable[row];
214:       if (proc != rank) { /* copy to the outgoing buf*/
215:         ctr[proc]++;
216:         *ptr[proc] = row;
217:         ptr[proc]++;
218:       }
219:     }
220:     /* Update the headers for the current IS */
221:     for (j=0; j<size; j++) { /* Can Optimise this loop too */
222:       if ((ctr_j = ctr[j])) {
223:         sbuf1_j        = sbuf1[j];
224:         k              = ++sbuf1_j[0];
225:         sbuf1_j[2*k]   = ctr_j;
226:         sbuf1_j[2*k-1] = i;
227:       }
228:     }
229:   }

231:   /*  Now  post the sends */
232:   PetscMalloc((nrqs+1)*sizeof(MPI_Request),&s_waits1);
233:   for (i=0; i<nrqs; ++i) {
234:     j = pa[i];
235:     MPI_Isend(sbuf1[j],w1[2*j],MPIU_INT,j,tag0,comm,s_waits1+i);
236:   }

238:   /* Post recieves to capture the row_data from other procs */
239:   PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits2);
240:   PetscMalloc((nrqs+1)*sizeof(PetscScalar*),&rbuf2);
241:   for (i=0; i<nrqs; i++) {
242:     j        = pa[i];
243:     count    = (w1[2*j] - (2*sbuf1[j][0] + 1))*N;
244:     PetscMalloc((count+1)*sizeof(PetscScalar),&rbuf2[i]);
245:     MPI_Irecv(rbuf2[i],count,MPIU_SCALAR,j,tag1,comm,r_waits2+i);
246:   }

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

251:   PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits2);
252:   PetscMalloc((nrqr+1)*sizeof(MPI_Status),&r_status1);
253:   PetscMalloc((nrqr+1)*sizeof(PetscScalar*),&sbuf2);
254: 
255:   {
256:     PetscScalar *sbuf2_i,*v_start;
257:     PetscInt         s_proc;
258:     for (i=0; i<nrqr; ++i) {
259:       MPI_Waitany(nrqr,r_waits1,&idex,r_status1+i);
260:       s_proc          = r_status1[i].MPI_SOURCE; /* send processor */
261:       rbuf1_i         = rbuf1[idex]; /* Actual message from s_proc */
262:       /* no of rows = end - start; since start is array idex[], 0idex, whel end
263:          is length of the buffer - which is 1idex */
264:       start           = 2*rbuf1_i[0] + 1;
265:       MPI_Get_count(r_status1+i,MPIU_INT,&end);
266:       /* allocate memory sufficinet to hold all the row values */
267:       PetscMalloc((end-start)*N*sizeof(PetscScalar),&sbuf2[idex]);
268:       sbuf2_i      = sbuf2[idex];
269:       /* Now pack the data */
270:       for (j=start; j<end; j++) {
271:         row = rbuf1_i[j] - rstart;
272:         v_start = a->v + row;
273:         for (k=0; k<N; k++) {
274:           sbuf2_i[0] = v_start[0];
275:           sbuf2_i++; v_start += C->rmap->n;
276:         }
277:       }
278:       /* Now send off the data */
279:       MPI_Isend(sbuf2[idex],(end-start)*N,MPIU_SCALAR,s_proc,tag1,comm,s_waits2+i);
280:     }
281:   }
282:   /* End Send-Recv of IS + row_numbers */
283:   PetscFree(r_status1);
284:   PetscFree(r_waits1);
285:   PetscMalloc((nrqs+1)*sizeof(MPI_Status),&s_status1);
286:   if (nrqs) {MPI_Waitall(nrqs,s_waits1,s_status1);}
287:   PetscFree(s_status1);
288:   PetscFree(s_waits1);

290:   /* Create the submatrices */
291:   if (scall == MAT_REUSE_MATRIX) {
292:     for (i=0; i<ismax; i++) {
293:       mat = (Mat_SeqDense *)(submats[i]->data);
294:       if ((submats[i]->rmap->n != nrow[i]) || (submats[i]->cmap->n != ncol[i])) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size");
295:       PetscMemzero(mat->v,submats[i]->rmap->n*submats[i]->cmap->n*sizeof(PetscScalar));
296:       submats[i]->factortype = C->factortype;
297:     }
298:   } else {
299:     for (i=0; i<ismax; i++) {
300:       MatCreate(PETSC_COMM_SELF,submats+i);
301:       MatSetSizes(submats[i],nrow[i],ncol[i],nrow[i],ncol[i]);
302:       MatSetType(submats[i],((PetscObject)A)->type_name);
303:       MatSeqDenseSetPreallocation(submats[i],PETSC_NULL);
304:     }
305:   }
306: 
307:   /* Assemble the matrices */
308:   {
309:     PetscInt         col;
310:     PetscScalar *imat_v,*mat_v,*imat_vi,*mat_vi;
311: 
312:     for (i=0; i<ismax; i++) {
313:       mat       = (Mat_SeqDense*)submats[i]->data;
314:       mat_v     = a->v;
315:       imat_v    = mat->v;
316:       irow_i    = irow[i];
317:       m         = nrow[i];
318:       for (j=0; j<m; j++) {
319:         row      = irow_i[j] ;
320:         proc     = rtable[row];
321:         if (proc == rank) {
322:           row      = row - rstart;
323:           mat_vi   = mat_v + row;
324:           imat_vi  = imat_v + j;
325:           for (k=0; k<ncol[i]; k++) {
326:             col = icol[i][k];
327:             imat_vi[k*m] = mat_vi[col*C->rmap->n];
328:           }
329:         }
330:       }
331:     }
332:   }

334:   /* Create row map-> This maps c->row to submat->row for each submat*/
335:   /* this is a very expensive operation wrt memory usage */
336:   PetscMalloc(ismax*sizeof(PetscInt*),&rmap);
337:   PetscMalloc(ismax*C->rmap->N*sizeof(PetscInt),&rmap[0]);
338:   PetscMemzero(rmap[0],ismax*C->rmap->N*sizeof(PetscInt));
339:   for (i=1; i<ismax; i++) { rmap[i] = rmap[i-1] + C->rmap->N;}
340:   for (i=0; i<ismax; i++) {
341:     rmap_i = rmap[i];
342:     irow_i = irow[i];
343:     jmax   = nrow[i];
344:     for (j=0; j<jmax; j++) {
345:       rmap_i[irow_i[j]] = j;
346:     }
347:   }
348: 
349:   /* Now Receive the row_values and assemble the rest of the matrix */
350:   PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status2);
351:   {
352:     PetscInt    is_max,tmp1,col,*sbuf1_i,is_sz;
353:     PetscScalar *rbuf2_i,*imat_v,*imat_vi;
354: 
355:     for (tmp1=0; tmp1<nrqs; tmp1++) { /* For each message */
356:       MPI_Waitany(nrqs,r_waits2,&i,r_status2+tmp1);
357:       /* Now dig out the corresponding sbuf1, which contains the IS data_structure */
358:       sbuf1_i = sbuf1[pa[i]];
359:       is_max  = sbuf1_i[0];
360:       ct1     = 2*is_max+1;
361:       rbuf2_i = rbuf2[i];
362:       for (j=1; j<=is_max; j++) { /* For each IS belonging to the message */
363:         is_no     = sbuf1_i[2*j-1];
364:         is_sz     = sbuf1_i[2*j];
365:         mat       = (Mat_SeqDense*)submats[is_no]->data;
366:         imat_v    = mat->v;
367:         rmap_i    = rmap[is_no];
368:         m         = nrow[is_no];
369:         for (k=0; k<is_sz; k++,rbuf2_i+=N) {  /* For each row */
370:           row      = sbuf1_i[ct1]; ct1++;
371:           row      = rmap_i[row];
372:           imat_vi  = imat_v + row;
373:           for (l=0; l<ncol[is_no]; l++) { /* For each col */
374:             col = icol[is_no][l];
375:             imat_vi[l*m] = rbuf2_i[col];
376:           }
377:         }
378:       }
379:     }
380:   }
381:   /* End Send-Recv of row_values */
382:   PetscFree(r_status2);
383:   PetscFree(r_waits2);
384:   PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status2);
385:   if (nrqr) {MPI_Waitall(nrqr,s_waits2,s_status2);}
386:   PetscFree(s_status2);
387:   PetscFree(s_waits2);

389:   /* Restore the indices */
390:   for (i=0; i<ismax; i++) {
391:     ISRestoreIndices(isrow[i],irow+i);
392:     ISRestoreIndices(iscol[i],icol+i);
393:   }

395:   /* Destroy allocated memory */
396:   PetscFree5(irow,icol,nrow,ncol,rtable);
397:   PetscFree3(w1,w3,w4);
398:   PetscFree(pa);

400:   for (i=0; i<nrqs; ++i) {
401:     PetscFree(rbuf2[i]);
402:   }
403:   PetscFree(rbuf2);
404:   PetscFree4(sbuf1,ptr,tmp,ctr);
405:   PetscFree(rbuf1[0]);
406:   PetscFree(rbuf1);

408:   for (i=0; i<nrqr; ++i) {
409:     PetscFree(sbuf2[i]);
410:   }

412:   PetscFree(sbuf2);
413:   PetscFree(rmap[0]);
414:   PetscFree(rmap);

416:   for (i=0; i<ismax; i++) {
417:     MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);
418:     MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);
419:   }

421:   return(0);
422: }

426: PetscErrorCode MatScale_MPIDense(Mat inA,PetscScalar alpha)
427: {
428:   Mat_MPIDense   *A = (Mat_MPIDense*)inA->data;
429:   Mat_SeqDense   *a = (Mat_SeqDense*)A->A->data;
430:   PetscScalar    oalpha = alpha;
432:   PetscBLASInt   one = 1,nz = PetscBLASIntCast(inA->rmap->n*inA->cmap->N);

435:   BLASscal_(&nz,&oalpha,a->v,&one);
436:   PetscLogFlops(nz);
437:   return(0);
438: }