Actual source code: mpisbaij.c

petsc-3.14.6 2021-03-30
Report Typos and Errors

  2: #include <../src/mat/impls/baij/mpi/mpibaij.h>
  3: #include <../src/mat/impls/sbaij/mpi/mpisbaij.h>
  4: #include <../src/mat/impls/sbaij/seq/sbaij.h>
  5: #include <petscblaslapack.h>

  7: #if defined(PETSC_HAVE_ELEMENTAL)
  8: PETSC_INTERN PetscErrorCode MatConvert_MPISBAIJ_Elemental(Mat,MatType,MatReuse,Mat*);
  9: #endif
 10: #if defined(PETSC_HAVE_SCALAPACK)
 11: PETSC_INTERN PetscErrorCode MatConvert_SBAIJ_ScaLAPACK(Mat,MatType,MatReuse,Mat*);
 12: #endif

 14: /* This could be moved to matimpl.h */
 15: static PetscErrorCode MatPreallocateWithMats_Private(Mat B, PetscInt nm, Mat X[], PetscBool symm[], PetscBool fill)
 16: {
 17:   Mat            preallocator;
 18:   PetscInt       r,rstart,rend;
 19:   PetscInt       bs,i,m,n,M,N;
 20:   PetscBool      cong = PETSC_TRUE;

 26:   for (i = 0; i < nm; i++) {
 28:     PetscLayoutCompare(B->rmap,X[i]->rmap,&cong);
 29:     if (!cong) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"Not for different layouts");
 30:   }
 32:   MatGetBlockSize(B,&bs);
 33:   MatGetSize(B,&M,&N);
 34:   MatGetLocalSize(B,&m,&n);
 35:   MatCreate(PetscObjectComm((PetscObject)B),&preallocator);
 36:   MatSetType(preallocator,MATPREALLOCATOR);
 37:   MatSetBlockSize(preallocator,bs);
 38:   MatSetSizes(preallocator,m,n,M,N);
 39:   MatSetUp(preallocator);
 40:   MatGetOwnershipRange(preallocator,&rstart,&rend);
 41:   for (r = rstart; r < rend; ++r) {
 42:     PetscInt          ncols;
 43:     const PetscInt    *row;
 44:     const PetscScalar *vals;

 46:     for (i = 0; i < nm; i++) {
 47:       MatGetRow(X[i],r,&ncols,&row,&vals);
 48:       MatSetValues(preallocator,1,&r,ncols,row,vals,INSERT_VALUES);
 49:       if (symm && symm[i]) {
 50:         MatSetValues(preallocator,ncols,row,1,&r,vals,INSERT_VALUES);
 51:       }
 52:       MatRestoreRow(X[i],r,&ncols,&row,&vals);
 53:     }
 54:   }
 55:   MatAssemblyBegin(preallocator,MAT_FINAL_ASSEMBLY);
 56:   MatAssemblyEnd(preallocator,MAT_FINAL_ASSEMBLY);
 57:   MatPreallocatorPreallocate(preallocator,fill,B);
 58:   MatDestroy(&preallocator);
 59:   return(0);
 60: }

 62: PETSC_INTERN PetscErrorCode MatConvert_MPISBAIJ_Basic(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
 63: {
 64:   Mat            B;
 66:   PetscInt       r;

 69:   if (reuse != MAT_REUSE_MATRIX) {
 70:     PetscBool symm = PETSC_TRUE,isdense;
 71:     PetscInt  bs;

 73:     MatCreate(PetscObjectComm((PetscObject)A),&B);
 74:     MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
 75:     MatSetType(B,newtype);
 76:     MatGetBlockSize(A,&bs);
 77:     MatSetBlockSize(B,bs);
 78:     PetscLayoutSetUp(B->rmap);
 79:     PetscLayoutSetUp(B->cmap);
 80:     PetscObjectTypeCompareAny((PetscObject)B,&isdense,MATSEQDENSE,MATMPIDENSE,MATSEQDENSECUDA,"");
 81:     if (!isdense) {
 82:       MatGetRowUpperTriangular(A);
 83:       MatPreallocateWithMats_Private(B,1,&A,&symm,PETSC_TRUE);
 84:       MatRestoreRowUpperTriangular(A);
 85:     } else {
 86:       MatSetUp(B);
 87:     }
 88:   } else {
 89:     B    = *newmat;
 90:     MatZeroEntries(B);
 91:   }

 93:   MatGetRowUpperTriangular(A);
 94:   for (r = A->rmap->rstart; r < A->rmap->rend; r++) {
 95:     PetscInt          ncols;
 96:     const PetscInt    *row;
 97:     const PetscScalar *vals;

 99:     MatGetRow(A,r,&ncols,&row,&vals);
100:     MatSetValues(B,1,&r,ncols,row,vals,INSERT_VALUES);
101: #if defined(PETSC_USE_COMPLEX)
102:     if (A->hermitian) {
103:       PetscInt i;
104:       for (i = 0; i < ncols; i++) {
105:         MatSetValue(B,row[i],r,PetscConj(vals[i]),INSERT_VALUES);
106:       }
107:     } else {
108:       MatSetValues(B,ncols,row,1,&r,vals,INSERT_VALUES);
109:     }
110: #else
111:     MatSetValues(B,ncols,row,1,&r,vals,INSERT_VALUES);
112: #endif
113:     MatRestoreRow(A,r,&ncols,&row,&vals);
114:   }
115:   MatRestoreRowUpperTriangular(A);
116:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
117:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

119:   if (reuse == MAT_INPLACE_MATRIX) {
120:     MatHeaderReplace(A,&B);
121:   } else {
122:     *newmat = B;
123:   }
124:   return(0);
125: }

127: PetscErrorCode  MatStoreValues_MPISBAIJ(Mat mat)
128: {
129:   Mat_MPISBAIJ   *aij = (Mat_MPISBAIJ*)mat->data;

133:   MatStoreValues(aij->A);
134:   MatStoreValues(aij->B);
135:   return(0);
136: }

138: PetscErrorCode  MatRetrieveValues_MPISBAIJ(Mat mat)
139: {
140:   Mat_MPISBAIJ   *aij = (Mat_MPISBAIJ*)mat->data;

144:   MatRetrieveValues(aij->A);
145:   MatRetrieveValues(aij->B);
146:   return(0);
147: }

149: #define  MatSetValues_SeqSBAIJ_A_Private(row,col,value,addv,orow,ocol)      \
150:   { \
151:     brow = row/bs;  \
152:     rp   = aj + ai[brow]; ap = aa + bs2*ai[brow]; \
153:     rmax = aimax[brow]; nrow = ailen[brow]; \
154:     bcol = col/bs; \
155:     ridx = row % bs; cidx = col % bs; \
156:     low  = 0; high = nrow; \
157:     while (high-low > 3) { \
158:       t = (low+high)/2; \
159:       if (rp[t] > bcol) high = t; \
160:       else              low  = t; \
161:     } \
162:     for (_i=low; _i<high; _i++) { \
163:       if (rp[_i] > bcol) break; \
164:       if (rp[_i] == bcol) { \
165:         bap = ap + bs2*_i + bs*cidx + ridx; \
166:         if (addv == ADD_VALUES) *bap += value;  \
167:         else                    *bap  = value;  \
168:         goto a_noinsert; \
169:       } \
170:     } \
171:     if (a->nonew == 1) goto a_noinsert; \
172:     if (a->nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero at global row/column (%D, %D) into matrix", orow, ocol); \
173:     MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,aimax,a->nonew,MatScalar); \
174:     N = nrow++ - 1;  \
175:     /* shift up all the later entries in this row */ \
176:     PetscArraymove(rp+_i+1,rp+_i,N-_i+1); \
177:     PetscArraymove(ap+bs2*(_i+1),ap+bs2*_i,bs2*(N-_i+1)); \
178:     PetscArrayzero(ap+bs2*_i,bs2);  \
179:     rp[_i]                      = bcol;  \
180:     ap[bs2*_i + bs*cidx + ridx] = value;  \
181:     A->nonzerostate++;\
182: a_noinsert:; \
183:     ailen[brow] = nrow; \
184:   }

186: #define  MatSetValues_SeqSBAIJ_B_Private(row,col,value,addv,orow,ocol) \
187:   { \
188:     brow = row/bs;  \
189:     rp   = bj + bi[brow]; ap = ba + bs2*bi[brow]; \
190:     rmax = bimax[brow]; nrow = bilen[brow]; \
191:     bcol = col/bs; \
192:     ridx = row % bs; cidx = col % bs; \
193:     low  = 0; high = nrow; \
194:     while (high-low > 3) { \
195:       t = (low+high)/2; \
196:       if (rp[t] > bcol) high = t; \
197:       else              low  = t; \
198:     } \
199:     for (_i=low; _i<high; _i++) { \
200:       if (rp[_i] > bcol) break; \
201:       if (rp[_i] == bcol) { \
202:         bap = ap + bs2*_i + bs*cidx + ridx; \
203:         if (addv == ADD_VALUES) *bap += value;  \
204:         else                    *bap  = value;  \
205:         goto b_noinsert; \
206:       } \
207:     } \
208:     if (b->nonew == 1) goto b_noinsert; \
209:     if (b->nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero at global row/column (%D, %D) into matrix", orow, ocol); \
210:     MatSeqXAIJReallocateAIJ(B,b->mbs,bs2,nrow,brow,bcol,rmax,ba,bi,bj,rp,ap,bimax,b->nonew,MatScalar); \
211:     N = nrow++ - 1;  \
212:     /* shift up all the later entries in this row */ \
213:     PetscArraymove(rp+_i+1,rp+_i,N-_i+1); \
214:     PetscArraymove(ap+bs2*(_i+1),ap+bs2*_i,bs2*(N-_i+1)); \
215:     PetscArrayzero(ap+bs2*_i,bs2); \
216:     rp[_i]                      = bcol;  \
217:     ap[bs2*_i + bs*cidx + ridx] = value;  \
218:     B->nonzerostate++;\
219: b_noinsert:; \
220:     bilen[brow] = nrow; \
221:   }

223: /* Only add/insert a(i,j) with i<=j (blocks).
224:    Any a(i,j) with i>j input by user is ingored or generates an error
225: */
226: PetscErrorCode MatSetValues_MPISBAIJ(Mat mat,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode addv)
227: {
228:   Mat_MPISBAIJ   *baij = (Mat_MPISBAIJ*)mat->data;
229:   MatScalar      value;
230:   PetscBool      roworiented = baij->roworiented;
232:   PetscInt       i,j,row,col;
233:   PetscInt       rstart_orig=mat->rmap->rstart;
234:   PetscInt       rend_orig  =mat->rmap->rend,cstart_orig=mat->cmap->rstart;
235:   PetscInt       cend_orig  =mat->cmap->rend,bs=mat->rmap->bs;

237:   /* Some Variables required in the macro */
238:   Mat          A     = baij->A;
239:   Mat_SeqSBAIJ *a    = (Mat_SeqSBAIJ*)(A)->data;
240:   PetscInt     *aimax=a->imax,*ai=a->i,*ailen=a->ilen,*aj=a->j;
241:   MatScalar    *aa   =a->a;

243:   Mat         B     = baij->B;
244:   Mat_SeqBAIJ *b    = (Mat_SeqBAIJ*)(B)->data;
245:   PetscInt    *bimax=b->imax,*bi=b->i,*bilen=b->ilen,*bj=b->j;
246:   MatScalar   *ba   =b->a;

248:   PetscInt  *rp,ii,nrow,_i,rmax,N,brow,bcol;
249:   PetscInt  low,high,t,ridx,cidx,bs2=a->bs2;
250:   MatScalar *ap,*bap;

252:   /* for stash */
253:   PetscInt  n_loc, *in_loc = NULL;
254:   MatScalar *v_loc = NULL;

257:   if (!baij->donotstash) {
258:     if (n > baij->n_loc) {
259:       PetscFree(baij->in_loc);
260:       PetscFree(baij->v_loc);
261:       PetscMalloc1(n,&baij->in_loc);
262:       PetscMalloc1(n,&baij->v_loc);

264:       baij->n_loc = n;
265:     }
266:     in_loc = baij->in_loc;
267:     v_loc  = baij->v_loc;
268:   }

270:   for (i=0; i<m; i++) {
271:     if (im[i] < 0) continue;
272:     if (PetscUnlikely(im[i] >= mat->rmap->N)) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",im[i],mat->rmap->N-1);
273:     if (im[i] >= rstart_orig && im[i] < rend_orig) { /* this processor entry */
274:       row = im[i] - rstart_orig;              /* local row index */
275:       for (j=0; j<n; j++) {
276:         if (im[i]/bs > in[j]/bs) {
277:           if (a->ignore_ltriangular) {
278:             continue;    /* ignore lower triangular blocks */
279:           } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Lower triangular value cannot be set for sbaij format. Ignoring these values, run with -mat_ignore_lower_triangular or call MatSetOption(mat,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE)");
280:         }
281:         if (in[j] >= cstart_orig && in[j] < cend_orig) {  /* diag entry (A) */
282:           col  = in[j] - cstart_orig;         /* local col index */
283:           brow = row/bs; bcol = col/bs;
284:           if (brow > bcol) continue;  /* ignore lower triangular blocks of A */
285:           if (roworiented) value = v[i*n+j];
286:           else             value = v[i+j*m];
287:           MatSetValues_SeqSBAIJ_A_Private(row,col,value,addv,im[i],in[j]);
288:           /* MatSetValues_SeqBAIJ(baij->A,1,&row,1,&col,&value,addv); */
289:         } else if (in[j] < 0) continue;
290:         else if (in[j] >= mat->cmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[j],mat->cmap->N-1);
291:         else {  /* off-diag entry (B) */
292:           if (mat->was_assembled) {
293:             if (!baij->colmap) {
294:               MatCreateColmap_MPIBAIJ_Private(mat);
295:             }
296: #if defined(PETSC_USE_CTABLE)
297:             PetscTableFind(baij->colmap,in[j]/bs + 1,&col);
298:             col  = col - 1;
299: #else
300:             col = baij->colmap[in[j]/bs] - 1;
301: #endif
302:             if (col < 0 && !((Mat_SeqSBAIJ*)(baij->A->data))->nonew) {
303:               MatDisAssemble_MPISBAIJ(mat);
304:               col  =  in[j];
305:               /* Reinitialize the variables required by MatSetValues_SeqBAIJ_B_Private() */
306:               B    = baij->B;
307:               b    = (Mat_SeqBAIJ*)(B)->data;
308:               bimax= b->imax;bi=b->i;bilen=b->ilen;bj=b->j;
309:               ba   = b->a;
310:             } else col += in[j]%bs;
311:           } else col = in[j];
312:           if (roworiented) value = v[i*n+j];
313:           else             value = v[i+j*m];
314:           MatSetValues_SeqSBAIJ_B_Private(row,col,value,addv,im[i],in[j]);
315:           /* MatSetValues_SeqBAIJ(baij->B,1,&row,1,&col,&value,addv); */
316:         }
317:       }
318:     } else {  /* off processor entry */
319:       if (mat->nooffprocentries) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Setting off process row %D even though MatSetOption(,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE) was set",im[i]);
320:       if (!baij->donotstash) {
321:         mat->assembled = PETSC_FALSE;
322:         n_loc          = 0;
323:         for (j=0; j<n; j++) {
324:           if (im[i]/bs > in[j]/bs) continue; /* ignore lower triangular blocks */
325:           in_loc[n_loc] = in[j];
326:           if (roworiented) {
327:             v_loc[n_loc] = v[i*n+j];
328:           } else {
329:             v_loc[n_loc] = v[j*m+i];
330:           }
331:           n_loc++;
332:         }
333:         MatStashValuesRow_Private(&mat->stash,im[i],n_loc,in_loc,v_loc,PETSC_FALSE);
334:       }
335:     }
336:   }
337:   return(0);
338: }

340: PETSC_STATIC_INLINE PetscErrorCode MatSetValuesBlocked_SeqSBAIJ_Inlined(Mat A,PetscInt row,PetscInt col,const PetscScalar v[],InsertMode is,PetscInt orow,PetscInt ocol)
341: {
342:   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
343:   PetscErrorCode    ierr;
344:   PetscInt          *rp,low,high,t,ii,jj,nrow,i,rmax,N;
345:   PetscInt          *imax      =a->imax,*ai=a->i,*ailen=a->ilen;
346:   PetscInt          *aj        =a->j,nonew=a->nonew,bs2=a->bs2,bs=A->rmap->bs;
347:   PetscBool         roworiented=a->roworiented;
348:   const PetscScalar *value     = v;
349:   MatScalar         *ap,*aa = a->a,*bap;

352:   if (col < row) {
353:     if (a->ignore_ltriangular) return(0); /* ignore lower triangular block */
354:     else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Lower triangular value cannot be set for sbaij format. Ignoring these values, run with -mat_ignore_lower_triangular or call MatSetOption(mat,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE)");
355:   }
356:   rp   = aj + ai[row];
357:   ap   = aa + bs2*ai[row];
358:   rmax = imax[row];
359:   nrow = ailen[row];
360:   value = v;
361:   low   = 0;
362:   high  = nrow;

364:   while (high-low > 7) {
365:     t = (low+high)/2;
366:     if (rp[t] > col) high = t;
367:     else             low  = t;
368:   }
369:   for (i=low; i<high; i++) {
370:     if (rp[i] > col) break;
371:     if (rp[i] == col) {
372:       bap = ap +  bs2*i;
373:       if (roworiented) {
374:         if (is == ADD_VALUES) {
375:           for (ii=0; ii<bs; ii++) {
376:             for (jj=ii; jj<bs2; jj+=bs) {
377:               bap[jj] += *value++;
378:             }
379:           }
380:         } else {
381:           for (ii=0; ii<bs; ii++) {
382:             for (jj=ii; jj<bs2; jj+=bs) {
383:               bap[jj] = *value++;
384:             }
385:           }
386:         }
387:       } else {
388:         if (is == ADD_VALUES) {
389:           for (ii=0; ii<bs; ii++) {
390:             for (jj=0; jj<bs; jj++) {
391:               *bap++ += *value++;
392:             }
393:           }
394:         } else {
395:           for (ii=0; ii<bs; ii++) {
396:             for (jj=0; jj<bs; jj++) {
397:               *bap++  = *value++;
398:             }
399:           }
400:         }
401:       }
402:       goto noinsert2;
403:     }
404:   }
405:   if (nonew == 1) goto noinsert2;
406:   if (nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new block index nonzero block (%D, %D) in the matrix", orow, ocol);
407:   MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar);
408:   N = nrow++ - 1; high++;
409:   /* shift up all the later entries in this row */
410:   PetscArraymove(rp+i+1,rp+i,N-i+1);
411:   PetscArraymove(ap+bs2*(i+1),ap+bs2*i,bs2*(N-i+1));
412:   rp[i] = col;
413:   bap   = ap +  bs2*i;
414:   if (roworiented) {
415:     for (ii=0; ii<bs; ii++) {
416:       for (jj=ii; jj<bs2; jj+=bs) {
417:         bap[jj] = *value++;
418:       }
419:     }
420:   } else {
421:     for (ii=0; ii<bs; ii++) {
422:       for (jj=0; jj<bs; jj++) {
423:         *bap++ = *value++;
424:       }
425:     }
426:   }
427:   noinsert2:;
428:   ailen[row] = nrow;
429:   return(0);
430: }

432: /*
433:    This routine is exactly duplicated in mpibaij.c
434: */
435: PETSC_STATIC_INLINE PetscErrorCode MatSetValuesBlocked_SeqBAIJ_Inlined(Mat A,PetscInt row,PetscInt col,const PetscScalar v[],InsertMode is,PetscInt orow,PetscInt ocol)
436: {
437:   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
438:   PetscInt          *rp,low,high,t,ii,jj,nrow,i,rmax,N;
439:   PetscInt          *imax=a->imax,*ai=a->i,*ailen=a->ilen;
440:   PetscErrorCode    ierr;
441:   PetscInt          *aj        =a->j,nonew=a->nonew,bs2=a->bs2,bs=A->rmap->bs;
442:   PetscBool         roworiented=a->roworiented;
443:   const PetscScalar *value     = v;
444:   MatScalar         *ap,*aa = a->a,*bap;

447:   rp   = aj + ai[row];
448:   ap   = aa + bs2*ai[row];
449:   rmax = imax[row];
450:   nrow = ailen[row];
451:   low  = 0;
452:   high = nrow;
453:   value = v;
454:   while (high-low > 7) {
455:     t = (low+high)/2;
456:     if (rp[t] > col) high = t;
457:     else             low  = t;
458:   }
459:   for (i=low; i<high; i++) {
460:     if (rp[i] > col) break;
461:     if (rp[i] == col) {
462:       bap = ap +  bs2*i;
463:       if (roworiented) {
464:         if (is == ADD_VALUES) {
465:           for (ii=0; ii<bs; ii++) {
466:             for (jj=ii; jj<bs2; jj+=bs) {
467:               bap[jj] += *value++;
468:             }
469:           }
470:         } else {
471:           for (ii=0; ii<bs; ii++) {
472:             for (jj=ii; jj<bs2; jj+=bs) {
473:               bap[jj] = *value++;
474:             }
475:           }
476:         }
477:       } else {
478:         if (is == ADD_VALUES) {
479:           for (ii=0; ii<bs; ii++,value+=bs) {
480:             for (jj=0; jj<bs; jj++) {
481:               bap[jj] += value[jj];
482:             }
483:             bap += bs;
484:           }
485:         } else {
486:           for (ii=0; ii<bs; ii++,value+=bs) {
487:             for (jj=0; jj<bs; jj++) {
488:               bap[jj]  = value[jj];
489:             }
490:             bap += bs;
491:           }
492:         }
493:       }
494:       goto noinsert2;
495:     }
496:   }
497:   if (nonew == 1) goto noinsert2;
498:   if (nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new global block indexed nonzero block (%D, %D) in the matrix", orow, ocol);
499:   MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar);
500:   N = nrow++ - 1; high++;
501:   /* shift up all the later entries in this row */
502:   PetscArraymove(rp+i+1,rp+i,N-i+1);
503:   PetscArraymove(ap+bs2*(i+1),ap+bs2*i,bs2*(N-i+1));
504:   rp[i] = col;
505:   bap   = ap +  bs2*i;
506:   if (roworiented) {
507:     for (ii=0; ii<bs; ii++) {
508:       for (jj=ii; jj<bs2; jj+=bs) {
509:         bap[jj] = *value++;
510:       }
511:     }
512:   } else {
513:     for (ii=0; ii<bs; ii++) {
514:       for (jj=0; jj<bs; jj++) {
515:         *bap++ = *value++;
516:       }
517:     }
518:   }
519:   noinsert2:;
520:   ailen[row] = nrow;
521:   return(0);
522: }

524: /*
525:     This routine could be optimized by removing the need for the block copy below and passing stride information
526:   to the above inline routines; similarly in MatSetValuesBlocked_MPIBAIJ()
527: */
528: PetscErrorCode MatSetValuesBlocked_MPISBAIJ(Mat mat,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const MatScalar v[],InsertMode addv)
529: {
530:   Mat_MPISBAIJ    *baij = (Mat_MPISBAIJ*)mat->data;
531:   const MatScalar *value;
532:   MatScalar       *barray     =baij->barray;
533:   PetscBool       roworiented = baij->roworiented,ignore_ltriangular = ((Mat_SeqSBAIJ*)baij->A->data)->ignore_ltriangular;
534:   PetscErrorCode  ierr;
535:   PetscInt        i,j,ii,jj,row,col,rstart=baij->rstartbs;
536:   PetscInt        rend=baij->rendbs,cstart=baij->cstartbs,stepval;
537:   PetscInt        cend=baij->cendbs,bs=mat->rmap->bs,bs2=baij->bs2;

540:   if (!barray) {
541:     PetscMalloc1(bs2,&barray);
542:     baij->barray = barray;
543:   }

545:   if (roworiented) {
546:     stepval = (n-1)*bs;
547:   } else {
548:     stepval = (m-1)*bs;
549:   }
550:   for (i=0; i<m; i++) {
551:     if (im[i] < 0) continue;
552:     if (PetscUnlikelyDebug(im[i] >= baij->Mbs)) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Block indexed row too large %D max %D",im[i],baij->Mbs-1);
553:     if (im[i] >= rstart && im[i] < rend) {
554:       row = im[i] - rstart;
555:       for (j=0; j<n; j++) {
556:         if (im[i] > in[j]) {
557:           if (ignore_ltriangular) continue; /* ignore lower triangular blocks */
558:           else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Lower triangular value cannot be set for sbaij format. Ignoring these values, run with -mat_ignore_lower_triangular or call MatSetOption(mat,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE)");
559:         }
560:         /* If NumCol = 1 then a copy is not required */
561:         if ((roworiented) && (n == 1)) {
562:           barray = (MatScalar*) v + i*bs2;
563:         } else if ((!roworiented) && (m == 1)) {
564:           barray = (MatScalar*) v + j*bs2;
565:         } else { /* Here a copy is required */
566:           if (roworiented) {
567:             value = v + i*(stepval+bs)*bs + j*bs;
568:           } else {
569:             value = v + j*(stepval+bs)*bs + i*bs;
570:           }
571:           for (ii=0; ii<bs; ii++,value+=stepval) {
572:             for (jj=0; jj<bs; jj++) {
573:               *barray++ = *value++;
574:             }
575:           }
576:           barray -=bs2;
577:         }

579:         if (in[j] >= cstart && in[j] < cend) {
580:           col  = in[j] - cstart;
581:           MatSetValuesBlocked_SeqSBAIJ_Inlined(baij->A,row,col,barray,addv,im[i],in[j]);
582:         } else if (in[j] < 0) continue;
583:         else if (PetscUnlikelyDebug(in[j] >= baij->Nbs)) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Block indexed column too large %D max %D",in[j],baij->Nbs-1);
584:         else {
585:           if (mat->was_assembled) {
586:             if (!baij->colmap) {
587:               MatCreateColmap_MPIBAIJ_Private(mat);
588:             }

590: #if defined(PETSC_USE_DEBUG)
591: #if defined(PETSC_USE_CTABLE)
592:             { PetscInt data;
593:               PetscTableFind(baij->colmap,in[j]+1,&data);
594:               if ((data - 1) % bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect colmap");
595:             }
596: #else
597:             if ((baij->colmap[in[j]] - 1) % bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect colmap");
598: #endif
599: #endif
600: #if defined(PETSC_USE_CTABLE)
601:             PetscTableFind(baij->colmap,in[j]+1,&col);
602:             col  = (col - 1)/bs;
603: #else
604:             col = (baij->colmap[in[j]] - 1)/bs;
605: #endif
606:             if (col < 0 && !((Mat_SeqBAIJ*)(baij->A->data))->nonew) {
607:               MatDisAssemble_MPISBAIJ(mat);
608:               col  = in[j];
609:             }
610:           } else col = in[j];
611:           MatSetValuesBlocked_SeqBAIJ_Inlined(baij->B,row,col,barray,addv,im[i],in[j]);
612:         }
613:       }
614:     } else {
615:       if (mat->nooffprocentries) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Setting off process block indexed row %D even though MatSetOption(,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE) was set",im[i]);
616:       if (!baij->donotstash) {
617:         if (roworiented) {
618:           MatStashValuesRowBlocked_Private(&mat->bstash,im[i],n,in,v,m,n,i);
619:         } else {
620:           MatStashValuesColBlocked_Private(&mat->bstash,im[i],n,in,v,m,n,i);
621:         }
622:       }
623:     }
624:   }
625:   return(0);
626: }

628: PetscErrorCode MatGetValues_MPISBAIJ(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],PetscScalar v[])
629: {
630:   Mat_MPISBAIJ   *baij = (Mat_MPISBAIJ*)mat->data;
632:   PetscInt       bs       = mat->rmap->bs,i,j,bsrstart = mat->rmap->rstart,bsrend = mat->rmap->rend;
633:   PetscInt       bscstart = mat->cmap->rstart,bscend = mat->cmap->rend,row,col,data;

636:   for (i=0; i<m; i++) {
637:     if (idxm[i] < 0) continue; /* SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",idxm[i]); */
638:     if (idxm[i] >= mat->rmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm[i],mat->rmap->N-1);
639:     if (idxm[i] >= bsrstart && idxm[i] < bsrend) {
640:       row = idxm[i] - bsrstart;
641:       for (j=0; j<n; j++) {
642:         if (idxn[j] < 0) continue; /* SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column %D",idxn[j]); */
643:         if (idxn[j] >= mat->cmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",idxn[j],mat->cmap->N-1);
644:         if (idxn[j] >= bscstart && idxn[j] < bscend) {
645:           col  = idxn[j] - bscstart;
646:           MatGetValues_SeqSBAIJ(baij->A,1,&row,1,&col,v+i*n+j);
647:         } else {
648:           if (!baij->colmap) {
649:             MatCreateColmap_MPIBAIJ_Private(mat);
650:           }
651: #if defined(PETSC_USE_CTABLE)
652:           PetscTableFind(baij->colmap,idxn[j]/bs+1,&data);
653:           data--;
654: #else
655:           data = baij->colmap[idxn[j]/bs]-1;
656: #endif
657:           if ((data < 0) || (baij->garray[data/bs] != idxn[j]/bs)) *(v+i*n+j) = 0.0;
658:           else {
659:             col  = data + idxn[j]%bs;
660:             MatGetValues_SeqBAIJ(baij->B,1,&row,1,&col,v+i*n+j);
661:           }
662:         }
663:       }
664:     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only local values currently supported");
665:   }
666:   return(0);
667: }

669: PetscErrorCode MatNorm_MPISBAIJ(Mat mat,NormType type,PetscReal *norm)
670: {
671:   Mat_MPISBAIJ   *baij = (Mat_MPISBAIJ*)mat->data;
673:   PetscReal      sum[2],*lnorm2;

676:   if (baij->size == 1) {
677:      MatNorm(baij->A,type,norm);
678:   } else {
679:     if (type == NORM_FROBENIUS) {
680:       PetscMalloc1(2,&lnorm2);
681:        MatNorm(baij->A,type,lnorm2);
682:       *lnorm2 = (*lnorm2)*(*lnorm2); lnorm2++;            /* squar power of norm(A) */
683:        MatNorm(baij->B,type,lnorm2);
684:       *lnorm2 = (*lnorm2)*(*lnorm2); lnorm2--;             /* squar power of norm(B) */
685:       MPIU_Allreduce(lnorm2,sum,2,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)mat));
686:       *norm   = PetscSqrtReal(sum[0] + 2*sum[1]);
687:       PetscFree(lnorm2);
688:     } else if (type == NORM_INFINITY || type == NORM_1) { /* max row/column sum */
689:       Mat_SeqSBAIJ *amat=(Mat_SeqSBAIJ*)baij->A->data;
690:       Mat_SeqBAIJ  *bmat=(Mat_SeqBAIJ*)baij->B->data;
691:       PetscReal    *rsum,*rsum2,vabs;
692:       PetscInt     *jj,*garray=baij->garray,rstart=baij->rstartbs,nz;
693:       PetscInt     brow,bcol,col,bs=baij->A->rmap->bs,row,grow,gcol,mbs=amat->mbs;
694:       MatScalar    *v;

696:       PetscMalloc2(mat->cmap->N,&rsum,mat->cmap->N,&rsum2);
697:       PetscArrayzero(rsum,mat->cmap->N);
698:       /* Amat */
699:       v = amat->a; jj = amat->j;
700:       for (brow=0; brow<mbs; brow++) {
701:         grow = bs*(rstart + brow);
702:         nz   = amat->i[brow+1] - amat->i[brow];
703:         for (bcol=0; bcol<nz; bcol++) {
704:           gcol = bs*(rstart + *jj); jj++;
705:           for (col=0; col<bs; col++) {
706:             for (row=0; row<bs; row++) {
707:               vabs            = PetscAbsScalar(*v); v++;
708:               rsum[gcol+col] += vabs;
709:               /* non-diagonal block */
710:               if (bcol > 0 && vabs > 0.0) rsum[grow+row] += vabs;
711:             }
712:           }
713:         }
714:         PetscLogFlops(nz*bs*bs);
715:       }
716:       /* Bmat */
717:       v = bmat->a; jj = bmat->j;
718:       for (brow=0; brow<mbs; brow++) {
719:         grow = bs*(rstart + brow);
720:         nz = bmat->i[brow+1] - bmat->i[brow];
721:         for (bcol=0; bcol<nz; bcol++) {
722:           gcol = bs*garray[*jj]; jj++;
723:           for (col=0; col<bs; col++) {
724:             for (row=0; row<bs; row++) {
725:               vabs            = PetscAbsScalar(*v); v++;
726:               rsum[gcol+col] += vabs;
727:               rsum[grow+row] += vabs;
728:             }
729:           }
730:         }
731:         PetscLogFlops(nz*bs*bs);
732:       }
733:       MPIU_Allreduce(rsum,rsum2,mat->cmap->N,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)mat));
734:       *norm = 0.0;
735:       for (col=0; col<mat->cmap->N; col++) {
736:         if (rsum2[col] > *norm) *norm = rsum2[col];
737:       }
738:       PetscFree2(rsum,rsum2);
739:     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for this norm yet");
740:   }
741:   return(0);
742: }

744: PetscErrorCode MatAssemblyBegin_MPISBAIJ(Mat mat,MatAssemblyType mode)
745: {
746:   Mat_MPISBAIJ   *baij = (Mat_MPISBAIJ*)mat->data;
748:   PetscInt       nstash,reallocs;

751:   if (baij->donotstash || mat->nooffprocentries) return(0);

753:   MatStashScatterBegin_Private(mat,&mat->stash,mat->rmap->range);
754:   MatStashScatterBegin_Private(mat,&mat->bstash,baij->rangebs);
755:   MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);
756:   PetscInfo2(mat,"Stash has %D entries,uses %D mallocs.\n",nstash,reallocs);
757:   MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);
758:   PetscInfo2(mat,"Block-Stash has %D entries, uses %D mallocs.\n",nstash,reallocs);
759:   return(0);
760: }

762: PetscErrorCode MatAssemblyEnd_MPISBAIJ(Mat mat,MatAssemblyType mode)
763: {
764:   Mat_MPISBAIJ   *baij=(Mat_MPISBAIJ*)mat->data;
765:   Mat_SeqSBAIJ   *a   =(Mat_SeqSBAIJ*)baij->A->data;
767:   PetscInt       i,j,rstart,ncols,flg,bs2=baij->bs2;
768:   PetscInt       *row,*col;
769:   PetscBool      other_disassembled;
770:   PetscMPIInt    n;
771:   PetscBool      r1,r2,r3;
772:   MatScalar      *val;

774:   /* do not use 'b=(Mat_SeqBAIJ*)baij->B->data' as B can be reset in disassembly */
776:   if (!baij->donotstash &&  !mat->nooffprocentries) {
777:     while (1) {
778:       MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);
779:       if (!flg) break;

781:       for (i=0; i<n;) {
782:         /* Now identify the consecutive vals belonging to the same row */
783:         for (j=i,rstart=row[j]; j<n; j++) {
784:           if (row[j] != rstart) break;
785:         }
786:         if (j < n) ncols = j-i;
787:         else       ncols = n-i;
788:         /* Now assemble all these values with a single function call */
789:         MatSetValues_MPISBAIJ(mat,1,row+i,ncols,col+i,val+i,mat->insertmode);
790:         i    = j;
791:       }
792:     }
793:     MatStashScatterEnd_Private(&mat->stash);
794:     /* Now process the block-stash. Since the values are stashed column-oriented,
795:        set the roworiented flag to column oriented, and after MatSetValues()
796:        restore the original flags */
797:     r1 = baij->roworiented;
798:     r2 = a->roworiented;
799:     r3 = ((Mat_SeqBAIJ*)baij->B->data)->roworiented;

801:     baij->roworiented = PETSC_FALSE;
802:     a->roworiented    = PETSC_FALSE;

804:     ((Mat_SeqBAIJ*)baij->B->data)->roworiented = PETSC_FALSE; /* b->roworinted */
805:     while (1) {
806:       MatStashScatterGetMesg_Private(&mat->bstash,&n,&row,&col,&val,&flg);
807:       if (!flg) break;

809:       for (i=0; i<n;) {
810:         /* Now identify the consecutive vals belonging to the same row */
811:         for (j=i,rstart=row[j]; j<n; j++) {
812:           if (row[j] != rstart) break;
813:         }
814:         if (j < n) ncols = j-i;
815:         else       ncols = n-i;
816:         MatSetValuesBlocked_MPISBAIJ(mat,1,row+i,ncols,col+i,val+i*bs2,mat->insertmode);
817:         i    = j;
818:       }
819:     }
820:     MatStashScatterEnd_Private(&mat->bstash);

822:     baij->roworiented = r1;
823:     a->roworiented    = r2;

825:     ((Mat_SeqBAIJ*)baij->B->data)->roworiented = r3; /* b->roworinted */
826:   }

828:   MatAssemblyBegin(baij->A,mode);
829:   MatAssemblyEnd(baij->A,mode);

831:   /* determine if any processor has disassembled, if so we must
832:      also disassemble ourselfs, in order that we may reassemble. */
833:   /*
834:      if nonzero structure of submatrix B cannot change then we know that
835:      no processor disassembled thus we can skip this stuff
836:   */
837:   if (!((Mat_SeqBAIJ*)baij->B->data)->nonew) {
838:     MPIU_Allreduce(&mat->was_assembled,&other_disassembled,1,MPIU_BOOL,MPI_PROD,PetscObjectComm((PetscObject)mat));
839:     if (mat->was_assembled && !other_disassembled) {
840:       MatDisAssemble_MPISBAIJ(mat);
841:     }
842:   }

844:   if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
845:     MatSetUpMultiply_MPISBAIJ(mat); /* setup Mvctx and sMvctx */
846:   }
847:   MatAssemblyBegin(baij->B,mode);
848:   MatAssemblyEnd(baij->B,mode);

850:   PetscFree2(baij->rowvalues,baij->rowindices);

852:   baij->rowvalues = NULL;

854:   /* if no new nonzero locations are allowed in matrix then only set the matrix state the first time through */
855:   if ((!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) || !((Mat_SeqBAIJ*)(baij->A->data))->nonew) {
856:     PetscObjectState state = baij->A->nonzerostate + baij->B->nonzerostate;
857:     MPIU_Allreduce(&state,&mat->nonzerostate,1,MPIU_INT64,MPI_SUM,PetscObjectComm((PetscObject)mat));
858:   }
859:   return(0);
860: }

862: extern PetscErrorCode MatSetValues_MPIBAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
863: #include <petscdraw.h>
864: static PetscErrorCode MatView_MPISBAIJ_ASCIIorDraworSocket(Mat mat,PetscViewer viewer)
865: {
866:   Mat_MPISBAIJ      *baij = (Mat_MPISBAIJ*)mat->data;
867:   PetscErrorCode    ierr;
868:   PetscInt          bs   = mat->rmap->bs;
869:   PetscMPIInt       rank = baij->rank;
870:   PetscBool         iascii,isdraw;
871:   PetscViewer       sviewer;
872:   PetscViewerFormat format;

875:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
876:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
877:   if (iascii) {
878:     PetscViewerGetFormat(viewer,&format);
879:     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
880:       MatInfo info;
881:       MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&rank);
882:       MatGetInfo(mat,MAT_LOCAL,&info);
883:       PetscViewerASCIIPushSynchronized(viewer);
884:       PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Local rows %D nz %D nz alloced %D bs %D mem %g\n",rank,mat->rmap->n,(PetscInt)info.nz_used,(PetscInt)info.nz_allocated,mat->rmap->bs,(double)info.memory);
885:       MatGetInfo(baij->A,MAT_LOCAL,&info);
886:       PetscViewerASCIISynchronizedPrintf(viewer,"[%d] on-diagonal part: nz %D \n",rank,(PetscInt)info.nz_used);
887:       MatGetInfo(baij->B,MAT_LOCAL,&info);
888:       PetscViewerASCIISynchronizedPrintf(viewer,"[%d] off-diagonal part: nz %D \n",rank,(PetscInt)info.nz_used);
889:       PetscViewerFlush(viewer);
890:       PetscViewerASCIIPopSynchronized(viewer);
891:       PetscViewerASCIIPrintf(viewer,"Information on VecScatter used in matrix-vector product: \n");
892:       VecScatterView(baij->Mvctx,viewer);
893:       return(0);
894:     } else if (format == PETSC_VIEWER_ASCII_INFO) {
895:       PetscViewerASCIIPrintf(viewer,"  block size is %D\n",bs);
896:       return(0);
897:     } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
898:       return(0);
899:     }
900:   }

902:   if (isdraw) {
903:     PetscDraw draw;
904:     PetscBool isnull;
905:     PetscViewerDrawGetDraw(viewer,0,&draw);
906:     PetscDrawIsNull(draw,&isnull);
907:     if (isnull) return(0);
908:   }

910:   {
911:     /* assemble the entire matrix onto first processor. */
912:     Mat          A;
913:     Mat_SeqSBAIJ *Aloc;
914:     Mat_SeqBAIJ  *Bloc;
915:     PetscInt     M = mat->rmap->N,N = mat->cmap->N,*ai,*aj,col,i,j,k,*rvals,mbs = baij->mbs;
916:     MatScalar    *a;
917:     const char   *matname;

919:     /* Should this be the same type as mat? */
920:     MatCreate(PetscObjectComm((PetscObject)mat),&A);
921:     if (!rank) {
922:       MatSetSizes(A,M,N,M,N);
923:     } else {
924:       MatSetSizes(A,0,0,M,N);
925:     }
926:     MatSetType(A,MATMPISBAIJ);
927:     MatMPISBAIJSetPreallocation(A,mat->rmap->bs,0,NULL,0,NULL);
928:     MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);
929:     PetscLogObjectParent((PetscObject)mat,(PetscObject)A);

931:     /* copy over the A part */
932:     Aloc = (Mat_SeqSBAIJ*)baij->A->data;
933:     ai   = Aloc->i; aj = Aloc->j; a = Aloc->a;
934:     PetscMalloc1(bs,&rvals);

936:     for (i=0; i<mbs; i++) {
937:       rvals[0] = bs*(baij->rstartbs + i);
938:       for (j=1; j<bs; j++) rvals[j] = rvals[j-1] + 1;
939:       for (j=ai[i]; j<ai[i+1]; j++) {
940:         col = (baij->cstartbs+aj[j])*bs;
941:         for (k=0; k<bs; k++) {
942:           MatSetValues_MPISBAIJ(A,bs,rvals,1,&col,a,INSERT_VALUES);
943:           col++;
944:           a += bs;
945:         }
946:       }
947:     }
948:     /* copy over the B part */
949:     Bloc = (Mat_SeqBAIJ*)baij->B->data;
950:     ai   = Bloc->i; aj = Bloc->j; a = Bloc->a;
951:     for (i=0; i<mbs; i++) {

953:       rvals[0] = bs*(baij->rstartbs + i);
954:       for (j=1; j<bs; j++) rvals[j] = rvals[j-1] + 1;
955:       for (j=ai[i]; j<ai[i+1]; j++) {
956:         col = baij->garray[aj[j]]*bs;
957:         for (k=0; k<bs; k++) {
958:           MatSetValues_MPIBAIJ(A,bs,rvals,1,&col,a,INSERT_VALUES);
959:           col++;
960:           a += bs;
961:         }
962:       }
963:     }
964:     PetscFree(rvals);
965:     MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
966:     MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
967:     /*
968:        Everyone has to call to draw the matrix since the graphics waits are
969:        synchronized across all processors that share the PetscDraw object
970:     */
971:     PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);
972:     PetscObjectGetName((PetscObject)mat,&matname);
973:     if (!rank) {
974:       PetscObjectSetName((PetscObject)((Mat_MPISBAIJ*)(A->data))->A,matname);
975:       MatView_SeqSBAIJ(((Mat_MPISBAIJ*)(A->data))->A,sviewer);
976:     }
977:     PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);
978:     PetscViewerFlush(viewer);
979:     MatDestroy(&A);
980:   }
981:   return(0);
982: }

984: /* Used for both MPIBAIJ and MPISBAIJ matrices */
985: #define MatView_MPISBAIJ_Binary MatView_MPIBAIJ_Binary

987: PetscErrorCode MatView_MPISBAIJ(Mat mat,PetscViewer viewer)
988: {
990:   PetscBool      iascii,isdraw,issocket,isbinary;

993:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
994:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
995:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSOCKET,&issocket);
996:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
997:   if (iascii || isdraw || issocket) {
998:     MatView_MPISBAIJ_ASCIIorDraworSocket(mat,viewer);
999:   } else if (isbinary) {
1000:     MatView_MPISBAIJ_Binary(mat,viewer);
1001:   }
1002:   return(0);
1003: }

1005: PetscErrorCode MatDestroy_MPISBAIJ(Mat mat)
1006: {
1007:   Mat_MPISBAIJ   *baij = (Mat_MPISBAIJ*)mat->data;

1011: #if defined(PETSC_USE_LOG)
1012:   PetscLogObjectState((PetscObject)mat,"Rows=%D,Cols=%D",mat->rmap->N,mat->cmap->N);
1013: #endif
1014:   MatStashDestroy_Private(&mat->stash);
1015:   MatStashDestroy_Private(&mat->bstash);
1016:   MatDestroy(&baij->A);
1017:   MatDestroy(&baij->B);
1018: #if defined(PETSC_USE_CTABLE)
1019:   PetscTableDestroy(&baij->colmap);
1020: #else
1021:   PetscFree(baij->colmap);
1022: #endif
1023:   PetscFree(baij->garray);
1024:   VecDestroy(&baij->lvec);
1025:   VecScatterDestroy(&baij->Mvctx);
1026:   VecDestroy(&baij->slvec0);
1027:   VecDestroy(&baij->slvec0b);
1028:   VecDestroy(&baij->slvec1);
1029:   VecDestroy(&baij->slvec1a);
1030:   VecDestroy(&baij->slvec1b);
1031:   VecScatterDestroy(&baij->sMvctx);
1032:   PetscFree2(baij->rowvalues,baij->rowindices);
1033:   PetscFree(baij->barray);
1034:   PetscFree(baij->hd);
1035:   VecDestroy(&baij->diag);
1036:   VecDestroy(&baij->bb1);
1037:   VecDestroy(&baij->xx1);
1038: #if defined(PETSC_USE_REAL_MAT_SINGLE)
1039:   PetscFree(baij->setvaluescopy);
1040: #endif
1041:   PetscFree(baij->in_loc);
1042:   PetscFree(baij->v_loc);
1043:   PetscFree(baij->rangebs);
1044:   PetscFree(mat->data);

1046:   PetscObjectChangeTypeName((PetscObject)mat,NULL);
1047:   PetscObjectComposeFunction((PetscObject)mat,"MatStoreValues_C",NULL);
1048:   PetscObjectComposeFunction((PetscObject)mat,"MatRetrieveValues_C",NULL);
1049:   PetscObjectComposeFunction((PetscObject)mat,"MatMPISBAIJSetPreallocation_C",NULL);
1050:   PetscObjectComposeFunction((PetscObject)mat,"MatMPISBAIJSetPreallocationCSR_C",NULL);
1051: #if defined(PETSC_HAVE_ELEMENTAL)
1052:   PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpisbaij_elemental_C",NULL);
1053: #endif
1054: #if defined(PETSC_HAVE_SCALAPACK)
1055:   PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpisbaij_scalapack_C",NULL);
1056: #endif
1057:   PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpisbaij_mpiaij_C",NULL);
1058:   PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpisbaij_mpibaij_C",NULL);
1059:   return(0);
1060: }

1062: PetscErrorCode MatMult_MPISBAIJ_Hermitian(Mat A,Vec xx,Vec yy)
1063: {
1064:   Mat_MPISBAIJ      *a = (Mat_MPISBAIJ*)A->data;
1065:   PetscErrorCode    ierr;
1066:   PetscInt          mbs=a->mbs,bs=A->rmap->bs;
1067:   PetscScalar       *from;
1068:   const PetscScalar *x;

1071:   /* diagonal part */
1072:   (*a->A->ops->mult)(a->A,xx,a->slvec1a);
1073:   VecSet(a->slvec1b,0.0);

1075:   /* subdiagonal part */
1076:   if (!a->B->ops->multhermitiantranspose) SETERRQ1(PetscObjectComm((PetscObject)a->B),PETSC_ERR_SUP,"Not for type %s\n",((PetscObject)a->B)->type_name);
1077:   (*a->B->ops->multhermitiantranspose)(a->B,xx,a->slvec0b);

1079:   /* copy x into the vec slvec0 */
1080:   VecGetArray(a->slvec0,&from);
1081:   VecGetArrayRead(xx,&x);

1083:   PetscArraycpy(from,x,bs*mbs);
1084:   VecRestoreArray(a->slvec0,&from);
1085:   VecRestoreArrayRead(xx,&x);

1087:   VecScatterBegin(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
1088:   VecScatterEnd(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
1089:   /* supperdiagonal part */
1090:   (*a->B->ops->multadd)(a->B,a->slvec1b,a->slvec1a,yy);
1091:   return(0);
1092: }

1094: PetscErrorCode MatMult_MPISBAIJ(Mat A,Vec xx,Vec yy)
1095: {
1096:   Mat_MPISBAIJ      *a = (Mat_MPISBAIJ*)A->data;
1097:   PetscErrorCode    ierr;
1098:   PetscInt          mbs=a->mbs,bs=A->rmap->bs;
1099:   PetscScalar       *from;
1100:   const PetscScalar *x;

1103:   /* diagonal part */
1104:   (*a->A->ops->mult)(a->A,xx,a->slvec1a);
1105:   VecSet(a->slvec1b,0.0);

1107:   /* subdiagonal part */
1108:   (*a->B->ops->multtranspose)(a->B,xx,a->slvec0b);

1110:   /* copy x into the vec slvec0 */
1111:   VecGetArray(a->slvec0,&from);
1112:   VecGetArrayRead(xx,&x);

1114:   PetscArraycpy(from,x,bs*mbs);
1115:   VecRestoreArray(a->slvec0,&from);
1116:   VecRestoreArrayRead(xx,&x);

1118:   VecScatterBegin(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
1119:   VecScatterEnd(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
1120:   /* supperdiagonal part */
1121:   (*a->B->ops->multadd)(a->B,a->slvec1b,a->slvec1a,yy);
1122:   return(0);
1123: }

1125: PetscErrorCode MatMultAdd_MPISBAIJ_Hermitian(Mat A,Vec xx,Vec yy,Vec zz)
1126: {
1127:   Mat_MPISBAIJ      *a = (Mat_MPISBAIJ*)A->data;
1128:   PetscErrorCode    ierr;
1129:   PetscInt          mbs=a->mbs,bs=A->rmap->bs;
1130:   PetscScalar       *from,zero=0.0;
1131:   const PetscScalar *x;

1134:   /* diagonal part */
1135:   (*a->A->ops->multadd)(a->A,xx,yy,a->slvec1a);
1136:   VecSet(a->slvec1b,zero);

1138:   /* subdiagonal part */
1139:   if (!a->B->ops->multhermitiantranspose) SETERRQ1(PetscObjectComm((PetscObject)a->B),PETSC_ERR_SUP,"Not for type %s\n",((PetscObject)a->B)->type_name);
1140:   (*a->B->ops->multhermitiantranspose)(a->B,xx,a->slvec0b);

1142:   /* copy x into the vec slvec0 */
1143:   VecGetArray(a->slvec0,&from);
1144:   VecGetArrayRead(xx,&x);
1145:   PetscArraycpy(from,x,bs*mbs);
1146:   VecRestoreArray(a->slvec0,&from);

1148:   VecScatterBegin(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
1149:   VecRestoreArrayRead(xx,&x);
1150:   VecScatterEnd(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);

1152:   /* supperdiagonal part */
1153:   (*a->B->ops->multadd)(a->B,a->slvec1b,a->slvec1a,zz);
1154:   return(0);
1155: }

1157: PetscErrorCode MatMultAdd_MPISBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
1158: {
1159:   Mat_MPISBAIJ      *a = (Mat_MPISBAIJ*)A->data;
1160:   PetscErrorCode    ierr;
1161:   PetscInt          mbs=a->mbs,bs=A->rmap->bs;
1162:   PetscScalar       *from,zero=0.0;
1163:   const PetscScalar *x;

1166:   /* diagonal part */
1167:   (*a->A->ops->multadd)(a->A,xx,yy,a->slvec1a);
1168:   VecSet(a->slvec1b,zero);

1170:   /* subdiagonal part */
1171:   (*a->B->ops->multtranspose)(a->B,xx,a->slvec0b);

1173:   /* copy x into the vec slvec0 */
1174:   VecGetArray(a->slvec0,&from);
1175:   VecGetArrayRead(xx,&x);
1176:   PetscArraycpy(from,x,bs*mbs);
1177:   VecRestoreArray(a->slvec0,&from);

1179:   VecScatterBegin(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
1180:   VecRestoreArrayRead(xx,&x);
1181:   VecScatterEnd(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);

1183:   /* supperdiagonal part */
1184:   (*a->B->ops->multadd)(a->B,a->slvec1b,a->slvec1a,zz);
1185:   return(0);
1186: }

1188: /*
1189:   This only works correctly for square matrices where the subblock A->A is the
1190:    diagonal block
1191: */
1192: PetscErrorCode MatGetDiagonal_MPISBAIJ(Mat A,Vec v)
1193: {
1194:   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;

1198:   /* if (a->rmap->N != a->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Supports only square matrix where A->A is diag block"); */
1199:   MatGetDiagonal(a->A,v);
1200:   return(0);
1201: }

1203: PetscErrorCode MatScale_MPISBAIJ(Mat A,PetscScalar aa)
1204: {
1205:   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;

1209:   MatScale(a->A,aa);
1210:   MatScale(a->B,aa);
1211:   return(0);
1212: }

1214: PetscErrorCode MatGetRow_MPISBAIJ(Mat matin,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
1215: {
1216:   Mat_MPISBAIJ   *mat = (Mat_MPISBAIJ*)matin->data;
1217:   PetscScalar    *vworkA,*vworkB,**pvA,**pvB,*v_p;
1219:   PetscInt       bs = matin->rmap->bs,bs2 = mat->bs2,i,*cworkA,*cworkB,**pcA,**pcB;
1220:   PetscInt       nztot,nzA,nzB,lrow,brstart = matin->rmap->rstart,brend = matin->rmap->rend;
1221:   PetscInt       *cmap,*idx_p,cstart = mat->rstartbs;

1224:   if (mat->getrowactive) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Already active");
1225:   mat->getrowactive = PETSC_TRUE;

1227:   if (!mat->rowvalues && (idx || v)) {
1228:     /*
1229:         allocate enough space to hold information from the longest row.
1230:     */
1231:     Mat_SeqSBAIJ *Aa = (Mat_SeqSBAIJ*)mat->A->data;
1232:     Mat_SeqBAIJ  *Ba = (Mat_SeqBAIJ*)mat->B->data;
1233:     PetscInt     max = 1,mbs = mat->mbs,tmp;
1234:     for (i=0; i<mbs; i++) {
1235:       tmp = Aa->i[i+1] - Aa->i[i] + Ba->i[i+1] - Ba->i[i]; /* row length */
1236:       if (max < tmp) max = tmp;
1237:     }
1238:     PetscMalloc2(max*bs2,&mat->rowvalues,max*bs2,&mat->rowindices);
1239:   }

1241:   if (row < brstart || row >= brend) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only local rows");
1242:   lrow = row - brstart;  /* local row index */

1244:   pvA = &vworkA; pcA = &cworkA; pvB = &vworkB; pcB = &cworkB;
1245:   if (!v)   {pvA = NULL; pvB = NULL;}
1246:   if (!idx) {pcA = NULL; if (!v) pcB = NULL;}
1247:   (*mat->A->ops->getrow)(mat->A,lrow,&nzA,pcA,pvA);
1248:   (*mat->B->ops->getrow)(mat->B,lrow,&nzB,pcB,pvB);
1249:   nztot = nzA + nzB;

1251:   cmap = mat->garray;
1252:   if (v  || idx) {
1253:     if (nztot) {
1254:       /* Sort by increasing column numbers, assuming A and B already sorted */
1255:       PetscInt imark = -1;
1256:       if (v) {
1257:         *v = v_p = mat->rowvalues;
1258:         for (i=0; i<nzB; i++) {
1259:           if (cmap[cworkB[i]/bs] < cstart) v_p[i] = vworkB[i];
1260:           else break;
1261:         }
1262:         imark = i;
1263:         for (i=0; i<nzA; i++)     v_p[imark+i] = vworkA[i];
1264:         for (i=imark; i<nzB; i++) v_p[nzA+i]   = vworkB[i];
1265:       }
1266:       if (idx) {
1267:         *idx = idx_p = mat->rowindices;
1268:         if (imark > -1) {
1269:           for (i=0; i<imark; i++) {
1270:             idx_p[i] = cmap[cworkB[i]/bs]*bs + cworkB[i]%bs;
1271:           }
1272:         } else {
1273:           for (i=0; i<nzB; i++) {
1274:             if (cmap[cworkB[i]/bs] < cstart) idx_p[i] = cmap[cworkB[i]/bs]*bs + cworkB[i]%bs;
1275:             else break;
1276:           }
1277:           imark = i;
1278:         }
1279:         for (i=0; i<nzA; i++)     idx_p[imark+i] = cstart*bs + cworkA[i];
1280:         for (i=imark; i<nzB; i++) idx_p[nzA+i]   = cmap[cworkB[i]/bs]*bs + cworkB[i]%bs ;
1281:       }
1282:     } else {
1283:       if (idx) *idx = NULL;
1284:       if (v)   *v   = NULL;
1285:     }
1286:   }
1287:   *nz  = nztot;
1288:   (*mat->A->ops->restorerow)(mat->A,lrow,&nzA,pcA,pvA);
1289:   (*mat->B->ops->restorerow)(mat->B,lrow,&nzB,pcB,pvB);
1290:   return(0);
1291: }

1293: PetscErrorCode MatRestoreRow_MPISBAIJ(Mat mat,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
1294: {
1295:   Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;

1298:   if (!baij->getrowactive) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"MatGetRow() must be called first");
1299:   baij->getrowactive = PETSC_FALSE;
1300:   return(0);
1301: }

1303: PetscErrorCode MatGetRowUpperTriangular_MPISBAIJ(Mat A)
1304: {
1305:   Mat_MPISBAIJ *a  = (Mat_MPISBAIJ*)A->data;
1306:   Mat_SeqSBAIJ *aA = (Mat_SeqSBAIJ*)a->A->data;

1309:   aA->getrow_utriangular = PETSC_TRUE;
1310:   return(0);
1311: }
1312: PetscErrorCode MatRestoreRowUpperTriangular_MPISBAIJ(Mat A)
1313: {
1314:   Mat_MPISBAIJ *a  = (Mat_MPISBAIJ*)A->data;
1315:   Mat_SeqSBAIJ *aA = (Mat_SeqSBAIJ*)a->A->data;

1318:   aA->getrow_utriangular = PETSC_FALSE;
1319:   return(0);
1320: }

1322: PetscErrorCode MatRealPart_MPISBAIJ(Mat A)
1323: {
1324:   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;

1328:   MatRealPart(a->A);
1329:   MatRealPart(a->B);
1330:   return(0);
1331: }

1333: PetscErrorCode MatImaginaryPart_MPISBAIJ(Mat A)
1334: {
1335:   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;

1339:   MatImaginaryPart(a->A);
1340:   MatImaginaryPart(a->B);
1341:   return(0);
1342: }

1344: /* Check if isrow is a subset of iscol_local, called by MatCreateSubMatrix_MPISBAIJ()
1345:    Input: isrow       - distributed(parallel),
1346:           iscol_local - locally owned (seq)
1347: */
1348: PetscErrorCode ISEqual_private(IS isrow,IS iscol_local,PetscBool  *flg)
1349: {
1351:   PetscInt       sz1,sz2,*a1,*a2,i,j,k,nmatch;
1352:   const PetscInt *ptr1,*ptr2;

1355:   ISGetLocalSize(isrow,&sz1);
1356:   ISGetLocalSize(iscol_local,&sz2);
1357:   if (sz1 > sz2) {
1358:     *flg = PETSC_FALSE;
1359:     return(0);
1360:   }

1362:   ISGetIndices(isrow,&ptr1);
1363:   ISGetIndices(iscol_local,&ptr2);

1365:   PetscMalloc1(sz1,&a1);
1366:   PetscMalloc1(sz2,&a2);
1367:   PetscArraycpy(a1,ptr1,sz1);
1368:   PetscArraycpy(a2,ptr2,sz2);
1369:   PetscSortInt(sz1,a1);
1370:   PetscSortInt(sz2,a2);

1372:   nmatch=0;
1373:   k     = 0;
1374:   for (i=0; i<sz1; i++){
1375:     for (j=k; j<sz2; j++){
1376:       if (a1[i] == a2[j]) {
1377:         k = j; nmatch++;
1378:         break;
1379:       }
1380:     }
1381:   }
1382:   ISRestoreIndices(isrow,&ptr1);
1383:   ISRestoreIndices(iscol_local,&ptr2);
1384:   PetscFree(a1);
1385:   PetscFree(a2);
1386:   if (nmatch < sz1) {
1387:     *flg = PETSC_FALSE;
1388:   } else {
1389:     *flg = PETSC_TRUE;
1390:   }
1391:   return(0);
1392: }

1394: PetscErrorCode MatCreateSubMatrix_MPISBAIJ(Mat mat,IS isrow,IS iscol,MatReuse call,Mat *newmat)
1395: {
1397:   IS             iscol_local;
1398:   PetscInt       csize;
1399:   PetscBool      isequal;

1402:   ISGetLocalSize(iscol,&csize);
1403:   if (call == MAT_REUSE_MATRIX) {
1404:     PetscObjectQuery((PetscObject)*newmat,"ISAllGather",(PetscObject*)&iscol_local);
1405:     if (!iscol_local) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Submatrix passed in was not used before, cannot reuse");
1406:   } else {
1407:     ISAllGather(iscol,&iscol_local);
1408:     ISEqual_private(isrow,iscol_local,&isequal);
1409:     if (!isequal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"For symmetric format, iscol must equal isrow");
1410:   }

1412:   /* now call MatCreateSubMatrix_MPIBAIJ() */
1413:   MatCreateSubMatrix_MPIBAIJ_Private(mat,isrow,iscol_local,csize,call,newmat);
1414:   if (call == MAT_INITIAL_MATRIX) {
1415:     PetscObjectCompose((PetscObject)*newmat,"ISAllGather",(PetscObject)iscol_local);
1416:     ISDestroy(&iscol_local);
1417:   }
1418:   return(0);
1419: }

1421: PetscErrorCode MatZeroEntries_MPISBAIJ(Mat A)
1422: {
1423:   Mat_MPISBAIJ   *l = (Mat_MPISBAIJ*)A->data;

1427:   MatZeroEntries(l->A);
1428:   MatZeroEntries(l->B);
1429:   return(0);
1430: }

1432: PetscErrorCode MatGetInfo_MPISBAIJ(Mat matin,MatInfoType flag,MatInfo *info)
1433: {
1434:   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)matin->data;
1435:   Mat            A  = a->A,B = a->B;
1437:   PetscLogDouble isend[5],irecv[5];

1440:   info->block_size = (PetscReal)matin->rmap->bs;

1442:   MatGetInfo(A,MAT_LOCAL,info);

1444:   isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
1445:   isend[3] = info->memory;  isend[4] = info->mallocs;

1447:   MatGetInfo(B,MAT_LOCAL,info);

1449:   isend[0] += info->nz_used; isend[1] += info->nz_allocated; isend[2] += info->nz_unneeded;
1450:   isend[3] += info->memory;  isend[4] += info->mallocs;
1451:   if (flag == MAT_LOCAL) {
1452:     info->nz_used      = isend[0];
1453:     info->nz_allocated = isend[1];
1454:     info->nz_unneeded  = isend[2];
1455:     info->memory       = isend[3];
1456:     info->mallocs      = isend[4];
1457:   } else if (flag == MAT_GLOBAL_MAX) {
1458:     MPIU_Allreduce(isend,irecv,5,MPIU_PETSCLOGDOUBLE,MPI_MAX,PetscObjectComm((PetscObject)matin));

1460:     info->nz_used      = irecv[0];
1461:     info->nz_allocated = irecv[1];
1462:     info->nz_unneeded  = irecv[2];
1463:     info->memory       = irecv[3];
1464:     info->mallocs      = irecv[4];
1465:   } else if (flag == MAT_GLOBAL_SUM) {
1466:     MPIU_Allreduce(isend,irecv,5,MPIU_PETSCLOGDOUBLE,MPI_SUM,PetscObjectComm((PetscObject)matin));

1468:     info->nz_used      = irecv[0];
1469:     info->nz_allocated = irecv[1];
1470:     info->nz_unneeded  = irecv[2];
1471:     info->memory       = irecv[3];
1472:     info->mallocs      = irecv[4];
1473:   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown MatInfoType argument %d",(int)flag);
1474:   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
1475:   info->fill_ratio_needed = 0;
1476:   info->factor_mallocs    = 0;
1477:   return(0);
1478: }

1480: PetscErrorCode MatSetOption_MPISBAIJ(Mat A,MatOption op,PetscBool flg)
1481: {
1482:   Mat_MPISBAIJ   *a  = (Mat_MPISBAIJ*)A->data;
1483:   Mat_SeqSBAIJ   *aA = (Mat_SeqSBAIJ*)a->A->data;

1487:   switch (op) {
1488:   case MAT_NEW_NONZERO_LOCATIONS:
1489:   case MAT_NEW_NONZERO_ALLOCATION_ERR:
1490:   case MAT_UNUSED_NONZERO_LOCATION_ERR:
1491:   case MAT_KEEP_NONZERO_PATTERN:
1492:   case MAT_SUBMAT_SINGLEIS:
1493:   case MAT_NEW_NONZERO_LOCATION_ERR:
1494:     MatCheckPreallocated(A,1);
1495:     MatSetOption(a->A,op,flg);
1496:     MatSetOption(a->B,op,flg);
1497:     break;
1498:   case MAT_ROW_ORIENTED:
1499:     MatCheckPreallocated(A,1);
1500:     a->roworiented = flg;

1502:     MatSetOption(a->A,op,flg);
1503:     MatSetOption(a->B,op,flg);
1504:     break;
1505:   case MAT_NEW_DIAGONALS:
1506:   case MAT_SORTED_FULL:
1507:     PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
1508:     break;
1509:   case MAT_IGNORE_OFF_PROC_ENTRIES:
1510:     a->donotstash = flg;
1511:     break;
1512:   case MAT_USE_HASH_TABLE:
1513:     a->ht_flag = flg;
1514:     break;
1515:   case MAT_HERMITIAN:
1516:     MatCheckPreallocated(A,1);
1517:     MatSetOption(a->A,op,flg);
1518: #if defined(PETSC_USE_COMPLEX)
1519:     if (flg) { /* need different mat-vec ops */
1520:       A->ops->mult             = MatMult_MPISBAIJ_Hermitian;
1521:       A->ops->multadd          = MatMultAdd_MPISBAIJ_Hermitian;
1522:       A->ops->multtranspose    = NULL;
1523:       A->ops->multtransposeadd = NULL;
1524:       A->symmetric = PETSC_FALSE;
1525:     }
1526: #endif
1527:     break;
1528:   case MAT_SPD:
1529:   case MAT_SYMMETRIC:
1530:     MatCheckPreallocated(A,1);
1531:     MatSetOption(a->A,op,flg);
1532: #if defined(PETSC_USE_COMPLEX)
1533:     if (flg) { /* restore to use default mat-vec ops */
1534:       A->ops->mult             = MatMult_MPISBAIJ;
1535:       A->ops->multadd          = MatMultAdd_MPISBAIJ;
1536:       A->ops->multtranspose    = MatMult_MPISBAIJ;
1537:       A->ops->multtransposeadd = MatMultAdd_MPISBAIJ;
1538:     }
1539: #endif
1540:     break;
1541:   case MAT_STRUCTURALLY_SYMMETRIC:
1542:     MatCheckPreallocated(A,1);
1543:     MatSetOption(a->A,op,flg);
1544:     break;
1545:   case MAT_SYMMETRY_ETERNAL:
1546:     if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix must be symmetric");
1547:     PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
1548:     break;
1549:   case MAT_IGNORE_LOWER_TRIANGULAR:
1550:     aA->ignore_ltriangular = flg;
1551:     break;
1552:   case MAT_ERROR_LOWER_TRIANGULAR:
1553:     aA->ignore_ltriangular = flg;
1554:     break;
1555:   case MAT_GETROW_UPPERTRIANGULAR:
1556:     aA->getrow_utriangular = flg;
1557:     break;
1558:   default:
1559:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %d",op);
1560:   }
1561:   return(0);
1562: }

1564: PetscErrorCode MatTranspose_MPISBAIJ(Mat A,MatReuse reuse,Mat *B)
1565: {

1569:   if (reuse == MAT_INITIAL_MATRIX) {
1570:     MatDuplicate(A,MAT_COPY_VALUES,B);
1571:   }  else if (reuse == MAT_REUSE_MATRIX) {
1572:     MatCopy(A,*B,SAME_NONZERO_PATTERN);
1573:   }
1574:   return(0);
1575: }

1577: PetscErrorCode MatDiagonalScale_MPISBAIJ(Mat mat,Vec ll,Vec rr)
1578: {
1579:   Mat_MPISBAIJ   *baij = (Mat_MPISBAIJ*)mat->data;
1580:   Mat            a     = baij->A, b=baij->B;
1582:   PetscInt       nv,m,n;
1583:   PetscBool      flg;

1586:   if (ll != rr) {
1587:     VecEqual(ll,rr,&flg);
1588:     if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"For symmetric format, left and right scaling vectors must be same\n");
1589:   }
1590:   if (!ll) return(0);

1592:   MatGetLocalSize(mat,&m,&n);
1593:   if (m != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"For symmetric format, local size %d %d must be same",m,n);

1595:   VecGetLocalSize(rr,&nv);
1596:   if (nv!=n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left and right vector non-conforming local size");

1598:   VecScatterBegin(baij->Mvctx,rr,baij->lvec,INSERT_VALUES,SCATTER_FORWARD);

1600:   /* left diagonalscale the off-diagonal part */
1601:   (*b->ops->diagonalscale)(b,ll,NULL);

1603:   /* scale the diagonal part */
1604:   (*a->ops->diagonalscale)(a,ll,rr);

1606:   /* right diagonalscale the off-diagonal part */
1607:   VecScatterEnd(baij->Mvctx,rr,baij->lvec,INSERT_VALUES,SCATTER_FORWARD);
1608:   (*b->ops->diagonalscale)(b,NULL,baij->lvec);
1609:   return(0);
1610: }

1612: PetscErrorCode MatSetUnfactored_MPISBAIJ(Mat A)
1613: {
1614:   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;

1618:   MatSetUnfactored(a->A);
1619:   return(0);
1620: }

1622: static PetscErrorCode MatDuplicate_MPISBAIJ(Mat,MatDuplicateOption,Mat*);

1624: PetscErrorCode MatEqual_MPISBAIJ(Mat A,Mat B,PetscBool  *flag)
1625: {
1626:   Mat_MPISBAIJ   *matB = (Mat_MPISBAIJ*)B->data,*matA = (Mat_MPISBAIJ*)A->data;
1627:   Mat            a,b,c,d;
1628:   PetscBool      flg;

1632:   a = matA->A; b = matA->B;
1633:   c = matB->A; d = matB->B;

1635:   MatEqual(a,c,&flg);
1636:   if (flg) {
1637:     MatEqual(b,d,&flg);
1638:   }
1639:   MPIU_Allreduce(&flg,flag,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));
1640:   return(0);
1641: }

1643: PetscErrorCode MatCopy_MPISBAIJ(Mat A,Mat B,MatStructure str)
1644: {
1646:   PetscBool      isbaij;

1649:   PetscObjectTypeCompareAny((PetscObject)B,&isbaij,MATSEQSBAIJ,MATMPISBAIJ,"");
1650:   if (!isbaij) SETERRQ1(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"Not for matrix type %s",((PetscObject)B)->type_name);
1651:   /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
1652:   if ((str != SAME_NONZERO_PATTERN) || (A->ops->copy != B->ops->copy)) {
1653:     MatGetRowUpperTriangular(A);
1654:     MatCopy_Basic(A,B,str);
1655:     MatRestoreRowUpperTriangular(A);
1656:   } else {
1657:     Mat_MPISBAIJ *a = (Mat_MPISBAIJ*)A->data;
1658:     Mat_MPISBAIJ *b = (Mat_MPISBAIJ*)B->data;

1660:     MatCopy(a->A,b->A,str);
1661:     MatCopy(a->B,b->B,str);
1662:   }
1663:   PetscObjectStateIncrease((PetscObject)B);
1664:   return(0);
1665: }

1667: PetscErrorCode MatSetUp_MPISBAIJ(Mat A)
1668: {

1672:   MatMPISBAIJSetPreallocation(A,A->rmap->bs,PETSC_DEFAULT,NULL,PETSC_DEFAULT,NULL);
1673:   return(0);
1674: }

1676: PetscErrorCode MatAXPY_MPISBAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str)
1677: {
1679:   Mat_MPISBAIJ   *xx=(Mat_MPISBAIJ*)X->data,*yy=(Mat_MPISBAIJ*)Y->data;
1680:   PetscBLASInt   bnz,one=1;
1681:   Mat_SeqSBAIJ   *xa,*ya;
1682:   Mat_SeqBAIJ    *xb,*yb;

1685:   if (str == SAME_NONZERO_PATTERN) {
1686:     PetscScalar alpha = a;
1687:     xa   = (Mat_SeqSBAIJ*)xx->A->data;
1688:     ya   = (Mat_SeqSBAIJ*)yy->A->data;
1689:     PetscBLASIntCast(xa->nz,&bnz);
1690:     PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&bnz,&alpha,xa->a,&one,ya->a,&one));
1691:     xb   = (Mat_SeqBAIJ*)xx->B->data;
1692:     yb   = (Mat_SeqBAIJ*)yy->B->data;
1693:     PetscBLASIntCast(xb->nz,&bnz);
1694:     PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&bnz,&alpha,xb->a,&one,yb->a,&one));
1695:     PetscObjectStateIncrease((PetscObject)Y);
1696:   } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
1697:     MatSetOption(X,MAT_GETROW_UPPERTRIANGULAR,PETSC_TRUE);
1698:     MatAXPY_Basic(Y,a,X,str);
1699:     MatSetOption(X,MAT_GETROW_UPPERTRIANGULAR,PETSC_FALSE);
1700:   } else {
1701:     Mat      B;
1702:     PetscInt *nnz_d,*nnz_o,bs=Y->rmap->bs;
1703:     if (bs != X->rmap->bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrices must have same block size");
1704:     MatGetRowUpperTriangular(X);
1705:     MatGetRowUpperTriangular(Y);
1706:     PetscMalloc1(yy->A->rmap->N,&nnz_d);
1707:     PetscMalloc1(yy->B->rmap->N,&nnz_o);
1708:     MatCreate(PetscObjectComm((PetscObject)Y),&B);
1709:     PetscObjectSetName((PetscObject)B,((PetscObject)Y)->name);
1710:     MatSetSizes(B,Y->rmap->n,Y->cmap->n,Y->rmap->N,Y->cmap->N);
1711:     MatSetBlockSizesFromMats(B,Y,Y);
1712:     MatSetType(B,MATMPISBAIJ);
1713:     MatAXPYGetPreallocation_SeqSBAIJ(yy->A,xx->A,nnz_d);
1714:     MatAXPYGetPreallocation_MPIBAIJ(yy->B,yy->garray,xx->B,xx->garray,nnz_o);
1715:     MatMPISBAIJSetPreallocation(B,bs,0,nnz_d,0,nnz_o);
1716:     MatAXPY_BasicWithPreallocation(B,Y,a,X,str);
1717:     MatHeaderReplace(Y,&B);
1718:     PetscFree(nnz_d);
1719:     PetscFree(nnz_o);
1720:     MatRestoreRowUpperTriangular(X);
1721:     MatRestoreRowUpperTriangular(Y);
1722:   }
1723:   return(0);
1724: }

1726: PetscErrorCode MatCreateSubMatrices_MPISBAIJ(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1727: {
1729:   PetscInt       i;
1730:   PetscBool      flg;

1733:   MatCreateSubMatrices_MPIBAIJ(A,n,irow,icol,scall,B); /* B[] are sbaij matrices */
1734:   for (i=0; i<n; i++) {
1735:     ISEqual(irow[i],icol[i],&flg);
1736:     if (!flg) {
1737:       MatSeqSBAIJZeroOps_Private(*B[i]);
1738:     }
1739:   }
1740:   return(0);
1741: }

1743: PetscErrorCode MatShift_MPISBAIJ(Mat Y,PetscScalar a)
1744: {
1746:   Mat_MPISBAIJ    *maij = (Mat_MPISBAIJ*)Y->data;
1747:   Mat_SeqSBAIJ    *aij = (Mat_SeqSBAIJ*)maij->A->data;

1750:   if (!Y->preallocated) {
1751:     MatMPISBAIJSetPreallocation(Y,Y->rmap->bs,1,NULL,0,NULL);
1752:   } else if (!aij->nz) {
1753:     PetscInt nonew = aij->nonew;
1754:     MatSeqSBAIJSetPreallocation(maij->A,Y->rmap->bs,1,NULL);
1755:     aij->nonew = nonew;
1756:   }
1757:   MatShift_Basic(Y,a);
1758:   return(0);
1759: }

1761: PetscErrorCode MatMissingDiagonal_MPISBAIJ(Mat A,PetscBool  *missing,PetscInt *d)
1762: {
1763:   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;

1767:   if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only works for square matrices");
1768:   MatMissingDiagonal(a->A,missing,d);
1769:   if (d) {
1770:     PetscInt rstart;
1771:     MatGetOwnershipRange(A,&rstart,NULL);
1772:     *d += rstart/A->rmap->bs;

1774:   }
1775:   return(0);
1776: }

1778: PetscErrorCode  MatGetDiagonalBlock_MPISBAIJ(Mat A,Mat *a)
1779: {
1781:   *a = ((Mat_MPISBAIJ*)A->data)->A;
1782:   return(0);
1783: }

1785: /* -------------------------------------------------------------------*/
1786: static struct _MatOps MatOps_Values = {MatSetValues_MPISBAIJ,
1787:                                        MatGetRow_MPISBAIJ,
1788:                                        MatRestoreRow_MPISBAIJ,
1789:                                        MatMult_MPISBAIJ,
1790:                                /*  4*/ MatMultAdd_MPISBAIJ,
1791:                                        MatMult_MPISBAIJ,       /* transpose versions are same as non-transpose */
1792:                                        MatMultAdd_MPISBAIJ,
1793:                                        NULL,
1794:                                        NULL,
1795:                                        NULL,
1796:                                /* 10*/ NULL,
1797:                                        NULL,
1798:                                        NULL,
1799:                                        MatSOR_MPISBAIJ,
1800:                                        MatTranspose_MPISBAIJ,
1801:                                /* 15*/ MatGetInfo_MPISBAIJ,
1802:                                        MatEqual_MPISBAIJ,
1803:                                        MatGetDiagonal_MPISBAIJ,
1804:                                        MatDiagonalScale_MPISBAIJ,
1805:                                        MatNorm_MPISBAIJ,
1806:                                /* 20*/ MatAssemblyBegin_MPISBAIJ,
1807:                                        MatAssemblyEnd_MPISBAIJ,
1808:                                        MatSetOption_MPISBAIJ,
1809:                                        MatZeroEntries_MPISBAIJ,
1810:                                /* 24*/ NULL,
1811:                                        NULL,
1812:                                        NULL,
1813:                                        NULL,
1814:                                        NULL,
1815:                                /* 29*/ MatSetUp_MPISBAIJ,
1816:                                        NULL,
1817:                                        NULL,
1818:                                        MatGetDiagonalBlock_MPISBAIJ,
1819:                                        NULL,
1820:                                /* 34*/ MatDuplicate_MPISBAIJ,
1821:                                        NULL,
1822:                                        NULL,
1823:                                        NULL,
1824:                                        NULL,
1825:                                /* 39*/ MatAXPY_MPISBAIJ,
1826:                                        MatCreateSubMatrices_MPISBAIJ,
1827:                                        MatIncreaseOverlap_MPISBAIJ,
1828:                                        MatGetValues_MPISBAIJ,
1829:                                        MatCopy_MPISBAIJ,
1830:                                /* 44*/ NULL,
1831:                                        MatScale_MPISBAIJ,
1832:                                        MatShift_MPISBAIJ,
1833:                                        NULL,
1834:                                        NULL,
1835:                                /* 49*/ NULL,
1836:                                        NULL,
1837:                                        NULL,
1838:                                        NULL,
1839:                                        NULL,
1840:                                /* 54*/ NULL,
1841:                                        NULL,
1842:                                        MatSetUnfactored_MPISBAIJ,
1843:                                        NULL,
1844:                                        MatSetValuesBlocked_MPISBAIJ,
1845:                                /* 59*/ MatCreateSubMatrix_MPISBAIJ,
1846:                                        NULL,
1847:                                        NULL,
1848:                                        NULL,
1849:                                        NULL,
1850:                                /* 64*/ NULL,
1851:                                        NULL,
1852:                                        NULL,
1853:                                        NULL,
1854:                                        NULL,
1855:                                /* 69*/ MatGetRowMaxAbs_MPISBAIJ,
1856:                                        NULL,
1857:                                        MatConvert_MPISBAIJ_Basic,
1858:                                        NULL,
1859:                                        NULL,
1860:                                /* 74*/ NULL,
1861:                                        NULL,
1862:                                        NULL,
1863:                                        NULL,
1864:                                        NULL,
1865:                                /* 79*/ NULL,
1866:                                        NULL,
1867:                                        NULL,
1868:                                        NULL,
1869:                                        MatLoad_MPISBAIJ,
1870:                                /* 84*/ NULL,
1871:                                        NULL,
1872:                                        NULL,
1873:                                        NULL,
1874:                                        NULL,
1875:                                /* 89*/ NULL,
1876:                                        NULL,
1877:                                        NULL,
1878:                                        NULL,
1879:                                        NULL,
1880:                                /* 94*/ NULL,
1881:                                        NULL,
1882:                                        NULL,
1883:                                        NULL,
1884:                                        NULL,
1885:                                /* 99*/ NULL,
1886:                                        NULL,
1887:                                        NULL,
1888:                                        NULL,
1889:                                        NULL,
1890:                                /*104*/ NULL,
1891:                                        MatRealPart_MPISBAIJ,
1892:                                        MatImaginaryPart_MPISBAIJ,
1893:                                        MatGetRowUpperTriangular_MPISBAIJ,
1894:                                        MatRestoreRowUpperTriangular_MPISBAIJ,
1895:                                /*109*/ NULL,
1896:                                        NULL,
1897:                                        NULL,
1898:                                        NULL,
1899:                                        MatMissingDiagonal_MPISBAIJ,
1900:                                /*114*/ NULL,
1901:                                        NULL,
1902:                                        NULL,
1903:                                        NULL,
1904:                                        NULL,
1905:                                /*119*/ NULL,
1906:                                        NULL,
1907:                                        NULL,
1908:                                        NULL,
1909:                                        NULL,
1910:                                /*124*/ NULL,
1911:                                        NULL,
1912:                                        NULL,
1913:                                        NULL,
1914:                                        NULL,
1915:                                /*129*/ NULL,
1916:                                        NULL,
1917:                                        NULL,
1918:                                        NULL,
1919:                                        NULL,
1920:                                /*134*/ NULL,
1921:                                        NULL,
1922:                                        NULL,
1923:                                        NULL,
1924:                                        NULL,
1925:                                /*139*/ MatSetBlockSizes_Default,
1926:                                        NULL,
1927:                                        NULL,
1928:                                        NULL,
1929:                                        NULL,
1930:                                 /*144*/MatCreateMPIMatConcatenateSeqMat_MPISBAIJ
1931: };

1933: PetscErrorCode  MatMPISBAIJSetPreallocation_MPISBAIJ(Mat B,PetscInt bs,PetscInt d_nz,const PetscInt *d_nnz,PetscInt o_nz,const PetscInt *o_nnz)
1934: {
1935:   Mat_MPISBAIJ   *b = (Mat_MPISBAIJ*)B->data;
1937:   PetscInt       i,mbs,Mbs;
1938:   PetscMPIInt    size;

1941:   MatSetBlockSize(B,PetscAbs(bs));
1942:   PetscLayoutSetUp(B->rmap);
1943:   PetscLayoutSetUp(B->cmap);
1944:   PetscLayoutGetBlockSize(B->rmap,&bs);
1945:   if (B->rmap->N > B->cmap->N) SETERRQ2(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"MPISBAIJ matrix cannot have more rows %D than columns %D",B->rmap->N,B->cmap->N);
1946:   if (B->rmap->n > B->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"MPISBAIJ matrix cannot have more local rows %D than columns %D",B->rmap->n,B->cmap->n);

1948:   mbs = B->rmap->n/bs;
1949:   Mbs = B->rmap->N/bs;
1950:   if (mbs*bs != B->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"No of local rows %D must be divisible by blocksize %D",B->rmap->N,bs);

1952:   B->rmap->bs = bs;
1953:   b->bs2      = bs*bs;
1954:   b->mbs      = mbs;
1955:   b->Mbs      = Mbs;
1956:   b->nbs      = B->cmap->n/bs;
1957:   b->Nbs      = B->cmap->N/bs;

1959:   for (i=0; i<=b->size; i++) {
1960:     b->rangebs[i] = B->rmap->range[i]/bs;
1961:   }
1962:   b->rstartbs = B->rmap->rstart/bs;
1963:   b->rendbs   = B->rmap->rend/bs;

1965:   b->cstartbs = B->cmap->rstart/bs;
1966:   b->cendbs   = B->cmap->rend/bs;

1968: #if defined(PETSC_USE_CTABLE)
1969:   PetscTableDestroy(&b->colmap);
1970: #else
1971:   PetscFree(b->colmap);
1972: #endif
1973:   PetscFree(b->garray);
1974:   VecDestroy(&b->lvec);
1975:   VecScatterDestroy(&b->Mvctx);
1976:   VecDestroy(&b->slvec0);
1977:   VecDestroy(&b->slvec0b);
1978:   VecDestroy(&b->slvec1);
1979:   VecDestroy(&b->slvec1a);
1980:   VecDestroy(&b->slvec1b);
1981:   VecScatterDestroy(&b->sMvctx);

1983:   /* Because the B will have been resized we simply destroy it and create a new one each time */
1984:   MPI_Comm_size(PetscObjectComm((PetscObject)B),&size);
1985:   MatDestroy(&b->B);
1986:   MatCreate(PETSC_COMM_SELF,&b->B);
1987:   MatSetSizes(b->B,B->rmap->n,size > 1 ? B->cmap->N : 0,B->rmap->n,size > 1 ? B->cmap->N : 0);
1988:   MatSetType(b->B,MATSEQBAIJ);
1989:   PetscLogObjectParent((PetscObject)B,(PetscObject)b->B);

1991:   if (!B->preallocated) {
1992:     MatCreate(PETSC_COMM_SELF,&b->A);
1993:     MatSetSizes(b->A,B->rmap->n,B->cmap->n,B->rmap->n,B->cmap->n);
1994:     MatSetType(b->A,MATSEQSBAIJ);
1995:     PetscLogObjectParent((PetscObject)B,(PetscObject)b->A);
1996:     MatStashCreate_Private(PetscObjectComm((PetscObject)B),bs,&B->bstash);
1997:   }

1999:   MatSeqSBAIJSetPreallocation(b->A,bs,d_nz,d_nnz);
2000:   MatSeqBAIJSetPreallocation(b->B,bs,o_nz,o_nnz);

2002:   B->preallocated  = PETSC_TRUE;
2003:   B->was_assembled = PETSC_FALSE;
2004:   B->assembled     = PETSC_FALSE;
2005:   return(0);
2006: }

2008: PetscErrorCode MatMPISBAIJSetPreallocationCSR_MPISBAIJ(Mat B,PetscInt bs,const PetscInt ii[],const PetscInt jj[],const PetscScalar V[])
2009: {
2010:   PetscInt       m,rstart,cend;
2011:   PetscInt       i,j,d,nz,bd, nz_max=0,*d_nnz=NULL,*o_nnz=NULL;
2012:   const PetscInt *JJ    =NULL;
2013:   PetscScalar    *values=NULL;
2014:   PetscBool      roworiented = ((Mat_MPISBAIJ*)B->data)->roworiented;
2016:   PetscBool      nooffprocentries;

2019:   if (bs < 1) SETERRQ1(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_OUTOFRANGE,"Invalid block size specified, must be positive but it is %D",bs);
2020:   PetscLayoutSetBlockSize(B->rmap,bs);
2021:   PetscLayoutSetBlockSize(B->cmap,bs);
2022:   PetscLayoutSetUp(B->rmap);
2023:   PetscLayoutSetUp(B->cmap);
2024:   PetscLayoutGetBlockSize(B->rmap,&bs);
2025:   m      = B->rmap->n/bs;
2026:   rstart = B->rmap->rstart/bs;
2027:   cend   = B->cmap->rend/bs;

2029:   if (ii[0]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"ii[0] must be 0 but it is %D",ii[0]);
2030:   PetscMalloc2(m,&d_nnz,m,&o_nnz);
2031:   for (i=0; i<m; i++) {
2032:     nz = ii[i+1] - ii[i];
2033:     if (nz < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local row %D has a negative number of columns %D",i,nz);
2034:     /* count the ones on the diagonal and above, split into diagonal and off diagonal portions. */
2035:     JJ     = jj + ii[i];
2036:     bd     = 0;
2037:     for (j=0; j<nz; j++) {
2038:       if (*JJ >= i + rstart) break;
2039:       JJ++;
2040:       bd++;
2041:     }
2042:     d  = 0;
2043:     for (; j<nz; j++) {
2044:       if (*JJ++ >= cend) break;
2045:       d++;
2046:     }
2047:     d_nnz[i] = d;
2048:     o_nnz[i] = nz - d - bd;
2049:     nz       = nz - bd;
2050:     nz_max = PetscMax(nz_max,nz);
2051:   }
2052:   MatMPISBAIJSetPreallocation(B,bs,0,d_nnz,0,o_nnz);
2053:   MatSetOption(B,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);
2054:   PetscFree2(d_nnz,o_nnz);

2056:   values = (PetscScalar*)V;
2057:   if (!values) {
2058:     PetscCalloc1(bs*bs*nz_max,&values);
2059:   }
2060:   for (i=0; i<m; i++) {
2061:     PetscInt          row    = i + rstart;
2062:     PetscInt          ncols  = ii[i+1] - ii[i];
2063:     const PetscInt    *icols = jj + ii[i];
2064:     if (bs == 1 || !roworiented) {         /* block ordering matches the non-nested layout of MatSetValues so we can insert entire rows */
2065:       const PetscScalar *svals = values + (V ? (bs*bs*ii[i]) : 0);
2066:       MatSetValuesBlocked_MPISBAIJ(B,1,&row,ncols,icols,svals,INSERT_VALUES);
2067:     } else {                    /* block ordering does not match so we can only insert one block at a time. */
2068:       PetscInt j;
2069:       for (j=0; j<ncols; j++) {
2070:         const PetscScalar *svals = values + (V ? (bs*bs*(ii[i]+j)) : 0);
2071:         MatSetValuesBlocked_MPISBAIJ(B,1,&row,1,&icols[j],svals,INSERT_VALUES);
2072:       }
2073:     }
2074:   }

2076:   if (!V) { PetscFree(values); }
2077:   nooffprocentries    = B->nooffprocentries;
2078:   B->nooffprocentries = PETSC_TRUE;
2079:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
2080:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
2081:   B->nooffprocentries = nooffprocentries;

2083:   MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
2084:   return(0);
2085: }

2087: /*MC
2088:    MATMPISBAIJ - MATMPISBAIJ = "mpisbaij" - A matrix type to be used for distributed symmetric sparse block matrices,
2089:    based on block compressed sparse row format.  Only the upper triangular portion of the "diagonal" portion of
2090:    the matrix is stored.

2092:    For complex numbers by default this matrix is symmetric, NOT Hermitian symmetric. To make it Hermitian symmetric you
2093:    can call MatSetOption(Mat, MAT_HERMITIAN);

2095:    Options Database Keys:
2096: . -mat_type mpisbaij - sets the matrix type to "mpisbaij" during a call to MatSetFromOptions()

2098:    Notes:
2099:      The number of rows in the matrix must be less than or equal to the number of columns. Similarly the number of rows in the
2100:      diagonal portion of the matrix of each process has to less than or equal the number of columns.

2102:    Level: beginner

2104: .seealso: MatCreateBAIJ(), MATSEQSBAIJ, MatType
2105: M*/

2107: PETSC_EXTERN PetscErrorCode MatCreate_MPISBAIJ(Mat B)
2108: {
2109:   Mat_MPISBAIJ   *b;
2111:   PetscBool      flg = PETSC_FALSE;

2114:   PetscNewLog(B,&b);
2115:   B->data = (void*)b;
2116:   PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));

2118:   B->ops->destroy = MatDestroy_MPISBAIJ;
2119:   B->ops->view    = MatView_MPISBAIJ;
2120:   B->assembled    = PETSC_FALSE;
2121:   B->insertmode   = NOT_SET_VALUES;

2123:   MPI_Comm_rank(PetscObjectComm((PetscObject)B),&b->rank);
2124:   MPI_Comm_size(PetscObjectComm((PetscObject)B),&b->size);

2126:   /* build local table of row and column ownerships */
2127:   PetscMalloc1(b->size+2,&b->rangebs);

2129:   /* build cache for off array entries formed */
2130:   MatStashCreate_Private(PetscObjectComm((PetscObject)B),1,&B->stash);

2132:   b->donotstash  = PETSC_FALSE;
2133:   b->colmap      = NULL;
2134:   b->garray      = NULL;
2135:   b->roworiented = PETSC_TRUE;

2137:   /* stuff used in block assembly */
2138:   b->barray = NULL;

2140:   /* stuff used for matrix vector multiply */
2141:   b->lvec    = NULL;
2142:   b->Mvctx   = NULL;
2143:   b->slvec0  = NULL;
2144:   b->slvec0b = NULL;
2145:   b->slvec1  = NULL;
2146:   b->slvec1a = NULL;
2147:   b->slvec1b = NULL;
2148:   b->sMvctx  = NULL;

2150:   /* stuff for MatGetRow() */
2151:   b->rowindices   = NULL;
2152:   b->rowvalues    = NULL;
2153:   b->getrowactive = PETSC_FALSE;

2155:   /* hash table stuff */
2156:   b->ht           = NULL;
2157:   b->hd           = NULL;
2158:   b->ht_size      = 0;
2159:   b->ht_flag      = PETSC_FALSE;
2160:   b->ht_fact      = 0;
2161:   b->ht_total_ct  = 0;
2162:   b->ht_insert_ct = 0;

2164:   /* stuff for MatCreateSubMatrices_MPIBAIJ_local() */
2165:   b->ijonly = PETSC_FALSE;

2167:   b->in_loc = NULL;
2168:   b->v_loc  = NULL;
2169:   b->n_loc  = 0;

2171:   PetscObjectComposeFunction((PetscObject)B,"MatStoreValues_C",MatStoreValues_MPISBAIJ);
2172:   PetscObjectComposeFunction((PetscObject)B,"MatRetrieveValues_C",MatRetrieveValues_MPISBAIJ);
2173:   PetscObjectComposeFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",MatMPISBAIJSetPreallocation_MPISBAIJ);
2174:   PetscObjectComposeFunction((PetscObject)B,"MatMPISBAIJSetPreallocationCSR_C",MatMPISBAIJSetPreallocationCSR_MPISBAIJ);
2175: #if defined(PETSC_HAVE_ELEMENTAL)
2176:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpisbaij_elemental_C",MatConvert_MPISBAIJ_Elemental);
2177: #endif
2178: #if defined(PETSC_HAVE_SCALAPACK)
2179:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpisbaij_scalapack_C",MatConvert_SBAIJ_ScaLAPACK);
2180: #endif
2181:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpisbaij_mpiaij_C",MatConvert_MPISBAIJ_Basic);
2182:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpisbaij_mpibaij_C",MatConvert_MPISBAIJ_Basic);

2184:   B->symmetric                  = PETSC_TRUE;
2185:   B->structurally_symmetric     = PETSC_TRUE;
2186:   B->symmetric_set              = PETSC_TRUE;
2187:   B->structurally_symmetric_set = PETSC_TRUE;
2188:   B->symmetric_eternal          = PETSC_TRUE;
2189: #if defined(PETSC_USE_COMPLEX)
2190:   B->hermitian                  = PETSC_FALSE;
2191:   B->hermitian_set              = PETSC_FALSE;
2192: #else
2193:   B->hermitian                  = PETSC_TRUE;
2194:   B->hermitian_set              = PETSC_TRUE;
2195: #endif

2197:   PetscObjectChangeTypeName((PetscObject)B,MATMPISBAIJ);
2198:   PetscOptionsBegin(PetscObjectComm((PetscObject)B),NULL,"Options for loading MPISBAIJ matrix 1","Mat");
2199:   PetscOptionsBool("-mat_use_hash_table","Use hash table to save memory in constructing matrix","MatSetOption",flg,&flg,NULL);
2200:   if (flg) {
2201:     PetscReal fact = 1.39;
2202:     MatSetOption(B,MAT_USE_HASH_TABLE,PETSC_TRUE);
2203:     PetscOptionsReal("-mat_use_hash_table","Use hash table factor","MatMPIBAIJSetHashTableFactor",fact,&fact,NULL);
2204:     if (fact <= 1.0) fact = 1.39;
2205:     MatMPIBAIJSetHashTableFactor(B,fact);
2206:     PetscInfo1(B,"Hash table Factor used %5.2f\n",fact);
2207:   }
2208:   PetscOptionsEnd();
2209:   return(0);
2210: }

2212: /*MC
2213:    MATSBAIJ - MATSBAIJ = "sbaij" - A matrix type to be used for symmetric block sparse matrices.

2215:    This matrix type is identical to MATSEQSBAIJ when constructed with a single process communicator,
2216:    and MATMPISBAIJ otherwise.

2218:    Options Database Keys:
2219: . -mat_type sbaij - sets the matrix type to "sbaij" during a call to MatSetFromOptions()

2221:   Level: beginner

2223: .seealso: MatCreateMPISBAIJ, MATSEQSBAIJ, MATMPISBAIJ
2224: M*/

2226: /*@C
2227:    MatMPISBAIJSetPreallocation - For good matrix assembly performance
2228:    the user should preallocate the matrix storage by setting the parameters
2229:    d_nz (or d_nnz) and o_nz (or o_nnz).  By setting these parameters accurately,
2230:    performance can be increased by more than a factor of 50.

2232:    Collective on Mat

2234:    Input Parameters:
2235: +  B - the matrix
2236: .  bs   - size of block, the blocks are ALWAYS square. One can use MatSetBlockSizes() to set a different row and column blocksize but the row
2237:           blocksize always defines the size of the blocks. The column blocksize sets the blocksize of the vectors obtained with MatCreateVecs()
2238: .  d_nz  - number of block nonzeros per block row in diagonal portion of local
2239:            submatrix  (same for all local rows)
2240: .  d_nnz - array containing the number of block nonzeros in the various block rows
2241:            in the upper triangular and diagonal part of the in diagonal portion of the local
2242:            (possibly different for each block row) or NULL.  If you plan to factor the matrix you must leave room
2243:            for the diagonal entry and set a value even if it is zero.
2244: .  o_nz  - number of block nonzeros per block row in the off-diagonal portion of local
2245:            submatrix (same for all local rows).
2246: -  o_nnz - array containing the number of nonzeros in the various block rows of the
2247:            off-diagonal portion of the local submatrix that is right of the diagonal
2248:            (possibly different for each block row) or NULL.


2251:    Options Database Keys:
2252: +   -mat_no_unroll - uses code that does not unroll the loops in the
2253:                      block calculations (much slower)
2254: -   -mat_block_size - size of the blocks to use

2256:    Notes:

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

2261:    If the *_nnz parameter is given then the *_nz parameter is ignored

2263:    Storage Information:
2264:    For a square global matrix we define each processor's diagonal portion
2265:    to be its local rows and the corresponding columns (a square submatrix);
2266:    each processor's off-diagonal portion encompasses the remainder of the
2267:    local matrix (a rectangular submatrix).

2269:    The user can specify preallocated storage for the diagonal part of
2270:    the local submatrix with either d_nz or d_nnz (not both).  Set
2271:    d_nz=PETSC_DEFAULT and d_nnz=NULL for PETSc to control dynamic
2272:    memory allocation.  Likewise, specify preallocated storage for the
2273:    off-diagonal part of the local submatrix with o_nz or o_nnz (not both).

2275:    You can call MatGetInfo() to get information on how effective the preallocation was;
2276:    for example the fields mallocs,nz_allocated,nz_used,nz_unneeded;
2277:    You can also run with the option -info and look for messages with the string
2278:    malloc in them to see if additional memory allocation was needed.

2280:    Consider a processor that owns rows 3, 4 and 5 of a parallel matrix. In
2281:    the figure below we depict these three local rows and all columns (0-11).

2283: .vb
2284:            0 1 2 3 4 5 6 7 8 9 10 11
2285:           --------------------------
2286:    row 3  |. . . d d d o o o o  o  o
2287:    row 4  |. . . d d d o o o o  o  o
2288:    row 5  |. . . d d d o o o o  o  o
2289:           --------------------------
2290: .ve

2292:    Thus, any entries in the d locations are stored in the d (diagonal)
2293:    submatrix, and any entries in the o locations are stored in the
2294:    o (off-diagonal) submatrix.  Note that the d matrix is stored in
2295:    MatSeqSBAIJ format and the o submatrix in MATSEQBAIJ format.

2297:    Now d_nz should indicate the number of block nonzeros per row in the upper triangular
2298:    plus the diagonal part of the d matrix,
2299:    and o_nz should indicate the number of block nonzeros per row in the o matrix

2301:    In general, for PDE problems in which most nonzeros are near the diagonal,
2302:    one expects d_nz >> o_nz.   For large problems you MUST preallocate memory
2303:    or you will get TERRIBLE performance; see the users' manual chapter on
2304:    matrices.

2306:    Level: intermediate

2308: .seealso: MatCreate(), MatCreateSeqSBAIJ(), MatSetValues(), MatCreateBAIJ(), PetscSplitOwnership()
2309: @*/
2310: PetscErrorCode  MatMPISBAIJSetPreallocation(Mat B,PetscInt bs,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
2311: {

2318:   PetscTryMethod(B,"MatMPISBAIJSetPreallocation_C",(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,bs,d_nz,d_nnz,o_nz,o_nnz));
2319:   return(0);
2320: }

2322: /*@C
2323:    MatCreateSBAIJ - Creates a sparse parallel matrix in symmetric block AIJ format
2324:    (block compressed row).  For good matrix assembly performance
2325:    the user should preallocate the matrix storage by setting the parameters
2326:    d_nz (or d_nnz) and o_nz (or o_nnz).  By setting these parameters accurately,
2327:    performance can be increased by more than a factor of 50.

2329:    Collective

2331:    Input Parameters:
2332: +  comm - MPI communicator
2333: .  bs   - size of block, the blocks are ALWAYS square. One can use MatSetBlockSizes() to set a different row and column blocksize but the row
2334:           blocksize always defines the size of the blocks. The column blocksize sets the blocksize of the vectors obtained with MatCreateVecs()
2335: .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
2336:            This value should be the same as the local size used in creating the
2337:            y vector for the matrix-vector product y = Ax.
2338: .  n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
2339:            This value should be the same as the local size used in creating the
2340:            x vector for the matrix-vector product y = Ax.
2341: .  M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
2342: .  N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
2343: .  d_nz  - number of block nonzeros per block row in diagonal portion of local
2344:            submatrix  (same for all local rows)
2345: .  d_nnz - array containing the number of block nonzeros in the various block rows
2346:            in the upper triangular portion of the in diagonal portion of the local
2347:            (possibly different for each block block row) or NULL.
2348:            If you plan to factor the matrix you must leave room for the diagonal entry and
2349:            set its value even if it is zero.
2350: .  o_nz  - number of block nonzeros per block row in the off-diagonal portion of local
2351:            submatrix (same for all local rows).
2352: -  o_nnz - array containing the number of nonzeros in the various block rows of the
2353:            off-diagonal portion of the local submatrix (possibly different for
2354:            each block row) or NULL.

2356:    Output Parameter:
2357: .  A - the matrix

2359:    Options Database Keys:
2360: +   -mat_no_unroll - uses code that does not unroll the loops in the
2361:                      block calculations (much slower)
2362: .   -mat_block_size - size of the blocks to use
2363: -   -mat_mpi - use the parallel matrix data structures even on one processor
2364:                (defaults to using SeqBAIJ format on one processor)

2366:    It is recommended that one use the MatCreate(), MatSetType() and/or MatSetFromOptions(),
2367:    MatXXXXSetPreallocation() paradigm instead of this routine directly.
2368:    [MatXXXXSetPreallocation() is, for example, MatSeqAIJSetPreallocation]

2370:    Notes:
2371:    The number of rows and columns must be divisible by blocksize.
2372:    This matrix type does not support complex Hermitian operation.

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

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

2380:    If the *_nnz parameter is given then the *_nz parameter is ignored

2382:    Storage Information:
2383:    For a square global matrix we define each processor's diagonal portion
2384:    to be its local rows and the corresponding columns (a square submatrix);
2385:    each processor's off-diagonal portion encompasses the remainder of the
2386:    local matrix (a rectangular submatrix).

2388:    The user can specify preallocated storage for the diagonal part of
2389:    the local submatrix with either d_nz or d_nnz (not both).  Set
2390:    d_nz=PETSC_DEFAULT and d_nnz=NULL for PETSc to control dynamic
2391:    memory allocation.  Likewise, specify preallocated storage for the
2392:    off-diagonal part of the local submatrix with o_nz or o_nnz (not both).

2394:    Consider a processor that owns rows 3, 4 and 5 of a parallel matrix. In
2395:    the figure below we depict these three local rows and all columns (0-11).

2397: .vb
2398:            0 1 2 3 4 5 6 7 8 9 10 11
2399:           --------------------------
2400:    row 3  |. . . d d d o o o o  o  o
2401:    row 4  |. . . d d d o o o o  o  o
2402:    row 5  |. . . d d d o o o o  o  o
2403:           --------------------------
2404: .ve

2406:    Thus, any entries in the d locations are stored in the d (diagonal)
2407:    submatrix, and any entries in the o locations are stored in the
2408:    o (off-diagonal) submatrix.  Note that the d matrix is stored in
2409:    MatSeqSBAIJ format and the o submatrix in MATSEQBAIJ format.

2411:    Now d_nz should indicate the number of block nonzeros per row in the upper triangular
2412:    plus the diagonal part of the d matrix,
2413:    and o_nz should indicate the number of block nonzeros per row in the o matrix.
2414:    In general, for PDE problems in which most nonzeros are near the diagonal,
2415:    one expects d_nz >> o_nz.   For large problems you MUST preallocate memory
2416:    or you will get TERRIBLE performance; see the users' manual chapter on
2417:    matrices.

2419:    Level: intermediate

2421: .seealso: MatCreate(), MatCreateSeqSBAIJ(), MatSetValues(), MatCreateBAIJ()
2422: @*/

2424: PetscErrorCode  MatCreateSBAIJ(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[],Mat *A)
2425: {
2427:   PetscMPIInt    size;

2430:   MatCreate(comm,A);
2431:   MatSetSizes(*A,m,n,M,N);
2432:   MPI_Comm_size(comm,&size);
2433:   if (size > 1) {
2434:     MatSetType(*A,MATMPISBAIJ);
2435:     MatMPISBAIJSetPreallocation(*A,bs,d_nz,d_nnz,o_nz,o_nnz);
2436:   } else {
2437:     MatSetType(*A,MATSEQSBAIJ);
2438:     MatSeqSBAIJSetPreallocation(*A,bs,d_nz,d_nnz);
2439:   }
2440:   return(0);
2441: }


2444: static PetscErrorCode MatDuplicate_MPISBAIJ(Mat matin,MatDuplicateOption cpvalues,Mat *newmat)
2445: {
2446:   Mat            mat;
2447:   Mat_MPISBAIJ   *a,*oldmat = (Mat_MPISBAIJ*)matin->data;
2449:   PetscInt       len=0,nt,bs=matin->rmap->bs,mbs=oldmat->mbs;
2450:   PetscScalar    *array;

2453:   *newmat = NULL;

2455:   MatCreate(PetscObjectComm((PetscObject)matin),&mat);
2456:   MatSetSizes(mat,matin->rmap->n,matin->cmap->n,matin->rmap->N,matin->cmap->N);
2457:   MatSetType(mat,((PetscObject)matin)->type_name);
2458:   PetscLayoutReference(matin->rmap,&mat->rmap);
2459:   PetscLayoutReference(matin->cmap,&mat->cmap);

2461:   mat->factortype   = matin->factortype;
2462:   mat->preallocated = PETSC_TRUE;
2463:   mat->assembled    = PETSC_TRUE;
2464:   mat->insertmode   = NOT_SET_VALUES;

2466:   a      = (Mat_MPISBAIJ*)mat->data;
2467:   a->bs2 = oldmat->bs2;
2468:   a->mbs = oldmat->mbs;
2469:   a->nbs = oldmat->nbs;
2470:   a->Mbs = oldmat->Mbs;
2471:   a->Nbs = oldmat->Nbs;

2473:   a->size         = oldmat->size;
2474:   a->rank         = oldmat->rank;
2475:   a->donotstash   = oldmat->donotstash;
2476:   a->roworiented  = oldmat->roworiented;
2477:   a->rowindices   = NULL;
2478:   a->rowvalues    = NULL;
2479:   a->getrowactive = PETSC_FALSE;
2480:   a->barray       = NULL;
2481:   a->rstartbs     = oldmat->rstartbs;
2482:   a->rendbs       = oldmat->rendbs;
2483:   a->cstartbs     = oldmat->cstartbs;
2484:   a->cendbs       = oldmat->cendbs;

2486:   /* hash table stuff */
2487:   a->ht           = NULL;
2488:   a->hd           = NULL;
2489:   a->ht_size      = 0;
2490:   a->ht_flag      = oldmat->ht_flag;
2491:   a->ht_fact      = oldmat->ht_fact;
2492:   a->ht_total_ct  = 0;
2493:   a->ht_insert_ct = 0;

2495:   PetscArraycpy(a->rangebs,oldmat->rangebs,a->size+2);
2496:   if (oldmat->colmap) {
2497: #if defined(PETSC_USE_CTABLE)
2498:     PetscTableCreateCopy(oldmat->colmap,&a->colmap);
2499: #else
2500:     PetscMalloc1(a->Nbs,&a->colmap);
2501:     PetscLogObjectMemory((PetscObject)mat,(a->Nbs)*sizeof(PetscInt));
2502:     PetscArraycpy(a->colmap,oldmat->colmap,a->Nbs);
2503: #endif
2504:   } else a->colmap = NULL;

2506:   if (oldmat->garray && (len = ((Mat_SeqBAIJ*)(oldmat->B->data))->nbs)) {
2507:     PetscMalloc1(len,&a->garray);
2508:     PetscLogObjectMemory((PetscObject)mat,len*sizeof(PetscInt));
2509:     PetscArraycpy(a->garray,oldmat->garray,len);
2510:   } else a->garray = NULL;

2512:   MatStashCreate_Private(PetscObjectComm((PetscObject)matin),matin->rmap->bs,&mat->bstash);
2513:   VecDuplicate(oldmat->lvec,&a->lvec);
2514:   PetscLogObjectParent((PetscObject)mat,(PetscObject)a->lvec);
2515:   VecScatterCopy(oldmat->Mvctx,&a->Mvctx);
2516:   PetscLogObjectParent((PetscObject)mat,(PetscObject)a->Mvctx);

2518:   VecDuplicate(oldmat->slvec0,&a->slvec0);
2519:   PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec0);
2520:   VecDuplicate(oldmat->slvec1,&a->slvec1);
2521:   PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec1);

2523:   VecGetLocalSize(a->slvec1,&nt);
2524:   VecGetArray(a->slvec1,&array);
2525:   VecCreateSeqWithArray(PETSC_COMM_SELF,1,bs*mbs,array,&a->slvec1a);
2526:   VecCreateSeqWithArray(PETSC_COMM_SELF,1,nt-bs*mbs,array+bs*mbs,&a->slvec1b);
2527:   VecRestoreArray(a->slvec1,&array);
2528:   VecGetArray(a->slvec0,&array);
2529:   VecCreateSeqWithArray(PETSC_COMM_SELF,1,nt-bs*mbs,array+bs*mbs,&a->slvec0b);
2530:   VecRestoreArray(a->slvec0,&array);
2531:   PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec0);
2532:   PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec1);
2533:   PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec0b);
2534:   PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec1a);
2535:   PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec1b);

2537:   /*  VecScatterCopy(oldmat->sMvctx,&a->sMvctx); - not written yet, replaced by the lazy trick: */
2538:   PetscObjectReference((PetscObject)oldmat->sMvctx);
2539:   a->sMvctx = oldmat->sMvctx;
2540:   PetscLogObjectParent((PetscObject)mat,(PetscObject)a->sMvctx);

2542:   MatDuplicate(oldmat->A,cpvalues,&a->A);
2543:   PetscLogObjectParent((PetscObject)mat,(PetscObject)a->A);
2544:   MatDuplicate(oldmat->B,cpvalues,&a->B);
2545:   PetscLogObjectParent((PetscObject)mat,(PetscObject)a->B);
2546:   PetscFunctionListDuplicate(((PetscObject)matin)->qlist,&((PetscObject)mat)->qlist);
2547:   *newmat = mat;
2548:   return(0);
2549: }

2551: /* Used for both MPIBAIJ and MPISBAIJ matrices */
2552: #define MatLoad_MPISBAIJ_Binary MatLoad_MPIBAIJ_Binary

2554: PetscErrorCode MatLoad_MPISBAIJ(Mat mat,PetscViewer viewer)
2555: {
2557:   PetscBool      isbinary;

2560:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
2561:   if (!isbinary) SETERRQ2(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Viewer type %s not yet supported for reading %s matrices",((PetscObject)viewer)->type_name,((PetscObject)mat)->type_name);
2562:   MatLoad_MPISBAIJ_Binary(mat,viewer);
2563:   return(0);
2564: }

2566: /*XXXXX@
2567:    MatMPISBAIJSetHashTableFactor - Sets the factor required to compute the size of the HashTable.

2569:    Input Parameters:
2570: .  mat  - the matrix
2571: .  fact - factor

2573:    Not Collective on Mat, each process can have a different hash factor

2575:    Level: advanced

2577:   Notes:
2578:    This can also be set by the command line option: -mat_use_hash_table fact

2580: .seealso: MatSetOption()
2581: @XXXXX*/


2584: PetscErrorCode MatGetRowMaxAbs_MPISBAIJ(Mat A,Vec v,PetscInt idx[])
2585: {
2586:   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;
2587:   Mat_SeqBAIJ    *b = (Mat_SeqBAIJ*)(a->B)->data;
2588:   PetscReal      atmp;
2589:   PetscReal      *work,*svalues,*rvalues;
2591:   PetscInt       i,bs,mbs,*bi,*bj,brow,j,ncols,krow,kcol,col,row,Mbs,bcol;
2592:   PetscMPIInt    rank,size;
2593:   PetscInt       *rowners_bs,dest,count,source;
2594:   PetscScalar    *va;
2595:   MatScalar      *ba;
2596:   MPI_Status     stat;

2599:   if (idx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Send email to petsc-maint@mcs.anl.gov");
2600:   MatGetRowMaxAbs(a->A,v,NULL);
2601:   VecGetArray(v,&va);

2603:   MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);
2604:   MPI_Comm_rank(PetscObjectComm((PetscObject)A),&rank);

2606:   bs  = A->rmap->bs;
2607:   mbs = a->mbs;
2608:   Mbs = a->Mbs;
2609:   ba  = b->a;
2610:   bi  = b->i;
2611:   bj  = b->j;

2613:   /* find ownerships */
2614:   rowners_bs = A->rmap->range;

2616:   /* each proc creates an array to be distributed */
2617:   PetscCalloc1(bs*Mbs,&work);

2619:   /* row_max for B */
2620:   if (rank != size-1) {
2621:     for (i=0; i<mbs; i++) {
2622:       ncols = bi[1] - bi[0]; bi++;
2623:       brow  = bs*i;
2624:       for (j=0; j<ncols; j++) {
2625:         bcol = bs*(*bj);
2626:         for (kcol=0; kcol<bs; kcol++) {
2627:           col  = bcol + kcol;                /* local col index */
2628:           col += rowners_bs[rank+1];      /* global col index */
2629:           for (krow=0; krow<bs; krow++) {
2630:             atmp = PetscAbsScalar(*ba); ba++;
2631:             row  = brow + krow;   /* local row index */
2632:             if (PetscRealPart(va[row]) < atmp) va[row] = atmp;
2633:             if (work[col] < atmp) work[col] = atmp;
2634:           }
2635:         }
2636:         bj++;
2637:       }
2638:     }

2640:     /* send values to its owners */
2641:     for (dest=rank+1; dest<size; dest++) {
2642:       svalues = work + rowners_bs[dest];
2643:       count   = rowners_bs[dest+1]-rowners_bs[dest];
2644:       MPI_Send(svalues,count,MPIU_REAL,dest,rank,PetscObjectComm((PetscObject)A));
2645:     }
2646:   }

2648:   /* receive values */
2649:   if (rank) {
2650:     rvalues = work;
2651:     count   = rowners_bs[rank+1]-rowners_bs[rank];
2652:     for (source=0; source<rank; source++) {
2653:       MPI_Recv(rvalues,count,MPIU_REAL,MPI_ANY_SOURCE,MPI_ANY_TAG,PetscObjectComm((PetscObject)A),&stat);
2654:       /* process values */
2655:       for (i=0; i<count; i++) {
2656:         if (PetscRealPart(va[i]) < rvalues[i]) va[i] = rvalues[i];
2657:       }
2658:     }
2659:   }

2661:   VecRestoreArray(v,&va);
2662:   PetscFree(work);
2663:   return(0);
2664: }

2666: PetscErrorCode MatSOR_MPISBAIJ(Mat matin,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
2667: {
2668:   Mat_MPISBAIJ      *mat = (Mat_MPISBAIJ*)matin->data;
2669:   PetscErrorCode    ierr;
2670:   PetscInt          mbs=mat->mbs,bs=matin->rmap->bs;
2671:   PetscScalar       *x,*ptr,*from;
2672:   Vec               bb1;
2673:   const PetscScalar *b;

2676:   if (its <= 0 || lits <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
2677:   if (bs > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented");

2679:   if (flag == SOR_APPLY_UPPER) {
2680:     (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx);
2681:     return(0);
2682:   }

2684:   if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP) {
2685:     if (flag & SOR_ZERO_INITIAL_GUESS) {
2686:       (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,lits,xx);
2687:       its--;
2688:     }

2690:     VecDuplicate(bb,&bb1);
2691:     while (its--) {

2693:       /* lower triangular part: slvec0b = - B^T*xx */
2694:       (*mat->B->ops->multtranspose)(mat->B,xx,mat->slvec0b);

2696:       /* copy xx into slvec0a */
2697:       VecGetArray(mat->slvec0,&ptr);
2698:       VecGetArray(xx,&x);
2699:       PetscArraycpy(ptr,x,bs*mbs);
2700:       VecRestoreArray(mat->slvec0,&ptr);

2702:       VecScale(mat->slvec0,-1.0);

2704:       /* copy bb into slvec1a */
2705:       VecGetArray(mat->slvec1,&ptr);
2706:       VecGetArrayRead(bb,&b);
2707:       PetscArraycpy(ptr,b,bs*mbs);
2708:       VecRestoreArray(mat->slvec1,&ptr);

2710:       /* set slvec1b = 0 */
2711:       VecSet(mat->slvec1b,0.0);

2713:       VecScatterBegin(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD);
2714:       VecRestoreArray(xx,&x);
2715:       VecRestoreArrayRead(bb,&b);
2716:       VecScatterEnd(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD);

2718:       /* upper triangular part: bb1 = bb1 - B*x */
2719:       (*mat->B->ops->multadd)(mat->B,mat->slvec1b,mat->slvec1a,bb1);

2721:       /* local diagonal sweep */
2722:       (*mat->A->ops->sor)(mat->A,bb1,omega,SOR_SYMMETRIC_SWEEP,fshift,lits,lits,xx);
2723:     }
2724:     VecDestroy(&bb1);
2725:   } else if ((flag & SOR_LOCAL_FORWARD_SWEEP) && (its == 1) && (flag & SOR_ZERO_INITIAL_GUESS)) {
2726:     (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx);
2727:   } else if ((flag & SOR_LOCAL_BACKWARD_SWEEP) && (its == 1) && (flag & SOR_ZERO_INITIAL_GUESS)) {
2728:     (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx);
2729:   } else if (flag & SOR_EISENSTAT) {
2730:     Vec               xx1;
2731:     PetscBool         hasop;
2732:     const PetscScalar *diag;
2733:     PetscScalar       *sl,scale = (omega - 2.0)/omega;
2734:     PetscInt          i,n;

2736:     if (!mat->xx1) {
2737:       VecDuplicate(bb,&mat->xx1);
2738:       VecDuplicate(bb,&mat->bb1);
2739:     }
2740:     xx1 = mat->xx1;
2741:     bb1 = mat->bb1;

2743:     (*mat->A->ops->sor)(mat->A,bb,omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_BACKWARD_SWEEP),fshift,lits,1,xx);

2745:     if (!mat->diag) {
2746:       /* this is wrong for same matrix with new nonzero values */
2747:       MatCreateVecs(matin,&mat->diag,NULL);
2748:       MatGetDiagonal(matin,mat->diag);
2749:     }
2750:     MatHasOperation(matin,MATOP_MULT_DIAGONAL_BLOCK,&hasop);

2752:     if (hasop) {
2753:       MatMultDiagonalBlock(matin,xx,bb1);
2754:       VecAYPX(mat->slvec1a,scale,bb);
2755:     } else {
2756:       /*
2757:           These two lines are replaced by code that may be a bit faster for a good compiler
2758:       VecPointwiseMult(mat->slvec1a,mat->diag,xx);
2759:       VecAYPX(mat->slvec1a,scale,bb);
2760:       */
2761:       VecGetArray(mat->slvec1a,&sl);
2762:       VecGetArrayRead(mat->diag,&diag);
2763:       VecGetArrayRead(bb,&b);
2764:       VecGetArray(xx,&x);
2765:       VecGetLocalSize(xx,&n);
2766:       if (omega == 1.0) {
2767:         for (i=0; i<n; i++) sl[i] = b[i] - diag[i]*x[i];
2768:         PetscLogFlops(2.0*n);
2769:       } else {
2770:         for (i=0; i<n; i++) sl[i] = b[i] + scale*diag[i]*x[i];
2771:         PetscLogFlops(3.0*n);
2772:       }
2773:       VecRestoreArray(mat->slvec1a,&sl);
2774:       VecRestoreArrayRead(mat->diag,&diag);
2775:       VecRestoreArrayRead(bb,&b);
2776:       VecRestoreArray(xx,&x);
2777:     }

2779:     /* multiply off-diagonal portion of matrix */
2780:     VecSet(mat->slvec1b,0.0);
2781:     (*mat->B->ops->multtranspose)(mat->B,xx,mat->slvec0b);
2782:     VecGetArray(mat->slvec0,&from);
2783:     VecGetArray(xx,&x);
2784:     PetscArraycpy(from,x,bs*mbs);
2785:     VecRestoreArray(mat->slvec0,&from);
2786:     VecRestoreArray(xx,&x);
2787:     VecScatterBegin(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD);
2788:     VecScatterEnd(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD);
2789:     (*mat->B->ops->multadd)(mat->B,mat->slvec1b,mat->slvec1a,mat->slvec1a);

2791:     /* local sweep */
2792:     (*mat->A->ops->sor)(mat->A,mat->slvec1a,omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_FORWARD_SWEEP),fshift,lits,1,xx1);
2793:     VecAXPY(xx,1.0,xx1);
2794:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatSORType is not supported for SBAIJ matrix format");
2795:   return(0);
2796: }

2798: /*@
2799:      MatCreateMPISBAIJWithArrays - creates a MPI SBAIJ matrix using arrays that contain in standard
2800:          CSR format the local rows.

2802:    Collective

2804:    Input Parameters:
2805: +  comm - MPI communicator
2806: .  bs - the block size, only a block size of 1 is supported
2807: .  m - number of local rows (Cannot be PETSC_DECIDE)
2808: .  n - This value should be the same as the local size used in creating the
2809:        x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have
2810:        calculated if N is given) For square matrices n is almost always m.
2811: .  M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
2812: .  N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
2813: .   i - row indices; that is i[0] = 0, i[row] = i[row-1] + number of block elements in that row block row of the matrix
2814: .   j - column indices
2815: -   a - matrix values

2817:    Output Parameter:
2818: .   mat - the matrix

2820:    Level: intermediate

2822:    Notes:
2823:        The i, j, and a arrays ARE copied by this routine into the internal format used by PETSc;
2824:      thus you CANNOT change the matrix entries by changing the values of a[] after you have
2825:      called this routine. Use MatCreateMPIAIJWithSplitArrays() to avoid needing to copy the arrays.

2827:        The i and j indices are 0 based, and i indices are indices corresponding to the local j array.

2829: .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatMPIAIJSetPreallocation(), MatMPIAIJSetPreallocationCSR(),
2830:           MPIAIJ, MatCreateAIJ(), MatCreateMPIAIJWithSplitArrays()
2831: @*/
2832: PetscErrorCode  MatCreateMPISBAIJWithArrays(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt i[],const PetscInt j[],const PetscScalar a[],Mat *mat)
2833: {


2838:   if (i[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"i (row indices) must start with 0");
2839:   if (m < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"local number of rows (m) cannot be PETSC_DECIDE, or negative");
2840:   MatCreate(comm,mat);
2841:   MatSetSizes(*mat,m,n,M,N);
2842:   MatSetType(*mat,MATMPISBAIJ);
2843:   MatMPISBAIJSetPreallocationCSR(*mat,bs,i,j,a);
2844:   return(0);
2845: }


2848: /*@C
2849:    MatMPISBAIJSetPreallocationCSR - Creates a sparse parallel matrix in SBAIJ format using the given nonzero structure and (optional) numerical values

2851:    Collective

2853:    Input Parameters:
2854: +  B - the matrix
2855: .  bs - the block size
2856: .  i - the indices into j for the start of each local row (starts with zero)
2857: .  j - the column indices for each local row (starts with zero) these must be sorted for each row
2858: -  v - optional values in the matrix

2860:    Level: advanced

2862:    Notes:
2863:    Though this routine has Preallocation() in the name it also sets the exact nonzero locations of the matrix entries
2864:    and usually the numerical values as well

2866:    Any entries below the diagonal are ignored

2868: .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatMPIBAIJSetPreallocation(), MatCreateAIJ(), MPIAIJ
2869: @*/
2870: PetscErrorCode  MatMPISBAIJSetPreallocationCSR(Mat B,PetscInt bs,const PetscInt i[],const PetscInt j[], const PetscScalar v[])
2871: {

2875:   PetscTryMethod(B,"MatMPISBAIJSetPreallocationCSR_C",(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]),(B,bs,i,j,v));
2876:   return(0);
2877: }

2879: PetscErrorCode MatCreateMPIMatConcatenateSeqMat_MPISBAIJ(MPI_Comm comm,Mat inmat,PetscInt n,MatReuse scall,Mat *outmat)
2880: {
2882:   PetscInt       m,N,i,rstart,nnz,Ii,bs,cbs;
2883:   PetscInt       *indx;
2884:   PetscScalar    *values;

2887:   MatGetSize(inmat,&m,&N);
2888:   if (scall == MAT_INITIAL_MATRIX) { /* symbolic phase */
2889:     Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)inmat->data;
2890:     PetscInt       *dnz,*onz,mbs,Nbs,nbs;
2891:     PetscInt       *bindx,rmax=a->rmax,j;
2892:     PetscMPIInt    rank,size;

2894:     MatGetBlockSizes(inmat,&bs,&cbs);
2895:     mbs = m/bs; Nbs = N/cbs;
2896:     if (n == PETSC_DECIDE) {
2897:       PetscSplitOwnershipBlock(comm,cbs,&n,&N);
2898:     }
2899:     nbs = n/cbs;

2901:     PetscMalloc1(rmax,&bindx);
2902:     MatPreallocateInitialize(comm,mbs,nbs,dnz,onz); /* inline function, output __end and __rstart are used below */

2904:     MPI_Comm_rank(comm,&rank);
2905:     MPI_Comm_rank(comm,&size);
2906:     if (rank == size-1) {
2907:       /* Check sum(nbs) = Nbs */
2908:       if (__end != Nbs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Sum of local block columns %D != global block columns %D",__end,Nbs);
2909:     }

2911:     rstart = __rstart; /* block rstart of *outmat; see inline function MatPreallocateInitialize */
2912:     MatSetOption(inmat,MAT_GETROW_UPPERTRIANGULAR,PETSC_TRUE);
2913:     for (i=0; i<mbs; i++) {
2914:       MatGetRow_SeqSBAIJ(inmat,i*bs,&nnz,&indx,NULL); /* non-blocked nnz and indx */
2915:       nnz  = nnz/bs;
2916:       for (j=0; j<nnz; j++) bindx[j] = indx[j*bs]/bs;
2917:       MatPreallocateSet(i+rstart,nnz,bindx,dnz,onz);
2918:       MatRestoreRow_SeqSBAIJ(inmat,i*bs,&nnz,&indx,NULL);
2919:     }
2920:     MatSetOption(inmat,MAT_GETROW_UPPERTRIANGULAR,PETSC_FALSE);
2921:     PetscFree(bindx);

2923:     MatCreate(comm,outmat);
2924:     MatSetSizes(*outmat,m,n,PETSC_DETERMINE,PETSC_DETERMINE);
2925:     MatSetBlockSizes(*outmat,bs,cbs);
2926:     MatSetType(*outmat,MATSBAIJ);
2927:     MatSeqSBAIJSetPreallocation(*outmat,bs,0,dnz);
2928:     MatMPISBAIJSetPreallocation(*outmat,bs,0,dnz,0,onz);
2929:     MatPreallocateFinalize(dnz,onz);
2930:   }

2932:   /* numeric phase */
2933:   MatGetBlockSizes(inmat,&bs,&cbs);
2934:   MatGetOwnershipRange(*outmat,&rstart,NULL);

2936:   MatSetOption(inmat,MAT_GETROW_UPPERTRIANGULAR,PETSC_TRUE);
2937:   for (i=0; i<m; i++) {
2938:     MatGetRow_SeqSBAIJ(inmat,i,&nnz,&indx,&values);
2939:     Ii   = i + rstart;
2940:     MatSetValues(*outmat,1,&Ii,nnz,indx,values,INSERT_VALUES);
2941:     MatRestoreRow_SeqSBAIJ(inmat,i,&nnz,&indx,&values);
2942:   }
2943:   MatSetOption(inmat,MAT_GETROW_UPPERTRIANGULAR,PETSC_FALSE);
2944:   MatAssemblyBegin(*outmat,MAT_FINAL_ASSEMBLY);
2945:   MatAssemblyEnd(*outmat,MAT_FINAL_ASSEMBLY);
2946:   return(0);
2947: }