Actual source code: mpisbaij.c

petsc-3.4.5 2014-06-29
  2: #include <../src/mat/impls/baij/mpi/mpibaij.h>    /*I "petscmat.h" I*/
  3: #include <../src/mat/impls/sbaij/mpi/mpisbaij.h>
  4: #include <../src/mat/impls/sbaij/seq/sbaij.h>
  5: #include <petscblaslapack.h>

  7: extern PetscErrorCode MatSetUpMultiply_MPISBAIJ(Mat);
  8: extern PetscErrorCode MatSetUpMultiply_MPISBAIJ_2comm(Mat);
  9: extern PetscErrorCode MatDisAssemble_MPISBAIJ(Mat);
 10: extern PetscErrorCode MatIncreaseOverlap_MPISBAIJ(Mat,PetscInt,IS[],PetscInt);
 11: extern PetscErrorCode MatGetValues_SeqSBAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar []);
 12: extern PetscErrorCode MatGetValues_SeqBAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar []);
 13: extern PetscErrorCode MatSetValues_SeqSBAIJ(Mat,PetscInt,const PetscInt [],PetscInt,const PetscInt [],const PetscScalar [],InsertMode);
 14: extern PetscErrorCode MatSetValuesBlocked_SeqSBAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
 15: extern PetscErrorCode MatSetValuesBlocked_SeqBAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
 16: extern PetscErrorCode MatGetRow_SeqSBAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
 17: extern PetscErrorCode MatRestoreRow_SeqSBAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
 18: extern PetscErrorCode MatZeroRows_SeqSBAIJ(Mat,IS,PetscScalar*,Vec,Vec);
 19: extern PetscErrorCode MatZeroRows_SeqBAIJ(Mat,IS,PetscScalar*,Vec,Vec);
 20: extern PetscErrorCode MatGetRowMaxAbs_MPISBAIJ(Mat,Vec,PetscInt[]);
 21: extern PetscErrorCode MatSOR_MPISBAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);

 25: PetscErrorCode  MatStoreValues_MPISBAIJ(Mat mat)
 26: {
 27:   Mat_MPISBAIJ   *aij = (Mat_MPISBAIJ*)mat->data;

 31:   MatStoreValues(aij->A);
 32:   MatStoreValues(aij->B);
 33:   return(0);
 34: }

 38: PetscErrorCode  MatRetrieveValues_MPISBAIJ(Mat mat)
 39: {
 40:   Mat_MPISBAIJ   *aij = (Mat_MPISBAIJ*)mat->data;

 44:   MatRetrieveValues(aij->A);
 45:   MatRetrieveValues(aij->B);
 46:   return(0);
 47: }

 49: #define  MatSetValues_SeqSBAIJ_A_Private(row,col,value,addv) \
 50:   { \
 51:  \
 52:     brow = row/bs;  \
 53:     rp   = aj + ai[brow]; ap = aa + bs2*ai[brow]; \
 54:     rmax = aimax[brow]; nrow = ailen[brow]; \
 55:     bcol = col/bs; \
 56:     ridx = row % bs; cidx = col % bs; \
 57:     low  = 0; high = nrow; \
 58:     while (high-low > 3) { \
 59:       t = (low+high)/2; \
 60:       if (rp[t] > bcol) high = t; \
 61:       else              low  = t; \
 62:     } \
 63:     for (_i=low; _i<high; _i++) { \
 64:       if (rp[_i] > bcol) break; \
 65:       if (rp[_i] == bcol) { \
 66:         bap = ap + bs2*_i + bs*cidx + ridx; \
 67:         if (addv == ADD_VALUES) *bap += value;  \
 68:         else                    *bap  = value;  \
 69:         goto a_noinsert; \
 70:       } \
 71:     } \
 72:     if (a->nonew == 1) goto a_noinsert; \
 73:     if (a->nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) into matrix", row, col); \
 74:     MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,aimax,a->nonew,MatScalar); \
 75:     N = nrow++ - 1;  \
 76:     /* shift up all the later entries in this row */ \
 77:     for (ii=N; ii>=_i; ii--) { \
 78:       rp[ii+1] = rp[ii]; \
 79:       PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar)); \
 80:     } \
 81:     if (N>=_i) { PetscMemzero(ap+bs2*_i,bs2*sizeof(MatScalar)); }  \
 82:     rp[_i]                      = bcol;  \
 83:     ap[bs2*_i + bs*cidx + ridx] = value;  \
 84: a_noinsert:; \
 85:     ailen[brow] = nrow; \
 86:   }

 88: #define  MatSetValues_SeqSBAIJ_B_Private(row,col,value,addv) \
 89:   { \
 90:     brow = row/bs;  \
 91:     rp   = bj + bi[brow]; ap = ba + bs2*bi[brow]; \
 92:     rmax = bimax[brow]; nrow = bilen[brow]; \
 93:     bcol = col/bs; \
 94:     ridx = row % bs; cidx = col % bs; \
 95:     low  = 0; high = nrow; \
 96:     while (high-low > 3) { \
 97:       t = (low+high)/2; \
 98:       if (rp[t] > bcol) high = t; \
 99:       else              low  = t; \
100:     } \
101:     for (_i=low; _i<high; _i++) { \
102:       if (rp[_i] > bcol) break; \
103:       if (rp[_i] == bcol) { \
104:         bap = ap + bs2*_i + bs*cidx + ridx; \
105:         if (addv == ADD_VALUES) *bap += value;  \
106:         else                    *bap  = value;  \
107:         goto b_noinsert; \
108:       } \
109:     } \
110:     if (b->nonew == 1) goto b_noinsert; \
111:     if (b->nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) into matrix", row, col); \
112:     MatSeqXAIJReallocateAIJ(B,b->mbs,bs2,nrow,brow,bcol,rmax,ba,bi,bj,rp,ap,bimax,b->nonew,MatScalar); \
113:     N = nrow++ - 1;  \
114:     /* shift up all the later entries in this row */ \
115:     for (ii=N; ii>=_i; ii--) { \
116:       rp[ii+1] = rp[ii]; \
117:       PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar)); \
118:     } \
119:     if (N>=_i) { PetscMemzero(ap+bs2*_i,bs2*sizeof(MatScalar));}  \
120:     rp[_i]                      = bcol;  \
121:     ap[bs2*_i + bs*cidx + ridx] = value;  \
122: b_noinsert:; \
123:     bilen[brow] = nrow; \
124:   }

126: /* Only add/insert a(i,j) with i<=j (blocks).
127:    Any a(i,j) with i>j input by user is ingored.
128: */
131: PetscErrorCode MatSetValues_MPISBAIJ(Mat mat,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode addv)
132: {
133:   Mat_MPISBAIJ   *baij = (Mat_MPISBAIJ*)mat->data;
134:   MatScalar      value;
135:   PetscBool      roworiented = baij->roworiented;
137:   PetscInt       i,j,row,col;
138:   PetscInt       rstart_orig=mat->rmap->rstart;
139:   PetscInt       rend_orig  =mat->rmap->rend,cstart_orig=mat->cmap->rstart;
140:   PetscInt       cend_orig  =mat->cmap->rend,bs=mat->rmap->bs;

142:   /* Some Variables required in the macro */
143:   Mat          A     = baij->A;
144:   Mat_SeqSBAIJ *a    = (Mat_SeqSBAIJ*)(A)->data;
145:   PetscInt     *aimax=a->imax,*ai=a->i,*ailen=a->ilen,*aj=a->j;
146:   MatScalar    *aa   =a->a;

148:   Mat         B     = baij->B;
149:   Mat_SeqBAIJ *b    = (Mat_SeqBAIJ*)(B)->data;
150:   PetscInt    *bimax=b->imax,*bi=b->i,*bilen=b->ilen,*bj=b->j;
151:   MatScalar   *ba   =b->a;

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

157:   /* for stash */
158:   PetscInt  n_loc, *in_loc = NULL;
159:   MatScalar *v_loc = NULL;

163:   if (!baij->donotstash) {
164:     if (n > baij->n_loc) {
165:       PetscFree(baij->in_loc);
166:       PetscFree(baij->v_loc);
167:       PetscMalloc(n*sizeof(PetscInt),&baij->in_loc);
168:       PetscMalloc(n*sizeof(MatScalar),&baij->v_loc);

170:       baij->n_loc = n;
171:     }
172:     in_loc = baij->in_loc;
173:     v_loc  = baij->v_loc;
174:   }

176:   for (i=0; i<m; i++) {
177:     if (im[i] < 0) continue;
178: #if defined(PETSC_USE_DEBUG)
179:     if (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);
180: #endif
181:     if (im[i] >= rstart_orig && im[i] < rend_orig) { /* this processor entry */
182:       row = im[i] - rstart_orig;              /* local row index */
183:       for (j=0; j<n; j++) {
184:         if (im[i]/bs > in[j]/bs) {
185:           if (a->ignore_ltriangular) {
186:             continue;    /* ignore lower triangular blocks */
187:           } 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)");
188:         }
189:         if (in[j] >= cstart_orig && in[j] < cend_orig) {  /* diag entry (A) */
190:           col  = in[j] - cstart_orig;         /* local col index */
191:           brow = row/bs; bcol = col/bs;
192:           if (brow > bcol) continue;  /* ignore lower triangular blocks of A */
193:           if (roworiented) value = v[i*n+j];
194:           else             value = v[i+j*m];
195:           MatSetValues_SeqSBAIJ_A_Private(row,col,value,addv);
196:           /* MatSetValues_SeqBAIJ(baij->A,1,&row,1,&col,&value,addv); */
197:         } else if (in[j] < 0) continue;
198: #if defined(PETSC_USE_DEBUG)
199:         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);
200: #endif
201:         else {  /* off-diag entry (B) */
202:           if (mat->was_assembled) {
203:             if (!baij->colmap) {
204:               MatCreateColmap_MPIBAIJ_Private(mat);
205:             }
206: #if defined(PETSC_USE_CTABLE)
207:             PetscTableFind(baij->colmap,in[j]/bs + 1,&col);
208:             col  = col - 1;
209: #else
210:             col = baij->colmap[in[j]/bs] - 1;
211: #endif
212:             if (col < 0 && !((Mat_SeqSBAIJ*)(baij->A->data))->nonew) {
213:               MatDisAssemble_MPISBAIJ(mat);
214:               col  =  in[j];
215:               /* Reinitialize the variables required by MatSetValues_SeqBAIJ_B_Private() */
216:               B    = baij->B;
217:               b    = (Mat_SeqBAIJ*)(B)->data;
218:               bimax= b->imax;bi=b->i;bilen=b->ilen;bj=b->j;
219:               ba   = b->a;
220:             } else col += in[j]%bs;
221:           } else col = in[j];
222:           if (roworiented) value = v[i*n+j];
223:           else             value = v[i+j*m];
224:           MatSetValues_SeqSBAIJ_B_Private(row,col,value,addv);
225:           /* MatSetValues_SeqBAIJ(baij->B,1,&row,1,&col,&value,addv); */
226:         }
227:       }
228:     } else {  /* off processor entry */
229:       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]);
230:       if (!baij->donotstash) {
231:         mat->assembled = PETSC_FALSE;
232:         n_loc          = 0;
233:         for (j=0; j<n; j++) {
234:           if (im[i]/bs > in[j]/bs) continue; /* ignore lower triangular blocks */
235:           in_loc[n_loc] = in[j];
236:           if (roworiented) {
237:             v_loc[n_loc] = v[i*n+j];
238:           } else {
239:             v_loc[n_loc] = v[j*m+i];
240:           }
241:           n_loc++;
242:         }
243:         MatStashValuesRow_Private(&mat->stash,im[i],n_loc,in_loc,v_loc,PETSC_FALSE);
244:       }
245:     }
246:   }
247:   return(0);
248: }

252: PetscErrorCode MatSetValuesBlocked_MPISBAIJ(Mat mat,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const MatScalar v[],InsertMode addv)
253: {
254:   Mat_MPISBAIJ    *baij = (Mat_MPISBAIJ*)mat->data;
255:   const MatScalar *value;
256:   MatScalar       *barray     =baij->barray;
257:   PetscBool       roworiented = baij->roworiented,ignore_ltriangular = ((Mat_SeqSBAIJ*)baij->A->data)->ignore_ltriangular;
258:   PetscErrorCode  ierr;
259:   PetscInt        i,j,ii,jj,row,col,rstart=baij->rstartbs;
260:   PetscInt        rend=baij->rendbs,cstart=baij->rstartbs,stepval;
261:   PetscInt        cend=baij->rendbs,bs=mat->rmap->bs,bs2=baij->bs2;

264:   if (!barray) {
265:     PetscMalloc(bs2*sizeof(MatScalar),&barray);
266:     baij->barray = barray;
267:   }

269:   if (roworiented) {
270:     stepval = (n-1)*bs;
271:   } else {
272:     stepval = (m-1)*bs;
273:   }
274:   for (i=0; i<m; i++) {
275:     if (im[i] < 0) continue;
276: #if defined(PETSC_USE_DEBUG)
277:     if (im[i] >= baij->Mbs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large, row %D max %D",im[i],baij->Mbs-1);
278: #endif
279:     if (im[i] >= rstart && im[i] < rend) {
280:       row = im[i] - rstart;
281:       for (j=0; j<n; j++) {
282:         if (im[i] > in[j]) {
283:           if (ignore_ltriangular) continue; /* ignore lower triangular blocks */
284:           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)");
285:         }
286:         /* If NumCol = 1 then a copy is not required */
287:         if ((roworiented) && (n == 1)) {
288:           barray = (MatScalar*) v + i*bs2;
289:         } else if ((!roworiented) && (m == 1)) {
290:           barray = (MatScalar*) v + j*bs2;
291:         } else { /* Here a copy is required */
292:           if (roworiented) {
293:             value = v + i*(stepval+bs)*bs + j*bs;
294:           } else {
295:             value = v + j*(stepval+bs)*bs + i*bs;
296:           }
297:           for (ii=0; ii<bs; ii++,value+=stepval) {
298:             for (jj=0; jj<bs; jj++) {
299:               *barray++ = *value++;
300:             }
301:           }
302:           barray -=bs2;
303:         }

305:         if (in[j] >= cstart && in[j] < cend) {
306:           col  = in[j] - cstart;
307:           MatSetValuesBlocked_SeqSBAIJ(baij->A,1,&row,1,&col,barray,addv);
308:         } else if (in[j] < 0) continue;
309: #if defined(PETSC_USE_DEBUG)
310:         else if (in[j] >= baij->Nbs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large, col %D max %D",in[j],baij->Nbs-1);
311: #endif
312:         else {
313:           if (mat->was_assembled) {
314:             if (!baij->colmap) {
315:               MatCreateColmap_MPIBAIJ_Private(mat);
316:             }

318: #if defined(PETSC_USE_DEBUG)
319: #if defined(PETSC_USE_CTABLE)
320:             { PetscInt data;
321:               PetscTableFind(baij->colmap,in[j]+1,&data);
322:               if ((data - 1) % bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect colmap");
323:             }
324: #else
325:             if ((baij->colmap[in[j]] - 1) % bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect colmap");
326: #endif
327: #endif
328: #if defined(PETSC_USE_CTABLE)
329:             PetscTableFind(baij->colmap,in[j]+1,&col);
330:             col  = (col - 1)/bs;
331: #else
332:             col = (baij->colmap[in[j]] - 1)/bs;
333: #endif
334:             if (col < 0 && !((Mat_SeqBAIJ*)(baij->A->data))->nonew) {
335:               MatDisAssemble_MPISBAIJ(mat);
336:               col  = in[j];
337:             }
338:           } else col = in[j];
339:           MatSetValuesBlocked_SeqBAIJ(baij->B,1,&row,1,&col,barray,addv);
340:         }
341:       }
342:     } else {
343:       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]);
344:       if (!baij->donotstash) {
345:         if (roworiented) {
346:           MatStashValuesRowBlocked_Private(&mat->bstash,im[i],n,in,v,m,n,i);
347:         } else {
348:           MatStashValuesColBlocked_Private(&mat->bstash,im[i],n,in,v,m,n,i);
349:         }
350:       }
351:     }
352:   }
353:   return(0);
354: }

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

366:   for (i=0; i<m; i++) {
367:     if (idxm[i] < 0) continue; /* SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",idxm[i]); */
368:     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);
369:     if (idxm[i] >= bsrstart && idxm[i] < bsrend) {
370:       row = idxm[i] - bsrstart;
371:       for (j=0; j<n; j++) {
372:         if (idxn[j] < 0) continue; /* SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column %D",idxn[j]); */
373:         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);
374:         if (idxn[j] >= bscstart && idxn[j] < bscend) {
375:           col  = idxn[j] - bscstart;
376:           MatGetValues_SeqSBAIJ(baij->A,1,&row,1,&col,v+i*n+j);
377:         } else {
378:           if (!baij->colmap) {
379:             MatCreateColmap_MPIBAIJ_Private(mat);
380:           }
381: #if defined(PETSC_USE_CTABLE)
382:           PetscTableFind(baij->colmap,idxn[j]/bs+1,&data);
383:           data--;
384: #else
385:           data = baij->colmap[idxn[j]/bs]-1;
386: #endif
387:           if ((data < 0) || (baij->garray[data/bs] != idxn[j]/bs)) *(v+i*n+j) = 0.0;
388:           else {
389:             col  = data + idxn[j]%bs;
390:             MatGetValues_SeqBAIJ(baij->B,1,&row,1,&col,v+i*n+j);
391:           }
392:         }
393:       }
394:     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only local values currently supported");
395:   }
396:   return(0);
397: }

401: PetscErrorCode MatNorm_MPISBAIJ(Mat mat,NormType type,PetscReal *norm)
402: {
403:   Mat_MPISBAIJ   *baij = (Mat_MPISBAIJ*)mat->data;
405:   PetscReal      sum[2],*lnorm2;

408:   if (baij->size == 1) {
409:      MatNorm(baij->A,type,norm);
410:   } else {
411:     if (type == NORM_FROBENIUS) {
412:       PetscMalloc(2*sizeof(PetscReal),&lnorm2);
413:        MatNorm(baij->A,type,lnorm2);
414:       *lnorm2 = (*lnorm2)*(*lnorm2); lnorm2++;            /* squar power of norm(A) */
415:        MatNorm(baij->B,type,lnorm2);
416:       *lnorm2 = (*lnorm2)*(*lnorm2); lnorm2--;             /* squar power of norm(B) */
417:       MPI_Allreduce(lnorm2,sum,2,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)mat));
418:       *norm   = PetscSqrtReal(sum[0] + 2*sum[1]);
419:       PetscFree(lnorm2);
420:     } else if (type == NORM_INFINITY || type == NORM_1) { /* max row/column sum */
421:       Mat_SeqSBAIJ *amat=(Mat_SeqSBAIJ*)baij->A->data;
422:       Mat_SeqBAIJ  *bmat=(Mat_SeqBAIJ*)baij->B->data;
423:       PetscReal    *rsum,*rsum2,vabs;
424:       PetscInt     *jj,*garray=baij->garray,rstart=baij->rstartbs,nz;
425:       PetscInt     brow,bcol,col,bs=baij->A->rmap->bs,row,grow,gcol,mbs=amat->mbs;
426:       MatScalar    *v;

428:       PetscMalloc2(mat->cmap->N,PetscReal,&rsum,mat->cmap->N,PetscReal,&rsum2);
429:       PetscMemzero(rsum,mat->cmap->N*sizeof(PetscReal));
430:       /* Amat */
431:       v = amat->a; jj = amat->j;
432:       for (brow=0; brow<mbs; brow++) {
433:         grow = bs*(rstart + brow);
434:         nz   = amat->i[brow+1] - amat->i[brow];
435:         for (bcol=0; bcol<nz; bcol++) {
436:           gcol = bs*(rstart + *jj); jj++;
437:           for (col=0; col<bs; col++) {
438:             for (row=0; row<bs; row++) {
439:               vabs            = PetscAbsScalar(*v); v++;
440:               rsum[gcol+col] += vabs;
441:               /* non-diagonal block */
442:               if (bcol > 0 && vabs > 0.0) rsum[grow+row] += vabs;
443:             }
444:           }
445:         }
446:       }
447:       /* Bmat */
448:       v = bmat->a; jj = bmat->j;
449:       for (brow=0; brow<mbs; brow++) {
450:         grow = bs*(rstart + brow);
451:         nz = bmat->i[brow+1] - bmat->i[brow];
452:         for (bcol=0; bcol<nz; bcol++) {
453:           gcol = bs*garray[*jj]; jj++;
454:           for (col=0; col<bs; col++) {
455:             for (row=0; row<bs; row++) {
456:               vabs            = PetscAbsScalar(*v); v++;
457:               rsum[gcol+col] += vabs;
458:               rsum[grow+row] += vabs;
459:             }
460:           }
461:         }
462:       }
463:       MPI_Allreduce(rsum,rsum2,mat->cmap->N,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)mat));
464:       *norm = 0.0;
465:       for (col=0; col<mat->cmap->N; col++) {
466:         if (rsum2[col] > *norm) *norm = rsum2[col];
467:       }
468:       PetscFree2(rsum,rsum2);
469:     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for this norm yet");
470:   }
471:   return(0);
472: }

476: PetscErrorCode MatAssemblyBegin_MPISBAIJ(Mat mat,MatAssemblyType mode)
477: {
478:   Mat_MPISBAIJ   *baij = (Mat_MPISBAIJ*)mat->data;
480:   PetscInt       nstash,reallocs;
481:   InsertMode     addv;

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

486:   /* make sure all processors are either in INSERTMODE or ADDMODE */
487:   MPI_Allreduce((PetscEnum*)&mat->insertmode,(PetscEnum*)&addv,1,MPIU_ENUM,MPI_BOR,PetscObjectComm((PetscObject)mat));
488:   if (addv == (ADD_VALUES|INSERT_VALUES)) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Some processors inserted others added");
489:   mat->insertmode = addv; /* in case this processor had no cache */

491:   MatStashScatterBegin_Private(mat,&mat->stash,mat->rmap->range);
492:   MatStashScatterBegin_Private(mat,&mat->bstash,baij->rangebs);
493:   MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);
494:   PetscInfo2(mat,"Stash has %D entries,uses %D mallocs.\n",nstash,reallocs);
495:   MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);
496:   PetscInfo2(mat,"Block-Stash has %D entries, uses %D mallocs.\n",nstash,reallocs);
497:   return(0);
498: }

502: PetscErrorCode MatAssemblyEnd_MPISBAIJ(Mat mat,MatAssemblyType mode)
503: {
504:   Mat_MPISBAIJ   *baij=(Mat_MPISBAIJ*)mat->data;
505:   Mat_SeqSBAIJ   *a   =(Mat_SeqSBAIJ*)baij->A->data;
507:   PetscInt       i,j,rstart,ncols,flg,bs2=baij->bs2;
508:   PetscInt       *row,*col;
509:   PetscBool      other_disassembled;
510:   PetscMPIInt    n;
511:   PetscBool      r1,r2,r3;
512:   MatScalar      *val;
513:   InsertMode     addv = mat->insertmode;

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

522:       for (i=0; i<n;) {
523:         /* Now identify the consecutive vals belonging to the same row */
524:         for (j=i,rstart=row[j]; j<n; j++) {
525:           if (row[j] != rstart) break;
526:         }
527:         if (j < n) ncols = j-i;
528:         else       ncols = n-i;
529:         /* Now assemble all these values with a single function call */
530:         MatSetValues_MPISBAIJ(mat,1,row+i,ncols,col+i,val+i,addv);
531:         i    = j;
532:       }
533:     }
534:     MatStashScatterEnd_Private(&mat->stash);
535:     /* Now process the block-stash. Since the values are stashed column-oriented,
536:        set the roworiented flag to column oriented, and after MatSetValues()
537:        restore the original flags */
538:     r1 = baij->roworiented;
539:     r2 = a->roworiented;
540:     r3 = ((Mat_SeqBAIJ*)baij->B->data)->roworiented;

542:     baij->roworiented = PETSC_FALSE;
543:     a->roworiented    = PETSC_FALSE;

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

550:       for (i=0; i<n;) {
551:         /* Now identify the consecutive vals belonging to the same row */
552:         for (j=i,rstart=row[j]; j<n; j++) {
553:           if (row[j] != rstart) break;
554:         }
555:         if (j < n) ncols = j-i;
556:         else       ncols = n-i;
557:         MatSetValuesBlocked_MPISBAIJ(mat,1,row+i,ncols,col+i,val+i*bs2,addv);
558:         i    = j;
559:       }
560:     }
561:     MatStashScatterEnd_Private(&mat->bstash);

563:     baij->roworiented = r1;
564:     a->roworiented    = r2;

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

569:   MatAssemblyBegin(baij->A,mode);
570:   MatAssemblyEnd(baij->A,mode);

572:   /* determine if any processor has disassembled, if so we must
573:      also disassemble ourselfs, in order that we may reassemble. */
574:   /*
575:      if nonzero structure of submatrix B cannot change then we know that
576:      no processor disassembled thus we can skip this stuff
577:   */
578:   if (!((Mat_SeqBAIJ*)baij->B->data)->nonew) {
579:     MPI_Allreduce(&mat->was_assembled,&other_disassembled,1,MPIU_BOOL,MPI_PROD,PetscObjectComm((PetscObject)mat));
580:     if (mat->was_assembled && !other_disassembled) {
581:       MatDisAssemble_MPISBAIJ(mat);
582:     }
583:   }

585:   if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
586:     MatSetUpMultiply_MPISBAIJ(mat); /* setup Mvctx and sMvctx */
587:   }
588:   MatSetOption(baij->B,MAT_CHECK_COMPRESSED_ROW,PETSC_TRUE);
589:   MatAssemblyBegin(baij->B,mode);
590:   MatAssemblyEnd(baij->B,mode);

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

594:   baij->rowvalues = 0;
595:   return(0);
596: }

598: extern PetscErrorCode MatSetValues_MPIBAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
599: #include <petscdraw.h>
602: static PetscErrorCode MatView_MPISBAIJ_ASCIIorDraworSocket(Mat mat,PetscViewer viewer)
603: {
604:   Mat_MPISBAIJ      *baij = (Mat_MPISBAIJ*)mat->data;
605:   PetscErrorCode    ierr;
606:   PetscInt          bs   = mat->rmap->bs;
607:   PetscMPIInt       size = baij->size,rank = baij->rank;
608:   PetscBool         iascii,isdraw;
609:   PetscViewer       sviewer;
610:   PetscViewerFormat format;

613:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
614:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
615:   if (iascii) {
616:     PetscViewerGetFormat(viewer,&format);
617:     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
618:       MatInfo info;
619:       MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&rank);
620:       MatGetInfo(mat,MAT_LOCAL,&info);
621:       PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);
622:       PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Local rows %D nz %D nz alloced %D bs %D mem %D\n",rank,mat->rmap->n,(PetscInt)info.nz_used,(PetscInt)info.nz_allocated,mat->rmap->bs,(PetscInt)info.memory);
623:       MatGetInfo(baij->A,MAT_LOCAL,&info);
624:       PetscViewerASCIISynchronizedPrintf(viewer,"[%d] on-diagonal part: nz %D \n",rank,(PetscInt)info.nz_used);
625:       MatGetInfo(baij->B,MAT_LOCAL,&info);
626:       PetscViewerASCIISynchronizedPrintf(viewer,"[%d] off-diagonal part: nz %D \n",rank,(PetscInt)info.nz_used);
627:       PetscViewerFlush(viewer);
628:       PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);
629:       PetscViewerASCIIPrintf(viewer,"Information on VecScatter used in matrix-vector product: \n");
630:       VecScatterView(baij->Mvctx,viewer);
631:       return(0);
632:     } else if (format == PETSC_VIEWER_ASCII_INFO) {
633:       PetscViewerASCIIPrintf(viewer,"  block size is %D\n",bs);
634:       return(0);
635:     } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
636:       return(0);
637:     }
638:   }

640:   if (isdraw) {
641:     PetscDraw draw;
642:     PetscBool isnull;
643:     PetscViewerDrawGetDraw(viewer,0,&draw);
644:     PetscDrawIsNull(draw,&isnull); if (isnull) return(0);
645:   }

647:   if (size == 1) {
648:     PetscObjectSetName((PetscObject)baij->A,((PetscObject)mat)->name);
649:     MatView(baij->A,viewer);
650:   } else {
651:     /* assemble the entire matrix onto first processor. */
652:     Mat          A;
653:     Mat_SeqSBAIJ *Aloc;
654:     Mat_SeqBAIJ  *Bloc;
655:     PetscInt     M = mat->rmap->N,N = mat->cmap->N,*ai,*aj,col,i,j,k,*rvals,mbs = baij->mbs;
656:     MatScalar    *a;

658:     /* Should this be the same type as mat? */
659:     MatCreate(PetscObjectComm((PetscObject)mat),&A);
660:     if (!rank) {
661:       MatSetSizes(A,M,N,M,N);
662:     } else {
663:       MatSetSizes(A,0,0,M,N);
664:     }
665:     MatSetType(A,MATMPISBAIJ);
666:     MatMPISBAIJSetPreallocation(A,mat->rmap->bs,0,NULL,0,NULL);
667:     MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);
668:     PetscLogObjectParent(mat,A);

670:     /* copy over the A part */
671:     Aloc = (Mat_SeqSBAIJ*)baij->A->data;
672:     ai   = Aloc->i; aj = Aloc->j; a = Aloc->a;
673:     PetscMalloc(bs*sizeof(PetscInt),&rvals);

675:     for (i=0; i<mbs; i++) {
676:       rvals[0] = bs*(baij->rstartbs + i);
677:       for (j=1; j<bs; j++) rvals[j] = rvals[j-1] + 1;
678:       for (j=ai[i]; j<ai[i+1]; j++) {
679:         col = (baij->cstartbs+aj[j])*bs;
680:         for (k=0; k<bs; k++) {
681:           MatSetValues_MPISBAIJ(A,bs,rvals,1,&col,a,INSERT_VALUES);
682:           col++;
683:           a += bs;
684:         }
685:       }
686:     }
687:     /* copy over the B part */
688:     Bloc = (Mat_SeqBAIJ*)baij->B->data;
689:     ai   = Bloc->i; aj = Bloc->j; a = Bloc->a;
690:     for (i=0; i<mbs; i++) {

692:       rvals[0] = bs*(baij->rstartbs + i);
693:       for (j=1; j<bs; j++) rvals[j] = rvals[j-1] + 1;
694:       for (j=ai[i]; j<ai[i+1]; j++) {
695:         col = baij->garray[aj[j]]*bs;
696:         for (k=0; k<bs; k++) {
697:           MatSetValues_MPIBAIJ(A,bs,rvals,1,&col,a,INSERT_VALUES);
698:           col++;
699:           a += bs;
700:         }
701:       }
702:     }
703:     PetscFree(rvals);
704:     MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
705:     MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
706:     /*
707:        Everyone has to call to draw the matrix since the graphics waits are
708:        synchronized across all processors that share the PetscDraw object
709:     */
710:     PetscViewerGetSingleton(viewer,&sviewer);
711:     if (!rank) {
712:       PetscObjectSetName((PetscObject)((Mat_MPISBAIJ*)(A->data))->A,((PetscObject)mat)->name);
713:       /* Set the type name to MATMPISBAIJ so that the correct type can be printed out by PetscObjectPrintClassNamePrefixType() in MatView_SeqSBAIJ_ASCII()*/
714:       PetscStrcpy(((PetscObject)((Mat_MPISBAIJ*)(A->data))->A)->type_name,MATMPISBAIJ);
715:       MatView(((Mat_MPISBAIJ*)(A->data))->A,sviewer);
716:     }
717:     PetscViewerRestoreSingleton(viewer,&sviewer);
718:     MatDestroy(&A);
719:   }
720:   return(0);
721: }

725: PetscErrorCode MatView_MPISBAIJ(Mat mat,PetscViewer viewer)
726: {
728:   PetscBool      iascii,isdraw,issocket,isbinary;

731:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
732:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
733:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSOCKET,&issocket);
734:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
735:   if (iascii || isdraw || issocket || isbinary) {
736:     MatView_MPISBAIJ_ASCIIorDraworSocket(mat,viewer);
737:   }
738:   return(0);
739: }

743: PetscErrorCode MatDestroy_MPISBAIJ(Mat mat)
744: {
745:   Mat_MPISBAIJ   *baij = (Mat_MPISBAIJ*)mat->data;

749: #if defined(PETSC_USE_LOG)
750:   PetscLogObjectState((PetscObject)mat,"Rows=%D,Cols=%D",mat->rmap->N,mat->cmap->N);
751: #endif
752:   MatStashDestroy_Private(&mat->stash);
753:   MatStashDestroy_Private(&mat->bstash);
754:   MatDestroy(&baij->A);
755:   MatDestroy(&baij->B);
756: #if defined(PETSC_USE_CTABLE)
757:   PetscTableDestroy(&baij->colmap);
758: #else
759:   PetscFree(baij->colmap);
760: #endif
761:   PetscFree(baij->garray);
762:   VecDestroy(&baij->lvec);
763:   VecScatterDestroy(&baij->Mvctx);
764:   VecDestroy(&baij->slvec0);
765:   VecDestroy(&baij->slvec0b);
766:   VecDestroy(&baij->slvec1);
767:   VecDestroy(&baij->slvec1a);
768:   VecDestroy(&baij->slvec1b);
769:   VecScatterDestroy(&baij->sMvctx);
770:   PetscFree2(baij->rowvalues,baij->rowindices);
771:   PetscFree(baij->barray);
772:   PetscFree(baij->hd);
773:   VecDestroy(&baij->diag);
774:   VecDestroy(&baij->bb1);
775:   VecDestroy(&baij->xx1);
776: #if defined(PETSC_USE_REAL_MAT_SINGLE)
777:   PetscFree(baij->setvaluescopy);
778: #endif
779:   PetscFree(baij->in_loc);
780:   PetscFree(baij->v_loc);
781:   PetscFree(baij->rangebs);
782:   PetscFree(mat->data);

784:   PetscObjectChangeTypeName((PetscObject)mat,0);
785:   PetscObjectComposeFunction((PetscObject)mat,"MatStoreValues_C",NULL);
786:   PetscObjectComposeFunction((PetscObject)mat,"MatRetrieveValues_C",NULL);
787:   PetscObjectComposeFunction((PetscObject)mat,"MatGetDiagonalBlock_C",NULL);
788:   PetscObjectComposeFunction((PetscObject)mat,"MatMPISBAIJSetPreallocation_C",NULL);
789:   PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpisbaij_mpisbstrm_C",NULL);
790:   return(0);
791: }

795: PetscErrorCode MatMult_MPISBAIJ_Hermitian(Mat A,Vec xx,Vec yy)
796: {
797:   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;
799:   PetscInt       nt,mbs=a->mbs,bs=A->rmap->bs;
800:   PetscScalar    *x,*from;

803:   VecGetLocalSize(xx,&nt);
804:   if (nt != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Incompatible partition of A and xx");

806:   /* diagonal part */
807:   (*a->A->ops->mult)(a->A,xx,a->slvec1a);
808:   VecSet(a->slvec1b,0.0);

810:   /* subdiagonal part */
811:   (*a->B->ops->multhermitiantranspose)(a->B,xx,a->slvec0b);

813:   /* copy x into the vec slvec0 */
814:   VecGetArray(a->slvec0,&from);
815:   VecGetArray(xx,&x);

817:   PetscMemcpy(from,x,bs*mbs*sizeof(MatScalar));
818:   VecRestoreArray(a->slvec0,&from);
819:   VecRestoreArray(xx,&x);

821:   VecScatterBegin(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
822:   VecScatterEnd(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
823:   /* supperdiagonal part */
824:   (*a->B->ops->multadd)(a->B,a->slvec1b,a->slvec1a,yy);
825:   return(0);
826: }

830: PetscErrorCode MatMult_MPISBAIJ(Mat A,Vec xx,Vec yy)
831: {
832:   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;
834:   PetscInt       nt,mbs=a->mbs,bs=A->rmap->bs;
835:   PetscScalar    *x,*from;

838:   VecGetLocalSize(xx,&nt);
839:   if (nt != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Incompatible partition of A and xx");

841:   /* diagonal part */
842:   (*a->A->ops->mult)(a->A,xx,a->slvec1a);
843:   VecSet(a->slvec1b,0.0);

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

848:   /* copy x into the vec slvec0 */
849:   VecGetArray(a->slvec0,&from);
850:   VecGetArray(xx,&x);

852:   PetscMemcpy(from,x,bs*mbs*sizeof(MatScalar));
853:   VecRestoreArray(a->slvec0,&from);
854:   VecRestoreArray(xx,&x);

856:   VecScatterBegin(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
857:   VecScatterEnd(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
858:   /* supperdiagonal part */
859:   (*a->B->ops->multadd)(a->B,a->slvec1b,a->slvec1a,yy);
860:   return(0);
861: }

865: PetscErrorCode MatMult_MPISBAIJ_2comm(Mat A,Vec xx,Vec yy)
866: {
867:   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;
869:   PetscInt       nt;

872:   VecGetLocalSize(xx,&nt);
873:   if (nt != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Incompatible partition of A and xx");

875:   VecGetLocalSize(yy,&nt);
876:   if (nt != A->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Incompatible parition of A and yy");

878:   VecScatterBegin(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);
879:   /* do diagonal part */
880:   (*a->A->ops->mult)(a->A,xx,yy);
881:   /* do supperdiagonal part */
882:   VecScatterEnd(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);
883:   (*a->B->ops->multadd)(a->B,a->lvec,yy,yy);
884:   /* do subdiagonal part */
885:   (*a->B->ops->multtranspose)(a->B,xx,a->lvec);
886:   VecScatterBegin(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);
887:   VecScatterEnd(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);
888:   return(0);
889: }

893: PetscErrorCode MatMultAdd_MPISBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
894: {
895:   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;
897:   PetscInt       mbs=a->mbs,bs=A->rmap->bs;
898:   PetscScalar    *x,*from,zero=0.0;

901:   /*
902:   PetscSynchronizedPrintf(PetscObjectComm((PetscObject)A)," MatMultAdd is called ...\n");
903:   PetscSynchronizedFlush(PetscObjectComm((PetscObject)A));
904:   */
905:   /* diagonal part */
906:   (*a->A->ops->multadd)(a->A,xx,yy,a->slvec1a);
907:   VecSet(a->slvec1b,zero);

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

912:   /* copy x into the vec slvec0 */
913:   VecGetArray(a->slvec0,&from);
914:   VecGetArray(xx,&x);
915:   PetscMemcpy(from,x,bs*mbs*sizeof(MatScalar));
916:   VecRestoreArray(a->slvec0,&from);

918:   VecScatterBegin(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
919:   VecRestoreArray(xx,&x);
920:   VecScatterEnd(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);

922:   /* supperdiagonal part */
923:   (*a->B->ops->multadd)(a->B,a->slvec1b,a->slvec1a,zz);
924:   return(0);
925: }

929: PetscErrorCode MatMultAdd_MPISBAIJ_2comm(Mat A,Vec xx,Vec yy,Vec zz)
930: {
931:   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;

935:   VecScatterBegin(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);
936:   /* do diagonal part */
937:   (*a->A->ops->multadd)(a->A,xx,yy,zz);
938:   /* do supperdiagonal part */
939:   VecScatterEnd(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);
940:   (*a->B->ops->multadd)(a->B,a->lvec,zz,zz);

942:   /* do subdiagonal part */
943:   (*a->B->ops->multtranspose)(a->B,xx,a->lvec);
944:   VecScatterBegin(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);
945:   VecScatterEnd(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);
946:   return(0);
947: }

949: /*
950:   This only works correctly for square matrices where the subblock A->A is the
951:    diagonal block
952: */
955: PetscErrorCode MatGetDiagonal_MPISBAIJ(Mat A,Vec v)
956: {
957:   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;

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

968: PetscErrorCode MatScale_MPISBAIJ(Mat A,PetscScalar aa)
969: {
970:   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;

974:   MatScale(a->A,aa);
975:   MatScale(a->B,aa);
976:   return(0);
977: }

981: PetscErrorCode MatGetRow_MPISBAIJ(Mat matin,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
982: {
983:   Mat_MPISBAIJ   *mat = (Mat_MPISBAIJ*)matin->data;
984:   PetscScalar    *vworkA,*vworkB,**pvA,**pvB,*v_p;
986:   PetscInt       bs = matin->rmap->bs,bs2 = mat->bs2,i,*cworkA,*cworkB,**pcA,**pcB;
987:   PetscInt       nztot,nzA,nzB,lrow,brstart = matin->rmap->rstart,brend = matin->rmap->rend;
988:   PetscInt       *cmap,*idx_p,cstart = mat->rstartbs;

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

994:   if (!mat->rowvalues && (idx || v)) {
995:     /*
996:         allocate enough space to hold information from the longest row.
997:     */
998:     Mat_SeqSBAIJ *Aa = (Mat_SeqSBAIJ*)mat->A->data;
999:     Mat_SeqBAIJ  *Ba = (Mat_SeqBAIJ*)mat->B->data;
1000:     PetscInt     max = 1,mbs = mat->mbs,tmp;
1001:     for (i=0; i<mbs; i++) {
1002:       tmp = Aa->i[i+1] - Aa->i[i] + Ba->i[i+1] - Ba->i[i]; /* row length */
1003:       if (max < tmp) max = tmp;
1004:     }
1005:     PetscMalloc2(max*bs2,PetscScalar,&mat->rowvalues,max*bs2,PetscInt,&mat->rowindices);
1006:   }

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

1011:   pvA = &vworkA; pcA = &cworkA; pvB = &vworkB; pcB = &cworkB;
1012:   if (!v)   {pvA = 0; pvB = 0;}
1013:   if (!idx) {pcA = 0; if (!v) pcB = 0;}
1014:   (*mat->A->ops->getrow)(mat->A,lrow,&nzA,pcA,pvA);
1015:   (*mat->B->ops->getrow)(mat->B,lrow,&nzB,pcB,pvB);
1016:   nztot = nzA + nzB;

1018:   cmap = mat->garray;
1019:   if (v  || idx) {
1020:     if (nztot) {
1021:       /* Sort by increasing column numbers, assuming A and B already sorted */
1022:       PetscInt imark = -1;
1023:       if (v) {
1024:         *v = v_p = mat->rowvalues;
1025:         for (i=0; i<nzB; i++) {
1026:           if (cmap[cworkB[i]/bs] < cstart) v_p[i] = vworkB[i];
1027:           else break;
1028:         }
1029:         imark = i;
1030:         for (i=0; i<nzA; i++)     v_p[imark+i] = vworkA[i];
1031:         for (i=imark; i<nzB; i++) v_p[nzA+i]   = vworkB[i];
1032:       }
1033:       if (idx) {
1034:         *idx = idx_p = mat->rowindices;
1035:         if (imark > -1) {
1036:           for (i=0; i<imark; i++) {
1037:             idx_p[i] = cmap[cworkB[i]/bs]*bs + cworkB[i]%bs;
1038:           }
1039:         } else {
1040:           for (i=0; i<nzB; i++) {
1041:             if (cmap[cworkB[i]/bs] < cstart) idx_p[i] = cmap[cworkB[i]/bs]*bs + cworkB[i]%bs;
1042:             else break;
1043:           }
1044:           imark = i;
1045:         }
1046:         for (i=0; i<nzA; i++)     idx_p[imark+i] = cstart*bs + cworkA[i];
1047:         for (i=imark; i<nzB; i++) idx_p[nzA+i]   = cmap[cworkB[i]/bs]*bs + cworkB[i]%bs ;
1048:       }
1049:     } else {
1050:       if (idx) *idx = 0;
1051:       if (v)   *v   = 0;
1052:     }
1053:   }
1054:   *nz  = nztot;
1055:   (*mat->A->ops->restorerow)(mat->A,lrow,&nzA,pcA,pvA);
1056:   (*mat->B->ops->restorerow)(mat->B,lrow,&nzB,pcB,pvB);
1057:   return(0);
1058: }

1062: PetscErrorCode MatRestoreRow_MPISBAIJ(Mat mat,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
1063: {
1064:   Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;

1067:   if (!baij->getrowactive) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"MatGetRow() must be called first");
1068:   baij->getrowactive = PETSC_FALSE;
1069:   return(0);
1070: }

1074: PetscErrorCode MatGetRowUpperTriangular_MPISBAIJ(Mat A)
1075: {
1076:   Mat_MPISBAIJ *a  = (Mat_MPISBAIJ*)A->data;
1077:   Mat_SeqSBAIJ *aA = (Mat_SeqSBAIJ*)a->A->data;

1080:   aA->getrow_utriangular = PETSC_TRUE;
1081:   return(0);
1082: }
1085: PetscErrorCode MatRestoreRowUpperTriangular_MPISBAIJ(Mat A)
1086: {
1087:   Mat_MPISBAIJ *a  = (Mat_MPISBAIJ*)A->data;
1088:   Mat_SeqSBAIJ *aA = (Mat_SeqSBAIJ*)a->A->data;

1091:   aA->getrow_utriangular = PETSC_FALSE;
1092:   return(0);
1093: }

1097: PetscErrorCode MatRealPart_MPISBAIJ(Mat A)
1098: {
1099:   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;

1103:   MatRealPart(a->A);
1104:   MatRealPart(a->B);
1105:   return(0);
1106: }

1110: PetscErrorCode MatImaginaryPart_MPISBAIJ(Mat A)
1111: {
1112:   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;

1116:   MatImaginaryPart(a->A);
1117:   MatImaginaryPart(a->B);
1118:   return(0);
1119: }

1123: PetscErrorCode MatZeroEntries_MPISBAIJ(Mat A)
1124: {
1125:   Mat_MPISBAIJ   *l = (Mat_MPISBAIJ*)A->data;

1129:   MatZeroEntries(l->A);
1130:   MatZeroEntries(l->B);
1131:   return(0);
1132: }

1136: PetscErrorCode MatGetInfo_MPISBAIJ(Mat matin,MatInfoType flag,MatInfo *info)
1137: {
1138:   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)matin->data;
1139:   Mat            A  = a->A,B = a->B;
1141:   PetscReal      isend[5],irecv[5];

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

1146:   MatGetInfo(A,MAT_LOCAL,info);

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

1151:   MatGetInfo(B,MAT_LOCAL,info);

1153:   isend[0] += info->nz_used; isend[1] += info->nz_allocated; isend[2] += info->nz_unneeded;
1154:   isend[3] += info->memory;  isend[4] += info->mallocs;
1155:   if (flag == MAT_LOCAL) {
1156:     info->nz_used      = isend[0];
1157:     info->nz_allocated = isend[1];
1158:     info->nz_unneeded  = isend[2];
1159:     info->memory       = isend[3];
1160:     info->mallocs      = isend[4];
1161:   } else if (flag == MAT_GLOBAL_MAX) {
1162:     MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)matin));

1164:     info->nz_used      = irecv[0];
1165:     info->nz_allocated = irecv[1];
1166:     info->nz_unneeded  = irecv[2];
1167:     info->memory       = irecv[3];
1168:     info->mallocs      = irecv[4];
1169:   } else if (flag == MAT_GLOBAL_SUM) {
1170:     MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)matin));

1172:     info->nz_used      = irecv[0];
1173:     info->nz_allocated = irecv[1];
1174:     info->nz_unneeded  = irecv[2];
1175:     info->memory       = irecv[3];
1176:     info->mallocs      = irecv[4];
1177:   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown MatInfoType argument %d",(int)flag);
1178:   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
1179:   info->fill_ratio_needed = 0;
1180:   info->factor_mallocs    = 0;
1181:   return(0);
1182: }

1186: PetscErrorCode MatSetOption_MPISBAIJ(Mat A,MatOption op,PetscBool flg)
1187: {
1188:   Mat_MPISBAIJ   *a  = (Mat_MPISBAIJ*)A->data;
1189:   Mat_SeqSBAIJ   *aA = (Mat_SeqSBAIJ*)a->A->data;

1193:   switch (op) {
1194:   case MAT_NEW_NONZERO_LOCATIONS:
1195:   case MAT_NEW_NONZERO_ALLOCATION_ERR:
1196:   case MAT_UNUSED_NONZERO_LOCATION_ERR:
1197:   case MAT_KEEP_NONZERO_PATTERN:
1198:   case MAT_NEW_NONZERO_LOCATION_ERR:
1199:     MatSetOption(a->A,op,flg);
1200:     MatSetOption(a->B,op,flg);
1201:     break;
1202:   case MAT_ROW_ORIENTED:
1203:     a->roworiented = flg;

1205:     MatSetOption(a->A,op,flg);
1206:     MatSetOption(a->B,op,flg);
1207:     break;
1208:   case MAT_NEW_DIAGONALS:
1209:     PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
1210:     break;
1211:   case MAT_IGNORE_OFF_PROC_ENTRIES:
1212:     a->donotstash = flg;
1213:     break;
1214:   case MAT_USE_HASH_TABLE:
1215:     a->ht_flag = flg;
1216:     break;
1217:   case MAT_HERMITIAN:
1218:     if (!A->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call MatAssemblyEnd() first");
1219:     MatSetOption(a->A,op,flg);

1221:     A->ops->mult = MatMult_MPISBAIJ_Hermitian;
1222:     break;
1223:   case MAT_SPD:
1224:     A->spd_set = PETSC_TRUE;
1225:     A->spd     = flg;
1226:     if (flg) {
1227:       A->symmetric                  = PETSC_TRUE;
1228:       A->structurally_symmetric     = PETSC_TRUE;
1229:       A->symmetric_set              = PETSC_TRUE;
1230:       A->structurally_symmetric_set = PETSC_TRUE;
1231:     }
1232:     break;
1233:   case MAT_SYMMETRIC:
1234:     MatSetOption(a->A,op,flg);
1235:     break;
1236:   case MAT_STRUCTURALLY_SYMMETRIC:
1237:     MatSetOption(a->A,op,flg);
1238:     break;
1239:   case MAT_SYMMETRY_ETERNAL:
1240:     if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix must be symmetric");
1241:     PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
1242:     break;
1243:   case MAT_IGNORE_LOWER_TRIANGULAR:
1244:     aA->ignore_ltriangular = flg;
1245:     break;
1246:   case MAT_ERROR_LOWER_TRIANGULAR:
1247:     aA->ignore_ltriangular = flg;
1248:     break;
1249:   case MAT_GETROW_UPPERTRIANGULAR:
1250:     aA->getrow_utriangular = flg;
1251:     break;
1252:   default:
1253:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %d",op);
1254:   }
1255:   return(0);
1256: }

1260: PetscErrorCode MatTranspose_MPISBAIJ(Mat A,MatReuse reuse,Mat *B)
1261: {

1265:   if (MAT_INITIAL_MATRIX || *B != A) {
1266:     MatDuplicate(A,MAT_COPY_VALUES,B);
1267:   }
1268:   return(0);
1269: }

1273: PetscErrorCode MatDiagonalScale_MPISBAIJ(Mat mat,Vec ll,Vec rr)
1274: {
1275:   Mat_MPISBAIJ   *baij = (Mat_MPISBAIJ*)mat->data;
1276:   Mat            a     = baij->A, b=baij->B;
1278:   PetscInt       nv,m,n;
1279:   PetscBool      flg;

1282:   if (ll != rr) {
1283:     VecEqual(ll,rr,&flg);
1284:     if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"For symmetric format, left and right scaling vectors must be same\n");
1285:   }
1286:   if (!ll) return(0);

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

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

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

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

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

1302:   /* right diagonalscale the off-diagonal part */
1303:   VecScatterEnd(baij->Mvctx,rr,baij->lvec,INSERT_VALUES,SCATTER_FORWARD);
1304:   (*b->ops->diagonalscale)(b,NULL,baij->lvec);
1305:   return(0);
1306: }

1310: PetscErrorCode MatSetUnfactored_MPISBAIJ(Mat A)
1311: {
1312:   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;

1316:   MatSetUnfactored(a->A);
1317:   return(0);
1318: }

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

1324: PetscErrorCode MatEqual_MPISBAIJ(Mat A,Mat B,PetscBool  *flag)
1325: {
1326:   Mat_MPISBAIJ   *matB = (Mat_MPISBAIJ*)B->data,*matA = (Mat_MPISBAIJ*)A->data;
1327:   Mat            a,b,c,d;
1328:   PetscBool      flg;

1332:   a = matA->A; b = matA->B;
1333:   c = matB->A; d = matB->B;

1335:   MatEqual(a,c,&flg);
1336:   if (flg) {
1337:     MatEqual(b,d,&flg);
1338:   }
1339:   MPI_Allreduce(&flg,flag,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));
1340:   return(0);
1341: }

1345: PetscErrorCode MatCopy_MPISBAIJ(Mat A,Mat B,MatStructure str)
1346: {
1348:   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;
1349:   Mat_MPISBAIJ   *b = (Mat_MPISBAIJ*)B->data;

1352:   /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
1353:   if ((str != SAME_NONZERO_PATTERN) || (A->ops->copy != B->ops->copy)) {
1354:     MatGetRowUpperTriangular(A);
1355:     MatCopy_Basic(A,B,str);
1356:     MatRestoreRowUpperTriangular(A);
1357:   } else {
1358:     MatCopy(a->A,b->A,str);
1359:     MatCopy(a->B,b->B,str);
1360:   }
1361:   return(0);
1362: }

1366: PetscErrorCode MatSetUp_MPISBAIJ(Mat A)
1367: {

1371:   MatMPISBAIJSetPreallocation(A,A->rmap->bs,PETSC_DEFAULT,0,PETSC_DEFAULT,0);
1372:   return(0);
1373: }

1377: PetscErrorCode MatAXPY_MPISBAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str)
1378: {
1380:   Mat_MPISBAIJ   *xx=(Mat_MPISBAIJ*)X->data,*yy=(Mat_MPISBAIJ*)Y->data;
1381:   PetscBLASInt   bnz,one=1;
1382:   Mat_SeqSBAIJ   *xa,*ya;
1383:   Mat_SeqBAIJ    *xb,*yb;

1386:   if (str == SAME_NONZERO_PATTERN) {
1387:     PetscScalar alpha = a;
1388:     xa   = (Mat_SeqSBAIJ*)xx->A->data;
1389:     ya   = (Mat_SeqSBAIJ*)yy->A->data;
1390:     PetscBLASIntCast(xa->nz,&bnz);
1391:     PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&bnz,&alpha,xa->a,&one,ya->a,&one));
1392:     xb   = (Mat_SeqBAIJ*)xx->B->data;
1393:     yb   = (Mat_SeqBAIJ*)yy->B->data;
1394:     PetscBLASIntCast(xb->nz,&bnz);
1395:     PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&bnz,&alpha,xb->a,&one,yb->a,&one));
1396:   } else {
1397:     MatGetRowUpperTriangular(X);
1398:     MatAXPY_Basic(Y,a,X,str);
1399:     MatRestoreRowUpperTriangular(X);
1400:   }
1401:   return(0);
1402: }

1406: PetscErrorCode MatGetSubMatrices_MPISBAIJ(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1407: {
1409:   PetscInt       i;
1410:   PetscBool      flg;

1413:   for (i=0; i<n; i++) {
1414:     ISEqual(irow[i],icol[i],&flg);
1415:     if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Can only get symmetric submatrix for MPISBAIJ matrices");
1416:   }
1417:   MatGetSubMatrices_MPIBAIJ(A,n,irow,icol,scall,B);
1418:   return(0);
1419: }


1422: /* -------------------------------------------------------------------*/
1423: static struct _MatOps MatOps_Values = {MatSetValues_MPISBAIJ,
1424:                                        MatGetRow_MPISBAIJ,
1425:                                        MatRestoreRow_MPISBAIJ,
1426:                                        MatMult_MPISBAIJ,
1427:                                /*  4*/ MatMultAdd_MPISBAIJ,
1428:                                        MatMult_MPISBAIJ,       /* transpose versions are same as non-transpose */
1429:                                        MatMultAdd_MPISBAIJ,
1430:                                        0,
1431:                                        0,
1432:                                        0,
1433:                                /* 10*/ 0,
1434:                                        0,
1435:                                        0,
1436:                                        MatSOR_MPISBAIJ,
1437:                                        MatTranspose_MPISBAIJ,
1438:                                /* 15*/ MatGetInfo_MPISBAIJ,
1439:                                        MatEqual_MPISBAIJ,
1440:                                        MatGetDiagonal_MPISBAIJ,
1441:                                        MatDiagonalScale_MPISBAIJ,
1442:                                        MatNorm_MPISBAIJ,
1443:                                /* 20*/ MatAssemblyBegin_MPISBAIJ,
1444:                                        MatAssemblyEnd_MPISBAIJ,
1445:                                        MatSetOption_MPISBAIJ,
1446:                                        MatZeroEntries_MPISBAIJ,
1447:                                /* 24*/ 0,
1448:                                        0,
1449:                                        0,
1450:                                        0,
1451:                                        0,
1452:                                /* 29*/ MatSetUp_MPISBAIJ,
1453:                                        0,
1454:                                        0,
1455:                                        0,
1456:                                        0,
1457:                                /* 34*/ MatDuplicate_MPISBAIJ,
1458:                                        0,
1459:                                        0,
1460:                                        0,
1461:                                        0,
1462:                                /* 39*/ MatAXPY_MPISBAIJ,
1463:                                        MatGetSubMatrices_MPISBAIJ,
1464:                                        MatIncreaseOverlap_MPISBAIJ,
1465:                                        MatGetValues_MPISBAIJ,
1466:                                        MatCopy_MPISBAIJ,
1467:                                /* 44*/ 0,
1468:                                        MatScale_MPISBAIJ,
1469:                                        0,
1470:                                        0,
1471:                                        0,
1472:                                /* 49*/ 0,
1473:                                        0,
1474:                                        0,
1475:                                        0,
1476:                                        0,
1477:                                /* 54*/ 0,
1478:                                        0,
1479:                                        MatSetUnfactored_MPISBAIJ,
1480:                                        0,
1481:                                        MatSetValuesBlocked_MPISBAIJ,
1482:                                /* 59*/ 0,
1483:                                        0,
1484:                                        0,
1485:                                        0,
1486:                                        0,
1487:                                /* 64*/ 0,
1488:                                        0,
1489:                                        0,
1490:                                        0,
1491:                                        0,
1492:                                /* 69*/ MatGetRowMaxAbs_MPISBAIJ,
1493:                                        0,
1494:                                        0,
1495:                                        0,
1496:                                        0,
1497:                                /* 74*/ 0,
1498:                                        0,
1499:                                        0,
1500:                                        0,
1501:                                        0,
1502:                                /* 79*/ 0,
1503:                                        0,
1504:                                        0,
1505:                                        0,
1506:                                        MatLoad_MPISBAIJ,
1507:                                /* 84*/ 0,
1508:                                        0,
1509:                                        0,
1510:                                        0,
1511:                                        0,
1512:                                /* 89*/ 0,
1513:                                        0,
1514:                                        0,
1515:                                        0,
1516:                                        0,
1517:                                /* 94*/ 0,
1518:                                        0,
1519:                                        0,
1520:                                        0,
1521:                                        0,
1522:                                /* 99*/ 0,
1523:                                        0,
1524:                                        0,
1525:                                        0,
1526:                                        0,
1527:                                /*104*/ 0,
1528:                                        MatRealPart_MPISBAIJ,
1529:                                        MatImaginaryPart_MPISBAIJ,
1530:                                        MatGetRowUpperTriangular_MPISBAIJ,
1531:                                        MatRestoreRowUpperTriangular_MPISBAIJ,
1532:                                /*109*/ 0,
1533:                                        0,
1534:                                        0,
1535:                                        0,
1536:                                        0,
1537:                                /*114*/ 0,
1538:                                        0,
1539:                                        0,
1540:                                        0,
1541:                                        0,
1542:                                /*119*/ 0,
1543:                                        0,
1544:                                        0,
1545:                                        0,
1546:                                        0,
1547:                                /*124*/ 0,
1548:                                        0,
1549:                                        0,
1550:                                        0,
1551:                                        0,
1552:                                /*129*/ 0,
1553:                                        0,
1554:                                        0,
1555:                                        0,
1556:                                        0,
1557:                                /*134*/ 0,
1558:                                        0,
1559:                                        0,
1560:                                        0,
1561:                                        0,
1562:                                /*139*/ 0,
1563:                                        0
1564: };


1569: PetscErrorCode  MatGetDiagonalBlock_MPISBAIJ(Mat A,Mat *a)
1570: {
1572:   *a = ((Mat_MPISBAIJ*)A->data)->A;
1573:   return(0);
1574: }

1578: PetscErrorCode  MatMPISBAIJSetPreallocation_MPISBAIJ(Mat B,PetscInt bs,PetscInt d_nz,const PetscInt *d_nnz,PetscInt o_nz,const PetscInt *o_nnz)
1579: {
1580:   Mat_MPISBAIJ   *b;
1582:   PetscInt       i,mbs,Mbs;

1585:   PetscLayoutSetBlockSize(B->rmap,bs);
1586:   PetscLayoutSetBlockSize(B->cmap,bs);
1587:   PetscLayoutSetUp(B->rmap);
1588:   PetscLayoutSetUp(B->cmap);
1589:   PetscLayoutGetBlockSize(B->rmap,&bs);

1591:   b   = (Mat_MPISBAIJ*)B->data;
1592:   mbs = B->rmap->n/bs;
1593:   Mbs = B->rmap->N/bs;
1594:   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);

1596:   B->rmap->bs = bs;
1597:   b->bs2      = bs*bs;
1598:   b->mbs      = mbs;
1599:   b->nbs      = mbs;
1600:   b->Mbs      = Mbs;
1601:   b->Nbs      = Mbs;

1603:   for (i=0; i<=b->size; i++) {
1604:     b->rangebs[i] = B->rmap->range[i]/bs;
1605:   }
1606:   b->rstartbs = B->rmap->rstart/bs;
1607:   b->rendbs   = B->rmap->rend/bs;

1609:   b->cstartbs = B->cmap->rstart/bs;
1610:   b->cendbs   = B->cmap->rend/bs;

1612:   if (!B->preallocated) {
1613:     MatCreate(PETSC_COMM_SELF,&b->A);
1614:     MatSetSizes(b->A,B->rmap->n,B->cmap->n,B->rmap->n,B->cmap->n);
1615:     MatSetType(b->A,MATSEQSBAIJ);
1616:     PetscLogObjectParent(B,b->A);
1617:     MatCreate(PETSC_COMM_SELF,&b->B);
1618:     MatSetSizes(b->B,B->rmap->n,B->cmap->N,B->rmap->n,B->cmap->N);
1619:     MatSetType(b->B,MATSEQBAIJ);
1620:     PetscLogObjectParent(B,b->B);
1621:     MatStashCreate_Private(PetscObjectComm((PetscObject)B),bs,&B->bstash);
1622:   }

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

1627:   B->preallocated = PETSC_TRUE;
1628:   return(0);
1629: }

1633: PetscErrorCode MatMPISBAIJSetPreallocationCSR_MPISBAIJ(Mat B,PetscInt bs,const PetscInt ii[],const PetscInt jj[],const PetscScalar V[])
1634: {
1635:   PetscInt       m,rstart,cstart,cend;
1636:   PetscInt       i,j,d,nz,nz_max=0,*d_nnz=0,*o_nnz=0;
1637:   const PetscInt *JJ    =0;
1638:   PetscScalar    *values=0;

1642:   if (bs < 1) SETERRQ1(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_OUTOFRANGE,"Invalid block size specified, must be positive but it is %D",bs);
1643:   PetscLayoutSetBlockSize(B->rmap,bs);
1644:   PetscLayoutSetBlockSize(B->cmap,bs);
1645:   PetscLayoutSetUp(B->rmap);
1646:   PetscLayoutSetUp(B->cmap);
1647:   PetscLayoutGetBlockSize(B->rmap,&bs);
1648:   m      = B->rmap->n/bs;
1649:   rstart = B->rmap->rstart/bs;
1650:   cstart = B->cmap->rstart/bs;
1651:   cend   = B->cmap->rend/bs;

1653:   if (ii[0]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"ii[0] must be 0 but it is %D",ii[0]);
1654:   PetscMalloc2(m,PetscInt,&d_nnz,m,PetscInt,&o_nnz);
1655:   for (i=0; i<m; i++) {
1656:     nz = ii[i+1] - ii[i];
1657:     if (nz < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local row %D has a negative number of columns %D",i,nz);
1658:     nz_max = PetscMax(nz_max,nz);
1659:     JJ     = jj + ii[i];
1660:     for (j=0; j<nz; j++) {
1661:       if (*JJ >= cstart) break;
1662:       JJ++;
1663:     }
1664:     d = 0;
1665:     for (; j<nz; j++) {
1666:       if (*JJ++ >= cend) break;
1667:       d++;
1668:     }
1669:     d_nnz[i] = d;
1670:     o_nnz[i] = nz - d;
1671:   }
1672:   MatMPISBAIJSetPreallocation(B,bs,0,d_nnz,0,o_nnz);
1673:   PetscFree2(d_nnz,o_nnz);

1675:   values = (PetscScalar*)V;
1676:   if (!values) {
1677:     PetscMalloc(bs*bs*nz_max*sizeof(PetscScalar),&values);
1678:     PetscMemzero(values,bs*bs*nz_max*sizeof(PetscScalar));
1679:   }
1680:   for (i=0; i<m; i++) {
1681:     PetscInt          row    = i + rstart;
1682:     PetscInt          ncols  = ii[i+1] - ii[i];
1683:     const PetscInt    *icols = jj + ii[i];
1684:     const PetscScalar *svals = values + (V ? (bs*bs*ii[i]) : 0);
1685:     MatSetValuesBlocked_MPISBAIJ(B,1,&row,ncols,icols,svals,INSERT_VALUES);
1686:   }

1688:   if (!V) { PetscFree(values); }
1689:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
1690:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
1691:   MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
1692:   return(0);
1693: }

1695: #if defined(PETSC_HAVE_MUMPS)
1696: PETSC_EXTERN PetscErrorCode MatGetFactor_sbaij_mumps(Mat,MatFactorType,Mat*);
1697: #endif
1698: #if defined(PETSC_HAVE_PASTIX)
1699: PETSC_EXTERN PetscErrorCode MatGetFactor_mpisbaij_pastix(Mat,MatFactorType,Mat*);
1700: #endif

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

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

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

1713:   Level: beginner

1715: .seealso: MatCreateMPISBAIJ
1716: M*/

1718: PETSC_EXTERN PetscErrorCode MatConvert_MPISBAIJ_MPISBSTRM(Mat,MatType,MatReuse,Mat*);

1722: PETSC_EXTERN PetscErrorCode MatCreate_MPISBAIJ(Mat B)
1723: {
1724:   Mat_MPISBAIJ   *b;
1726:   PetscBool      flg;

1729:   PetscNewLog(B,Mat_MPISBAIJ,&b);
1730:   B->data = (void*)b;
1731:   PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));

1733:   B->ops->destroy = MatDestroy_MPISBAIJ;
1734:   B->ops->view    = MatView_MPISBAIJ;
1735:   B->assembled    = PETSC_FALSE;
1736:   B->insertmode   = NOT_SET_VALUES;

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

1741:   /* build local table of row and column ownerships */
1742:   PetscMalloc((b->size+2)*sizeof(PetscInt),&b->rangebs);

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

1747:   b->donotstash  = PETSC_FALSE;
1748:   b->colmap      = NULL;
1749:   b->garray      = NULL;
1750:   b->roworiented = PETSC_TRUE;

1752:   /* stuff used in block assembly */
1753:   b->barray = 0;

1755:   /* stuff used for matrix vector multiply */
1756:   b->lvec    = 0;
1757:   b->Mvctx   = 0;
1758:   b->slvec0  = 0;
1759:   b->slvec0b = 0;
1760:   b->slvec1  = 0;
1761:   b->slvec1a = 0;
1762:   b->slvec1b = 0;
1763:   b->sMvctx  = 0;

1765:   /* stuff for MatGetRow() */
1766:   b->rowindices   = 0;
1767:   b->rowvalues    = 0;
1768:   b->getrowactive = PETSC_FALSE;

1770:   /* hash table stuff */
1771:   b->ht           = 0;
1772:   b->hd           = 0;
1773:   b->ht_size      = 0;
1774:   b->ht_flag      = PETSC_FALSE;
1775:   b->ht_fact      = 0;
1776:   b->ht_total_ct  = 0;
1777:   b->ht_insert_ct = 0;

1779:   /* stuff for MatGetSubMatrices_MPIBAIJ_local() */
1780:   b->ijonly = PETSC_FALSE;

1782:   b->in_loc = 0;
1783:   b->v_loc  = 0;
1784:   b->n_loc  = 0;
1785:   PetscOptionsBegin(PetscObjectComm((PetscObject)B),NULL,"Options for loading MPISBAIJ matrix 1","Mat");
1786:   PetscOptionsBool("-mat_use_hash_table","Use hash table to save memory in constructing matrix","MatSetOption",PETSC_FALSE,&flg,NULL);
1787:   if (flg) {
1788:     PetscReal fact = 1.39;
1789:     MatSetOption(B,MAT_USE_HASH_TABLE,PETSC_TRUE);
1790:     PetscOptionsReal("-mat_use_hash_table","Use hash table factor","MatMPIBAIJSetHashTableFactor",fact,&fact,NULL);
1791:     if (fact <= 1.0) fact = 1.39;
1792:     MatMPIBAIJSetHashTableFactor(B,fact);
1793:     PetscInfo1(B,"Hash table Factor used %5.2f\n",fact);
1794:   }
1795:   PetscOptionsEnd();

1797: #if defined(PETSC_HAVE_PASTIX)
1798:   PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_pastix_C",MatGetFactor_mpisbaij_pastix);
1799: #endif
1800: #if defined(PETSC_HAVE_MUMPS)
1801:   PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_mumps_C",MatGetFactor_sbaij_mumps);
1802: #endif
1803:   PetscObjectComposeFunction((PetscObject)B,"MatStoreValues_C",MatStoreValues_MPISBAIJ);
1804:   PetscObjectComposeFunction((PetscObject)B,"MatRetrieveValues_C",MatRetrieveValues_MPISBAIJ);
1805:   PetscObjectComposeFunction((PetscObject)B,"MatGetDiagonalBlock_C",MatGetDiagonalBlock_MPISBAIJ);
1806:   PetscObjectComposeFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",MatMPISBAIJSetPreallocation_MPISBAIJ);
1807:   PetscObjectComposeFunction((PetscObject)B,"MatMPISBAIJSetPreallocationCSR_C",MatMPISBAIJSetPreallocationCSR_MPISBAIJ);
1808:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpisbaij_mpisbstrm_C",MatConvert_MPISBAIJ_MPISBSTRM);

1810:   B->symmetric                  = PETSC_TRUE;
1811:   B->structurally_symmetric     = PETSC_TRUE;
1812:   B->symmetric_set              = PETSC_TRUE;
1813:   B->structurally_symmetric_set = PETSC_TRUE;

1815:   PetscObjectChangeTypeName((PetscObject)B,MATMPISBAIJ);
1816:   return(0);
1817: }

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

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

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

1828:   Level: beginner

1830: .seealso: MatCreateMPISBAIJ,MATSEQSBAIJ,MATMPISBAIJ
1831: M*/

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

1841:    Collective on Mat

1843:    Input Parameters:
1844: +  A - the matrix
1845: .  bs   - size of blockk
1846: .  d_nz  - number of block nonzeros per block row in diagonal portion of local
1847:            submatrix  (same for all local rows)
1848: .  d_nnz - array containing the number of block nonzeros in the various block rows
1849:            in the upper triangular and diagonal part of the in diagonal portion of the local
1850:            (possibly different for each block row) or NULL.  If you plan to factor the matrix you must leave room
1851:            for the diagonal entry and set a value even if it is zero.
1852: .  o_nz  - number of block nonzeros per block row in the off-diagonal portion of local
1853:            submatrix (same for all local rows).
1854: -  o_nnz - array containing the number of nonzeros in the various block rows of the
1855:            off-diagonal portion of the local submatrix that is right of the diagonal
1856:            (possibly different for each block row) or NULL.


1859:    Options Database Keys:
1860: .   -mat_no_unroll - uses code that does not unroll the loops in the
1861:                      block calculations (much slower)
1862: .   -mat_block_size - size of the blocks to use

1864:    Notes:

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

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

1871:    Storage Information:
1872:    For a square global matrix we define each processor's diagonal portion
1873:    to be its local rows and the corresponding columns (a square submatrix);
1874:    each processor's off-diagonal portion encompasses the remainder of the
1875:    local matrix (a rectangular submatrix).

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

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

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

1891: .vb
1892:            0 1 2 3 4 5 6 7 8 9 10 11
1893:           --------------------------
1894:    row 3  |. . . d d d o o o o  o  o
1895:    row 4  |. . . d d d o o o o  o  o
1896:    row 5  |. . . d d d o o o o  o  o
1897:           --------------------------
1898: .ve

1900:    Thus, any entries in the d locations are stored in the d (diagonal)
1901:    submatrix, and any entries in the o locations are stored in the
1902:    o (off-diagonal) submatrix.  Note that the d matrix is stored in
1903:    MatSeqSBAIJ format and the o submatrix in MATSEQBAIJ format.

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

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

1914:    Level: intermediate

1916: .keywords: matrix, block, aij, compressed row, sparse, parallel

1918: .seealso: MatCreate(), MatCreateSeqSBAIJ(), MatSetValues(), MatCreateBAIJ(), PetscSplitOwnership()
1919: @*/
1920: PetscErrorCode  MatMPISBAIJSetPreallocation(Mat B,PetscInt bs,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
1921: {

1928:   PetscTryMethod(B,"MatMPISBAIJSetPreallocation_C",(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,bs,d_nz,d_nnz,o_nz,o_nnz));
1929:   return(0);
1930: }

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

1941:    Collective on MPI_Comm

1943:    Input Parameters:
1944: +  comm - MPI communicator
1945: .  bs   - size of blockk
1946: .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
1947:            This value should be the same as the local size used in creating the
1948:            y vector for the matrix-vector product y = Ax.
1949: .  n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
1950:            This value should be the same as the local size used in creating the
1951:            x vector for the matrix-vector product y = Ax.
1952: .  M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
1953: .  N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
1954: .  d_nz  - number of block nonzeros per block row in diagonal portion of local
1955:            submatrix  (same for all local rows)
1956: .  d_nnz - array containing the number of block nonzeros in the various block rows
1957:            in the upper triangular portion of the in diagonal portion of the local
1958:            (possibly different for each block block row) or NULL.
1959:            If you plan to factor the matrix you must leave room for the diagonal entry and
1960:            set its value even if it is zero.
1961: .  o_nz  - number of block nonzeros per block row in the off-diagonal portion of local
1962:            submatrix (same for all local rows).
1963: -  o_nnz - array containing the number of nonzeros in the various block rows of the
1964:            off-diagonal portion of the local submatrix (possibly different for
1965:            each block row) or NULL.

1967:    Output Parameter:
1968: .  A - the matrix

1970:    Options Database Keys:
1971: .   -mat_no_unroll - uses code that does not unroll the loops in the
1972:                      block calculations (much slower)
1973: .   -mat_block_size - size of the blocks to use
1974: .   -mat_mpi - use the parallel matrix data structures even on one processor
1975:                (defaults to using SeqBAIJ format on one processor)

1977:    It is recommended that one use the MatCreate(), MatSetType() and/or MatSetFromOptions(),
1978:    MatXXXXSetPreallocation() paradgm instead of this routine directly.
1979:    [MatXXXXSetPreallocation() is, for example, MatSeqAIJSetPreallocation]

1981:    Notes:
1982:    The number of rows and columns must be divisible by blocksize.
1983:    This matrix type does not support complex Hermitian operation.

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

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

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

1993:    Storage Information:
1994:    For a square global matrix we define each processor's diagonal portion
1995:    to be its local rows and the corresponding columns (a square submatrix);
1996:    each processor's off-diagonal portion encompasses the remainder of the
1997:    local matrix (a rectangular submatrix).

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

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

2008: .vb
2009:            0 1 2 3 4 5 6 7 8 9 10 11
2010:           --------------------------
2011:    row 3  |. . . d d d o o o o  o  o
2012:    row 4  |. . . d d d o o o o  o  o
2013:    row 5  |. . . d d d o o o o  o  o
2014:           --------------------------
2015: .ve

2017:    Thus, any entries in the d locations are stored in the d (diagonal)
2018:    submatrix, and any entries in the o locations are stored in the
2019:    o (off-diagonal) submatrix.  Note that the d matrix is stored in
2020:    MatSeqSBAIJ format and the o submatrix in MATSEQBAIJ format.

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

2030:    Level: intermediate

2032: .keywords: matrix, block, aij, compressed row, sparse, parallel

2034: .seealso: MatCreate(), MatCreateSeqSBAIJ(), MatSetValues(), MatCreateBAIJ()
2035: @*/

2037: 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)
2038: {
2040:   PetscMPIInt    size;

2043:   MatCreate(comm,A);
2044:   MatSetSizes(*A,m,n,M,N);
2045:   MPI_Comm_size(comm,&size);
2046:   if (size > 1) {
2047:     MatSetType(*A,MATMPISBAIJ);
2048:     MatMPISBAIJSetPreallocation(*A,bs,d_nz,d_nnz,o_nz,o_nnz);
2049:   } else {
2050:     MatSetType(*A,MATSEQSBAIJ);
2051:     MatSeqSBAIJSetPreallocation(*A,bs,d_nz,d_nnz);
2052:   }
2053:   return(0);
2054: }


2059: static PetscErrorCode MatDuplicate_MPISBAIJ(Mat matin,MatDuplicateOption cpvalues,Mat *newmat)
2060: {
2061:   Mat            mat;
2062:   Mat_MPISBAIJ   *a,*oldmat = (Mat_MPISBAIJ*)matin->data;
2064:   PetscInt       len=0,nt,bs=matin->rmap->bs,mbs=oldmat->mbs;
2065:   PetscScalar    *array;

2068:   *newmat = 0;

2070:   MatCreate(PetscObjectComm((PetscObject)matin),&mat);
2071:   MatSetSizes(mat,matin->rmap->n,matin->cmap->n,matin->rmap->N,matin->cmap->N);
2072:   MatSetType(mat,((PetscObject)matin)->type_name);
2073:   PetscMemcpy(mat->ops,matin->ops,sizeof(struct _MatOps));
2074:   PetscLayoutReference(matin->rmap,&mat->rmap);
2075:   PetscLayoutReference(matin->cmap,&mat->cmap);

2077:   mat->factortype   = matin->factortype;
2078:   mat->preallocated = PETSC_TRUE;
2079:   mat->assembled    = PETSC_TRUE;
2080:   mat->insertmode   = NOT_SET_VALUES;

2082:   a      = (Mat_MPISBAIJ*)mat->data;
2083:   a->bs2 = oldmat->bs2;
2084:   a->mbs = oldmat->mbs;
2085:   a->nbs = oldmat->nbs;
2086:   a->Mbs = oldmat->Mbs;
2087:   a->Nbs = oldmat->Nbs;


2090:   a->size         = oldmat->size;
2091:   a->rank         = oldmat->rank;
2092:   a->donotstash   = oldmat->donotstash;
2093:   a->roworiented  = oldmat->roworiented;
2094:   a->rowindices   = 0;
2095:   a->rowvalues    = 0;
2096:   a->getrowactive = PETSC_FALSE;
2097:   a->barray       = 0;
2098:   a->rstartbs     = oldmat->rstartbs;
2099:   a->rendbs       = oldmat->rendbs;
2100:   a->cstartbs     = oldmat->cstartbs;
2101:   a->cendbs       = oldmat->cendbs;

2103:   /* hash table stuff */
2104:   a->ht           = 0;
2105:   a->hd           = 0;
2106:   a->ht_size      = 0;
2107:   a->ht_flag      = oldmat->ht_flag;
2108:   a->ht_fact      = oldmat->ht_fact;
2109:   a->ht_total_ct  = 0;
2110:   a->ht_insert_ct = 0;

2112:   PetscMemcpy(a->rangebs,oldmat->rangebs,(a->size+2)*sizeof(PetscInt));
2113:   if (oldmat->colmap) {
2114: #if defined(PETSC_USE_CTABLE)
2115:     PetscTableCreateCopy(oldmat->colmap,&a->colmap);
2116: #else
2117:     PetscMalloc((a->Nbs)*sizeof(PetscInt),&a->colmap);
2118:     PetscLogObjectMemory(mat,(a->Nbs)*sizeof(PetscInt));
2119:     PetscMemcpy(a->colmap,oldmat->colmap,(a->Nbs)*sizeof(PetscInt));
2120: #endif
2121:   } else a->colmap = 0;

2123:   if (oldmat->garray && (len = ((Mat_SeqBAIJ*)(oldmat->B->data))->nbs)) {
2124:     PetscMalloc(len*sizeof(PetscInt),&a->garray);
2125:     PetscLogObjectMemory(mat,len*sizeof(PetscInt));
2126:     PetscMemcpy(a->garray,oldmat->garray,len*sizeof(PetscInt));
2127:   } else a->garray = 0;

2129:   MatStashCreate_Private(PetscObjectComm((PetscObject)matin),matin->rmap->bs,&mat->bstash);
2130:   VecDuplicate(oldmat->lvec,&a->lvec);
2131:   PetscLogObjectParent(mat,a->lvec);
2132:   VecScatterCopy(oldmat->Mvctx,&a->Mvctx);
2133:   PetscLogObjectParent(mat,a->Mvctx);

2135:    VecDuplicate(oldmat->slvec0,&a->slvec0);
2136:   PetscLogObjectParent(mat,a->slvec0);
2137:    VecDuplicate(oldmat->slvec1,&a->slvec1);
2138:   PetscLogObjectParent(mat,a->slvec1);

2140:   VecGetLocalSize(a->slvec1,&nt);
2141:   VecGetArray(a->slvec1,&array);
2142:   VecCreateSeqWithArray(PETSC_COMM_SELF,1,bs*mbs,array,&a->slvec1a);
2143:   VecCreateSeqWithArray(PETSC_COMM_SELF,1,nt-bs*mbs,array+bs*mbs,&a->slvec1b);
2144:   VecRestoreArray(a->slvec1,&array);
2145:   VecGetArray(a->slvec0,&array);
2146:   VecCreateSeqWithArray(PETSC_COMM_SELF,1,nt-bs*mbs,array+bs*mbs,&a->slvec0b);
2147:   VecRestoreArray(a->slvec0,&array);
2148:   PetscLogObjectParent(mat,a->slvec0);
2149:   PetscLogObjectParent(mat,a->slvec1);
2150:   PetscLogObjectParent(mat,a->slvec0b);
2151:   PetscLogObjectParent(mat,a->slvec1a);
2152:   PetscLogObjectParent(mat,a->slvec1b);

2154:   /*  VecScatterCopy(oldmat->sMvctx,&a->sMvctx); - not written yet, replaced by the lazy trick: */
2155:   PetscObjectReference((PetscObject)oldmat->sMvctx);
2156:   a->sMvctx = oldmat->sMvctx;
2157:   PetscLogObjectParent(mat,a->sMvctx);

2159:    MatDuplicate(oldmat->A,cpvalues,&a->A);
2160:   PetscLogObjectParent(mat,a->A);
2161:    MatDuplicate(oldmat->B,cpvalues,&a->B);
2162:   PetscLogObjectParent(mat,a->B);
2163:   PetscFunctionListDuplicate(((PetscObject)matin)->qlist,&((PetscObject)mat)->qlist);
2164:   *newmat = mat;
2165:   return(0);
2166: }

2170: PetscErrorCode MatLoad_MPISBAIJ(Mat newmat,PetscViewer viewer)
2171: {
2173:   PetscInt       i,nz,j,rstart,rend;
2174:   PetscScalar    *vals,*buf;
2175:   MPI_Comm       comm;
2176:   MPI_Status     status;
2177:   PetscMPIInt    rank,size,tag = ((PetscObject)viewer)->tag,*sndcounts = 0,*browners,maxnz,*rowners,mmbs;
2178:   PetscInt       header[4],*rowlengths = 0,M,N,m,*cols,*locrowlens;
2179:   PetscInt       *procsnz = 0,jj,*mycols,*ibuf;
2180:   PetscInt       bs       =1,Mbs,mbs,extra_rows;
2181:   PetscInt       *dlens,*odlens,*mask,*masked1,*masked2,rowcount,odcount;
2182:   PetscInt       dcount,kmax,k,nzcount,tmp,sizesset=1,grows,gcols;
2183:   int            fd;

2186:   PetscObjectGetComm((PetscObject)viewer,&comm);
2187:   PetscOptionsBegin(comm,NULL,"Options for loading MPISBAIJ matrix 2","Mat");
2188:   PetscOptionsInt("-matload_block_size","Set the blocksize used to store the matrix","MatLoad",bs,&bs,NULL);
2189:   PetscOptionsEnd();

2191:   MPI_Comm_size(comm,&size);
2192:   MPI_Comm_rank(comm,&rank);
2193:   if (!rank) {
2194:     PetscViewerBinaryGetDescriptor(viewer,&fd);
2195:     PetscBinaryRead(fd,(char*)header,4,PETSC_INT);
2196:     if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"not matrix object");
2197:     if (header[3] < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as MPISBAIJ");
2198:   }

2200:   if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) sizesset = 0;

2202:   MPI_Bcast(header+1,3,MPIU_INT,0,comm);
2203:   M    = header[1];
2204:   N    = header[2];

2206:   /* If global rows/cols are set to PETSC_DECIDE, set it to the sizes given in the file */
2207:   if (sizesset && newmat->rmap->N < 0) newmat->rmap->N = M;
2208:   if (sizesset && newmat->cmap->N < 0) newmat->cmap->N = N;

2210:   /* If global sizes are set, check if they are consistent with that given in the file */
2211:   if (sizesset) {
2212:     MatGetSize(newmat,&grows,&gcols);
2213:   }
2214:   if (sizesset && newmat->rmap->N != grows) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Inconsistent # of rows:Matrix in file has (%d) and input matrix has (%d)",M,grows);
2215:   if (sizesset && newmat->cmap->N != gcols) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Inconsistent # of cols:Matrix in file has (%d) and input matrix has (%d)",N,gcols);

2217:   if (M != N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Can only do square matrices");

2219:   /*
2220:      This code adds extra rows to make sure the number of rows is
2221:      divisible by the blocksize
2222:   */
2223:   Mbs        = M/bs;
2224:   extra_rows = bs - M + bs*(Mbs);
2225:   if (extra_rows == bs) extra_rows = 0;
2226:   else                  Mbs++;
2227:   if (extra_rows &&!rank) {
2228:     PetscInfo(viewer,"Padding loaded matrix to match blocksize\n");
2229:   }

2231:   /* determine ownership of all rows */
2232:   if (newmat->rmap->n < 0) { /* PETSC_DECIDE */
2233:     mbs = Mbs/size + ((Mbs % size) > rank);
2234:     m   = mbs*bs;
2235:   } else { /* User Set */
2236:     m   = newmat->rmap->n;
2237:     mbs = m/bs;
2238:   }
2239:   PetscMalloc2(size+1,PetscMPIInt,&rowners,size+1,PetscMPIInt,&browners);
2240:   PetscMPIIntCast(mbs,&mmbs);
2241:   MPI_Allgather(&mmbs,1,MPI_INT,rowners+1,1,MPI_INT,comm);
2242:   rowners[0] = 0;
2243:   for (i=2; i<=size; i++) rowners[i] += rowners[i-1];
2244:   for (i=0; i<=size; i++) browners[i] = rowners[i]*bs;
2245:   rstart = rowners[rank];
2246:   rend   = rowners[rank+1];

2248:   /* distribute row lengths to all processors */
2249:   PetscMalloc((rend-rstart)*bs*sizeof(PetscInt),&locrowlens);
2250:   if (!rank) {
2251:     PetscMalloc((M+extra_rows)*sizeof(PetscInt),&rowlengths);
2252:     PetscBinaryRead(fd,rowlengths,M,PETSC_INT);
2253:     for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;
2254:     PetscMalloc(size*sizeof(PetscMPIInt),&sndcounts);
2255:     for (i=0; i<size; i++) sndcounts[i] = browners[i+1] - browners[i];
2256:     MPI_Scatterv(rowlengths,sndcounts,browners,MPIU_INT,locrowlens,(rend-rstart)*bs,MPIU_INT,0,comm);
2257:     PetscFree(sndcounts);
2258:   } else {
2259:     MPI_Scatterv(0,0,0,MPIU_INT,locrowlens,(rend-rstart)*bs,MPIU_INT,0,comm);
2260:   }

2262:   if (!rank) {   /* procs[0] */
2263:     /* calculate the number of nonzeros on each processor */
2264:     PetscMalloc(size*sizeof(PetscInt),&procsnz);
2265:     PetscMemzero(procsnz,size*sizeof(PetscInt));
2266:     for (i=0; i<size; i++) {
2267:       for (j=rowners[i]*bs; j< rowners[i+1]*bs; j++) {
2268:         procsnz[i] += rowlengths[j];
2269:       }
2270:     }
2271:     PetscFree(rowlengths);

2273:     /* determine max buffer needed and allocate it */
2274:     maxnz = 0;
2275:     for (i=0; i<size; i++) {
2276:       maxnz = PetscMax(maxnz,procsnz[i]);
2277:     }
2278:     PetscMalloc(maxnz*sizeof(PetscInt),&cols);

2280:     /* read in my part of the matrix column indices  */
2281:     nz     = procsnz[0];
2282:     PetscMalloc(nz*sizeof(PetscInt),&ibuf);
2283:     mycols = ibuf;
2284:     if (size == 1) nz -= extra_rows;
2285:     PetscBinaryRead(fd,mycols,nz,PETSC_INT);
2286:     if (size == 1) {
2287:       for (i=0; i< extra_rows; i++) mycols[nz+i] = M+i;
2288:     }

2290:     /* read in every ones (except the last) and ship off */
2291:     for (i=1; i<size-1; i++) {
2292:       nz   = procsnz[i];
2293:       PetscBinaryRead(fd,cols,nz,PETSC_INT);
2294:       MPI_Send(cols,nz,MPIU_INT,i,tag,comm);
2295:     }
2296:     /* read in the stuff for the last proc */
2297:     if (size != 1) {
2298:       nz   = procsnz[size-1] - extra_rows;  /* the extra rows are not on the disk */
2299:       PetscBinaryRead(fd,cols,nz,PETSC_INT);
2300:       for (i=0; i<extra_rows; i++) cols[nz+i] = M+i;
2301:       MPI_Send(cols,nz+extra_rows,MPIU_INT,size-1,tag,comm);
2302:     }
2303:     PetscFree(cols);
2304:   } else {  /* procs[i], i>0 */
2305:     /* determine buffer space needed for message */
2306:     nz = 0;
2307:     for (i=0; i<m; i++) nz += locrowlens[i];
2308:     PetscMalloc(nz*sizeof(PetscInt),&ibuf);
2309:     mycols = ibuf;
2310:     /* receive message of column indices*/
2311:     MPI_Recv(mycols,nz,MPIU_INT,0,tag,comm,&status);
2312:     MPI_Get_count(&status,MPIU_INT,&maxnz);
2313:     if (maxnz != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
2314:   }

2316:   /* loop over local rows, determining number of off diagonal entries */
2317:   PetscMalloc2(rend-rstart,PetscInt,&dlens,rend-rstart,PetscInt,&odlens);
2318:   PetscMalloc3(Mbs,PetscInt,&mask,Mbs,PetscInt,&masked1,Mbs,PetscInt,&masked2);
2319:   PetscMemzero(mask,Mbs*sizeof(PetscInt));
2320:   PetscMemzero(masked1,Mbs*sizeof(PetscInt));
2321:   PetscMemzero(masked2,Mbs*sizeof(PetscInt));
2322:   rowcount = 0;
2323:   nzcount  = 0;
2324:   for (i=0; i<mbs; i++) {
2325:     dcount  = 0;
2326:     odcount = 0;
2327:     for (j=0; j<bs; j++) {
2328:       kmax = locrowlens[rowcount];
2329:       for (k=0; k<kmax; k++) {
2330:         tmp = mycols[nzcount++]/bs; /* block col. index */
2331:         if (!mask[tmp]) {
2332:           mask[tmp] = 1;
2333:           if (tmp < rstart || tmp >= rend) masked2[odcount++] = tmp; /* entry in off-diag portion */
2334:           else masked1[dcount++] = tmp; /* entry in diag portion */
2335:         }
2336:       }
2337:       rowcount++;
2338:     }

2340:     dlens[i]  = dcount;  /* d_nzz[i] */
2341:     odlens[i] = odcount; /* o_nzz[i] */

2343:     /* zero out the mask elements we set */
2344:     for (j=0; j<dcount; j++) mask[masked1[j]] = 0;
2345:     for (j=0; j<odcount; j++) mask[masked2[j]] = 0;
2346:   }
2347:   if (!sizesset) {
2348:     MatSetSizes(newmat,m,m,M+extra_rows,N+extra_rows);
2349:   }
2350:   MatMPISBAIJSetPreallocation(newmat,bs,0,dlens,0,odlens);
2351:   MatSetOption(newmat,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);

2353:   if (!rank) {
2354:     PetscMalloc(maxnz*sizeof(PetscScalar),&buf);
2355:     /* read in my part of the matrix numerical values  */
2356:     nz     = procsnz[0];
2357:     vals   = buf;
2358:     mycols = ibuf;
2359:     if (size == 1) nz -= extra_rows;
2360:     PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
2361:     if (size == 1) {
2362:       for (i=0; i< extra_rows; i++) vals[nz+i] = 1.0;
2363:     }

2365:     /* insert into matrix */
2366:     jj = rstart*bs;
2367:     for (i=0; i<m; i++) {
2368:       MatSetValues(newmat,1,&jj,locrowlens[i],mycols,vals,INSERT_VALUES);
2369:       mycols += locrowlens[i];
2370:       vals   += locrowlens[i];
2371:       jj++;
2372:     }

2374:     /* read in other processors (except the last one) and ship out */
2375:     for (i=1; i<size-1; i++) {
2376:       nz   = procsnz[i];
2377:       vals = buf;
2378:       PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
2379:       MPI_Send(vals,nz,MPIU_SCALAR,i,((PetscObject)newmat)->tag,comm);
2380:     }
2381:     /* the last proc */
2382:     if (size != 1) {
2383:       nz   = procsnz[i] - extra_rows;
2384:       vals = buf;
2385:       PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
2386:       for (i=0; i<extra_rows; i++) vals[nz+i] = 1.0;
2387:       MPI_Send(vals,nz+extra_rows,MPIU_SCALAR,size-1,((PetscObject)newmat)->tag,comm);
2388:     }
2389:     PetscFree(procsnz);

2391:   } else {
2392:     /* receive numeric values */
2393:     PetscMalloc(nz*sizeof(PetscScalar),&buf);

2395:     /* receive message of values*/
2396:     vals   = buf;
2397:     mycols = ibuf;
2398:     MPI_Recv(vals,nz,MPIU_SCALAR,0,((PetscObject)newmat)->tag,comm,&status);
2399:     MPI_Get_count(&status,MPIU_SCALAR,&maxnz);
2400:     if (maxnz != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");

2402:     /* insert into matrix */
2403:     jj = rstart*bs;
2404:     for (i=0; i<m; i++) {
2405:       MatSetValues_MPISBAIJ(newmat,1,&jj,locrowlens[i],mycols,vals,INSERT_VALUES);
2406:       mycols += locrowlens[i];
2407:       vals   += locrowlens[i];
2408:       jj++;
2409:     }
2410:   }

2412:   PetscFree(locrowlens);
2413:   PetscFree(buf);
2414:   PetscFree(ibuf);
2415:   PetscFree2(rowners,browners);
2416:   PetscFree2(dlens,odlens);
2417:   PetscFree3(mask,masked1,masked2);
2418:   MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);
2419:   MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);
2420:   return(0);
2421: }

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

2428:    Input Parameters:
2429: .  mat  - the matrix
2430: .  fact - factor

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

2434:    Level: advanced

2436:   Notes:
2437:    This can also be set by the command line option: -mat_use_hash_table fact

2439: .keywords: matrix, hashtable, factor, HT

2441: .seealso: MatSetOption()
2442: @XXXXX*/


2447: PetscErrorCode MatGetRowMaxAbs_MPISBAIJ(Mat A,Vec v,PetscInt idx[])
2448: {
2449:   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;
2450:   Mat_SeqBAIJ    *b = (Mat_SeqBAIJ*)(a->B)->data;
2451:   PetscReal      atmp;
2452:   PetscReal      *work,*svalues,*rvalues;
2454:   PetscInt       i,bs,mbs,*bi,*bj,brow,j,ncols,krow,kcol,col,row,Mbs,bcol;
2455:   PetscMPIInt    rank,size;
2456:   PetscInt       *rowners_bs,dest,count,source;
2457:   PetscScalar    *va;
2458:   MatScalar      *ba;
2459:   MPI_Status     stat;

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

2466:   MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);
2467:   MPI_Comm_rank(PetscObjectComm((PetscObject)A),&rank);

2469:   bs  = A->rmap->bs;
2470:   mbs = a->mbs;
2471:   Mbs = a->Mbs;
2472:   ba  = b->a;
2473:   bi  = b->i;
2474:   bj  = b->j;

2476:   /* find ownerships */
2477:   rowners_bs = A->rmap->range;

2479:   /* each proc creates an array to be distributed */
2480:   PetscMalloc(bs*Mbs*sizeof(PetscReal),&work);
2481:   PetscMemzero(work,bs*Mbs*sizeof(PetscReal));

2483:   /* row_max for B */
2484:   if (rank != size-1) {
2485:     for (i=0; i<mbs; i++) {
2486:       ncols = bi[1] - bi[0]; bi++;
2487:       brow  = bs*i;
2488:       for (j=0; j<ncols; j++) {
2489:         bcol = bs*(*bj);
2490:         for (kcol=0; kcol<bs; kcol++) {
2491:           col  = bcol + kcol;                /* local col index */
2492:           col += rowners_bs[rank+1];      /* global col index */
2493:           for (krow=0; krow<bs; krow++) {
2494:             atmp = PetscAbsScalar(*ba); ba++;
2495:             row  = brow + krow;   /* local row index */
2496:             if (PetscRealPart(va[row]) < atmp) va[row] = atmp;
2497:             if (work[col] < atmp) work[col] = atmp;
2498:           }
2499:         }
2500:         bj++;
2501:       }
2502:     }

2504:     /* send values to its owners */
2505:     for (dest=rank+1; dest<size; dest++) {
2506:       svalues = work + rowners_bs[dest];
2507:       count   = rowners_bs[dest+1]-rowners_bs[dest];
2508:       MPI_Send(svalues,count,MPIU_REAL,dest,rank,PetscObjectComm((PetscObject)A));
2509:     }
2510:   }

2512:   /* receive values */
2513:   if (rank) {
2514:     rvalues = work;
2515:     count   = rowners_bs[rank+1]-rowners_bs[rank];
2516:     for (source=0; source<rank; source++) {
2517:       MPI_Recv(rvalues,count,MPIU_REAL,MPI_ANY_SOURCE,MPI_ANY_TAG,PetscObjectComm((PetscObject)A),&stat);
2518:       /* process values */
2519:       for (i=0; i<count; i++) {
2520:         if (PetscRealPart(va[i]) < rvalues[i]) va[i] = rvalues[i];
2521:       }
2522:     }
2523:   }

2525:   VecRestoreArray(v,&va);
2526:   PetscFree(work);
2527:   return(0);
2528: }

2532: PetscErrorCode MatSOR_MPISBAIJ(Mat matin,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
2533: {
2534:   Mat_MPISBAIJ      *mat = (Mat_MPISBAIJ*)matin->data;
2535:   PetscErrorCode    ierr;
2536:   PetscInt          mbs=mat->mbs,bs=matin->rmap->bs;
2537:   PetscScalar       *x,*ptr,*from;
2538:   Vec               bb1;
2539:   const PetscScalar *b;

2542:   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);
2543:   if (bs > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented");

2545:   if (flag == SOR_APPLY_UPPER) {
2546:     (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx);
2547:     return(0);
2548:   }

2550:   if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP) {
2551:     if (flag & SOR_ZERO_INITIAL_GUESS) {
2552:       (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,lits,xx);
2553:       its--;
2554:     }

2556:     VecDuplicate(bb,&bb1);
2557:     while (its--) {

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

2562:       /* copy xx into slvec0a */
2563:       VecGetArray(mat->slvec0,&ptr);
2564:       VecGetArray(xx,&x);
2565:       PetscMemcpy(ptr,x,bs*mbs*sizeof(MatScalar));
2566:       VecRestoreArray(mat->slvec0,&ptr);

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

2570:       /* copy bb into slvec1a */
2571:       VecGetArray(mat->slvec1,&ptr);
2572:       VecGetArrayRead(bb,&b);
2573:       PetscMemcpy(ptr,b,bs*mbs*sizeof(MatScalar));
2574:       VecRestoreArray(mat->slvec1,&ptr);

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

2579:       VecScatterBegin(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD);
2580:       VecRestoreArray(xx,&x);
2581:       VecRestoreArrayRead(bb,&b);
2582:       VecScatterEnd(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD);

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

2587:       /* local diagonal sweep */
2588:       (*mat->A->ops->sor)(mat->A,bb1,omega,SOR_SYMMETRIC_SWEEP,fshift,lits,lits,xx);
2589:     }
2590:     VecDestroy(&bb1);
2591:   } else if ((flag & SOR_LOCAL_FORWARD_SWEEP) && (its == 1) && (flag & SOR_ZERO_INITIAL_GUESS)) {
2592:     (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx);
2593:   } else if ((flag & SOR_LOCAL_BACKWARD_SWEEP) && (its == 1) && (flag & SOR_ZERO_INITIAL_GUESS)) {
2594:     (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx);
2595:   } else if (flag & SOR_EISENSTAT) {
2596:     Vec               xx1;
2597:     PetscBool         hasop;
2598:     const PetscScalar *diag;
2599:     PetscScalar       *sl,scale = (omega - 2.0)/omega;
2600:     PetscInt          i,n;

2602:     if (!mat->xx1) {
2603:       VecDuplicate(bb,&mat->xx1);
2604:       VecDuplicate(bb,&mat->bb1);
2605:     }
2606:     xx1 = mat->xx1;
2607:     bb1 = mat->bb1;

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

2611:     if (!mat->diag) {
2612:       /* this is wrong for same matrix with new nonzero values */
2613:       MatGetVecs(matin,&mat->diag,NULL);
2614:       MatGetDiagonal(matin,mat->diag);
2615:     }
2616:     MatHasOperation(matin,MATOP_MULT_DIAGONAL_BLOCK,&hasop);

2618:     if (hasop) {
2619:       MatMultDiagonalBlock(matin,xx,bb1);
2620:       VecAYPX(mat->slvec1a,scale,bb);
2621:     } else {
2622:       /*
2623:           These two lines are replaced by code that may be a bit faster for a good compiler
2624:       VecPointwiseMult(mat->slvec1a,mat->diag,xx);
2625:       VecAYPX(mat->slvec1a,scale,bb);
2626:       */
2627:       VecGetArray(mat->slvec1a,&sl);
2628:       VecGetArrayRead(mat->diag,&diag);
2629:       VecGetArrayRead(bb,&b);
2630:       VecGetArray(xx,&x);
2631:       VecGetLocalSize(xx,&n);
2632:       if (omega == 1.0) {
2633:         for (i=0; i<n; i++) sl[i] = b[i] - diag[i]*x[i];
2634:         PetscLogFlops(2.0*n);
2635:       } else {
2636:         for (i=0; i<n; i++) sl[i] = b[i] + scale*diag[i]*x[i];
2637:         PetscLogFlops(3.0*n);
2638:       }
2639:       VecRestoreArray(mat->slvec1a,&sl);
2640:       VecRestoreArrayRead(mat->diag,&diag);
2641:       VecRestoreArrayRead(bb,&b);
2642:       VecRestoreArray(xx,&x);
2643:     }

2645:     /* multiply off-diagonal portion of matrix */
2646:     VecSet(mat->slvec1b,0.0);
2647:     (*mat->B->ops->multtranspose)(mat->B,xx,mat->slvec0b);
2648:     VecGetArray(mat->slvec0,&from);
2649:     VecGetArray(xx,&x);
2650:     PetscMemcpy(from,x,bs*mbs*sizeof(MatScalar));
2651:     VecRestoreArray(mat->slvec0,&from);
2652:     VecRestoreArray(xx,&x);
2653:     VecScatterBegin(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD);
2654:     VecScatterEnd(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD);
2655:     (*mat->B->ops->multadd)(mat->B,mat->slvec1b,mat->slvec1a,mat->slvec1a);

2657:     /* local sweep */
2658:     (*mat->A->ops->sor)(mat->A,mat->slvec1a,omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_FORWARD_SWEEP),fshift,lits,1,xx1);
2659:     VecAXPY(xx,1.0,xx1);
2660:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatSORType is not supported for SBAIJ matrix format");
2661:   return(0);
2662: }

2666: PetscErrorCode MatSOR_MPISBAIJ_2comm(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;
2670:   Vec            lvec1,bb1;

2673:   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);
2674:   if (matin->rmap->bs > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented");

2676:   if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP) {
2677:     if (flag & SOR_ZERO_INITIAL_GUESS) {
2678:       (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,lits,xx);
2679:       its--;
2680:     }

2682:     VecDuplicate(mat->lvec,&lvec1);
2683:     VecDuplicate(bb,&bb1);
2684:     while (its--) {
2685:       VecScatterBegin(mat->Mvctx,xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD);

2687:       /* lower diagonal part: bb1 = bb - B^T*xx */
2688:       (*mat->B->ops->multtranspose)(mat->B,xx,lvec1);
2689:       VecScale(lvec1,-1.0);

2691:       VecScatterEnd(mat->Mvctx,xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD);
2692:       VecCopy(bb,bb1);
2693:       VecScatterBegin(mat->Mvctx,lvec1,bb1,ADD_VALUES,SCATTER_REVERSE);

2695:       /* upper diagonal part: bb1 = bb1 - B*x */
2696:       VecScale(mat->lvec,-1.0);
2697:       (*mat->B->ops->multadd)(mat->B,mat->lvec,bb1,bb1);

2699:       VecScatterEnd(mat->Mvctx,lvec1,bb1,ADD_VALUES,SCATTER_REVERSE);

2701:       /* diagonal sweep */
2702:       (*mat->A->ops->sor)(mat->A,bb1,omega,SOR_SYMMETRIC_SWEEP,fshift,lits,lits,xx);
2703:     }
2704:     VecDestroy(&lvec1);
2705:     VecDestroy(&bb1);
2706:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatSORType is not supported for SBAIJ matrix format");
2707:   return(0);
2708: }

2712: /*@
2713:      MatCreateMPISBAIJWithArrays - creates a MPI SBAIJ matrix using arrays that contain in standard
2714:          CSR format the local rows.

2716:    Collective on MPI_Comm

2718:    Input Parameters:
2719: +  comm - MPI communicator
2720: .  bs - the block size, only a block size of 1 is supported
2721: .  m - number of local rows (Cannot be PETSC_DECIDE)
2722: .  n - This value should be the same as the local size used in creating the
2723:        x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have
2724:        calculated if N is given) For square matrices n is almost always m.
2725: .  M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
2726: .  N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
2727: .   i - row indices
2728: .   j - column indices
2729: -   a - matrix values

2731:    Output Parameter:
2732: .   mat - the matrix

2734:    Level: intermediate

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

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

2743: .keywords: matrix, aij, compressed row, sparse, parallel

2745: .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatMPIAIJSetPreallocation(), MatMPIAIJSetPreallocationCSR(),
2746:           MPIAIJ, MatCreateAIJ(), MatCreateMPIAIJWithSplitArrays()
2747: @*/
2748: 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)
2749: {


2754:   if (i[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"i (row indices) must start with 0");
2755:   if (m < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"local number of rows (m) cannot be PETSC_DECIDE, or negative");
2756:   MatCreate(comm,mat);
2757:   MatSetSizes(*mat,m,n,M,N);
2758:   MatSetType(*mat,MATMPISBAIJ);
2759:   MatMPISBAIJSetPreallocationCSR(*mat,bs,i,j,a);
2760:   return(0);
2761: }


2766: /*@C
2767:    MatMPISBAIJSetPreallocationCSR - Allocates memory for a sparse parallel matrix in BAIJ format
2768:    (the default parallel PETSc format).

2770:    Collective on MPI_Comm

2772:    Input Parameters:
2773: +  A - the matrix
2774: .  bs - the block size
2775: .  i - the indices into j for the start of each local row (starts with zero)
2776: .  j - the column indices for each local row (starts with zero) these must be sorted for each row
2777: -  v - optional values in the matrix

2779:    Level: developer

2781: .keywords: matrix, aij, compressed row, sparse, parallel

2783: .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatMPIBAIJSetPreallocation(), MatCreateAIJ(), MPIAIJ
2784: @*/
2785: PetscErrorCode  MatMPISBAIJSetPreallocationCSR(Mat B,PetscInt bs,const PetscInt i[],const PetscInt j[], const PetscScalar v[])
2786: {

2790:   PetscTryMethod(B,"MatMPISBAIJSetPreallocationCSR_C",(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]),(B,bs,i,j,v));
2791:   return(0);
2792: }