Actual source code: mpirowbs.c

  1: #define PETSCMAT_DLL

 3:  #include src/mat/impls/rowbs/mpi/mpirowbs.h

  5: #define CHUNCKSIZE_LOCAL   10

  9: static PetscErrorCode MatFreeRowbs_Private(Mat A,int n,int *i,PetscScalar *v)
 10: {

 14:   if (v) {
 15: #if defined(PETSC_USE_LOG)
 16:     int len = -n*(sizeof(int)+sizeof(PetscScalar));
 17: #endif
 18:     PetscFree(v);
 19:     PetscLogObjectMemory(A,len);
 20:   }
 21:   return(0);
 22: }

 26: static PetscErrorCode MatMallocRowbs_Private(Mat A,int n,int **i,PetscScalar **v)
 27: {
 29:   int len;

 32:   if (!n) {
 33:     *i = 0; *v = 0;
 34:   } else {
 35:     len = n*(sizeof(int) + sizeof(PetscScalar));
 36:     PetscMalloc(len,v);
 37:     PetscLogObjectMemory(A,len);
 38:     *i = (int*)(*v + n);
 39:   }
 40:   return(0);
 41: }

 45: PetscErrorCode MatScale_MPIRowbs(Mat inA,PetscScalar alpha)
 46: {
 47:   Mat_MPIRowbs   *a = (Mat_MPIRowbs*)inA->data;
 48:   BSspmat        *A = a->A;
 49:   BSsprow        *vs;
 50:   PetscScalar    *ap;
 51:   int            i,m = inA->rmap.n,nrow,j;

 55:   for (i=0; i<m; i++) {
 56:     vs   = A->rows[i];
 57:     nrow = vs->length;
 58:     ap   = vs->nz;
 59:     for (j=0; j<nrow; j++) {
 60:       ap[j] *= alpha;
 61:     }
 62:   }
 63:   PetscLogFlops(a->nz);
 64:   return(0);
 65: }

 67: /* ----------------------------------------------------------------- */
 70: static PetscErrorCode MatCreateMPIRowbs_local(Mat A,int nz,const int nnz[])
 71: {
 72:   Mat_MPIRowbs *bsif = (Mat_MPIRowbs*)A->data;
 74:   int   i,len,m = A->rmap.n,*tnnz;
 75:   BSspmat      *bsmat;
 76:   BSsprow      *vs;

 79:   PetscMalloc((m+1)*sizeof(int),&tnnz);
 80:   if (!nnz) {
 81:     if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
 82:     if (nz <= 0)             nz = 1;
 83:     for (i=0; i<m; i++) tnnz[i] = nz;
 84:     nz      = nz*m;
 85:   } else {
 86:     nz = 0;
 87:     for (i=0; i<m; i++) {
 88:       if (nnz[i] <= 0) tnnz[i] = 1;
 89:       else             tnnz[i] = nnz[i];
 90:       nz += tnnz[i];
 91:     }
 92:   }

 94:   /* Allocate BlockSolve matrix context */
 95:   PetscNew(BSspmat,&bsif->A);
 96:   bsmat = bsif->A;
 97:   BSset_mat_icc_storage(bsmat,PETSC_FALSE);
 98:   BSset_mat_symmetric(bsmat,PETSC_FALSE);
 99:   len                    = m*(sizeof(BSsprow*)+ sizeof(BSsprow)) + 1;
100:   PetscMalloc(len,&bsmat->rows);
101:   bsmat->num_rows        = m;
102:   bsmat->global_num_rows = A->rmap.N;
103:   bsmat->map             = bsif->bsmap;
104:   vs                     = (BSsprow*)(bsmat->rows + m);
105:   for (i=0; i<m; i++) {
106:     bsmat->rows[i]  = vs;
107:     bsif->imax[i]   = tnnz[i];
108:     vs->diag_ind    = -1;
109:     MatMallocRowbs_Private(A,tnnz[i],&(vs->col),&(vs->nz));
110:     /* put zero on diagonal */
111:     /*vs->length            = 1;
112:     vs->col[0]      = i + bsif->rstart;
113:     vs->nz[0]       = 0.0;*/
114:     vs->length = 0;
115:     vs++;
116:   }
117:   PetscLogObjectMemory(A,sizeof(BSspmat) + len);
118:   bsif->nz               = 0;
119:   bsif->maxnz            = nz;
120:   bsif->sorted           = 0;
121:   bsif->roworiented      = PETSC_TRUE;
122:   bsif->nonew            = 0;
123:   bsif->bs_color_single  = 0;

125:   PetscFree(tnnz);
126:   return(0);
127: }

131: static PetscErrorCode MatSetValues_MPIRowbs_local(Mat AA,int m,const int im[],int n,const int in[],const PetscScalar v[],InsertMode addv)
132: {
133:   Mat_MPIRowbs *mat = (Mat_MPIRowbs*)AA->data;
134:   BSspmat      *A = mat->A;
135:   BSsprow      *vs;
137:   int          *rp,k,a,b,t,ii,row,nrow,i,col,l,rmax;
138:   int          *imax = mat->imax,nonew = mat->nonew,sorted = mat->sorted;
139:   PetscScalar  *ap,value;

142:   for (k=0; k<m; k++) { /* loop over added rows */
143:     row = im[k];
144:     if (row < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %d",row);
145:     if (row >= AA->rmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,AA->rmap.n-1);
146:     vs   = A->rows[row];
147:     ap   = vs->nz; rp = vs->col;
148:     rmax = imax[row]; nrow = vs->length;
149:     a    = 0;
150:     for (l=0; l<n; l++) { /* loop over added columns */
151:       if (in[l] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative col: %d",in[l]);
152:       if (in[l] >= AA->cmap.N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],AA->cmap.N-1);
153:       col = in[l]; value = *v++;
154:       if (!sorted) a = 0; b = nrow;
155:       while (b-a > 5) {
156:         t = (b+a)/2;
157:         if (rp[t] > col) b = t;
158:         else             a = t;
159:       }
160:       for (i=a; i<b; i++) {
161:         if (rp[i] > col) break;
162:         if (rp[i] == col) {
163:           if (addv == ADD_VALUES) ap[i] += value;
164:           else                    ap[i] = value;
165:           goto noinsert;
166:         }
167:       }
168:       if (nonew) goto noinsert;
169:       if (nrow >= rmax) {
170:         /* there is no extra room in row, therefore enlarge */
171:         int    *itemp,*iout,*iin = vs->col;
172:         PetscScalar *vout,*vin = vs->nz,*vtemp;

174:         /* malloc new storage space */
175:         imax[row] += CHUNCKSIZE_LOCAL;
176:         MatMallocRowbs_Private(AA,imax[row],&itemp,&vtemp);
177:         vout = vtemp; iout = itemp;
178:         for (ii=0; ii<i; ii++) {
179:           vout[ii] = vin[ii];
180:           iout[ii] = iin[ii];
181:         }
182:         vout[i] = value;
183:         iout[i] = col;
184:         for (ii=i+1; ii<=nrow; ii++) {
185:           vout[ii] = vin[ii-1];
186:           iout[ii] = iin[ii-1];
187:         }
188:         /* free old row storage */
189:         if (rmax > 0) {
190:           MatFreeRowbs_Private(AA,rmax,vs->col,vs->nz);
191:         }
192:         vs->col           =  iout; vs->nz = vout;
193:         rmax              =  imax[row];
194:         mat->maxnz        += CHUNCKSIZE_LOCAL;
195:         mat->reallocs++;
196:       } else {
197:         /* shift higher columns over to make room for newie */
198:         for (ii=nrow-1; ii>=i; ii--) {
199:           rp[ii+1] = rp[ii];
200:           ap[ii+1] = ap[ii];
201:         }
202:         rp[i] = col;
203:         ap[i] = value;
204:       }
205:       nrow++;
206:       mat->nz++;
207:       AA->same_nonzero = PETSC_FALSE;
208:       noinsert:;
209:       a = i + 1;
210:     }
211:     vs->length = nrow;
212:   }
213:   return(0);
214: }


219: static PetscErrorCode MatAssemblyBegin_MPIRowbs_local(Mat A,MatAssemblyType mode)
220: {
222:   return(0);
223: }

227: static PetscErrorCode MatAssemblyEnd_MPIRowbs_local(Mat AA,MatAssemblyType mode)
228: {
229:   Mat_MPIRowbs *a = (Mat_MPIRowbs*)AA->data;
230:   BSspmat      *A = a->A;
231:   BSsprow      *vs;
232:   int          i,j,rstart = AA->rmap.rstart;

235:   if (mode == MAT_FLUSH_ASSEMBLY) return(0);

237:   /* Mark location of diagonal */
238:   for (i=0; i<AA->rmap.n; i++) {
239:     vs = A->rows[i];
240:     for (j=0; j<vs->length; j++) {
241:       if (vs->col[j] == i + rstart) {
242:         vs->diag_ind = j;
243:         break;
244:       }
245:     }
246:     if (vs->diag_ind == -1) {
247:       SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"no diagonal entry");
248:     }
249:   }
250:   return(0);
251: }

255: static PetscErrorCode MatZeroRows_MPIRowbs_local(Mat A,PetscInt N,const PetscInt rz[],PetscScalar diag)
256: {
257:   Mat_MPIRowbs   *a = (Mat_MPIRowbs*)A->data;
258:   BSspmat        *l = a->A;
260:   int            i,m = A->rmap.n - 1,col,base=A->rmap.rstart;

263:   if (a->keepzeroedrows) {
264:     for (i=0; i<N; i++) {
265:       if (rz[i] < 0 || rz[i] > m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"row out of range");
266:       PetscMemzero(l->rows[rz[i]]->nz,l->rows[rz[i]]->length*sizeof(PetscScalar));
267:       if (diag != 0.0) {
268:         col=rz[i]+base;
269:         MatSetValues_MPIRowbs_local(A,1,&rz[i],1,&col,&diag,INSERT_VALUES);
270:       }
271:     }
272:   } else {
273:     if (diag != 0.0) {
274:       for (i=0; i<N; i++) {
275:         if (rz[i] < 0 || rz[i] > m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Out of range");
276:         if (l->rows[rz[i]]->length > 0) { /* in case row was completely empty */
277:           l->rows[rz[i]]->length = 1;
278:           l->rows[rz[i]]->nz[0]  = diag;
279:           l->rows[rz[i]]->col[0] = A->rmap.rstart + rz[i];
280:         } else {
281:           col=rz[i]+base;
282:           MatSetValues_MPIRowbs_local(A,1,&rz[i],1,&col,&diag,INSERT_VALUES);
283:         }
284:       }
285:     } else {
286:       for (i=0; i<N; i++) {
287:         if (rz[i] < 0 || rz[i] > m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Out of range");
288:         l->rows[rz[i]]->length = 0;
289:       }
290:     }
291:     A->same_nonzero = PETSC_FALSE;
292:   }
293:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
294:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
295:   return(0);
296: }

300: static PetscErrorCode MatNorm_MPIRowbs_local(Mat A,NormType type,PetscReal *norm)
301: {
302:   Mat_MPIRowbs *mat = (Mat_MPIRowbs*)A->data;
303:   BSsprow      *vs,**rs;
304:   PetscScalar  *xv;
305:   PetscReal    sum = 0.0;
307:   int          *xi,nz,i,j;

310:   rs = mat->A->rows;
311:   if (type == NORM_FROBENIUS) {
312:     for (i=0; i<A->rmap.n; i++) {
313:       vs = *rs++;
314:       nz = vs->length;
315:       xv = vs->nz;
316:       while (nz--) {
317: #if defined(PETSC_USE_COMPLEX)
318:         sum += PetscRealPart(PetscConj(*xv)*(*xv)); xv++;
319: #else
320:         sum += (*xv)*(*xv); xv++;
321: #endif
322:       }
323:     }
324:     *norm = sqrt(sum);
325:   } else if (type == NORM_1) { /* max column norm */
326:     PetscReal *tmp;
327:     PetscMalloc(A->cmap.n*sizeof(PetscReal),&tmp);
328:     PetscMemzero(tmp,A->cmap.n*sizeof(PetscReal));
329:     *norm = 0.0;
330:     for (i=0; i<A->rmap.n; i++) {
331:       vs = *rs++;
332:       nz = vs->length;
333:       xi = vs->col;
334:       xv = vs->nz;
335:       while (nz--) {
336:         tmp[*xi] += PetscAbsScalar(*xv);
337:         xi++; xv++;
338:       }
339:     }
340:     for (j=0; j<A->rmap.n; j++) {
341:       if (tmp[j] > *norm) *norm = tmp[j];
342:     }
343:     PetscFree(tmp);
344:   } else if (type == NORM_INFINITY) { /* max row norm */
345:     *norm = 0.0;
346:     for (i=0; i<A->rmap.n; i++) {
347:       vs = *rs++;
348:       nz = vs->length;
349:       xv = vs->nz;
350:       sum = 0.0;
351:       while (nz--) {
352:         sum += PetscAbsScalar(*xv); xv++;
353:       }
354:       if (sum > *norm) *norm = sum;
355:     }
356:   } else {
357:     SETERRQ(PETSC_ERR_SUP,"No support for the two norm");
358:   }
359:   return(0);
360: }

362: /* ----------------------------------------------------------------- */

366: PetscErrorCode MatSetValues_MPIRowbs(Mat mat,int m,const int im[],int n,const int in[],const PetscScalar v[],InsertMode av)
367: {
368:   Mat_MPIRowbs *a = (Mat_MPIRowbs*)mat->data;
370:   int   i,j,row,col,rstart = mat->rmap.rstart,rend = mat->rmap.rend;
371:   PetscTruth   roworiented = a->roworiented;

374:   /* Note:  There's no need to "unscale" the matrix, since scaling is
375:      confined to a->pA, and we're working with a->A here */
376:   for (i=0; i<m; i++) {
377:     if (im[i] < 0) continue;
378:     if (im[i] >= mat->rmap.N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",im[i],mat->rmap.N-1);
379:     if (im[i] >= rstart && im[i] < rend) {
380:       row = im[i] - rstart;
381:       for (j=0; j<n; j++) {
382:         if (in[j] < 0) continue;
383:         if (in[j] >= mat->cmap.N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[j],mat->cmap.N-1);
384:         if (in[j] >= 0 && in[j] < mat->cmap.N){
385:           col = in[j];
386:           if (roworiented) {
387:             MatSetValues_MPIRowbs_local(mat,1,&row,1,&col,v+i*n+j,av);
388:           } else {
389:             MatSetValues_MPIRowbs_local(mat,1,&row,1,&col,v+i+j*m,av);
390:           }
391:         } else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid column");}
392:       }
393:     } else {
394:       if (!a->donotstash) {
395:         if (roworiented) {
396:           MatStashValuesRow_Private(&mat->stash,im[i],n,in,v+i*n);
397:         } else {
398:           MatStashValuesCol_Private(&mat->stash,im[i],n,in,v+i,m);
399:         }
400:       }
401:     }
402:   }
403:   return(0);
404: }

408: PetscErrorCode MatAssemblyBegin_MPIRowbs(Mat mat,MatAssemblyType mode)
409: {
410:   MPI_Comm      comm = mat->comm;
412:   int         nstash,reallocs;
413:   InsertMode    addv;

416:   /* Note:  There's no need to "unscale" the matrix, since scaling is
417:             confined to a->pA, and we're working with a->A here */

419:   /* make sure all processors are either in INSERTMODE or ADDMODE */
420:   MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,comm);
421:   if (addv == (ADD_VALUES|INSERT_VALUES)) {
422:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Some procs inserted; others added");
423:   }
424:   mat->insertmode = addv; /* in case this processor had no cache */

426:   MatStashScatterBegin_Private(&mat->stash,mat->rmap.range);
427:   MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);
428:   PetscInfo2(0,"Block-Stash has %d entries, uses %d mallocs.\n",nstash,reallocs);
429:   return(0);
430: }

434: static PetscErrorCode MatView_MPIRowbs_ASCII(Mat mat,PetscViewer viewer)
435: {
436:   Mat_MPIRowbs      *a = (Mat_MPIRowbs*)mat->data;
438:   int               i,j;
439:   PetscTruth        iascii;
440:   BSspmat           *A = a->A;
441:   BSsprow           **rs = A->rows;
442:   PetscViewerFormat format;

445:   PetscViewerGetFormat(viewer,&format);
446:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);

448:   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
449:     int ind_l,ind_g,clq_l,clq_g,color;
450:     ind_l = BSlocal_num_inodes(a->pA);CHKERRBS(0);
451:     ind_g = BSglobal_num_inodes(a->pA);CHKERRBS(0);
452:     clq_l = BSlocal_num_cliques(a->pA);CHKERRBS(0);
453:     clq_g = BSglobal_num_cliques(a->pA);CHKERRBS(0);
454:     color = BSnum_colors(a->pA);CHKERRBS(0);
455:     PetscViewerASCIIPrintf(viewer,"  %d global inode(s), %d global clique(s), %d color(s)\n",ind_g,clq_g,color);
456:     PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d local inode(s), %d local clique(s)\n",a->rank,ind_l,clq_l);
457:   } else  if (format == PETSC_VIEWER_ASCII_COMMON) {
458:     for (i=0; i<A->num_rows; i++) {
459:       PetscViewerASCIISynchronizedPrintf(viewer,"row %d:",i+mat->rmap.rstart);
460:       for (j=0; j<rs[i]->length; j++) {
461:         if (rs[i]->nz[j]) {PetscViewerASCIISynchronizedPrintf(viewer," %d %g ",rs[i]->col[j],rs[i]->nz[j]);}
462:       }
463:       PetscViewerASCIISynchronizedPrintf(viewer,"\n");
464:     }
465:   } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
466:     SETERRQ(PETSC_ERR_SUP,"Matlab format not supported");
467:   } else {
468:     PetscViewerASCIIUseTabs(viewer,PETSC_NO);
469:     for (i=0; i<A->num_rows; i++) {
470:       PetscViewerASCIISynchronizedPrintf(viewer,"row %d:",i+mat->rmap.rstart);
471:       for (j=0; j<rs[i]->length; j++) {
472:         PetscViewerASCIISynchronizedPrintf(viewer," %d %g ",rs[i]->col[j],rs[i]->nz[j]);
473:       }
474:       PetscViewerASCIISynchronizedPrintf(viewer,"\n");
475:     }
476:     PetscViewerASCIIUseTabs(viewer,PETSC_YES);
477:   }
478:   PetscViewerFlush(viewer);
479:   return(0);
480: }

484: static PetscErrorCode MatView_MPIRowbs_Binary(Mat mat,PetscViewer viewer)
485: {
486:   Mat_MPIRowbs   *a = (Mat_MPIRowbs*)mat->data;
488:   PetscMPIInt    rank,size;
489:   PetscInt       i,M,m,*sbuff,*rowlengths;
490:   PetscInt       *recvcts,*recvdisp,fd,*cols,maxnz,nz,j;
491:   BSspmat        *A = a->A;
492:   BSsprow        **rs = A->rows;
493:   MPI_Comm       comm = mat->comm;
494:   MPI_Status     status;
495:   PetscScalar    *vals;
496:   MatInfo        info;

499:   MPI_Comm_size(comm,&size);
500:   MPI_Comm_rank(comm,&rank);

502:   M = mat->rmap.N; m = mat->rmap.n;
503:   /* First gather together on the first processor the lengths of 
504:      each row, and write them out to the file */
505:   PetscMalloc(m*sizeof(int),&sbuff);
506:   for (i=0; i<A->num_rows; i++) {
507:     sbuff[i] = rs[i]->length;
508:   }
509:   MatGetInfo(mat,MAT_GLOBAL_SUM,&info);
510:   if (!rank) {
511:     PetscViewerBinaryGetDescriptor(viewer,&fd);
512:     PetscMalloc((4+M)*sizeof(int),&rowlengths);
513:     PetscMalloc(size*sizeof(int),&recvcts);
514:     recvdisp = mat->rmap.range;
515:     for (i=0; i<size; i++) {
516:       recvcts[i] = recvdisp[i+1] - recvdisp[i];
517:     }
518:     /* first four elements of rowlength are the header */
519:     rowlengths[0] = mat->cookie;
520:     rowlengths[1] = mat->rmap.N;
521:     rowlengths[2] = mat->cmap.N;
522:     rowlengths[3] = (int)info.nz_used;
523:     MPI_Gatherv(sbuff,m,MPI_INT,rowlengths+4,recvcts,recvdisp,MPI_INT,0,comm);
524:     PetscFree(sbuff);
525:     PetscBinaryWrite(fd,rowlengths,4+M,PETSC_INT,PETSC_FALSE);
526:     /* count the number of nonzeros on each processor */
527:     PetscMemzero(recvcts,size*sizeof(int));
528:     for (i=0; i<size; i++) {
529:       for (j=recvdisp[i]; j<recvdisp[i+1]; j++) {
530:         recvcts[i] += rowlengths[j+3];
531:       }
532:     }
533:     /* allocate buffer long enough to hold largest one */
534:     maxnz = 0;
535:     for (i=0; i<size; i++) {
536:       maxnz = PetscMax(maxnz,recvcts[i]);
537:     }
538:     PetscFree(rowlengths);
539:     PetscFree(recvcts);
540:     PetscMalloc(maxnz*sizeof(int),&cols);

542:     /* binary store column indices for 0th processor */
543:     nz = 0;
544:     for (i=0; i<A->num_rows; i++) {
545:       for (j=0; j<rs[i]->length; j++) {
546:         cols[nz++] = rs[i]->col[j];
547:       }
548:     }
549:     PetscBinaryWrite(fd,cols,nz,PETSC_INT,PETSC_FALSE);

551:     /* receive and store column indices for all other processors */
552:     for (i=1; i<size; i++) {
553:       /* should tell processor that I am now ready and to begin the send */
554:       MPI_Recv(cols,maxnz,MPI_INT,i,mat->tag,comm,&status);
555:       MPI_Get_count(&status,MPI_INT,&nz);
556:       PetscBinaryWrite(fd,cols,nz,PETSC_INT,PETSC_FALSE);
557:     }
558:     PetscFree(cols);
559:     PetscMalloc(maxnz*sizeof(PetscScalar),&vals);

561:     /* binary store values for 0th processor */
562:     nz = 0;
563:     for (i=0; i<A->num_rows; i++) {
564:       for (j=0; j<rs[i]->length; j++) {
565:         vals[nz++] = rs[i]->nz[j];
566:       }
567:     }
568:     PetscBinaryWrite(fd,vals,nz,PETSC_SCALAR,PETSC_FALSE);

570:     /* receive and store nonzeros for all other processors */
571:     for (i=1; i<size; i++) {
572:       /* should tell processor that I am now ready and to begin the send */
573:       MPI_Recv(vals,maxnz,MPIU_SCALAR,i,mat->tag,comm,&status);
574:       MPI_Get_count(&status,MPIU_SCALAR,&nz);
575:       PetscBinaryWrite(fd,vals,nz,PETSC_SCALAR,PETSC_FALSE);
576:     }
577:     PetscFree(vals);
578:   } else {
579:     MPI_Gatherv(sbuff,m,MPI_INT,0,0,0,MPI_INT,0,comm);
580:     PetscFree(sbuff);

582:     /* count local nonzeros */
583:     nz = 0;
584:     for (i=0; i<A->num_rows; i++) {
585:       for (j=0; j<rs[i]->length; j++) {
586:         nz++;
587:       }
588:     }
589:     /* copy into buffer column indices */
590:     PetscMalloc(nz*sizeof(int),&cols);
591:     nz = 0;
592:     for (i=0; i<A->num_rows; i++) {
593:       for (j=0; j<rs[i]->length; j++) {
594:         cols[nz++] = rs[i]->col[j];
595:       }
596:     }
597:     /* send */  /* should wait until processor zero tells me to go */
598:     MPI_Send(cols,nz,MPI_INT,0,mat->tag,comm);
599:     PetscFree(cols);

601:     /* copy into buffer column values */
602:     PetscMalloc(nz*sizeof(PetscScalar),&vals);
603:     nz   = 0;
604:     for (i=0; i<A->num_rows; i++) {
605:       for (j=0; j<rs[i]->length; j++) {
606:         vals[nz++] = rs[i]->nz[j];
607:       }
608:     }
609:     /* send */  /* should wait until processor zero tells me to go */
610:     MPI_Send(vals,nz,MPIU_SCALAR,0,mat->tag,comm);
611:     PetscFree(vals);
612:   }

614:   return(0);
615: }

619: PetscErrorCode MatView_MPIRowbs(Mat mat,PetscViewer viewer)
620: {
621:   Mat_MPIRowbs *bsif = (Mat_MPIRowbs*)mat->data;
623:   PetscTruth   iascii,isbinary;

626:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
627:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
628:   if (!bsif->blocksolveassembly) {
629:     MatAssemblyEnd_MPIRowbs_ForBlockSolve(mat);
630:   }
631:   if (iascii) {
632:     MatView_MPIRowbs_ASCII(mat,viewer);
633:   } else if (isbinary) {
634:     MatView_MPIRowbs_Binary(mat,viewer);
635:   } else {
636:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by MPIRowbs matrices",((PetscObject)viewer)->type_name);
637:   }
638:   return(0);
639: }
640: 
643: static PetscErrorCode MatAssemblyEnd_MPIRowbs_MakeSymmetric(Mat mat)
644: {
645:   Mat_MPIRowbs *a = (Mat_MPIRowbs*)mat->data;
646:   BSspmat      *A = a->A;
647:   BSsprow      *vs;
648:   int          size,rank,M,rstart,tag,i,j,*rtable,*w1,*w3,*w4,len,proc,nrqs;
649:   int          msz,*pa,bsz,nrqr,**rbuf1,**sbuf1,**ptr,*tmp,*ctr,col,idx,row;
651:   int          ctr_j,*sbuf1_j,k;
652:   PetscScalar  val=0.0;
653:   MPI_Comm     comm;
654:   MPI_Request  *s_waits1,*r_waits1;
655:   MPI_Status   *s_status,*r_status;

658:   comm   = mat->comm;
659:   tag    = mat->tag;
660:   size   = a->size;
661:   rank   = a->rank;
662:   M      = mat->rmap.N;
663:   rstart = mat->rmap.rstart;

665:   PetscMalloc(M*sizeof(int),&rtable);
666:   /* Create hash table for the mapping :row -> proc */
667:   for (i=0,j=0; i<size; i++) {
668:     len = mat->rmap.range[i+1];
669:     for (; j<len; j++) {
670:       rtable[j] = i;
671:     }
672:   }

674:   /* Evaluate communication - mesg to whom, length of mesg, and buffer space
675:      required. Based on this, buffers are allocated, and data copied into them. */
676:   PetscMalloc(size*4*sizeof(int),&w1);/*  mesg size */
677:   w3   = w1 + 2*size;       /* no of IS that needs to be sent to proc i */
678:   w4   = w3 + size;       /* temp work space used in determining w1,  w3 */
679:   PetscMemzero(w1,size*3*sizeof(int)); /* initialize work vector */

681:   for (i=0;  i<mat->rmap.n; i++) {
682:     PetscMemzero(w4,size*sizeof(int)); /* initialize work vector */
683:     vs = A->rows[i];
684:     for (j=0; j<vs->length; j++) {
685:       proc = rtable[vs->col[j]];
686:       w4[proc]++;
687:     }
688:     for (j=0; j<size; j++) {
689:       if (w4[j]) { w1[2*j] += w4[j]; w3[j]++;}
690:     }
691:   }
692: 
693:   nrqs       = 0;              /* number of outgoing messages */
694:   msz        = 0;              /* total mesg length (for all proc */
695:   w1[2*rank] = 0;              /* no mesg sent to itself */
696:   w3[rank]   = 0;
697:   for (i=0; i<size; i++) {
698:     if (w1[2*i])  {w1[2*i+1] = 1; nrqs++;} /* there exists a message to proc i */
699:   }
700:   /* pa - is list of processors to communicate with */
701:   PetscMalloc((nrqs+1)*sizeof(int),&pa);
702:   for (i=0,j=0; i<size; i++) {
703:     if (w1[2*i]) {pa[j] = i; j++;}
704:   }

706:   /* Each message would have a header = 1 + 2*(no of ROWS) + data */
707:   for (i=0; i<nrqs; i++) {
708:     j       = pa[i];
709:     w1[2*j] += w1[2*j+1] + 2*w3[j];
710:     msz     += w1[2*j];
711:   }
712: 
713:   /* Do a global reduction to determine how many messages to expect */
714:   PetscMaxSum(comm,w1,&bsz,&nrqr);

716:   /* Allocate memory for recv buffers . Prob none if nrqr = 0 ???? */
717:   len      = (nrqr+1)*sizeof(int*) + nrqr*bsz*sizeof(int);
718:   PetscMalloc(len,&rbuf1);
719:   rbuf1[0] = (int*)(rbuf1 + nrqr);
720:   for (i=1; i<nrqr; ++i) rbuf1[i] = rbuf1[i-1] + bsz;

722:   /* Post the receives */
723:   PetscMalloc((nrqr+1)*sizeof(MPI_Request),&r_waits1);
724:   for (i=0; i<nrqr; ++i){
725:     MPI_Irecv(rbuf1[i],bsz,MPI_INT,MPI_ANY_SOURCE,tag,comm,r_waits1+i);
726:   }
727: 
728:   /* Allocate Memory for outgoing messages */
729:   len   = 2*size*sizeof(int*) + (size+msz)*sizeof(int);
730:   PetscMalloc(len,&sbuf1);
731:   ptr   = sbuf1 + size;     /* Pointers to the data in outgoing buffers */
732:   PetscMemzero(sbuf1,2*size*sizeof(int*));
733:   tmp   = (int*)(sbuf1 + 2*size);
734:   ctr   = tmp + msz;

736:   {
737:     int *iptr = tmp,ict  = 0;
738:     for (i=0; i<nrqs; i++) {
739:       j        = pa[i];
740:       iptr    += ict;
741:       sbuf1[j] = iptr;
742:       ict      = w1[2*j];
743:     }
744:   }

746:   /* Form the outgoing messages */
747:   /* Clean up the header space */
748:   for (i=0; i<nrqs; i++) {
749:     j           = pa[i];
750:     sbuf1[j][0] = 0;
751:     PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(int));
752:     ptr[j]      = sbuf1[j] + 2*w3[j] + 1;
753:   }

755:   /* Parse the matrix and copy the data into sbuf1 */
756:   for (i=0; i<mat->rmap.n; i++) {
757:     PetscMemzero(ctr,size*sizeof(int));
758:     vs = A->rows[i];
759:     for (j=0; j<vs->length; j++) {
760:       col  = vs->col[j];
761:       proc = rtable[col];
762:       if (proc != rank) { /* copy to the outgoing buffer */
763:         ctr[proc]++;
764:           *ptr[proc] = col;
765:           ptr[proc]++;
766:       } else {
767:         row = col - rstart;
768:         col = i + rstart;
769:         MatSetValues_MPIRowbs_local(mat,1,&row,1,&col,&val,ADD_VALUES);
770:       }
771:     }
772:     /* Update the headers for the current row */
773:     for (j=0; j<size; j++) { /* Can Optimise this loop by using pa[] */
774:       if ((ctr_j = ctr[j])) {
775:         sbuf1_j        = sbuf1[j];
776:         k               = ++sbuf1_j[0];
777:         sbuf1_j[2*k]   = ctr_j;
778:         sbuf1_j[2*k-1] = i + rstart;
779:       }
780:     }
781:   }
782:    /* Check Validity of the outgoing messages */
783:   {
784:     int sum;
785:     for (i=0 ; i<nrqs ; i++) {
786:       j = pa[i];
787:       if (w3[j] != sbuf1[j][0]) {SETERRQ(PETSC_ERR_PLIB,"Blew it! Header[1] mismatch!\n"); }
788:     }

790:     for (i=0 ; i<nrqs ; i++) {
791:       j = pa[i];
792:       sum = 1;
793:       for (k = 1; k <= w3[j]; k++) sum += sbuf1[j][2*k]+2;
794:       if (sum != w1[2*j]) { SETERRQ(PETSC_ERR_PLIB,"Blew it! Header[2-n] mismatch!\n"); }
795:     }
796:   }
797: 
798:   /* Now post the sends */
799:   PetscMalloc((nrqs+1)*sizeof(MPI_Request),&s_waits1);
800:   for (i=0; i<nrqs; ++i) {
801:     j    = pa[i];
802:     MPI_Isend(sbuf1[j],w1[2*j],MPI_INT,j,tag,comm,s_waits1+i);
803:   }
804: 
805:   /* Receive messages*/
806:   PetscMalloc((nrqr+1)*sizeof(MPI_Status),&r_status);
807:   for (i=0; i<nrqr; ++i) {
808:     MPI_Waitany(nrqr,r_waits1,&idx,r_status+i);
809:     /* Process the Message */
810:     {
811:       int    *rbuf1_i,n_row,ct1;

813:       rbuf1_i = rbuf1[idx];
814:       n_row   = rbuf1_i[0];
815:       ct1     = 2*n_row+1;
816:       val     = 0.0;
817:       /* Optimise this later */
818:       for (j=1; j<=n_row; j++) {
819:         col = rbuf1_i[2*j-1];
820:         for (k=0; k<rbuf1_i[2*j]; k++,ct1++) {
821:           row = rbuf1_i[ct1] - rstart;
822:           MatSetValues_MPIRowbs_local(mat,1,&row,1,&col,&val,ADD_VALUES);
823:         }
824:       }
825:     }
826:   }

828:   PetscMalloc((nrqs+1)*sizeof(MPI_Status),&s_status);
829:   if (nrqs) {MPI_Waitall(nrqs,s_waits1,s_status);}

831:   PetscFree(rtable);
832:   PetscFree(w1);
833:   PetscFree(pa);
834:   PetscFree(rbuf1);
835:   PetscFree(sbuf1);
836:   PetscFree(r_waits1);
837:   PetscFree(s_waits1);
838:   PetscFree(r_status);
839:   PetscFree(s_status);
840:   return(0);
841: }

843: /*
844:      This does the BlockSolve portion of the matrix assembly.
845:    It is provided in a separate routine so that users can
846:    operate on the matrix (using MatScale(), MatShift() etc.) after 
847:    the matrix has been assembled but before BlockSolve has sucked it
848:    in and devoured it.
849: */
852: PetscErrorCode MatAssemblyEnd_MPIRowbs_ForBlockSolve(Mat mat)
853: {
854:   Mat_MPIRowbs *a = (Mat_MPIRowbs*)mat->data;
856:   int          ldim,low,high,i;
857:   PetscScalar  *diag;

860:   if ((mat->was_assembled) && (!mat->same_nonzero)) {  /* Free the old info */
861:     if (a->pA)       {BSfree_par_mat(a->pA);CHKERRBS(0);}
862:     if (a->comm_pA)  {BSfree_comm(a->comm_pA);CHKERRBS(0);}
863:   }

865:   if ((!mat->same_nonzero) || (!mat->was_assembled)) {
866:     /* Indicates bypassing cliques in coloring */
867:     if (a->bs_color_single) {
868:       BSctx_set_si(a->procinfo,100);
869:     }
870:     /* Form permuted matrix for efficient parallel execution */
871:     a->pA = BSmain_perm(a->procinfo,a->A);CHKERRBS(0);
872:     /* Set up the communication */
873:     a->comm_pA = BSsetup_forward(a->pA,a->procinfo);CHKERRBS(0);
874:   } else {
875:     /* Repermute the matrix */
876:     BSmain_reperm(a->procinfo,a->A,a->pA);CHKERRBS(0);
877:   }

879:   /* Symmetrically scale the matrix by the diagonal */
880:   BSscale_diag(a->pA,a->pA->diag,a->procinfo);CHKERRBS(0);

882:   /* Store inverse of square root of permuted diagonal scaling matrix */
883:   VecGetLocalSize(a->diag,&ldim);
884:   VecGetOwnershipRange(a->diag,&low,&high);
885:   VecGetArray(a->diag,&diag);
886:   for (i=0; i<ldim; i++) {
887:     if (a->pA->scale_diag[i] != 0.0) {
888:       diag[i] = 1.0/sqrt(PetscAbsScalar(a->pA->scale_diag[i]));
889:     } else {
890:       diag[i] = 1.0;
891:     }
892:   }
893:   VecRestoreArray(a->diag,&diag);
894:   a->assembled_icc_storage = a->A->icc_storage;
895:   a->blocksolveassembly = 1;
896:   mat->was_assembled    = PETSC_TRUE;
897:   mat->same_nonzero     = PETSC_TRUE;
898:   PetscInfo(mat,"Completed BlockSolve95 matrix assembly\n");
899:   return(0);
900: }

904: PetscErrorCode MatAssemblyEnd_MPIRowbs(Mat mat,MatAssemblyType mode)
905: {
906:   Mat_MPIRowbs *a = (Mat_MPIRowbs*)mat->data;
908:   int          i,n,row,col,*rows,*cols,rstart,nzcount,flg,j,ncols;
909:   PetscScalar  *vals,val;
910:   InsertMode   addv = mat->insertmode;

913:   while (1) {
914:     MatStashScatterGetMesg_Private(&mat->stash,&n,&rows,&cols,&vals,&flg);
915:     if (!flg) break;
916: 
917:     for (i=0; i<n;) {
918:       /* Now identify the consecutive vals belonging to the same row */
919:       for (j=i,rstart=rows[j]; j<n; j++) { if (rows[j] != rstart) break; }
920:       if (j < n) ncols = j-i;
921:       else       ncols = n-i;
922:       /* Now assemble all these values with a single function call */
923:       MatSetValues_MPIRowbs(mat,1,rows+i,ncols,cols+i,vals+i,addv);
924:       i = j;
925:     }
926:   }
927:   MatStashScatterEnd_Private(&mat->stash);

929:   rstart = mat->rmap.rstart;
930:   nzcount = a->nz; /* This is the number of nonzeros entered by the user */
931:   /* BlockSolve requires that the matrix is structurally symmetric */
932:   if (mode == MAT_FINAL_ASSEMBLY && !mat->structurally_symmetric) {
933:     MatAssemblyEnd_MPIRowbs_MakeSymmetric(mat);
934:   }
935: 
936:   /* BlockSolve requires that all the diagonal elements are set */
937:   val  = 0.0;
938:   for (i=0; i<mat->rmap.n; i++) {
939:     row = i; col = i + rstart;
940:     MatSetValues_MPIRowbs_local(mat,1,&row,1,&col,&val,ADD_VALUES);
941:   }
942: 
943:   MatAssemblyBegin_MPIRowbs_local(mat,mode);
944:   MatAssemblyEnd_MPIRowbs_local(mat,mode);
945: 
946:   a->blocksolveassembly = 0;
947:   PetscInfo4(mat,"Matrix size: %d X %d; storage space: %d unneeded,%d used\n",mat->rmap.n,mat->cmap.n,a->maxnz-a->nz,a->nz);
948:   PetscInfo2(mat,"User entered %d nonzeros, PETSc added %d\n",nzcount,a->nz-nzcount);
949:   PetscInfo1(mat,"Number of mallocs during MatSetValues is %d\n",a->reallocs);
950:   return(0);
951: }

955: PetscErrorCode MatZeroEntries_MPIRowbs(Mat mat)
956: {
957:   Mat_MPIRowbs *l = (Mat_MPIRowbs*)mat->data;
958:   BSspmat      *A = l->A;
959:   BSsprow      *vs;
960:   int          i,j;

963:   for (i=0; i <mat->rmap.n; i++) {
964:     vs = A->rows[i];
965:     for (j=0; j< vs->length; j++) vs->nz[j] = 0.0;
966:   }
967:   return(0);
968: }

970: /* the code does not do the diagonal entries correctly unless the 
971:    matrix is square and the column and row owerships are identical.
972:    This is a BUG.
973: */

977: PetscErrorCode MatZeroRows_MPIRowbs(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag)
978: {
979:   Mat_MPIRowbs   *l = (Mat_MPIRowbs*)A->data;
981:   int            i,*owners = A->rmap.range,size = l->size;
982:   int            *nprocs,j,idx,nsends;
983:   int            nmax,*svalues,*starts,*owner,nrecvs,rank = l->rank;
984:   int            *rvalues,tag = A->tag,count,base,slen,n,*source;
985:   int            *lens,imdex,*lrows,*values;
986:   MPI_Comm       comm = A->comm;
987:   MPI_Request    *send_waits,*recv_waits;
988:   MPI_Status     recv_status,*send_status;
989:   PetscTruth     found;

992:   /*  first count number of contributors to each processor */
993:   PetscMalloc(2*size*sizeof(int),&nprocs);
994:   PetscMemzero(nprocs,2*size*sizeof(int));
995:   PetscMalloc((N+1)*sizeof(int),&owner); /* see note*/
996:   for (i=0; i<N; i++) {
997:     idx   = rows[i];
998:     found = PETSC_FALSE;
999:     for (j=0; j<size; j++) {
1000:       if (idx >= owners[j] && idx < owners[j+1]) {
1001:         nprocs[2*j]++; nprocs[2*j+1] = 1; owner[i] = j; found = PETSC_TRUE; break;
1002:       }
1003:     }
1004:     if (!found) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row out of range");
1005:   }
1006:   nsends = 0;  for (i=0; i<size; i++) {nsends += nprocs[2*i+1];}

1008:   /* inform other processors of number of messages and max length*/
1009:   PetscMaxSum(comm,nprocs,&nmax,&nrecvs);

1011:   /* post receives:   */
1012:   PetscMalloc((nrecvs+1)*(nmax+1)*sizeof(int),&rvalues);
1013:   PetscMalloc((nrecvs+1)*sizeof(MPI_Request),&recv_waits);
1014:   for (i=0; i<nrecvs; i++) {
1015:     MPI_Irecv(rvalues+nmax*i,nmax,MPI_INT,MPI_ANY_SOURCE,tag,comm,recv_waits+i);
1016:   }

1018:   /* do sends:
1019:       1) starts[i] gives the starting index in svalues for stuff going to 
1020:          the ith processor
1021:   */
1022:   PetscMalloc((N+1)*sizeof(int),&svalues);
1023:   PetscMalloc((nsends+1)*sizeof(MPI_Request),&send_waits);
1024:   PetscMalloc((size+1)*sizeof(int),&starts);
1025:   starts[0] = 0;
1026:   for (i=1; i<size; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];}
1027:   for (i=0; i<N; i++) {
1028:     svalues[starts[owner[i]]++] = rows[i];
1029:   }

1031:   starts[0] = 0;
1032:   for (i=1; i<size+1; i++) { starts[i] = starts[i-1] + nprocs[2*i-2];}
1033:   count = 0;
1034:   for (i=0; i<size; i++) {
1035:     if (nprocs[2*i+1]) {
1036:       MPI_Isend(svalues+starts[i],nprocs[2*i],MPI_INT,i,tag,comm,send_waits+count++);
1037:     }
1038:   }
1039:   PetscFree(starts);

1041:   base = owners[rank];

1043:   /*  wait on receives */
1044:   PetscMalloc(2*(nrecvs+1)*sizeof(int),&lens);
1045:   source = lens + nrecvs;
1046:   count = nrecvs; slen = 0;
1047:   while (count) {
1048:     MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);
1049:     /* unpack receives into our local space */
1050:     MPI_Get_count(&recv_status,MPI_INT,&n);
1051:     source[imdex]  = recv_status.MPI_SOURCE;
1052:     lens[imdex]    = n;
1053:     slen           += n;
1054:     count--;
1055:   }
1056:   PetscFree(recv_waits);
1057: 
1058:   /* move the data into the send scatter */
1059:   PetscMalloc((slen+1)*sizeof(int),&lrows);
1060:   count = 0;
1061:   for (i=0; i<nrecvs; i++) {
1062:     values = rvalues + i*nmax;
1063:     for (j=0; j<lens[i]; j++) {
1064:       lrows[count++] = values[j] - base;
1065:     }
1066:   }
1067:   PetscFree(rvalues);
1068:   PetscFree(lens);
1069:   PetscFree(owner);
1070:   PetscFree(nprocs);
1071: 
1072:   /* actually zap the local rows */
1073:   MatZeroRows_MPIRowbs_local(A,slen,lrows,diag);
1074:   PetscFree(lrows);

1076:   /* wait on sends */
1077:   if (nsends) {
1078:     PetscMalloc(nsends*sizeof(MPI_Status),&send_status);
1079:     MPI_Waitall(nsends,send_waits,send_status);
1080:     PetscFree(send_status);
1081:   }
1082:   PetscFree(send_waits);
1083:   PetscFree(svalues);

1085:   return(0);
1086: }

1090: PetscErrorCode MatNorm_MPIRowbs(Mat mat,NormType type,PetscReal *norm)
1091: {
1092:   Mat_MPIRowbs *a = (Mat_MPIRowbs*)mat->data;
1093:   BSsprow      *vs,**rs;
1094:   PetscScalar  *xv;
1095:   PetscReal    sum = 0.0;
1097:   int          *xi,nz,i,j;

1100:   if (a->size == 1) {
1101:     MatNorm_MPIRowbs_local(mat,type,norm);
1102:   } else {
1103:     rs = a->A->rows;
1104:     if (type == NORM_FROBENIUS) {
1105:       for (i=0; i<mat->rmap.n; i++) {
1106:         vs = *rs++;
1107:         nz = vs->length;
1108:         xv = vs->nz;
1109:         while (nz--) {
1110: #if defined(PETSC_USE_COMPLEX)
1111:           sum += PetscRealPart(PetscConj(*xv)*(*xv)); xv++;
1112: #else
1113:           sum += (*xv)*(*xv); xv++;
1114: #endif
1115:         }
1116:       }
1117:       MPI_Allreduce(&sum,norm,1,MPIU_REAL,MPI_SUM,mat->comm);
1118:       *norm = sqrt(*norm);
1119:     } else if (type == NORM_1) { /* max column norm */
1120:       PetscReal *tmp,*tmp2;
1121:       PetscMalloc(mat->cmap.n*sizeof(PetscReal),&tmp);
1122:       PetscMalloc(mat->cmap.n*sizeof(PetscReal),&tmp2);
1123:       PetscMemzero(tmp,mat->cmap.n*sizeof(PetscReal));
1124:       *norm = 0.0;
1125:       for (i=0; i<mat->rmap.n; i++) {
1126:         vs = *rs++;
1127:         nz = vs->length;
1128:         xi = vs->col;
1129:         xv = vs->nz;
1130:         while (nz--) {
1131:           tmp[*xi] += PetscAbsScalar(*xv);
1132:           xi++; xv++;
1133:         }
1134:       }
1135:       MPI_Allreduce(tmp,tmp2,mat->cmap.N,MPIU_REAL,MPI_SUM,mat->comm);
1136:       for (j=0; j<mat->cmap.n; j++) {
1137:         if (tmp2[j] > *norm) *norm = tmp2[j];
1138:       }
1139:       PetscFree(tmp);
1140:       PetscFree(tmp2);
1141:     } else if (type == NORM_INFINITY) { /* max row norm */
1142:       PetscReal ntemp = 0.0;
1143:       for (i=0; i<mat->rmap.n; i++) {
1144:         vs = *rs++;
1145:         nz = vs->length;
1146:         xv = vs->nz;
1147:         sum = 0.0;
1148:         while (nz--) {
1149:           sum += PetscAbsScalar(*xv); xv++;
1150:         }
1151:         if (sum > ntemp) ntemp = sum;
1152:       }
1153:       MPI_Allreduce(&ntemp,norm,1,MPIU_REAL,MPI_MAX,mat->comm);
1154:     } else {
1155:       SETERRQ(PETSC_ERR_SUP,"No support for two norm");
1156:     }
1157:   }
1158:   return(0);
1159: }

1163: PetscErrorCode MatMult_MPIRowbs(Mat mat,Vec xx,Vec yy)
1164: {
1165:   Mat_MPIRowbs *bsif = (Mat_MPIRowbs*)mat->data;
1166:   BSprocinfo   *bspinfo = bsif->procinfo;
1167:   PetscScalar  *xxa,*xworka,*yya;

1171:   if (!bsif->blocksolveassembly) {
1172:     MatAssemblyEnd_MPIRowbs_ForBlockSolve(mat);
1173:   }

1175:   /* Permute and apply diagonal scaling:  [ xwork = D^{1/2} * x ] */
1176:   if (!bsif->vecs_permscale) {
1177:     VecGetArray(bsif->xwork,&xworka);
1178:     VecGetArray(xx,&xxa);
1179:     BSperm_dvec(xxa,xworka,bsif->pA->perm);CHKERRBS(0);
1180:     VecRestoreArray(bsif->xwork,&xworka);
1181:     VecRestoreArray(xx,&xxa);
1182:     VecPointwiseDivide(xx,bsif->xwork,bsif->diag);
1183:   }

1185:   VecGetArray(xx,&xxa);
1186:   VecGetArray(yy,&yya);
1187:   /* Do lower triangular multiplication:  [ y = L * xwork ] */
1188:   if (bspinfo->single) {
1189:     BSforward1(bsif->pA,xxa,yya,bsif->comm_pA,bspinfo);CHKERRBS(0);
1190:   }  else {
1191:     BSforward(bsif->pA,xxa,yya,bsif->comm_pA,bspinfo);CHKERRBS(0);
1192:   }
1193: 
1194:   /* Do upper triangular multiplication:  [ y = y + L^{T} * xwork ] */
1195:   if (mat->symmetric) {
1196:     if (bspinfo->single){
1197:       BSbackward1(bsif->pA,xxa,yya,bsif->comm_pA,bspinfo);CHKERRBS(0);
1198:     } else {
1199:       BSbackward(bsif->pA,xxa,yya,bsif->comm_pA,bspinfo);CHKERRBS(0);
1200:     }
1201:   }
1202:   /* not needed for ILU version since forward does it all */
1203:   VecRestoreArray(xx,&xxa);
1204:   VecRestoreArray(yy,&yya);

1206:   /* Apply diagonal scaling to vector:  [  y = D^{1/2} * y ] */
1207:   if (!bsif->vecs_permscale) {
1208:     VecGetArray(bsif->xwork,&xworka);
1209:     VecGetArray(xx,&xxa);
1210:     BSiperm_dvec(xworka,xxa,bsif->pA->perm);CHKERRBS(0);
1211:     VecRestoreArray(bsif->xwork,&xworka);
1212:     VecRestoreArray(xx,&xxa);
1213:     VecPointwiseDivide(bsif->xwork,yy,bsif->diag);
1214:     VecGetArray(bsif->xwork,&xworka);
1215:     VecGetArray(yy,&yya);
1216:     BSiperm_dvec(xworka,yya,bsif->pA->perm);CHKERRBS(0);
1217:     VecRestoreArray(bsif->xwork,&xworka);
1218:     VecRestoreArray(yy,&yya);
1219:   }
1220:   PetscLogFlops(2*bsif->nz - mat->cmap.n);

1222:   return(0);
1223: }

1227: PetscErrorCode MatMultAdd_MPIRowbs(Mat mat,Vec xx,Vec yy,Vec zz)
1228: {
1230:   PetscScalar  one = 1.0;

1233:   (*mat->ops->mult)(mat,xx,zz);
1234:   VecAXPY(zz,one,yy);
1235:   return(0);
1236: }

1240: PetscErrorCode MatGetInfo_MPIRowbs(Mat A,MatInfoType flag,MatInfo *info)
1241: {
1242:   Mat_MPIRowbs *mat = (Mat_MPIRowbs*)A->data;
1243:   PetscReal    isend[5],irecv[5];

1247:   info->rows_global    = (double)A->rmap.N;
1248:   info->columns_global = (double)A->cmap.N;
1249:   info->rows_local     = (double)A->cmap.n;
1250:   info->columns_local  = (double)A->rmap.n;
1251:   info->block_size     = 1.0;
1252:   info->mallocs        = (double)mat->reallocs;
1253:   isend[0] = mat->nz; isend[1] = mat->maxnz; isend[2] =  mat->maxnz -  mat->nz;
1254:   isend[3] = A->mem;  isend[4] = info->mallocs;

1256:   if (flag == MAT_LOCAL) {
1257:     info->nz_used      = isend[0];
1258:     info->nz_allocated = isend[1];
1259:     info->nz_unneeded  = isend[2];
1260:     info->memory       = isend[3];
1261:     info->mallocs      = isend[4];
1262:   } else if (flag == MAT_GLOBAL_MAX) {
1263:     MPI_Allreduce(isend,irecv,3,MPIU_REAL,MPI_MAX,A->comm);
1264:     info->nz_used      = irecv[0];
1265:     info->nz_allocated = irecv[1];
1266:     info->nz_unneeded  = irecv[2];
1267:     info->memory       = irecv[3];
1268:     info->mallocs      = irecv[4];
1269:   } else if (flag == MAT_GLOBAL_SUM) {
1270:     MPI_Allreduce(isend,irecv,3,MPIU_REAL,MPI_SUM,A->comm);
1271:     info->nz_used      = irecv[0];
1272:     info->nz_allocated = irecv[1];
1273:     info->nz_unneeded  = irecv[2];
1274:     info->memory       = irecv[3];
1275:     info->mallocs      = irecv[4];
1276:   }
1277:   return(0);
1278: }

1282: PetscErrorCode MatGetDiagonal_MPIRowbs(Mat mat,Vec v)
1283: {
1284:   Mat_MPIRowbs *a = (Mat_MPIRowbs*)mat->data;
1285:   BSsprow      **rs = a->A->rows;
1287:   int          i,n;
1288:   PetscScalar  *x,zero = 0.0;

1291:   if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1292:   if (!a->blocksolveassembly) {
1293:     MatAssemblyEnd_MPIRowbs_ForBlockSolve(mat);
1294:   }

1296:   VecSet(v,zero);
1297:   VecGetLocalSize(v,&n);
1298:   if (n != mat->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec");
1299:   VecGetArray(v,&x);
1300:   for (i=0; i<mat->rmap.n; i++) {
1301:     x[i] = rs[i]->nz[rs[i]->diag_ind];
1302:   }
1303:   VecRestoreArray(v,&x);
1304:   return(0);
1305: }

1309: PetscErrorCode MatDestroy_MPIRowbs(Mat mat)
1310: {
1311:   Mat_MPIRowbs *a = (Mat_MPIRowbs*)mat->data;
1312:   BSspmat      *A = a->A;
1313:   BSsprow      *vs;
1315:   int          i;

1318: #if defined(PETSC_USE_LOG)
1319:   PetscLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d",mat->rmap.N,mat->cmap.N);
1320: #endif
1321:   MatStashDestroy_Private(&mat->stash);
1322:   if (a->bsmap) {
1323:     PetscFree(a->bsmap->vlocal2global);
1324:     PetscFree(a->bsmap->vglobal2local);
1325:     if (a->bsmap->vglobal2proc)  (*a->bsmap->free_g2p)(a->bsmap->vglobal2proc);
1326:     PetscFree(a->bsmap);
1327:   }

1329:   if (A) {
1330:     for (i=0; i<mat->rmap.n; i++) {
1331:       vs = A->rows[i];
1332:       MatFreeRowbs_Private(mat,vs->length,vs->col,vs->nz);
1333:     }
1334:     /* Note: A->map = a->bsmap is freed above */
1335:     PetscFree(A->rows);
1336:     PetscFree(A);
1337:   }
1338:   if (a->procinfo) {BSfree_ctx(a->procinfo);CHKERRBS(0);}
1339:   if (a->diag)     {VecDestroy(a->diag);}
1340:   if (a->xwork)    {VecDestroy(a->xwork);}
1341:   if (a->pA)       {BSfree_par_mat(a->pA);CHKERRBS(0);}
1342:   if (a->fpA)      {BSfree_copy_par_mat(a->fpA);CHKERRBS(0);}
1343:   if (a->comm_pA)  {BSfree_comm(a->comm_pA);CHKERRBS(0);}
1344:   if (a->comm_fpA) {BSfree_comm(a->comm_fpA);CHKERRBS(0);}
1345:   PetscFree(a->imax);
1346:   MPI_Comm_free(&(a->comm_mpirowbs));
1347:   PetscFree(a);
1348:   PetscObjectComposeFunction((PetscObject)mat,"MatMPIRowbsSetPreallocation_C","",PETSC_NULL);
1349:   return(0);
1350: }

1354: PetscErrorCode MatSetOption_MPIRowbs(Mat A,MatOption op)
1355: {
1356:   Mat_MPIRowbs   *a = (Mat_MPIRowbs*)A->data;

1360:   switch (op) {
1361:   case MAT_ROW_ORIENTED:
1362:     a->roworiented = PETSC_TRUE;
1363:     break;
1364:   case MAT_COLUMN_ORIENTED:
1365:     a->roworiented = PETSC_FALSE;
1366:     break;
1367:   case MAT_COLUMNS_SORTED:
1368:     a->sorted      = 1;
1369:     break;
1370:   case MAT_COLUMNS_UNSORTED:
1371:     a->sorted      = 0;
1372:     break;
1373:   case MAT_NO_NEW_NONZERO_LOCATIONS:
1374:     a->nonew       = 1;
1375:     break;
1376:   case MAT_YES_NEW_NONZERO_LOCATIONS:
1377:     a->nonew       = 0;
1378:     break;
1379:   case MAT_DO_NOT_USE_INODES:
1380:     a->bs_color_single = 1;
1381:     break;
1382:   case MAT_YES_NEW_DIAGONALS:
1383:   case MAT_ROWS_SORTED:
1384:   case MAT_NEW_NONZERO_LOCATION_ERR:
1385:   case MAT_NEW_NONZERO_ALLOCATION_ERR:
1386:   case MAT_ROWS_UNSORTED:
1387:   case MAT_USE_HASH_TABLE:
1388:     PetscInfo1(A,"Option %d ignored\n",op);
1389:     break;
1390:   case MAT_IGNORE_OFF_PROC_ENTRIES:
1391:     a->donotstash = PETSC_TRUE;
1392:     break;
1393:   case MAT_NO_NEW_DIAGONALS:
1394:     SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
1395:     break;
1396:   case MAT_KEEP_ZEROED_ROWS:
1397:     a->keepzeroedrows    = PETSC_TRUE;
1398:     break;
1399:   case MAT_SYMMETRIC:
1400:     BSset_mat_symmetric(a->A,PETSC_TRUE);CHKERRBS(0);
1401:     break;
1402:   case MAT_STRUCTURALLY_SYMMETRIC:
1403:   case MAT_NOT_SYMMETRIC:
1404:   case MAT_NOT_STRUCTURALLY_SYMMETRIC:
1405:   case MAT_HERMITIAN:
1406:   case MAT_NOT_HERMITIAN:
1407:   case MAT_SYMMETRY_ETERNAL:
1408:   case MAT_NOT_SYMMETRY_ETERNAL:
1409:     break;
1410:   default:
1411:     SETERRQ1(PETSC_ERR_SUP,"unknown option %d",op);
1412:     break;
1413:   }
1414:   return(0);
1415: }

1419: PetscErrorCode MatGetRow_MPIRowbs(Mat AA,int row,int *nz,int **idx,PetscScalar **v)
1420: {
1421:   Mat_MPIRowbs *mat = (Mat_MPIRowbs*)AA->data;
1422:   BSspmat      *A = mat->A;
1423:   BSsprow      *rs;
1424: 
1426:   if (row < AA->rmap.rstart || row >= AA->rmap.rend) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Only local rows");

1428:   rs  = A->rows[row - AA->rmap.rstart];
1429:   *nz = rs->length;
1430:   if (v)   *v   = rs->nz;
1431:   if (idx) *idx = rs->col;
1432:   return(0);
1433: }

1437: PetscErrorCode MatRestoreRow_MPIRowbs(Mat A,int row,int *nz,int **idx,PetscScalar **v)
1438: {
1440:   return(0);
1441: }

1443: /* ------------------------------------------------------------------ */

1447: PetscErrorCode MatSetUpPreallocation_MPIRowbs(Mat A)
1448: {

1452:    MatMPIRowbsSetPreallocation(A,PETSC_DEFAULT,0);
1453:   return(0);
1454: }

1456: /* -------------------------------------------------------------------*/
1457: static struct _MatOps MatOps_Values = {MatSetValues_MPIRowbs,
1458:        MatGetRow_MPIRowbs,
1459:        MatRestoreRow_MPIRowbs,
1460:        MatMult_MPIRowbs,
1461: /* 4*/ MatMultAdd_MPIRowbs,
1462:        MatMult_MPIRowbs,
1463:        MatMultAdd_MPIRowbs,
1464:        MatSolve_MPIRowbs,
1465:        0,
1466:        0,
1467: /*10*/ 0,
1468:        0,
1469:        0,
1470:        0,
1471:        0,
1472: /*15*/ MatGetInfo_MPIRowbs,
1473:        0,
1474:        MatGetDiagonal_MPIRowbs,
1475:        0,
1476:        MatNorm_MPIRowbs,
1477: /*20*/ MatAssemblyBegin_MPIRowbs,
1478:        MatAssemblyEnd_MPIRowbs,
1479:        0,
1480:        MatSetOption_MPIRowbs,
1481:        MatZeroEntries_MPIRowbs,
1482: /*25*/ MatZeroRows_MPIRowbs,
1483:        0,
1484:        MatLUFactorNumeric_MPIRowbs,
1485:        0,
1486:        MatCholeskyFactorNumeric_MPIRowbs,
1487: /*30*/ MatSetUpPreallocation_MPIRowbs,
1488:        MatILUFactorSymbolic_MPIRowbs,
1489:        MatIncompleteCholeskyFactorSymbolic_MPIRowbs,
1490:        0,
1491:        0,
1492: /*35*/ 0,
1493:        MatForwardSolve_MPIRowbs,
1494:        MatBackwardSolve_MPIRowbs,
1495:        0,
1496:        0,
1497: /*40*/ 0,
1498:        MatGetSubMatrices_MPIRowbs,
1499:        0,
1500:        0,
1501:        0,
1502: /*45*/ 0,
1503:        MatScale_MPIRowbs,
1504:        0,
1505:        0,
1506:        0,
1507: /*50*/ 0,
1508:        0,
1509:        0,
1510:        0,
1511:        0,
1512: /*55*/ 0,
1513:        0,
1514:        0,
1515:        0,
1516:        0,
1517: /*60*/ MatGetSubMatrix_MPIRowbs,
1518:        MatDestroy_MPIRowbs,
1519:        MatView_MPIRowbs,
1520:        0,
1521:        MatUseScaledForm_MPIRowbs,
1522: /*65*/ MatScaleSystem_MPIRowbs,
1523:        MatUnScaleSystem_MPIRowbs,
1524:        0,
1525:        0,
1526:        0,
1527: /*70*/ 0,
1528:        0,
1529:        0,
1530:        0,
1531:        0,
1532: /*75*/ 0,
1533:        0,
1534:        0,
1535:        0,
1536:        0,
1537: /*80*/ 0,
1538:        0,
1539:        0,
1540:        0,
1541:        MatLoad_MPIRowbs,
1542: /*85*/ 0,
1543:        0,
1544:        0,
1545:        0,
1546:        0,
1547: /*90*/ 0,
1548:        0,
1549:        0,
1550:        0,
1551:        0,
1552: /*95*/ 0,
1553:        0,
1554:        0,
1555:        0};

1557: /* ------------------------------------------------------------------- */

1562: PetscErrorCode  MatMPIRowbsSetPreallocation_MPIRowbs(Mat mat,int nz,const int nnz[])
1563: {

1567:   mat->preallocated = PETSC_TRUE;
1568:   MatCreateMPIRowbs_local(mat,nz,nnz);
1569:   return(0);
1570: }

1573: /*MC
1574:    MATMPIROWBS - MATMPIROWBS = "mpirowbs" - A matrix type providing ILU and ICC for distributed sparse matrices for use
1575:    with the external package BlockSolve95.  If BlockSolve95 is installed (see the manual for instructions
1576:    on how to declare the existence of external packages), a matrix type can be constructed which invokes
1577:    BlockSolve95 preconditioners and solvers. 

1579:    Options Database Keys:
1580: . -mat_type mpirowbs - sets the matrix type to "mpirowbs" during a call to MatSetFromOptions()

1582:   Level: beginner

1584: .seealso: MatCreateMPIRowbs
1585: M*/

1590: PetscErrorCode  MatCreate_MPIRowbs(Mat A)
1591: {
1592:   Mat_MPIRowbs *a;
1593:   BSmapping    *bsmap;
1594:   BSoff_map    *bsoff;
1596:   int          *offset,m,M;
1597:   PetscTruth   flg1,flg3;
1598:   BSprocinfo   *bspinfo;
1599:   MPI_Comm     comm;
1600: 
1602:   comm = A->comm;

1604:   PetscNew(Mat_MPIRowbs,&a);
1605:   A->data               = (void*)a;
1606:   PetscMemcpy(A->ops,&MatOps_Values,sizeof(struct _MatOps));
1607:   A->factor             = 0;
1608:   A->mapping            = 0;
1609:   a->vecs_permscale     = PETSC_FALSE;
1610:   A->insertmode         = NOT_SET_VALUES;
1611:   a->blocksolveassembly = 0;
1612:   a->keepzeroedrows     = PETSC_FALSE;

1614:   MPI_Comm_rank(comm,&a->rank);
1615:   MPI_Comm_size(comm,&a->size);


1618:   PetscMapInitialize(comm,&A->rmap);
1619:   PetscMapInitialize(comm,&A->cmap);
1620:   m    = A->rmap.n;
1621:   M    = A->rmap.N;

1623:   PetscMalloc((A->rmap.n+1)*sizeof(int),&a->imax);
1624:   a->reallocs                      = 0;

1626:   /* build cache for off array entries formed */
1627:   MatStashCreate_Private(A->comm,1,&A->stash);
1628:   a->donotstash = PETSC_FALSE;

1630:   /* Initialize BlockSolve information */
1631:   a->A              = 0;
1632:   a->pA              = 0;
1633:   a->comm_pA  = 0;
1634:   a->fpA      = 0;
1635:   a->comm_fpA = 0;
1636:   a->alpha    = 1.0;
1637:   a->0;
1638:   a->failures = 0;
1639:   MPI_Comm_dup(A->comm,&(a->comm_mpirowbs));
1640:   VecCreateMPI(A->comm,A->rmap.n,A->rmap.N,&(a->diag));
1641:   VecDuplicate(a->diag,&(a->xwork));
1642:   PetscLogObjectParent(A,a->diag);  PetscLogObjectParent(A,a->xwork);
1643:   PetscLogObjectMemory(A,(A->rmap.n+1)*sizeof(PetscScalar));
1644:   bspinfo = BScreate_ctx();CHKERRBS(0);
1645:   a->procinfo = bspinfo;
1646:   BSctx_set_id(bspinfo,a->rank);CHKERRBS(0);
1647:   BSctx_set_np(bspinfo,a->size);CHKERRBS(0);
1648:   BSctx_set_ps(bspinfo,a->comm_mpirowbs);CHKERRBS(0);
1649:   BSctx_set_cs(bspinfo,INT_MAX);CHKERRBS(0);
1650:   BSctx_set_is(bspinfo,INT_MAX);CHKERRBS(0);
1651:   BSctx_set_ct(bspinfo,IDO);CHKERRBS(0);
1652: #if defined(PETSC_USE_DEBUG)
1653:   BSctx_set_err(bspinfo,1);CHKERRBS(0);  /* BS error checking */
1654: #endif
1655:   BSctx_set_rt(bspinfo,1);CHKERRBS(0);
1656: #if defined (PETSC_USE_INFO)
1657:   PetscOptionsHasName(PETSC_NULL,"-info",&flg1);
1658:   if (flg1) {
1659:     BSctx_set_pr(bspinfo,1);CHKERRBS(0);
1660:   }
1661: #endif
1662:   PetscOptionsBegin(A->comm,PETSC_NULL,"Options for MPIROWBS matrix","Mat");
1663:     PetscOptionsTruth("-pc_factor_factorpointwise","Do not optimize for inodes (slow)",PETSC_NULL,PETSC_FALSE,&flg1,PETSC_NULL);
1664:     PetscOptionsTruth("-mat_rowbs_no_inode","Do not optimize for inodes (slow)",PETSC_NULL,PETSC_FALSE,&flg3,PETSC_NULL);
1665:   PetscOptionsEnd();
1666:   if (flg1 || flg3) {
1667:     BSctx_set_si(bspinfo,1);CHKERRBS(0);
1668:   } else {
1669:     BSctx_set_si(bspinfo,0);CHKERRBS(0);
1670:   }
1671: #if defined(PETSC_USE_LOG)
1672:   MLOG_INIT();  /* Initialize logging */
1673: #endif

1675:   /* Compute global offsets */
1676:   offset = &A->rmap.rstart;

1678:   PetscNew(BSmapping,&a->bsmap);
1679:   PetscLogObjectMemory(A,sizeof(BSmapping));
1680:   bsmap = a->bsmap;
1681:   PetscMalloc(sizeof(int),&bsmap->vlocal2global);
1682:   *((int*)bsmap->vlocal2global) = (*offset);
1683:   bsmap->flocal2global                 = BSloc2glob;
1684:   bsmap->free_l2g                = 0;
1685:   PetscMalloc(sizeof(int),&bsmap->vglobal2local);
1686:   *((int*)bsmap->vglobal2local) = (*offset);
1687:   bsmap->fglobal2local                 = BSglob2loc;
1688:   bsmap->free_g2l                 = 0;
1689:   bsoff                          = BSmake_off_map(*offset,bspinfo,A->rmap.N);
1690:   bsmap->vglobal2proc                 = (void*)bsoff;
1691:   bsmap->fglobal2proc                 = BSglob2proc;
1692:   bsmap->free_g2p                = (void(*)(void*)) BSfree_off_map;
1693:   PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMPIRowbsSetPreallocation_C",
1694:                                     "MatMPIRowbsSetPreallocation_MPIRowbs",
1695:                                      MatMPIRowbsSetPreallocation_MPIRowbs);
1696:   PetscObjectChangeTypeName((PetscObject)A,MATMPIROWBS);
1697:   return(0);
1698: }

1703: /* @
1704:   MatMPIRowbsSetPreallocation - Sets the number of expected nonzeros 
1705:   per row in the matrix.

1707:   Input Parameter:
1708: +  mat - matrix
1709: .  nz - maximum expected for any row
1710: -  nzz - number expected in each row

1712:   Note:
1713:   This routine is valid only for matrices stored in the MATMPIROWBS
1714:   format.
1715: @ */
1716: PetscErrorCode  MatMPIRowbsSetPreallocation(Mat mat,int nz,const int nnz[])
1717: {
1718:   PetscErrorCode ierr,(*f)(Mat,int,const int[]);

1721:   PetscObjectQueryFunction((PetscObject)mat,"MatMPIRowbsSetPreallocation_C",(void (**)(void))&f);
1722:   if (f) {
1723:     (*f)(mat,nz,nnz);
1724:   }
1725:   return(0);
1726: }

1728: /* --------------- extra BlockSolve-specific routines -------------- */
1731: /* @
1732:   MatGetBSProcinfo - Gets the BlockSolve BSprocinfo context, which the
1733:   user can then manipulate to alter the default parameters.

1735:   Input Parameter:
1736:   mat - matrix

1738:   Output Parameter:
1739:   procinfo - processor information context

1741:   Note:
1742:   This routine is valid only for matrices stored in the MATMPIROWBS
1743:   format.
1744: @ */
1745: PetscErrorCode  MatGetBSProcinfo(Mat mat,BSprocinfo *procinfo)
1746: {
1747:   Mat_MPIRowbs *a = (Mat_MPIRowbs*)mat->data;
1748:   PetscTruth   ismpirowbs;

1752:   PetscTypeCompare((PetscObject)mat,MATMPIROWBS,&ismpirowbs);
1753:   if (!ismpirowbs) SETERRQ(PETSC_ERR_ARG_WRONG,"For MATMPIROWBS matrix type");
1754:   procinfo = a->procinfo;
1755:   return(0);
1756: }

1760: PetscErrorCode MatLoad_MPIRowbs(PetscViewer viewer,MatType type,Mat *newmat)
1761: {
1762:   Mat_MPIRowbs *a;
1763:   BSspmat      *A;
1764:   BSsprow      **rs;
1765:   Mat          mat;
1767:   int          i,nz,j,rstart,rend,fd,*ourlens,*sndcounts = 0,*procsnz;
1768:   int          header[4],rank,size,*rowlengths = 0,M,m,*rowners,maxnz,*cols;
1769:   PetscScalar  *vals;
1770:   MPI_Comm     comm = ((PetscObject)viewer)->comm;
1771:   MPI_Status   status;

1774:   MPI_Comm_size(comm,&size);
1775:   MPI_Comm_rank(comm,&rank);
1776:   if (!rank) {
1777:     PetscViewerBinaryGetDescriptor(viewer,&fd);
1778:     PetscBinaryRead(fd,(char *)header,4,PETSC_INT);
1779:     if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object");
1780:     if (header[3] < 0) {
1781:       SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format,cannot load as MPIRowbs");
1782:     }
1783:   }

1785:   MPI_Bcast(header+1,3,MPI_INT,0,comm);
1786:   M = header[1];

1788:   /* determine ownership of all rows */
1789:   m          = M/size + ((M % size) > rank);
1790:   PetscMalloc((size+2)*sizeof(int),&rowners);
1791:   MPI_Allgather(&m,1,MPI_INT,rowners+1,1,MPI_INT,comm);
1792:   rowners[0] = 0;
1793:   for (i=2; i<=size; i++) {
1794:     rowners[i] += rowners[i-1];
1795:   }
1796:   rstart = rowners[rank];
1797:   rend   = rowners[rank+1];

1799:   /* distribute row lengths to all processors */
1800:   PetscMalloc((rend-rstart)*sizeof(int),&ourlens);
1801:   if (!rank) {
1802:     PetscMalloc(M*sizeof(int),&rowlengths);
1803:     PetscBinaryRead(fd,rowlengths,M,PETSC_INT);
1804:     PetscMalloc(size*sizeof(int),&sndcounts);
1805:     for (i=0; i<size; i++) sndcounts[i] = rowners[i+1] - rowners[i];
1806:     MPI_Scatterv(rowlengths,sndcounts,rowners,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);
1807:     PetscFree(sndcounts);
1808:   } else {
1809:     MPI_Scatterv(0,0,0,MPI_INT,ourlens,rend-rstart,MPI_INT,0,comm);
1810:   }

1812:   /* create our matrix */
1813:   MatCreate(comm,newmat);
1814:   MatSetSizes(*newmat,m,m,M,M);
1815:   MatSetType(*newmat,type);
1816:   MatMPIRowbsSetPreallocation(*newmat,0,ourlens);
1817:   mat = *newmat;
1818:   PetscFree(ourlens);

1820:   a = (Mat_MPIRowbs*)mat->data;
1821:   A = a->A;
1822:   rs = A->rows;

1824:   if (!rank) {
1825:     /* calculate the number of nonzeros on each processor */
1826:     PetscMalloc(size*sizeof(int),&procsnz);
1827:     PetscMemzero(procsnz,size*sizeof(int));
1828:     for (i=0; i<size; i++) {
1829:       for (j=rowners[i]; j< rowners[i+1]; j++) {
1830:         procsnz[i] += rowlengths[j];
1831:       }
1832:     }
1833:     PetscFree(rowlengths);

1835:     /* determine max buffer needed and allocate it */
1836:     maxnz = 0;
1837:     for (i=0; i<size; i++) {
1838:       maxnz = PetscMax(maxnz,procsnz[i]);
1839:     }
1840:     PetscMalloc(maxnz*sizeof(int),&cols);

1842:     /* read in my part of the matrix column indices  */
1843:     nz = procsnz[0];
1844:     PetscBinaryRead(fd,cols,nz,PETSC_INT);
1845: 
1846:     /* insert it into my part of matrix */
1847:     nz = 0;
1848:     for (i=0; i<A->num_rows; i++) {
1849:       for (j=0; j<a->imax[i]; j++) {
1850:         rs[i]->col[j] = cols[nz++];
1851:       }
1852:       rs[i]->length = a->imax[i];
1853:     }
1854:     /* read in parts for all other processors */
1855:     for (i=1; i<size; i++) {
1856:       nz   = procsnz[i];
1857:       PetscBinaryRead(fd,cols,nz,PETSC_INT);
1858:       MPI_Send(cols,nz,MPI_INT,i,mat->tag,comm);
1859:     }
1860:     PetscFree(cols);
1861:     PetscMalloc(maxnz*sizeof(PetscScalar),&vals);

1863:     /* read in my part of the matrix numerical values  */
1864:     nz   = procsnz[0];
1865:     PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
1866: 
1867:     /* insert it into my part of matrix */
1868:     nz = 0;
1869:     for (i=0; i<A->num_rows; i++) {
1870:       for (j=0; j<a->imax[i]; j++) {
1871:         rs[i]->nz[j] = vals[nz++];
1872:       }
1873:     }
1874:     /* read in parts for all other processors */
1875:     for (i=1; i<size; i++) {
1876:       nz   = procsnz[i];
1877:       PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
1878:       MPI_Send(vals,nz,MPIU_SCALAR,i,mat->tag,comm);
1879:     }
1880:     PetscFree(vals);
1881:     PetscFree(procsnz);
1882:   } else {
1883:     /* determine buffer space needed for message */
1884:     nz = 0;
1885:     for (i=0; i<A->num_rows; i++) {
1886:       nz += a->imax[i];
1887:     }
1888:     PetscMalloc(nz*sizeof(int),&cols);

1890:     /* receive message of column indices*/
1891:     MPI_Recv(cols,nz,MPI_INT,0,mat->tag,comm,&status);
1892:     MPI_Get_count(&status,MPI_INT,&maxnz);
1893:     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong");

1895:     /* insert it into my part of matrix */
1896:     nz = 0;
1897:     for (i=0; i<A->num_rows; i++) {
1898:       for (j=0; j<a->imax[i]; j++) {
1899:         rs[i]->col[j] = cols[nz++];
1900:       }
1901:       rs[i]->length = a->imax[i];
1902:     }
1903:     PetscFree(cols);
1904:     PetscMalloc(nz*sizeof(PetscScalar),&vals);

1906:     /* receive message of values*/
1907:     MPI_Recv(vals,nz,MPIU_SCALAR,0,mat->tag,comm,&status);
1908:     MPI_Get_count(&status,MPIU_SCALAR,&maxnz);
1909:     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong");

1911:     /* insert it into my part of matrix */
1912:     nz = 0;
1913:     for (i=0; i<A->num_rows; i++) {
1914:       for (j=0; j<a->imax[i]; j++) {
1915:         rs[i]->nz[j] = vals[nz++];
1916:       }
1917:       rs[i]->length = a->imax[i];
1918:     }
1919:     PetscFree(vals);
1920:   }
1921:   PetscFree(rowners);
1922:   a->nz = a->maxnz;
1923:   MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
1924:   MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);
1925:   return(0);
1926: }

1928: /* 
1929:     Special destroy and view routines for factored matrices 
1930: */
1933: static PetscErrorCode MatDestroy_MPIRowbs_Factored(Mat mat)
1934: {
1936: #if defined(PETSC_USE_LOG)
1937:   PetscLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d",mat->rmap.N,mat->cmap.N);
1938: #endif
1939:   return(0);
1940: }

1944: static PetscErrorCode MatView_MPIRowbs_Factored(Mat mat,PetscViewer viewer)
1945: {

1949:   MatView((Mat) mat->data,viewer);
1950:   return(0);
1951: }

1955: PetscErrorCode MatIncompleteCholeskyFactorSymbolic_MPIRowbs(Mat mat,IS isrow,MatFactorInfo *info,Mat *newfact)
1956: {
1957:   /* Note:  f is not currently used in BlockSolve */
1958:   Mat          newmat;
1959:   Mat_MPIRowbs *mbs = (Mat_MPIRowbs*)mat->data;
1961:   PetscTruth   idn;

1964:   if (isrow) {
1965:     ISIdentity(isrow,&idn);
1966:     if (!idn) SETERRQ(PETSC_ERR_SUP,"Only identity row permutation supported");
1967:   }

1969:   if (!mat->symmetric) {
1970:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"To use incomplete Cholesky \n\
1971:         preconditioning with a MATMPIROWBS matrix you must declare it to be \n\
1972:         symmetric using the option MatSetOption(A,MAT_SYMMETRIC)");
1973:   }

1975:   /* If the icc_storage flag wasn't set before the last blocksolveassembly,          */
1976:   /* we must completely redo the assembly as a different storage format is required. */
1977:   if (mbs->blocksolveassembly && !mbs->assembled_icc_storage) {
1978:     mat->same_nonzero       = PETSC_FALSE;
1979:     mbs->blocksolveassembly = 0;
1980:   }

1982:   if (!mbs->blocksolveassembly) {
1983:     BSset_mat_icc_storage(mbs->A,PETSC_TRUE);CHKERRBS(0);
1984:     BSset_mat_symmetric(mbs->A,PETSC_TRUE);CHKERRBS(0);
1985:     MatAssemblyEnd_MPIRowbs_ForBlockSolve(mat);
1986:   }

1988:   /* Copy permuted matrix */
1989:   if (mbs->fpA) {BSfree_copy_par_mat(mbs->fpA);CHKERRBS(0);}
1990:   mbs->fpA = BScopy_par_mat(mbs->pA);CHKERRBS(0);

1992:   /* Set up the communication for factorization */
1993:   if (mbs->comm_fpA) {BSfree_comm(mbs->comm_fpA);CHKERRBS(0);}
1994:   mbs->comm_fpA = BSsetup_factor(mbs->fpA,mbs->procinfo);CHKERRBS(0);

1996:   /* 
1997:       Create a new Mat structure to hold the "factored" matrix, 
1998:     not this merely contains a pointer to the original matrix, since
1999:     the original matrix contains the factor information.
2000:   */
2001:   PetscHeaderCreate(newmat,_p_Mat,struct _MatOps,MAT_COOKIE,-1,"Mat",mat->comm,MatDestroy,MatView);
2002:   PetscLogObjectMemory(newmat,sizeof(struct _p_Mat));

2004:   newmat->data         = (void*)mat;
2005:   PetscMemcpy(newmat->ops,&MatOps_Values,sizeof(struct _MatOps));
2006:   newmat->ops->destroy = MatDestroy_MPIRowbs_Factored;
2007:   newmat->ops->view    = MatView_MPIRowbs_Factored;
2008:   newmat->factor       = 1;
2009:   newmat->preallocated = PETSC_TRUE;
2010:   PetscMapCopy(mat->comm,&mat->rmap,&newmat->rmap);
2011:   PetscMapCopy(mat->comm,&mat->cmap,&newmat->cmap);

2013:   PetscStrallocpy(MATMPIROWBS,&newmat->type_name);

2015:   *newfact = newmat;
2016:   return(0);
2017: }

2021: PetscErrorCode MatILUFactorSymbolic_MPIRowbs(Mat mat,IS isrow,IS iscol,MatFactorInfo* info,Mat *newfact)
2022: {
2023:   Mat          newmat;
2024:   Mat_MPIRowbs *mbs = (Mat_MPIRowbs*)mat->data;
2026:   PetscTruth   idn;

2029:   if (info->levels) SETERRQ(PETSC_ERR_SUP,"Blocksolve ILU only supports 0 fill");
2030:   if (isrow) {
2031:     ISIdentity(isrow,&idn);
2032:     if (!idn) SETERRQ(PETSC_ERR_SUP,"Only identity row permutation supported");
2033:   }
2034:   if (iscol) {
2035:     ISIdentity(iscol,&idn);
2036:     if (!idn) SETERRQ(PETSC_ERR_SUP,"Only identity column permutation supported");
2037:   }

2039:   if (!mbs->blocksolveassembly) {
2040:     MatAssemblyEnd_MPIRowbs_ForBlockSolve(mat);
2041:   }
2042: 
2043: /*   if (mat->symmetric) { */
2044: /*     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"To use ILU preconditioner with \n\ */
2045: /*         MatCreateMPIRowbs() matrix you CANNOT declare it to be a symmetric matrix\n\ */
2046: /*         using the option MatSetOption(A,MAT_SYMMETRIC)"); */
2047: /*   } */

2049:   /* Copy permuted matrix */
2050:   if (mbs->fpA) {BSfree_copy_par_mat(mbs->fpA);CHKERRBS(0);}
2051:   mbs->fpA = BScopy_par_mat(mbs->pA);CHKERRBS(0);

2053:   /* Set up the communication for factorization */
2054:   if (mbs->comm_fpA) {BSfree_comm(mbs->comm_fpA);CHKERRBS(0);}
2055:   mbs->comm_fpA = BSsetup_factor(mbs->fpA,mbs->procinfo);CHKERRBS(0);

2057:   /* 
2058:       Create a new Mat structure to hold the "factored" matrix,
2059:     not this merely contains a pointer to the original matrix, since
2060:     the original matrix contains the factor information.
2061:   */
2062:   PetscHeaderCreate(newmat,_p_Mat,struct _MatOps,MAT_COOKIE,-1,"Mat",mat->comm,MatDestroy,MatView);
2063:   PetscLogObjectMemory(newmat,sizeof(struct _p_Mat));

2065:   newmat->data         = (void*)mat;
2066:   PetscMemcpy(newmat->ops,&MatOps_Values,sizeof(struct _MatOps));
2067:   newmat->ops->destroy = MatDestroy_MPIRowbs_Factored;
2068:   newmat->ops->view    = MatView_MPIRowbs_Factored;
2069:   newmat->factor       = 1;
2070:   newmat->preallocated = PETSC_TRUE;

2072:   PetscMapCopy(mat->comm,&mat->rmap,&newmat->rmap);
2073:   PetscMapCopy(mat->comm,&mat->cmap,&newmat->cmap);

2075:   PetscStrallocpy(MATMPIROWBS,&newmat->type_name);

2077:   *newfact = newmat;
2078:   return(0);
2079: }

2083: /*@C
2084:    MatCreateMPIRowbs - Creates a sparse parallel matrix in the MATMPIROWBS
2085:    format.  This format is intended primarily as an interface for BlockSolve95.

2087:    Collective on MPI_Comm

2089:    Input Parameters:
2090: +  comm - MPI communicator
2091: .  m - number of local rows (or PETSC_DECIDE to have calculated)
2092: .  M - number of global rows (or PETSC_DECIDE to have calculated)
2093: .  nz - number of nonzeros per row (same for all local rows)
2094: -  nnz - number of nonzeros per row (possibly different for each row).

2096:    Output Parameter:
2097: .  newA - the matrix 

2099:    Notes:
2100:    If PETSC_DECIDE or  PETSC_DETERMINE is used for a particular argument on one processor
2101:    than it must be used on all processors that share the object for that argument.

2103:    The user MUST specify either the local or global matrix dimensions
2104:    (possibly both).

2106:    Specify the preallocated storage with either nz or nnz (not both).  Set 
2107:    nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 
2108:    allocation.

2110:    Notes:
2111:    By default, the matrix is assumed to be nonsymmetric; the user can
2112:    take advantage of special optimizations for symmetric matrices by calling
2113: $     MatSetOption(mat,MAT_SYMMETRIC)
2114: $     MatSetOption(mat,MAT_SYMMETRY_ETERNAL)
2115:    BEFORE calling the routine MatAssemblyBegin().

2117:    Internally, the MATMPIROWBS format inserts zero elements to the
2118:    matrix if necessary, so that nonsymmetric matrices are considered
2119:    to be symmetric in terms of their sparsity structure; this format
2120:    is required for use of the parallel communication routines within
2121:    BlockSolve95. In particular, if the matrix element A[i,j] exists,
2122:    then PETSc will internally allocate a 0 value for the element
2123:    A[j,i] during MatAssemblyEnd() if the user has not already set
2124:    a value for the matrix element A[j,i].

2126:    Options Database Keys:
2127: .  -mat_rowbs_no_inode - Do not use inodes.

2129:    Level: intermediate
2130:   
2131: .keywords: matrix, row, symmetric, sparse, parallel, BlockSolve

2133: .seealso: MatCreate(), MatSetValues()
2134: @*/
2135: PetscErrorCode  MatCreateMPIRowbs(MPI_Comm comm,int m,int M,int nz,const int nnz[],Mat *newA)
2136: {
2138: 
2140:   MatCreate(comm,newA);
2141:   MatSetSizes(*newA,m,m,M,M);
2142:   MatSetType(*newA,MATMPIROWBS);
2143:   MatMPIRowbsSetPreallocation(*newA,nz,nnz);
2144:   return(0);
2145: }


2148: /* -------------------------------------------------------------------------*/

2150:  #include src/mat/impls/aij/seq/aij.h
2151:  #include src/mat/impls/aij/mpi/mpiaij.h

2155: PetscErrorCode MatGetSubMatrices_MPIRowbs(Mat C,int ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submat[])
2156: {
2158:   int         nmax,nstages_local,nstages,i,pos,max_no;


2162:   /* Allocate memory to hold all the submatrices */
2163:   if (scall != MAT_REUSE_MATRIX) {
2164:     PetscMalloc((ismax+1)*sizeof(Mat),submat);
2165:   }
2166: 
2167:   /* Determine the number of stages through which submatrices are done */
2168:   nmax          = 20*1000000 / (C->cmap.N * sizeof(int));
2169:   if (!nmax) nmax = 1;
2170:   nstages_local = ismax/nmax + ((ismax % nmax)?1:0);

2172:   /* Make sure every processor loops through the nstages */
2173:   MPI_Allreduce(&nstages_local,&nstages,1,MPI_INT,MPI_MAX,C->comm);

2175:   for (i=0,pos=0; i<nstages; i++) {
2176:     if (pos+nmax <= ismax) max_no = nmax;
2177:     else if (pos == ismax) max_no = 0;
2178:     else                   max_no = ismax-pos;
2179:     MatGetSubMatrices_MPIRowbs_Local(C,max_no,isrow+pos,iscol+pos,scall,*submat+pos);
2180:     pos += max_no;
2181:   }
2182:   return(0);
2183: }
2184: /* -------------------------------------------------------------------------*/
2185: /* for now MatGetSubMatrices_MPIRowbs_Local get MPIAij submatrices of input
2186:    matrix and preservs zeroes from structural symetry
2187:  */
2190: PetscErrorCode MatGetSubMatrices_MPIRowbs_Local(Mat C,int ismax,const IS isrow[],const IS iscol[],MatReuse scall,Mat *submats)
2191: {
2192:   Mat_MPIRowbs  *c = (Mat_MPIRowbs *)(C->data);
2193:   BSspmat       *A = c->A;
2194:   Mat_SeqAIJ    *mat;
2196:   int         **irow,**icol,*nrow,*ncol,*w1,*w2,*w3,*w4,*rtable,start,end,size;
2197:   int         **sbuf1,**sbuf2,rank,m,i,j,k,l,ct1,ct2,**rbuf1,row,proc;
2198:   int         nrqs,msz,**ptr,idx,*req_size,*ctr,*pa,*tmp,tcol,nrqr;
2199:   int         **rbuf3,*req_source,**sbuf_aj,**rbuf2,max1,max2,**rmap;
2200:   int         **cmap,**lens,is_no,ncols,*cols,mat_i,*mat_j,tmp2,jmax,*irow_i;
2201:   int         len,ctr_j,*sbuf1_j,*sbuf_aj_i,*rbuf1_i,kmax,*cmap_i,*lens_i;
2202:   int         *rmap_i,tag0,tag1,tag2,tag3;
2203:   MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3;
2204:   MPI_Request *r_waits4,*s_waits3,*s_waits4;
2205:   MPI_Status  *r_status1,*r_status2,*s_status1,*s_status3,*s_status2;
2206:   MPI_Status  *r_status3,*r_status4,*s_status4;
2207:   MPI_Comm    comm;
2208:   FLOAT       **rbuf4,**sbuf_aa,*vals,*sbuf_aa_i;
2209:   PetscScalar *mat_a;
2210:   PetscTruth  sorted;
2211:   int         *onodes1,*olengths1;

2214:   comm   = C->comm;
2215:   tag0   = C->tag;
2216:   size   = c->size;
2217:   rank   = c->rank;
2218:   m      = C->rmap.N;
2219: 
2220:   /* Get some new tags to keep the communication clean */
2221:   PetscObjectGetNewTag((PetscObject)C,&tag1);
2222:   PetscObjectGetNewTag((PetscObject)C,&tag2);
2223:   PetscObjectGetNewTag((PetscObject)C,&tag3);

2225:     /* Check if the col indices are sorted */
2226:   for (i=0; i<ismax; i++) {
2227:     ISSorted(isrow[i],&sorted);
2228:     if (!sorted) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"ISrow is not sorted");
2229:     ISSorted(iscol[i],&sorted);
2230:     /*    if (!sorted) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IScol is not sorted"); */
2231:   }

2233:   len    = (2*ismax+1)*(sizeof(int*)+ sizeof(int)) + (m+1)*sizeof(int);
2234:   PetscMalloc(len,&irow);
2235:   icol   = irow + ismax;
2236:   nrow   = (int*)(icol + ismax);
2237:   ncol   = nrow + ismax;
2238:   rtable = ncol + ismax;

2240:   for (i=0; i<ismax; i++) {
2241:     ISGetIndices(isrow[i],&irow[i]);
2242:     ISGetIndices(iscol[i],&icol[i]);
2243:     ISGetLocalSize(isrow[i],&nrow[i]);
2244:     ISGetLocalSize(iscol[i],&ncol[i]);
2245:   }

2247:   /* Create hash table for the mapping :row -> proc*/
2248:   for (i=0,j=0; i<size; i++) {
2249:     jmax = C->rmap.range[i+1];
2250:     for (; j<jmax; j++) {
2251:       rtable[j] = i;
2252:     }
2253:   }

2255:   /* evaluate communication - mesg to who, length of mesg, and buffer space
2256:      required. Based on this, buffers are allocated, and data copied into them*/
2257:   PetscMalloc(size*4*sizeof(int),&w1); /* mesg size */
2258:   w2     = w1 + size;      /* if w2[i] marked, then a message to proc i*/
2259:   w3     = w2 + size;      /* no of IS that needs to be sent to proc i */
2260:   w4     = w3 + size;      /* temp work space used in determining w1, w2, w3 */
2261:   PetscMemzero(w1,size*3*sizeof(int)); /* initialize work vector*/
2262:   for (i=0; i<ismax; i++) {
2263:     PetscMemzero(w4,size*sizeof(int)); /* initialize work vector*/
2264:     jmax   = nrow[i];
2265:     irow_i = irow[i];
2266:     for (j=0; j<jmax; j++) {
2267:       row  = irow_i[j];
2268:       proc = rtable[row];
2269:       w4[proc]++;
2270:     }
2271:     for (j=0; j<size; j++) {
2272:       if (w4[j]) { w1[j] += w4[j];  w3[j]++;}
2273:     }
2274:   }
2275: 
2276:   nrqs     = 0;              /* no of outgoing messages */
2277:   msz      = 0;              /* total mesg length (for all procs) */
2278:   w1[rank] = 0;              /* no mesg sent to self */
2279:   w3[rank] = 0;
2280:   for (i=0; i<size; i++) {
2281:     if (w1[i])  { w2[i] = 1; nrqs++;} /* there exists a message to proc i */
2282:   }
2283:   PetscMalloc((nrqs+1)*sizeof(int),&pa); /*(proc -array)*/
2284:   for (i=0,j=0; i<size; i++) {
2285:     if (w1[i]) { pa[j] = i; j++; }
2286:   }

2288:   /* Each message would have a header = 1 + 2*(no of IS) + data */
2289:   for (i=0; i<nrqs; i++) {
2290:     j     = pa[i];
2291:     w1[j] += w2[j] + 2* w3[j];
2292:     msz   += w1[j];
2293:   }

2295:   /* Determine the number of messages to expect, their lengths, from from-ids */
2296:   PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);
2297:   PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);

2299:   /* Now post the Irecvs corresponding to these messages */
2300:   PetscPostIrecvInt(comm,tag0,nrqr,onodes1,olengths1,&rbuf1,&r_waits1);
2301: 
2302:   PetscFree(onodes1);
2303:   PetscFree(olengths1);
2304: 
2305:   /* Allocate Memory for outgoing messages */
2306:   len      = 2*size*sizeof(int*) + 2*msz*sizeof(int) + size*sizeof(int);
2307:   PetscMalloc(len,&sbuf1);
2308:   ptr      = sbuf1 + size;   /* Pointers to the data in outgoing buffers */
2309:   PetscMemzero(sbuf1,2*size*sizeof(int*));
2310:   /* allocate memory for outgoing data + buf to receive the first reply */
2311:   tmp      = (int*)(ptr + size);
2312:   ctr      = tmp + 2*msz;

2314:   {
2315:     int *iptr = tmp,ict = 0;
2316:     for (i=0; i<nrqs; i++) {
2317:       j         = pa[i];
2318:       iptr     += ict;
2319:       sbuf1[j]  = iptr;
2320:       ict       = w1[j];
2321:     }
2322:   }

2324:   /* Form the outgoing messages */
2325:   /* Initialize the header space */
2326:   for (i=0; i<nrqs; i++) {
2327:     j           = pa[i];
2328:     sbuf1[j][0] = 0;
2329:     PetscMemzero(sbuf1[j]+1,2*w3[j]*sizeof(int));
2330:     ptr[j]      = sbuf1[j] + 2*w3[j] + 1;
2331:   }
2332: 
2333:   /* Parse the isrow and copy data into outbuf */
2334:   for (i=0; i<ismax; i++) {
2335:     PetscMemzero(ctr,size*sizeof(int));
2336:     irow_i = irow[i];
2337:     jmax   = nrow[i];
2338:     for (j=0; j<jmax; j++) {  /* parse the indices of each IS */
2339:       row  = irow_i[j];
2340:       proc = rtable[row];
2341:       if (proc != rank) { /* copy to the outgoing buf*/
2342:         ctr[proc]++;
2343:         *ptr[proc] = row;
2344:         ptr[proc]++;
2345:       }
2346:     }
2347:     /* Update the headers for the current IS */
2348:     for (j=0; j<size; j++) { /* Can Optimise this loop too */
2349:       if ((ctr_j = ctr[j])) {
2350:         sbuf1_j        = sbuf1[j];
2351:         k              = ++sbuf1_j[0];
2352:         sbuf1_j[2*k]   = ctr_j;
2353:         sbuf1_j[2*k-1] = i;
2354:       }
2355:     }
2356:   }

2358:   /*  Now  post the sends */
2359:   PetscMalloc((nrqs+1)*sizeof(MPI_Request),&s_waits1);
2360:   for (i=0; i<nrqs; ++i) {
2361:     j    = pa[i];
2362:     MPI_Isend(sbuf1[j],w1[j],MPI_INT,j,tag0,comm,s_waits1+i);
2363:   }

2365:   /* Post Receives to capture the buffer size */
2366:   PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits2);
2367:   PetscMalloc((nrqs+1)*sizeof(int*),&rbuf2);
2368:   rbuf2[0] = tmp + msz;
2369:   for (i=1; i<nrqs; ++i) {
2370:     rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]];
2371:   }
2372:   for (i=0; i<nrqs; ++i) {
2373:     j    = pa[i];
2374:     MPI_Irecv(rbuf2[i],w1[j],MPI_INT,j,tag1,comm,r_waits2+i);
2375:   }

2377:   /* Send to other procs the buf size they should allocate */
2378: 

2380:   /* Receive messages*/
2381:   PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits2);
2382:   PetscMalloc((nrqr+1)*sizeof(MPI_Status),&r_status1);
2383:   len         = 2*nrqr*sizeof(int) + (nrqr+1)*sizeof(int*);
2384:   PetscMalloc(len,&sbuf2);
2385:   req_size    = (int*)(sbuf2 + nrqr);
2386:   req_source  = req_size + nrqr;
2387: 
2388:   {
2389:     BSsprow    **sAi = A->rows;
2390:     int        id,rstart = C->rmap.rstart;
2391:     int        *sbuf2_i;

2393:     for (i=0; i<nrqr; ++i) {
2394:       MPI_Waitany(nrqr,r_waits1,&idx,r_status1+i);
2395:       req_size[idx]   = 0;
2396:       rbuf1_i         = rbuf1[idx];
2397:       start           = 2*rbuf1_i[0] + 1;
2398:       MPI_Get_count(r_status1+i,MPI_INT,&end);
2399:       PetscMalloc((end+1)*sizeof(int),&sbuf2[idx]);
2400:       sbuf2_i         = sbuf2[idx];
2401:       for (j=start; j<end; j++) {
2402:         id               = rbuf1_i[j] - rstart;
2403:         ncols            = (sAi[id])->length;
2404:         sbuf2_i[j]       = ncols;
2405:         req_size[idx]   += ncols;
2406:       }
2407:       req_source[idx] = r_status1[i].MPI_SOURCE;
2408:       /* form the header */
2409:       sbuf2_i[0]   = req_size[idx];
2410:       for (j=1; j<start; j++) { sbuf2_i[j] = rbuf1_i[j]; }
2411:       MPI_Isend(sbuf2_i,end,MPI_INT,req_source[idx],tag1,comm,s_waits2+i);
2412:     }
2413:   }
2414:   PetscFree(r_status1);
2415:   PetscFree(r_waits1);

2417:   /*  recv buffer sizes */
2418:   /* Receive messages*/
2419: 
2420:   PetscMalloc((nrqs+1)*sizeof(int*),&rbuf3);
2421:   PetscMalloc((nrqs+1)*sizeof(FLOAT *),&rbuf4);
2422:   PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits3);
2423:   PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits4);
2424:   PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status2);

2426:   for (i=0; i<nrqs; ++i) {
2427:     MPI_Waitany(nrqs,r_waits2,&idx,r_status2+i);
2428:     PetscMalloc((rbuf2[idx][0]+1)*sizeof(int),&rbuf3[idx]);
2429:     PetscMalloc((rbuf2[idx][0]+1)*sizeof(FLOAT),&rbuf4[idx]);
2430:     MPI_Irecv(rbuf3[idx],rbuf2[idx][0],MPI_INT,r_status2[i].MPI_SOURCE,tag2,comm,r_waits3+idx);
2431:     MPI_Irecv(rbuf4[idx],rbuf2[idx][0],MPIU_SCALAR,r_status2[i].MPI_SOURCE,tag3,comm,r_waits4+idx);
2432:   }
2433:   PetscFree(r_status2);
2434:   PetscFree(r_waits2);
2435: 
2436:   /* Wait on sends1 and sends2 */
2437:   PetscMalloc((nrqs+1)*sizeof(MPI_Status),&s_status1);
2438:   PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status2);

2440:   if (nrqs) {MPI_Waitall(nrqs,s_waits1,s_status1);}
2441:   if (nrqr) {MPI_Waitall(nrqr,s_waits2,s_status2);}
2442:   PetscFree(s_status1);
2443:   PetscFree(s_status2);
2444:   PetscFree(s_waits1);
2445:   PetscFree(s_waits2);

2447:   /* Now allocate buffers for a->j, and send them off */
2448:   PetscMalloc((nrqr+1)*sizeof(int*),&sbuf_aj);
2449:   for (i=0,j=0; i<nrqr; i++) j += req_size[i];
2450:   PetscMalloc((j+1)*sizeof(int),&sbuf_aj[0]);
2451:   for (i=1; i<nrqr; i++)  sbuf_aj[i] = sbuf_aj[i-1] + req_size[i-1];
2452: 
2453:   PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits3);
2454:   {
2455:     BSsprow *brow;
2456:     int *Acol;
2457:     int rstart = C->rmap.rstart;

2459:     for (i=0; i<nrqr; i++) {
2460:       rbuf1_i   = rbuf1[i];
2461:       sbuf_aj_i = sbuf_aj[i];
2462:       ct1       = 2*rbuf1_i[0] + 1;
2463:       ct2       = 0;
2464:       for (j=1,max1=rbuf1_i[0]; j<=max1; j++) {
2465:         kmax = rbuf1[i][2*j];
2466:         for (k=0; k<kmax; k++,ct1++) {
2467:           brow   = A->rows[rbuf1_i[ct1] - rstart];
2468:           ncols  = brow->length;
2469:           Acol   = brow->col;
2470:           /* load the column indices for this row into cols*/
2471:           cols  = sbuf_aj_i + ct2;
2472:           PetscMemcpy(cols,Acol,ncols*sizeof(int));
2473:           /*for (l=0; l<ncols;l++) cols[l]=Acol[l]; */ /* How is it with
2474:                                                           mappings?? */
2475:           ct2 += ncols;
2476:         }
2477:       }
2478:       MPI_Isend(sbuf_aj_i,req_size[i],MPI_INT,req_source[i],tag2,comm,s_waits3+i);
2479:     }
2480:   }
2481:   PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status3);
2482:   PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status3);

2484:   /* Allocate buffers for a->a, and send them off */
2485:   PetscMalloc((nrqr+1)*sizeof(FLOAT*),&sbuf_aa);
2486:   for (i=0,j=0; i<nrqr; i++) j += req_size[i];
2487:   PetscMalloc((j+1)*sizeof(FLOAT),&sbuf_aa[0]);
2488:   for (i=1; i<nrqr; i++)  sbuf_aa[i] = sbuf_aa[i-1] + req_size[i-1];
2489: 
2490:   PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits4);
2491:   {
2492:     BSsprow *brow;
2493:     FLOAT *Aval;
2494:     int rstart = C->rmap.rstart;
2495: 
2496:     for (i=0; i<nrqr; i++) {
2497:       rbuf1_i   = rbuf1[i];
2498:       sbuf_aa_i = sbuf_aa[i];
2499:       ct1       = 2*rbuf1_i[0]+1;
2500:       ct2       = 0;
2501:       for (j=1,max1=rbuf1_i[0]; j<=max1; j++) {
2502:         kmax = rbuf1_i[2*j];
2503:         for (k=0; k<kmax; k++,ct1++) {
2504:           brow  = A->rows[rbuf1_i[ct1] - rstart];
2505:           ncols = brow->length;
2506:           Aval  = brow->nz;
2507:           /* load the column values for this row into vals*/
2508:           vals  = sbuf_aa_i+ct2;
2509:           PetscMemcpy(vals,Aval,ncols*sizeof(FLOAT));
2510:           ct2 += ncols;
2511:         }
2512:       }
2513:       MPI_Isend(sbuf_aa_i,req_size[i],MPIU_SCALAR,req_source[i],tag3,comm,s_waits4+i);
2514:     }
2515:   }
2516:   PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status4);
2517:   PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status4);
2518:   PetscFree(rbuf1);

2520:   /* Form the matrix */
2521:   /* create col map */
2522:   {
2523:     int *icol_i;
2524: 
2525:     len     = (1+ismax)*sizeof(int*)+ ismax*C->cmap.N*sizeof(int);
2526:     PetscMalloc(len,&cmap);
2527:     cmap[0] = (int*)(cmap + ismax);
2528:     PetscMemzero(cmap[0],(1+ismax*C->cmap.N)*sizeof(int));
2529:     for (i=1; i<ismax; i++) { cmap[i] = cmap[i-1] + C->cmap.N; }
2530:     for (i=0; i<ismax; i++) {
2531:       jmax   = ncol[i];
2532:       icol_i = icol[i];
2533:       cmap_i = cmap[i];
2534:       for (j=0; j<jmax; j++) {
2535:         cmap_i[icol_i[j]] = j+1;
2536:       }
2537:     }
2538:   }

2540:   /* Create lens which is required for MatCreate... */
2541:   for (i=0,j=0; i<ismax; i++) { j += nrow[i]; }
2542:   len     = (1+ismax)*sizeof(int*)+ j*sizeof(int);
2543:   PetscMalloc(len,&lens);
2544:   lens[0] = (int*)(lens + ismax);
2545:   PetscMemzero(lens[0],j*sizeof(int));
2546:   for (i=1; i<ismax; i++) { lens[i] = lens[i-1] + nrow[i-1]; }
2547: 
2548:   /* Update lens from local data */
2549:   { BSsprow *Arow;
2550:     for (i=0; i<ismax; i++) {
2551:       jmax   = nrow[i];
2552:       cmap_i = cmap[i];
2553:       irow_i = irow[i];
2554:       lens_i = lens[i];
2555:       for (j=0; j<jmax; j++) {
2556:         row  = irow_i[j];
2557:         proc = rtable[row];
2558:         if (proc == rank) {
2559:           Arow=A->rows[row-C->rmap.rstart];
2560:           ncols=Arow->length;
2561:           cols=Arow->col;
2562:           for (k=0; k<ncols; k++) {
2563:             if (cmap_i[cols[k]]) { lens_i[j]++;}
2564:           }
2565:         }
2566:       }
2567:     }
2568:   }
2569: 
2570:   /* Create row map*/
2571:   len     = (1+ismax)*sizeof(int*)+ ismax*C->rmap.N*sizeof(int);
2572:   PetscMalloc(len,&rmap);
2573:   rmap[0] = (int*)(rmap + ismax);
2574:   PetscMemzero(rmap[0],ismax*C->rmap.N*sizeof(int));
2575:   for (i=1; i<ismax; i++) { rmap[i] = rmap[i-1] + C->rmap.N;}
2576:   for (i=0; i<ismax; i++) {
2577:     rmap_i = rmap[i];
2578:     irow_i = irow[i];
2579:     jmax   = nrow[i];
2580:     for (j=0; j<jmax; j++) {
2581:       rmap_i[irow_i[j]] = j;
2582:     }
2583:   }
2584: 
2585:   /* Update lens from offproc data */
2586:   {
2587:     int *rbuf2_i,*rbuf3_i,*sbuf1_i;

2589:     for (tmp2=0; tmp2<nrqs; tmp2++) {
2590:       MPI_Waitany(nrqs,r_waits3,&i,r_status3+tmp2);
2591:       idx     = pa[i];
2592:       sbuf1_i = sbuf1[idx];
2593:       jmax    = sbuf1_i[0];
2594:       ct1     = 2*jmax+1;
2595:       ct2     = 0;
2596:       rbuf2_i = rbuf2[i];
2597:       rbuf3_i = rbuf3[i];
2598:       for (j=1; j<=jmax; j++) {
2599:         is_no   = sbuf1_i[2*j-1];
2600:         max1    = sbuf1_i[2*j];
2601:         lens_i  = lens[is_no];
2602:         cmap_i  = cmap[is_no];
2603:         rmap_i  = rmap[is_no];
2604:         for (k=0; k<max1; k++,ct1++) {
2605:           row  = rmap_i[sbuf1_i[ct1]]; /* the val in the new matrix to be */
2606:           max2 = rbuf2_i[ct1];
2607:           for (l=0; l<max2; l++,ct2++) {
2608:             if (cmap_i[rbuf3_i[ct2]]) {
2609:               lens_i[row]++;
2610:             }
2611:           }
2612:         }
2613:       }
2614:     }
2615:   }
2616:   PetscFree(r_status3);
2617:   PetscFree(r_waits3);
2618:   if (nrqr) {MPI_Waitall(nrqr,s_waits3,s_status3);}
2619:   PetscFree(s_status3);
2620:   PetscFree(s_waits3);

2622:   /* Create the submatrices */
2623:   if (scall == MAT_REUSE_MATRIX) {
2624:     PetscTruth same;
2625: 
2626:     /*
2627:         Assumes new rows are same length as the old rows,hence bug!
2628:     */
2629:     for (i=0; i<ismax; i++) {
2630:       PetscTypeCompare((PetscObject)(submats[i]),MATSEQAIJ,&same);
2631:       if (!same) {
2632:         SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong type");
2633:       }
2634:       mat = (Mat_SeqAIJ*)(submats[i]->data);
2635:       if ((submats[i]->rmap.n != nrow[i]) || (submats[i]->cmap.n != ncol[i])) {
2636:         SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size");
2637:       }
2638:       PetscMemcmp(mat->ilen,lens[i],submats[i]->rmap.n*sizeof(int),&same);
2639:       if (!same) {
2640:         SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros");
2641:       }
2642:       /* Initial matrix as if empty */
2643:       PetscMemzero(mat->ilen,submats[i]->rmap.n*sizeof(int));
2644:       submats[i]->factor = C->factor;
2645:     }
2646:   } else {
2647:     for (i=0; i<ismax; i++) {
2648:       /* Here we want to explicitly generate SeqAIJ matrices */
2649:       MatCreate(PETSC_COMM_SELF,submats+i);
2650:       MatSetSizes(submats[i],nrow[i],ncol[i],nrow[i],ncol[i]);
2651:       MatSetType(submats[i],MATSEQAIJ);
2652:       MatSeqAIJSetPreallocation(submats[i],0,lens[i]);
2653:     }
2654:   }

2656:   /* Assemble the matrices */
2657:   /* First assemble the local rows */
2658:   {
2659:     int    ilen_row,*imat_ilen,*imat_j,*imat_i,old_row;
2660:     PetscScalar *imat_a;
2661:     BSsprow *Arow;
2662: 
2663:     for (i=0; i<ismax; i++) {
2664:       mat       = (Mat_SeqAIJ*)submats[i]->data;
2665:       imat_ilen = mat->ilen;
2666:       imat_j    = mat->j;
2667:       imat_i    = mat->i;
2668:       imat_a    = mat->a;
2669:       cmap_i    = cmap[i];
2670:       rmap_i    = rmap[i];
2671:       irow_i    = irow[i];
2672:       jmax      = nrow[i];
2673:       for (j=0; j<jmax; j++) {
2674:         row      = irow_i[j];
2675:         proc     = rtable[row];
2676:         if (proc == rank) {
2677:           old_row  = row;
2678:           row      = rmap_i[row];
2679:           ilen_row = imat_ilen[row];
2680: 
2681:           Arow=A->rows[old_row-C->rmap.rstart];
2682:           ncols=Arow->length;
2683:           cols=Arow->col;
2684:           vals=Arow->nz;
2685: 
2686:           mat_i    = imat_i[row];
2687:           mat_a    = imat_a + mat_i;
2688:           mat_j    = imat_j + mat_i;
2689:           for (k=0; k<ncols; k++) {
2690:             if ((tcol = cmap_i[cols[k]])) {
2691:               *mat_j++ = tcol - 1;
2692:               *mat_a++ = (PetscScalar)vals[k];
2693:               ilen_row++;
2694:             }
2695:           }
2696:           imat_ilen[row] = ilen_row;
2697:         }
2698:       }
2699:     }
2700:   }

2702:   /*   Now assemble the off proc rows*/
2703:   {
2704:     int    *sbuf1_i,*rbuf2_i,*rbuf3_i,*imat_ilen,ilen;
2705:     int    *imat_j,*imat_i;
2706:     PetscScalar *imat_a;
2707:     FLOAT *rbuf4_i;
2708: 
2709:     for (tmp2=0; tmp2<nrqs; tmp2++) {
2710:       MPI_Waitany(nrqs,r_waits4,&i,r_status4+tmp2);
2711:       idx     = pa[i];
2712:       sbuf1_i = sbuf1[idx];
2713:       jmax    = sbuf1_i[0];
2714:       ct1     = 2*jmax + 1;
2715:       ct2     = 0;
2716:       rbuf2_i = rbuf2[i];
2717:       rbuf3_i = rbuf3[i];
2718:       rbuf4_i = rbuf4[i];
2719:       for (j=1; j<=jmax; j++) {
2720:         is_no     = sbuf1_i[2*j-1];
2721:         rmap_i    = rmap[is_no];
2722:         cmap_i    = cmap[is_no];
2723:         mat       = (Mat_SeqAIJ*)submats[is_no]->data;
2724:         imat_ilen = mat->ilen;
2725:         imat_j    = mat->j;
2726:         imat_i    = mat->i;
2727:         imat_a    = mat->a;
2728:         max1      = sbuf1_i[2*j];
2729:         for (k=0; k<max1; k++,ct1++) {
2730:           row   = sbuf1_i[ct1];
2731:           row   = rmap_i[row];
2732:           ilen  = imat_ilen[row];
2733:           mat_i = imat_i[row];
2734:           mat_a = imat_a + mat_i;
2735:           mat_j = imat_j + mat_i;
2736:           max2 = rbuf2_i[ct1];
2737:           for (l=0; l<max2; l++,ct2++) {
2738:             if ((tcol = cmap_i[rbuf3_i[ct2]])) {
2739:               *mat_j++ = tcol - 1;
2740:               *mat_a++ = (PetscScalar)rbuf4_i[ct2];
2741:               ilen++;
2742:             }
2743:           }
2744:           imat_ilen[row] = ilen;
2745:         }
2746:       }
2747:     }
2748:   }
2749:   PetscFree(r_status4);
2750:   PetscFree(r_waits4);
2751:   if (nrqr) {MPI_Waitall(nrqr,s_waits4,s_status4);}
2752:   PetscFree(s_waits4);
2753:   PetscFree(s_status4);

2755:   /* Restore the indices */
2756:   for (i=0; i<ismax; i++) {
2757:     ISRestoreIndices(isrow[i],irow+i);
2758:     ISRestoreIndices(iscol[i],icol+i);
2759:   }

2761:   /* Destroy allocated memory */
2762:   PetscFree(irow);
2763:   PetscFree(w1);
2764:   PetscFree(pa);

2766:   PetscFree(sbuf1);
2767:   PetscFree(rbuf2);
2768:   for (i=0; i<nrqr; ++i) {
2769:     PetscFree(sbuf2[i]);
2770:   }
2771:   for (i=0; i<nrqs; ++i) {
2772:     PetscFree(rbuf3[i]);
2773:     PetscFree(rbuf4[i]);
2774:   }

2776:   PetscFree(sbuf2);
2777:   PetscFree(rbuf3);
2778:   PetscFree(rbuf4);
2779:   PetscFree(sbuf_aj[0]);
2780:   PetscFree(sbuf_aj);
2781:   PetscFree(sbuf_aa[0]);
2782:   PetscFree(sbuf_aa);
2783: 
2784:   PetscFree(cmap);
2785:   PetscFree(rmap);
2786:   PetscFree(lens);

2788:   for (i=0; i<ismax; i++) {
2789:     MatAssemblyBegin(submats[i],MAT_FINAL_ASSEMBLY);
2790:     MatAssemblyEnd(submats[i],MAT_FINAL_ASSEMBLY);
2791:   }
2792:   return(0);
2793: }

2795: /*
2796:   can be optimized by send only non-zeroes in iscol IS  -
2797:   so prebuild submatrix on sending side including A,B partitioning
2798:   */
2801:  #include src/vec/is/impls/general/general.h
2802: PetscErrorCode MatGetSubMatrix_MPIRowbs(Mat C,IS isrow,IS iscol,int csize,MatReuse scall,Mat *submat)
2803: {
2804:   Mat_MPIRowbs  *c = (Mat_MPIRowbs*)C->data;
2805:   BSspmat       *A = c->A;
2806:   BSsprow *Arow;
2807:   Mat_SeqAIJ    *matA,*matB; /* on prac , off proc part of submat */
2808:   Mat_MPIAIJ    *mat;  /* submat->data */
2810:   int    *irow,*icol,nrow,ncol,*rtable,size,rank,tag0,tag1,tag2,tag3;
2811:   int    *w1,*w2,*pa,nrqs,nrqr,msz,row_t;
2812:   int    i,j,k,l,len,jmax,proc,idx;
2813:   int    **sbuf1,**sbuf2,**rbuf1,**rbuf2,*req_size,**sbuf3,**rbuf3;
2814:   FLOAT  **rbuf4,**sbuf4; /* FLOAT is from Block Solve 95 library */

2816:   int    *cmap,*rmap,nlocal,*o_nz,*d_nz,cstart,cend;
2817:   int    *req_source;
2818:   int    ncols_t;
2819: 
2820: 
2821:   MPI_Request *s_waits1,*r_waits1,*s_waits2,*r_waits2,*r_waits3;
2822:   MPI_Request *r_waits4,*s_waits3,*s_waits4;
2823: 
2824:   MPI_Status  *r_status1,*r_status2,*s_status1,*s_status3,*s_status2;
2825:   MPI_Status  *r_status3,*r_status4,*s_status4;
2826:   MPI_Comm    comm;


2830:   comm   = C->comm;
2831:   tag0   = C->tag;
2832:   size   = c->size;
2833:   rank   = c->rank;

2835:   if (size==1) {
2836:     if (scall == MAT_REUSE_MATRIX) {
2837:       ierr=MatGetSubMatrices(C,1,&isrow,&iscol,MAT_REUSE_MATRIX,&submat);
2838:       return(0);
2839:     } else {
2840:       Mat *newsubmat;
2841: 
2842:       ierr=MatGetSubMatrices(C,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&newsubmat);
2843:       *submat=*newsubmat;
2844:       ierr=PetscFree(newsubmat);
2845:       return(0);
2846:     }
2847:   }
2848: 
2849:   /* Get some new tags to keep the communication clean */
2850:   PetscObjectGetNewTag((PetscObject)C,&tag1);
2851:   PetscObjectGetNewTag((PetscObject)C,&tag2);
2852:   PetscObjectGetNewTag((PetscObject)C,&tag3);

2854:   /* Check if the col indices are sorted */
2855:   {PetscTruth sorted;
2856:   ISSorted(isrow,&sorted);
2857:   if (!sorted) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"ISrow is not sorted");
2858:   ISSorted(iscol,&sorted);
2859:   if (!sorted) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IScol is not sorted");
2860:   }
2861: 
2862:   ISGetIndices(isrow,&irow);
2863:   ISGetIndices(iscol,&icol);
2864:   ISGetLocalSize(isrow,&nrow);
2865:   ISGetLocalSize(iscol,&ncol);
2866: 
2867:   if (!isrow) SETERRQ(PETSC_ERR_ARG_SIZ,"Empty ISrow");
2868:   if (!iscol) SETERRQ(PETSC_ERR_ARG_SIZ,"Empty IScol");
2869: 
2870: 
2871:   len    = (C->rmap.N+1)*sizeof(int);
2872:   PetscMalloc(len,&rtable);
2873:   /* Create hash table for the mapping :row -> proc*/
2874:   for (i=0,j=0; i<size; i++) {
2875:     jmax = C->rmap.range[i+1];
2876:     for (; j<jmax; j++) {
2877:       rtable[j] = i;
2878:     }
2879:   }

2881:   /* evaluate communication - mesg to who, length of mesg, and buffer space
2882:      required. Based on this, buffers are allocated, and data copied into them*/
2883:   PetscMalloc(size*2*sizeof(int),&w1); /* mesg size */
2884:   w2     = w1 + size;      /* if w2[i] marked, then a message to proc i*/
2885:   PetscMemzero(w1,size*2*sizeof(int)); /* initialize work vector*/
2886:   for (j=0; j<nrow; j++) {
2887:     row_t  = irow[j];
2888:     proc   = rtable[row_t];
2889:     w1[proc]++;
2890:   }
2891:   nrqs     = 0;              /* no of outgoing messages */
2892:   msz      = 0;              /* total mesg length (for all procs) */
2893:   w1[rank] = 0;              /* no mesg sent to self */
2894:   for (i=0; i<size; i++) {
2895:     if (w1[i])  { w2[i] = 1; nrqs++;} /* there exists a message to proc i */
2896:   }
2897: 
2898:   PetscMalloc((nrqs+1)*sizeof(int),&pa); /*(proc -array)*/
2899:   for (i=0,j=0; i<size; i++) {
2900:     if (w1[i]) {
2901:       pa[j++] = i;
2902:       w1[i]++;  /* header for return data */
2903:       msz+=w1[i];
2904:     }
2905:   }
2906: 
2907:   {int  *onodes1,*olengths1;
2908:   /* Determine the number of messages to expect, their lengths, from from-ids */
2909:   PetscGatherNumberOfMessages(comm,w2,w1,&nrqr);
2910:   PetscGatherMessageLengths(comm,nrqs,nrqr,w1,&onodes1,&olengths1);
2911:   /* Now post the Irecvs corresponding to these messages */
2912:   PetscPostIrecvInt(comm,tag0,nrqr,onodes1,olengths1,&rbuf1,&r_waits1);
2913:   PetscFree(onodes1);
2914:   PetscFree(olengths1);
2915:   }
2916: 
2917: { int **ptr,*iptr,*tmp;
2918:   /* Allocate Memory for outgoing messages */
2919:   len      = 2*size*sizeof(int*) + msz*sizeof(int);
2920:   PetscMalloc(len,&sbuf1);
2921:   ptr      = sbuf1 + size;   /* Pointers to the data in outgoing buffers */
2922:   PetscMemzero(sbuf1,2*size*sizeof(int*));
2923:   /* allocate memory for outgoing data + buf to receive the first reply */
2924:   tmp      = (int*)(ptr + size);

2926:   for (i=0,iptr=tmp; i<nrqs; i++) {
2927:     j         = pa[i];
2928:     sbuf1[j]  = iptr;
2929:     iptr     += w1[j];
2930:   }

2932:   /* Form the outgoing messages */
2933:   for (i=0; i<nrqs; i++) {
2934:     j           = pa[i];
2935:     sbuf1[j][0] = 0;   /*header */
2936:     ptr[j]      = sbuf1[j] + 1;
2937:   }
2938: 
2939:   /* Parse the isrow and copy data into outbuf */
2940:   for (j=0; j<nrow; j++) {
2941:     row_t  = irow[j];
2942:     proc = rtable[row_t];
2943:     if (proc != rank) { /* copy to the outgoing buf*/
2944:       sbuf1[proc][0]++;
2945:       *ptr[proc] = row_t;
2946:       ptr[proc]++;
2947:     }
2948:   }
2949: } /* block */

2951:   /*  Now  post the sends */
2952: 
2953:   /* structure of sbuf1[i]/rbuf1[i] : 1 (num of rows) + nrow-local rows (nuberes
2954:    * of requested rows)*/

2956:   PetscMalloc((nrqs+1)*sizeof(MPI_Request),&s_waits1);
2957:   for (i=0; i<nrqs; ++i) {
2958:     j    = pa[i];
2959:     MPI_Isend(sbuf1[j],w1[j],MPI_INT,j,tag0,comm,s_waits1+i);
2960:   }

2962:   /* Post Receives to capture the buffer size */
2963:   PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits2);
2964:   PetscMalloc((nrqs+1)*sizeof(int*),&rbuf2);
2965:   PetscMalloc(msz*sizeof(int)+1,&(rbuf2[0]));
2966:   for (i=1; i<nrqs; ++i) {
2967:     rbuf2[i] = rbuf2[i-1]+w1[pa[i-1]];
2968:   }
2969:   for (i=0; i<nrqs; ++i) {
2970:     j    = pa[i];
2971:     MPI_Irecv(rbuf2[i],w1[j],MPI_INT,j,tag1,comm,r_waits2+i);
2972:   }

2974:   /* Send to other procs the buf size they should allocate */
2975:   /* structure of sbuf2[i]/rbuf2[i]: 1 (total size to allocate) + nrow-locrow
2976:    * (row sizes) */

2978:   /* Receive messages*/
2979:   PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits2);
2980:   PetscMalloc((nrqr+1)*sizeof(MPI_Status),&r_status1);
2981:   len         = 2*nrqr*sizeof(int) + (nrqr+1)*sizeof(int*);
2982:   PetscMalloc(len,&sbuf2);
2983:   req_size    = (int*)(sbuf2 + nrqr);
2984:   req_source  = req_size + nrqr;
2985: 
2986:   {
2987:     BSsprow    **sAi = A->rows;
2988:     int        id,rstart = C->rmap.rstart;
2989:     int        *sbuf2_i,*rbuf1_i,end;

2991:     for (i=0; i<nrqr; ++i) {
2992:       MPI_Waitany(nrqr,r_waits1,&idx,r_status1+i);
2993:       req_size[idx]   = 0;
2994:       rbuf1_i         = rbuf1[idx];
2995:       MPI_Get_count(r_status1+i,MPI_INT,&end);
2996:       PetscMalloc((end+1)*sizeof(int),&sbuf2[idx]);
2997:       sbuf2_i         = sbuf2[idx];
2998:       for (j=1; j<end; j++) {
2999:         id               = rbuf1_i[j] - rstart;
3000:         ncols_t          = (sAi[id])->length;
3001:         sbuf2_i[j]       = ncols_t;
3002:         req_size[idx]   += ncols_t;
3003:       }
3004:       req_source[idx] = r_status1[i].MPI_SOURCE;
3005:       /* form the header */
3006:       sbuf2_i[0]   = req_size[idx];
3007:       MPI_Isend(sbuf2_i,end,MPI_INT,req_source[idx],tag1,comm,s_waits2+i);
3008:     }
3009:   }
3010:   PetscFree(r_status1);
3011:   PetscFree(r_waits1);

3013:   /*  recv buffer sizes */
3014:   /* Receive messages*/
3015: 
3016:   PetscMalloc((nrqs+1)*sizeof(int*),&rbuf3);
3017:   PetscMalloc((nrqs+1)*sizeof(FLOAT*),&rbuf4);
3018:   PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits3);
3019:   PetscMalloc((nrqs+1)*sizeof(MPI_Request),&r_waits4);
3020:   PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status2);

3022:   for (i=0; i<nrqs; ++i) {
3023:     MPI_Waitany(nrqs,r_waits2,&idx,r_status2+i);
3024:     PetscMalloc((rbuf2[idx][0]+1)*sizeof(int),&rbuf3[idx]);
3025:     PetscMalloc((rbuf2[idx][0]+1)*sizeof(FLOAT),&rbuf4[idx]);
3026:     MPI_Irecv(rbuf3[idx],rbuf2[idx][0],MPI_INT,r_status2[i].MPI_SOURCE,tag2,comm,r_waits3+idx);
3027:     MPI_Irecv(rbuf4[idx],rbuf2[idx][0],MPIU_SCALAR,r_status2[i].MPI_SOURCE,tag3,comm,r_waits4+idx);
3028:   }
3029:   PetscFree(r_status2);
3030:   PetscFree(r_waits2);
3031: 
3032:   /* Wait on sends1 and sends2 */
3033:   PetscMalloc((nrqs+1)*sizeof(MPI_Status),&s_status1);
3034:   PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status2);

3036:   if (nrqs) {MPI_Waitall(nrqs,s_waits1,s_status1);}
3037:   if (nrqr) {MPI_Waitall(nrqr,s_waits2,s_status2);}
3038:   PetscFree(s_status1);
3039:   PetscFree(s_status2);
3040:   PetscFree(s_waits1);
3041:   PetscFree(s_waits2);

3043:   /* Now allocate buffers for a->j, and send them off */
3044:   /* structure of sbuf3[i]/rbuf3[i],sbuf4[i]/rbuf4[i]: reqsize[i] (cols resp.
3045:    * vals of all req. rows; row sizes was in rbuf2; vals are of FLOAT type */
3046: 
3047:   PetscMalloc((nrqr+1)*sizeof(int*),&sbuf3);
3048:   for (i=0,j=0; i<nrqr; i++) j += req_size[i];
3049:   PetscMalloc((j+1)*sizeof(int),&sbuf3[0]);
3050:   for (i=1; i<nrqr; i++)  sbuf3[i] = sbuf3[i-1] + req_size[i-1];
3051: 
3052:   PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits3);
3053:   {
3054:     int *Acol,*rbuf1_i,*sbuf3_i,rqrow,noutcols,kmax,*cols,ncols;
3055:     int rstart = C->rmap.rstart;

3057:     for (i=0; i<nrqr; i++) {
3058:       rbuf1_i   = rbuf1[i];
3059:       sbuf3_i   = sbuf3[i];
3060:       noutcols  = 0;
3061:       kmax = rbuf1_i[0];  /* num. of req. rows */
3062:       for (k=0,rqrow=1; k<kmax; k++,rqrow++) {
3063:         Arow    = A->rows[rbuf1_i[rqrow] - rstart];
3064:         ncols  = Arow->length;
3065:         Acol   = Arow->col;
3066:         /* load the column indices for this row into cols*/
3067:         cols  = sbuf3_i + noutcols;
3068:         PetscMemcpy(cols,Acol,ncols*sizeof(int));
3069:         /*for (l=0; l<ncols;l++) cols[l]=Acol[l]; */ /* How is it with mappings?? */
3070:         noutcols += ncols;
3071:       }
3072:       MPI_Isend(sbuf3_i,req_size[i],MPI_INT,req_source[i],tag2,comm,s_waits3+i);
3073:     }
3074:   }
3075:   PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status3);
3076:   PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status3);

3078:   /* Allocate buffers for a->a, and send them off */
3079:   /* can be optimized by conect with previous block */
3080:   PetscMalloc((nrqr+1)*sizeof(FLOAT*),&sbuf4);
3081:   for (i=0,j=0; i<nrqr; i++) j += req_size[i];
3082:   PetscMalloc((j+1)*sizeof(FLOAT),&sbuf4[0]);
3083:   for (i=1; i<nrqr; i++)  sbuf4[i] = sbuf4[i-1] + req_size[i-1];
3084: 
3085:   PetscMalloc((nrqr+1)*sizeof(MPI_Request),&s_waits4);
3086:   {
3087:     FLOAT *Aval,*vals,*sbuf4_i;
3088:     int rstart = C->rmap.rstart,*rbuf1_i,rqrow,noutvals,kmax,ncols;
3089: 
3090: 
3091:     for (i=0; i<nrqr; i++) {
3092:       rbuf1_i   = rbuf1[i];
3093:       sbuf4_i   = sbuf4[i];
3094:       rqrow     = 1;
3095:       noutvals  = 0;
3096:       kmax      = rbuf1_i[0];  /* num of req. rows */
3097:       for (k=0; k<kmax; k++,rqrow++) {
3098:         Arow    = A->rows[rbuf1_i[rqrow] - rstart];
3099:         ncols  = Arow->length;
3100:         Aval = Arow->nz;
3101:         /* load the column values for this row into vals*/
3102:         vals  = sbuf4_i+noutvals;
3103:         PetscMemcpy(vals,Aval,ncols*sizeof(FLOAT));
3104:         noutvals += ncols;
3105:       }
3106:       MPI_Isend(sbuf4_i,req_size[i],MPIU_SCALAR,req_source[i],tag3,comm,s_waits4+i);
3107:     }
3108:   }
3109:   PetscMalloc((nrqs+1)*sizeof(MPI_Status),&r_status4);
3110:   PetscMalloc((nrqr+1)*sizeof(MPI_Status),&s_status4);
3111:   PetscFree(rbuf1);

3113:   /* Form the matrix */

3115:   /* create col map */
3116:   len     = C->cmap.N*sizeof(int)+1;
3117:   PetscMalloc(len,&cmap);
3118:   PetscMemzero(cmap,C->cmap.N*sizeof(int));
3119:   for (j=0; j<ncol; j++) {
3120:       cmap[icol[j]] = j+1;
3121:   }
3122: 
3123:   /* Create row map / maybe I will need global rowmap but here is local rowmap*/
3124:   len     = C->rmap.N*sizeof(int)+1;
3125:   PetscMalloc(len,&rmap);
3126:   PetscMemzero(rmap,C->rmap.N*sizeof(int));
3127:   for (j=0; j<nrow; j++) {
3128:     rmap[irow[j]] = j;
3129:   }

3131:   /*
3132:      Determine the number of non-zeros in the diagonal and off-diagonal 
3133:      portions of the matrix in order to do correct preallocation
3134:    */

3136:   /* first get start and end of "diagonal" columns */
3137:   if (csize == PETSC_DECIDE) {
3138:     nlocal = ncol/size + ((ncol % size) > rank);
3139:   } else {
3140:     nlocal = csize;
3141:   }
3142:   {
3143:     int ncols,*cols,olen,dlen,thecol;
3144:     int *rbuf2_i,*rbuf3_i,*sbuf1_i,row,kmax,cidx;
3145: 
3146:     MPI_Scan(&nlocal,&cend,1,MPI_INT,MPI_SUM,comm);
3147:     cstart = cend - nlocal;
3148:     if (rank == size - 1 && cend != ncol) {
3149:       SETERRQ(PETSC_ERR_ARG_SIZ,"Local column sizes do not add up to total number of columns");
3150:     }

3152:     PetscMalloc((2*nrow+1)*sizeof(int),&d_nz);
3153:     o_nz = d_nz + nrow;
3154: 
3155:     /* Update lens from local data */
3156:     for (j=0; j<nrow; j++) {
3157:       row  = irow[j];
3158:       proc = rtable[row];
3159:       if (proc == rank) {
3160:         Arow=A->rows[row-C->rmap.rstart];
3161:         ncols=Arow->length;
3162:         cols=Arow->col;
3163:         olen=dlen=0;
3164:         for (k=0; k<ncols; k++) {
3165:           if ((thecol=cmap[cols[k]])) {
3166:             if (cstart<thecol && thecol<=cend) dlen++; /* thecol is from 1 */
3167:             else olen++;
3168:           }
3169:         }
3170:         o_nz[j]=olen;
3171:         d_nz[j]=dlen;
3172:       } else d_nz[j]=o_nz[j]=0;
3173:     }
3174:     /* Update lens from offproc data and done waits */
3175:     /* this will be much simplier after sending only appropriate columns */
3176:     for (j=0; j<nrqs;j++) {
3177:       MPI_Waitany(nrqs,r_waits3,&i,r_status3+j);
3178:       proc   = pa[i];
3179:       sbuf1_i = sbuf1[proc];
3180:       cidx    = 0;
3181:       rbuf2_i = rbuf2[i];
3182:       rbuf3_i = rbuf3[i];
3183:       kmax    = sbuf1_i[0]; /*num of rq. rows*/
3184:       for (k=1; k<=kmax; k++) {
3185:         row  = rmap[sbuf1_i[k]]; /* the val in the new matrix to be */
3186:         for (l=0; l<rbuf2_i[k]; l++,cidx++) {
3187:           if ((thecol=cmap[rbuf3_i[cidx]])) {
3188: 
3189:             if (cstart<thecol && thecol<=cend) d_nz[row]++; /* thecol is from 1 */
3190:             else o_nz[row]++;
3191:           }
3192:         }
3193:       }
3194:     }
3195:   }
3196:   PetscFree(r_status3);
3197:   PetscFree(r_waits3);
3198:   if (nrqr) {MPI_Waitall(nrqr,s_waits3,s_status3);}
3199:   PetscFree(s_status3);
3200:   PetscFree(s_waits3);

3202:   if (scall ==  MAT_INITIAL_MATRIX) {
3203:     MatCreate(comm,submat);
3204:     MatSetSizes(*submat,nrow,nlocal,PETSC_DECIDE,ncol);
3205:     MatSetType(*submat,C->type_name);
3206:     MatMPIAIJSetPreallocation(*submat,0,d_nz,0,o_nz);
3207:     mat=(Mat_MPIAIJ *)((*submat)->data);
3208:     matA=(Mat_SeqAIJ *)(mat->A->data);
3209:     matB=(Mat_SeqAIJ *)(mat->B->data);
3210: 
3211:   } else {
3212:     PetscTruth same;
3213:     /* folowing code can be optionaly dropped for debuged versions of users
3214:      * program, but I don't know PETSc option which can switch off such safety
3215:      * tests - in a same way counting of o_nz,d_nz can be droped for  REUSE
3216:      * matrix */
3217: 
3218:     PetscTypeCompare((PetscObject)(*submat),MATMPIAIJ,&same);
3219:     if (!same) {
3220:       SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong type");
3221:     }
3222:     if (((*submat)->rmap.n != nrow) || ((*submat)->cmap.N != ncol)) {
3223:         SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size");
3224:     }
3225:     mat=(Mat_MPIAIJ *)((*submat)->data);
3226:     matA=(Mat_SeqAIJ *)(mat->A->data);
3227:     matB=(Mat_SeqAIJ *)(mat->B->data);
3228:     PetscMemcmp(matA->ilen,d_nz,nrow*sizeof(int),&same);
3229:     if (!same) {
3230:       SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros");
3231:     }
3232:     PetscMemcmp(matB->ilen,o_nz,nrow*sizeof(int),&same);
3233:     if (!same) {
3234:       SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros");
3235:     }
3236:   /* Initial matrix as if empty */
3237:     PetscMemzero(matA->ilen,nrow*sizeof(int));
3238:     PetscMemzero(matB->ilen,nrow*sizeof(int));
3239:     /* Perhaps MatZeroEnteries may be better - look what it is exactly doing - I must
3240:      * delete all possibly nonactual inforamtion */
3241:     /*submats[i]->factor = C->factor; !!! ??? if factor will be same then I must
3242:      * copy some factor information - where are thay */
3243:     (*submat)->was_assembled=PETSC_FALSE;
3244:     (*submat)->assembled=PETSC_FALSE;
3245: 
3246:   }
3247:   PetscFree(d_nz);

3249:   /* Assemble the matrix */
3250:   /* First assemble from local rows */
3251:   {
3252:     int    i_row,oldrow,row,ncols,*cols,*matA_j,*matB_j,ilenA,ilenB,tcol;
3253:     FLOAT  *vals;
3254:     PetscScalar *matA_a,*matB_a;
3255: 
3256:     for (j=0; j<nrow; j++) {
3257:       oldrow = irow[j];
3258:       proc   = rtable[oldrow];
3259:       if (proc == rank) {
3260:         row  = rmap[oldrow];
3261: 
3262:         Arow  = A->rows[oldrow-C->rmap.rstart];
3263:         ncols = Arow->length;
3264:         cols  = Arow->col;
3265:         vals  = Arow->nz;
3266: 
3267:         i_row   = matA->i[row];
3268:         matA_a = matA->a + i_row;
3269:         matA_j = matA->j + i_row;
3270:         i_row   = matB->i[row];
3271:         matB_a = matB->a + i_row;
3272:         matB_j = matB->j + i_row;
3273:         for (k=0,ilenA=0,ilenB=0; k<ncols; k++) {
3274:           if ((tcol = cmap[cols[k]])) {
3275:             if (tcol<=cstart) {
3276:               *matB_j++ = tcol-1;
3277:               *matB_a++ = vals[k];
3278:               ilenB++;
3279:             } else if (tcol<=cend) {
3280:               *matA_j++ = (tcol-1)-cstart;
3281:               *matA_a++ = (PetscScalar)(vals[k]);
3282:               ilenA++;
3283:             } else {
3284:               *matB_j++ = tcol-1;
3285:               *matB_a++ = vals[k];
3286:               ilenB++;
3287:             }
3288:           }
3289:         }
3290:         matA->ilen[row]=ilenA;
3291:         matB->ilen[row]=ilenB;
3292: 
3293:       }
3294:     }
3295:   }

3297:   /*   Now assemble the off proc rows*/
3298:   {
3299:     int  *sbuf1_i,*rbuf2_i,*rbuf3_i,cidx,kmax,row,i_row;
3300:     int  *matA_j,*matB_j,lmax,tcol,ilenA,ilenB;
3301:     PetscScalar *matA_a,*matB_a;
3302:     FLOAT *rbuf4_i;

3304:     for (j=0; j<nrqs; j++) {
3305:       MPI_Waitany(nrqs,r_waits4,&i,r_status4+j);
3306:       proc   = pa[i];
3307:       sbuf1_i = sbuf1[proc];
3308: 
3309:       cidx    = 0;
3310:       rbuf2_i = rbuf2[i];
3311:       rbuf3_i = rbuf3[i];
3312:       rbuf4_i = rbuf4[i];
3313:       kmax    = sbuf1_i[0];
3314:       for (k=1; k<=kmax; k++) {
3315:         row = rmap[sbuf1_i[k]];
3316: 
3317:         i_row  = matA->i[row];
3318:         matA_a = matA->a + i_row;
3319:         matA_j = matA->j + i_row;
3320:         i_row  = matB->i[row];
3321:         matB_a = matB->a + i_row;
3322:         matB_j = matB->j + i_row;
3323: 
3324:         lmax = rbuf2_i[k];
3325:         for (l=0,ilenA=0,ilenB=0; l<lmax; l++,cidx++) {
3326:           if ((tcol = cmap[rbuf3_i[cidx]])) {
3327:             if (tcol<=cstart) {
3328:               *matB_j++ = tcol-1;
3329:               *matB_a++ = (PetscScalar)(rbuf4_i[cidx]);;
3330:               ilenB++;
3331:             } else if (tcol<=cend) {
3332:               *matA_j++ = (tcol-1)-cstart;
3333:               *matA_a++ = (PetscScalar)(rbuf4_i[cidx]);
3334:               ilenA++;
3335:             } else {
3336:               *matB_j++ = tcol-1;
3337:               *matB_a++ = (PetscScalar)(rbuf4_i[cidx]);
3338:               ilenB++;
3339:             }
3340:           }
3341:         }
3342:         matA->ilen[row]=ilenA;
3343:         matB->ilen[row]=ilenB;
3344:       }
3345:     }
3346:   }

3348:   PetscFree(r_status4);
3349:   PetscFree(r_waits4);
3350:   if (nrqr) {MPI_Waitall(nrqr,s_waits4,s_status4);}
3351:   PetscFree(s_waits4);
3352:   PetscFree(s_status4);

3354:   /* Restore the indices */
3355:   ISRestoreIndices(isrow,&irow);
3356:   ISRestoreIndices(iscol,&icol);

3358:   /* Destroy allocated memory */
3359:   PetscFree(rtable);
3360:   PetscFree(w1);
3361:   PetscFree(pa);

3363:   PetscFree(sbuf1);
3364:   PetscFree(rbuf2[0]);
3365:   PetscFree(rbuf2);
3366:   for (i=0; i<nrqr; ++i) {
3367:     PetscFree(sbuf2[i]);
3368:   }
3369:   for (i=0; i<nrqs; ++i) {
3370:     PetscFree(rbuf3[i]);
3371:     PetscFree(rbuf4[i]);
3372:   }

3374:   PetscFree(sbuf2);
3375:   PetscFree(rbuf3);
3376:   PetscFree(rbuf4);
3377:   PetscFree(sbuf3[0]);
3378:   PetscFree(sbuf3);
3379:   PetscFree(sbuf4[0]);
3380:   PetscFree(sbuf4);
3381: 
3382:   PetscFree(cmap);
3383:   PetscFree(rmap);


3386:   MatAssemblyBegin(*submat,MAT_FINAL_ASSEMBLY);
3387:   MatAssemblyEnd(*submat,MAT_FINAL_ASSEMBLY);


3390:   return(0);
3391: }