Actual source code: mpisbaij.c

petsc-3.5.1 2014-07-24
Report Typos and Errors
  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->nonzerostate++;\
 85: a_noinsert:; \
 86:     ailen[brow] = nrow; \
 87:   }

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

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

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

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

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

159:   /* for stash */
160:   PetscInt  n_loc, *in_loc = NULL;
161:   MatScalar *v_loc = NULL;

164:   if (!baij->donotstash) {
165:     if (n > baij->n_loc) {
166:       PetscFree(baij->in_loc);
167:       PetscFree(baij->v_loc);
168:       PetscMalloc1(n,&baij->in_loc);
169:       PetscMalloc1(n,&baij->v_loc);

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

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

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

265:   if (!barray) {
266:     PetscMalloc1(bs2,&barray);
267:     baij->barray = barray;
268:   }

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

306:         if (in[j] >= cstart && in[j] < cend) {
307:           col  = in[j] - cstart;
308:           MatSetValuesBlocked_SeqSBAIJ(baij->A,1,&row,1,&col,barray,addv);
309:         } else if (in[j] < 0) continue;
310: #if defined(PETSC_USE_DEBUG)
311:         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);
312: #endif
313:         else {
314:           if (mat->was_assembled) {
315:             if (!baij->colmap) {
316:               MatCreateColmap_MPIBAIJ_Private(mat);
317:             }

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

594:   baij->rowvalues = 0;

596:   /* if no new nonzero locations are allowed in matrix then only set the matrix state the first time through */
597:   if ((!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) || !((Mat_SeqBAIJ*)(baij->A->data))->nonew) {
598:     PetscObjectState state = baij->A->nonzerostate + baij->B->nonzerostate;
599:     MPI_Allreduce(&state,&mat->nonzerostate,1,MPIU_INT64,MPI_SUM,PetscObjectComm((PetscObject)mat));
600:   }
601:   return(0);
602: }

604: extern PetscErrorCode MatView_SeqSBAIJ_ASCII(Mat,PetscViewer);
605: extern PetscErrorCode MatSetValues_MPIBAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
606: #include <petscdraw.h>
609: static PetscErrorCode MatView_MPISBAIJ_ASCIIorDraworSocket(Mat mat,PetscViewer viewer)
610: {
611:   Mat_MPISBAIJ      *baij = (Mat_MPISBAIJ*)mat->data;
612:   PetscErrorCode    ierr;
613:   PetscInt          bs   = mat->rmap->bs;
614:   PetscMPIInt       rank = baij->rank;
615:   PetscBool         iascii,isdraw;
616:   PetscViewer       sviewer;
617:   PetscViewerFormat format;

620:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
621:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
622:   if (iascii) {
623:     PetscViewerGetFormat(viewer,&format);
624:     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
625:       MatInfo info;
626:       MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&rank);
627:       MatGetInfo(mat,MAT_LOCAL,&info);
628:       PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);
629:       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);
630:       MatGetInfo(baij->A,MAT_LOCAL,&info);
631:       PetscViewerASCIISynchronizedPrintf(viewer,"[%d] on-diagonal part: nz %D \n",rank,(PetscInt)info.nz_used);
632:       MatGetInfo(baij->B,MAT_LOCAL,&info);
633:       PetscViewerASCIISynchronizedPrintf(viewer,"[%d] off-diagonal part: nz %D \n",rank,(PetscInt)info.nz_used);
634:       PetscViewerFlush(viewer);
635:       PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);
636:       PetscViewerASCIIPrintf(viewer,"Information on VecScatter used in matrix-vector product: \n");
637:       VecScatterView(baij->Mvctx,viewer);
638:       return(0);
639:     } else if (format == PETSC_VIEWER_ASCII_INFO) {
640:       PetscViewerASCIIPrintf(viewer,"  block size is %D\n",bs);
641:       return(0);
642:     } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
643:       return(0);
644:     }
645:   }

647:   if (isdraw) {
648:     PetscDraw draw;
649:     PetscBool isnull;
650:     PetscViewerDrawGetDraw(viewer,0,&draw);
651:     PetscDrawIsNull(draw,&isnull); if (isnull) return(0);
652:   }

654:   {
655:     /* assemble the entire matrix onto first processor. */
656:     Mat          A;
657:     Mat_SeqSBAIJ *Aloc;
658:     Mat_SeqBAIJ  *Bloc;
659:     PetscInt     M = mat->rmap->N,N = mat->cmap->N,*ai,*aj,col,i,j,k,*rvals,mbs = baij->mbs;
660:     MatScalar    *a;

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

674:     /* copy over the A part */
675:     Aloc = (Mat_SeqSBAIJ*)baij->A->data;
676:     ai   = Aloc->i; aj = Aloc->j; a = Aloc->a;
677:     PetscMalloc1(bs,&rvals);

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

1147:   MatGetInfo(A,MAT_LOCAL,info);

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

1152:   MatGetInfo(B,MAT_LOCAL,info);

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

1367: PetscErrorCode MatSetUp_MPISBAIJ(Mat A)
1368: {

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

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

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

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

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


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


1572: PetscErrorCode  MatGetDiagonalBlock_MPISBAIJ(Mat A,Mat *a)
1573: {
1575:   *a = ((Mat_MPISBAIJ*)A->data)->A;
1576:   return(0);
1577: }

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

1588:   MatSetBlockSize(B,PetscAbs(bs));
1589:   PetscLayoutSetUp(B->rmap);
1590:   PetscLayoutSetUp(B->cmap);
1591:   PetscLayoutGetBlockSize(B->rmap,&bs);

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

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

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

1611:   b->cstartbs = B->cmap->rstart/bs;
1612:   b->cendbs   = B->cmap->rend/bs;

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

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

1629:   B->preallocated = PETSC_TRUE;
1630:   return(0);
1631: }

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

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

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

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

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

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

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

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

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

1715:   Level: beginner

1717: .seealso: MatCreateMPISBAIJ
1718: M*/

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

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

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

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

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

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

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

1749:   b->donotstash  = PETSC_FALSE;
1750:   b->colmap      = NULL;
1751:   b->garray      = NULL;
1752:   b->roworiented = PETSC_TRUE;

1754:   /* stuff used in block assembly */
1755:   b->barray = 0;

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

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

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

1781:   /* stuff for MatGetSubMatrices_MPIBAIJ_local() */
1782:   b->ijonly = PETSC_FALSE;

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

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

1812:   B->symmetric                  = PETSC_TRUE;
1813:   B->structurally_symmetric     = PETSC_TRUE;
1814:   B->symmetric_set              = PETSC_TRUE;
1815:   B->structurally_symmetric_set = PETSC_TRUE;

1817:   PetscObjectChangeTypeName((PetscObject)B,MATMPISBAIJ);
1818:   return(0);
1819: }

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

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

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

1830:   Level: beginner

1832: .seealso: MatCreateMPISBAIJ,MATSEQSBAIJ,MATMPISBAIJ
1833: M*/

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

1843:    Collective on Mat

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


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

1866:    Notes:

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

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

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

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

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

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

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

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

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

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

1916:    Level: intermediate

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

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

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

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

1943:    Collective on MPI_Comm

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

1969:    Output Parameter:
1970: .  A - the matrix

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

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

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

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

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

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

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

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

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

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

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

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

2032:    Level: intermediate

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

2036: .seealso: MatCreate(), MatCreateSeqSBAIJ(), MatSetValues(), MatCreateBAIJ()
2037: @*/

2039: 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)
2040: {
2042:   PetscMPIInt    size;

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


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

2070:   *newmat = 0;

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

2079:   mat->factortype   = matin->factortype;
2080:   mat->preallocated = PETSC_TRUE;
2081:   mat->assembled    = PETSC_TRUE;
2082:   mat->insertmode   = NOT_SET_VALUES;

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


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

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

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

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

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

2137:    VecDuplicate(oldmat->slvec0,&a->slvec0);
2138:   PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec0);
2139:    VecDuplicate(oldmat->slvec1,&a->slvec1);
2140:   PetscLogObjectParent((PetscObject)mat,(PetscObject)a->slvec1);

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

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

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

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

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

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

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

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

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

2212:   /* If global sizes are set, check if they are consistent with that given in the file */
2213:   if (sizesset) {
2214:     MatGetSize(newmat,&grows,&gcols);
2215:   }
2216:   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);
2217:   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);

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

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

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

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

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

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

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

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

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

2342:     dlens[i]  = dcount;  /* d_nzz[i] */
2343:     odlens[i] = odcount; /* o_nzz[i] */

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

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

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

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

2393:   } else {
2394:     /* receive numeric values */
2395:     PetscMalloc1(nz,&buf);

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

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

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

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

2430:    Input Parameters:
2431: .  mat  - the matrix
2432: .  fact - factor

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

2436:    Level: advanced

2438:   Notes:
2439:    This can also be set by the command line option: -mat_use_hash_table fact

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

2443: .seealso: MatSetOption()
2444: @XXXXX*/


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

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

2468:   MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);
2469:   MPI_Comm_rank(PetscObjectComm((PetscObject)A),&rank);

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

2478:   /* find ownerships */
2479:   rowners_bs = A->rmap->range;

2481:   /* each proc creates an array to be distributed */
2482:   PetscMalloc1(bs*Mbs,&work);
2483:   PetscMemzero(work,bs*Mbs*sizeof(PetscReal));

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

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

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

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

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

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

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

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

2558:     VecDuplicate(bb,&bb1);
2559:     while (its--) {

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

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

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

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

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

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

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

2589:       /* local diagonal sweep */
2590:       (*mat->A->ops->sor)(mat->A,bb1,omega,SOR_SYMMETRIC_SWEEP,fshift,lits,lits,xx);
2591:     }
2592:     VecDestroy(&bb1);
2593:   } else if ((flag & SOR_LOCAL_FORWARD_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_LOCAL_BACKWARD_SWEEP) && (its == 1) && (flag & SOR_ZERO_INITIAL_GUESS)) {
2596:     (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx);
2597:   } else if (flag & SOR_EISENSTAT) {
2598:     Vec               xx1;
2599:     PetscBool         hasop;
2600:     const PetscScalar *diag;
2601:     PetscScalar       *sl,scale = (omega - 2.0)/omega;
2602:     PetscInt          i,n;

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

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

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

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

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

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

2668: PetscErrorCode MatSOR_MPISBAIJ_2comm(Mat matin,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
2669: {
2670:   Mat_MPISBAIJ   *mat = (Mat_MPISBAIJ*)matin->data;
2672:   Vec            lvec1,bb1;

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

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

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

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

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

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

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

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

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

2718:    Collective on MPI_Comm

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

2733:    Output Parameter:
2734: .   mat - the matrix

2736:    Level: intermediate

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

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

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

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


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


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

2772:    Collective on MPI_Comm

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

2781:    Level: developer

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

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

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