Actual source code: mpisbaij.c

  1: #define PETSCMAT_DLL

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

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

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

 33:   MatStoreValues(aij->A);
 34:   MatStoreValues(aij->B);
 35:   return(0);
 36: }

 42: PetscErrorCode  MatRetrieveValues_MPISBAIJ(Mat mat)
 43: {
 44:   Mat_MPISBAIJ   *aij = (Mat_MPISBAIJ *)mat->data;

 48:   MatRetrieveValues(aij->A);
 49:   MatRetrieveValues(aij->B);
 50:   return(0);
 51: }


 55: #define CHUNKSIZE  10

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

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

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

150:   /* Some Variables required in the macro */
151:   Mat            A = baij->A;
152:   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)(A)->data;
153:   PetscInt       *aimax=a->imax,*ai=a->i,*ailen=a->ilen,*aj=a->j;
154:   MatScalar      *aa=a->a;

156:   Mat            B = baij->B;
157:   Mat_SeqBAIJ   *b = (Mat_SeqBAIJ*)(B)->data;
158:   PetscInt      *bimax=b->imax,*bi=b->i,*bilen=b->ilen,*bj=b->j;
159:   MatScalar     *ba=b->a;

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

165:   /* for stash */
166:   PetscInt      n_loc, *in_loc = PETSC_NULL;
167:   MatScalar     *v_loc = PETSC_NULL;

171:   if (!baij->donotstash){
172:     if (n > baij->n_loc) {
173:       PetscFree(baij->in_loc);
174:       PetscFree(baij->v_loc);
175:       PetscMalloc(n*sizeof(PetscInt),&baij->in_loc);
176:       PetscMalloc(n*sizeof(MatScalar),&baij->v_loc);
177:       baij->n_loc = n;
178:     }
179:     in_loc = baij->in_loc;
180:     v_loc  = baij->v_loc;
181:   }

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

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

269:   if(!barray) {
270:     PetscMalloc(bs2*sizeof(MatScalar),&barray);
271:     baij->barray = barray;
272:   }

274:   if (roworiented) {
275:     stepval = (n-1)*bs;
276:   } else {
277:     stepval = (m-1)*bs;
278:   }
279:   for (i=0; i<m; i++) {
280:     if (im[i] < 0) continue;
281: #if defined(PETSC_USE_DEBUG)
282:     if (im[i] >= baij->Mbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large, row %D max %D",im[i],baij->Mbs-1);
283: #endif
284:     if (im[i] >= rstart && im[i] < rend) {
285:       row = im[i] - rstart;
286:       for (j=0; j<n; j++) {
287:         if (im[i] > in[j]) {
288:           if (ignore_ltriangular) continue; /* ignore lower triangular blocks */
289:           else SETERRQ(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)");
290:         }
291:         /* If NumCol = 1 then a copy is not required */
292:         if ((roworiented) && (n == 1)) {
293:           barray = (MatScalar*) v + i*bs2;
294:         } else if((!roworiented) && (m == 1)) {
295:           barray = (MatScalar*) v + j*bs2;
296:         } else { /* Here a copy is required */
297:           if (roworiented) {
298:             value = v + i*(stepval+bs)*bs + j*bs;
299:           } else {
300:             value = v + j*(stepval+bs)*bs + i*bs;
301:           }
302:           for (ii=0; ii<bs; ii++,value+=stepval) {
303:             for (jj=0; jj<bs; jj++) {
304:               *barray++  = *value++;
305:             }
306:           }
307:           barray -=bs2;
308:         }
309: 
310:         if (in[j] >= cstart && in[j] < cend){
311:           col  = in[j] - cstart;
312:           MatSetValuesBlocked_SeqSBAIJ(baij->A,1,&row,1,&col,barray,addv);
313:         }
314:         else if (in[j] < 0) continue;
315: #if defined(PETSC_USE_DEBUG)
316:         else if (in[j] >= baij->Nbs) {SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large, col %D max %D",in[j],baij->Nbs-1);}
317: #endif
318:         else {
319:           if (mat->was_assembled) {
320:             if (!baij->colmap) {
321:               CreateColmap_MPIBAIJ_Private(mat);
322:             }

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

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

372:   for (i=0; i<m; i++) {
373:     if (idxm[i] < 0) continue; /* SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",idxm[i]); */
374:     if (idxm[i] >= mat->rmap->N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm[i],mat->rmap->N-1);
375:     if (idxm[i] >= bsrstart && idxm[i] < bsrend) {
376:       row = idxm[i] - bsrstart;
377:       for (j=0; j<n; j++) {
378:         if (idxn[j] < 0) continue; /* SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative column %D",idxn[j]); */
379:         if (idxn[j] >= mat->cmap->N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",idxn[j],mat->cmap->N-1);
380:         if (idxn[j] >= bscstart && idxn[j] < bscend){
381:           col = idxn[j] - bscstart;
382:           MatGetValues_SeqSBAIJ(baij->A,1,&row,1,&col,v+i*n+j);
383:         } else {
384:           if (!baij->colmap) {
385:             CreateColmap_MPIBAIJ_Private(mat);
386:           }
387: #if defined (PETSC_USE_CTABLE)
388:           PetscTableFind(baij->colmap,idxn[j]/bs+1,&data);
389:           data --;
390: #else
391:           data = baij->colmap[idxn[j]/bs]-1;
392: #endif
393:           if((data < 0) || (baij->garray[data/bs] != idxn[j]/bs)) *(v+i*n+j) = 0.0;
394:           else {
395:             col  = data + idxn[j]%bs;
396:             MatGetValues_SeqBAIJ(baij->B,1,&row,1,&col,v+i*n+j);
397:           }
398:         }
399:       }
400:     } else {
401:       SETERRQ(PETSC_ERR_SUP,"Only local values currently supported");
402:     }
403:   }
404:  return(0);
405: }

409: PetscErrorCode MatNorm_MPISBAIJ(Mat mat,NormType type,PetscReal *norm)
410: {
411:   Mat_MPISBAIJ   *baij = (Mat_MPISBAIJ*)mat->data;
413:   PetscReal      sum[2],*lnorm2;

416:   if (baij->size == 1) {
417:      MatNorm(baij->A,type,norm);
418:   } else {
419:     if (type == NORM_FROBENIUS) {
420:       PetscMalloc(2*sizeof(PetscReal),&lnorm2);
421:        MatNorm(baij->A,type,lnorm2);
422:       *lnorm2 = (*lnorm2)*(*lnorm2); lnorm2++;            /* squar power of norm(A) */
423:        MatNorm(baij->B,type,lnorm2);
424:       *lnorm2 = (*lnorm2)*(*lnorm2); lnorm2--;             /* squar power of norm(B) */
425:       MPI_Allreduce(lnorm2,&sum,2,MPIU_REAL,MPI_SUM,((PetscObject)mat)->comm);
426:       *norm = sqrt(sum[0] + 2*sum[1]);
427:       PetscFree(lnorm2);
428:     } else if (type == NORM_INFINITY || type == NORM_1) { /* max row/column sum */
429:       Mat_SeqSBAIJ *amat=(Mat_SeqSBAIJ*)baij->A->data;
430:       Mat_SeqBAIJ  *bmat=(Mat_SeqBAIJ*)baij->B->data;
431:       PetscReal    *rsum,*rsum2,vabs;
432:       PetscInt     *jj,*garray=baij->garray,rstart=baij->rstartbs,nz;
433:       PetscInt     brow,bcol,col,bs=baij->A->rmap->bs,row,grow,gcol,mbs=amat->mbs;
434:       MatScalar    *v;

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

486: PetscErrorCode MatAssemblyBegin_MPISBAIJ(Mat mat,MatAssemblyType mode)
487: {
488:   Mat_MPISBAIJ   *baij = (Mat_MPISBAIJ*)mat->data;
490:   PetscInt       nstash,reallocs;
491:   InsertMode     addv;

494:   if (baij->donotstash) {
495:     return(0);
496:   }

498:   /* make sure all processors are either in INSERTMODE or ADDMODE */
499:   MPI_Allreduce(&mat->insertmode,&addv,1,MPI_INT,MPI_BOR,((PetscObject)mat)->comm);
500:   if (addv == (ADD_VALUES|INSERT_VALUES)) {
501:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Some processors inserted others added");
502:   }
503:   mat->insertmode = addv; /* in case this processor had no cache */

505:   MatStashScatterBegin_Private(mat,&mat->stash,mat->rmap->range);
506:   MatStashScatterBegin_Private(mat,&mat->bstash,baij->rangebs);
507:   MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);
508:   PetscInfo2(mat,"Stash has %D entries,uses %D mallocs.\n",nstash,reallocs);
509:   MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);
510:   PetscInfo2(mat,"Block-Stash has %D entries, uses %D mallocs.\n",nstash,reallocs);
511:   return(0);
512: }

516: PetscErrorCode MatAssemblyEnd_MPISBAIJ(Mat mat,MatAssemblyType mode)
517: {
518:   Mat_MPISBAIJ   *baij=(Mat_MPISBAIJ*)mat->data;
519:   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)baij->A->data;
521:   PetscInt       i,j,rstart,ncols,flg,bs2=baij->bs2;
522:   PetscInt       *row,*col;
523:   PetscTruth     other_disassembled;
524:   PetscMPIInt    n;
525:   PetscTruth     r1,r2,r3;
526:   MatScalar      *val;
527:   InsertMode     addv = mat->insertmode;

529:   /* do not use 'b=(Mat_SeqBAIJ*)baij->B->data' as B can be reset in disassembly */

532:   if (!baij->donotstash) {
533:     while (1) {
534:       MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);
535:       if (!flg) break;

537:       for (i=0; i<n;) {
538:         /* Now identify the consecutive vals belonging to the same row */
539:         for (j=i,rstart=row[j]; j<n; j++) { if (row[j] != rstart) break; }
540:         if (j < n) ncols = j-i;
541:         else       ncols = n-i;
542:         /* Now assemble all these values with a single function call */
543:         MatSetValues_MPISBAIJ(mat,1,row+i,ncols,col+i,val+i,addv);
544:         i = j;
545:       }
546:     }
547:     MatStashScatterEnd_Private(&mat->stash);
548:     /* Now process the block-stash. Since the values are stashed column-oriented,
549:        set the roworiented flag to column oriented, and after MatSetValues() 
550:        restore the original flags */
551:     r1 = baij->roworiented;
552:     r2 = a->roworiented;
553:     r3 = ((Mat_SeqBAIJ*)baij->B->data)->roworiented;
554:     baij->roworiented = PETSC_FALSE;
555:     a->roworiented    = PETSC_FALSE;
556:     ((Mat_SeqBAIJ*)baij->B->data)->roworiented    = PETSC_FALSE; /* b->roworinted */
557:     while (1) {
558:       MatStashScatterGetMesg_Private(&mat->bstash,&n,&row,&col,&val,&flg);
559:       if (!flg) break;
560: 
561:       for (i=0; i<n;) {
562:         /* Now identify the consecutive vals belonging to the same row */
563:         for (j=i,rstart=row[j]; j<n; j++) { if (row[j] != rstart) break; }
564:         if (j < n) ncols = j-i;
565:         else       ncols = n-i;
566:         MatSetValuesBlocked_MPISBAIJ(mat,1,row+i,ncols,col+i,val+i*bs2,addv);
567:         i = j;
568:       }
569:     }
570:     MatStashScatterEnd_Private(&mat->bstash);
571:     baij->roworiented = r1;
572:     a->roworiented    = r2;
573:     ((Mat_SeqBAIJ*)baij->B->data)->roworiented    = r3; /* b->roworinted */
574:   }

576:   MatAssemblyBegin(baij->A,mode);
577:   MatAssemblyEnd(baij->A,mode);

579:   /* determine if any processor has disassembled, if so we must 
580:      also disassemble ourselfs, in order that we may reassemble. */
581:   /*
582:      if nonzero structure of submatrix B cannot change then we know that
583:      no processor disassembled thus we can skip this stuff
584:   */
585:   if (!((Mat_SeqBAIJ*)baij->B->data)->nonew)  {
586:     MPI_Allreduce(&mat->was_assembled,&other_disassembled,1,MPI_INT,MPI_PROD,((PetscObject)mat)->comm);
587:     if (mat->was_assembled && !other_disassembled) {
588:       DisAssemble_MPISBAIJ(mat);
589:     }
590:   }

592:   if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
593:     MatSetUpMultiply_MPISBAIJ(mat); /* setup Mvctx and sMvctx */
594:   }
595:   ((Mat_SeqBAIJ*)baij->B->data)->compressedrow.use = PETSC_TRUE; /* b->compressedrow.use */
596:   MatAssemblyBegin(baij->B,mode);
597:   MatAssemblyEnd(baij->B,mode);
598: 
599:   PetscFree2(baij->rowvalues,baij->rowindices);
600:   baij->rowvalues = 0;

602:   return(0);
603: }

608: static PetscErrorCode MatView_MPISBAIJ_ASCIIorDraworSocket(Mat mat,PetscViewer viewer)
609: {
610:   Mat_MPISBAIJ      *baij = (Mat_MPISBAIJ*)mat->data;
611:   PetscErrorCode    ierr;
612:   PetscInt          bs = mat->rmap->bs;
613:   PetscMPIInt       size = baij->size,rank = baij->rank;
614:   PetscTruth        iascii,isdraw;
615:   PetscViewer       sviewer;
616:   PetscViewerFormat format;

619:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
620:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
621:   if (iascii) {
622:     PetscViewerGetFormat(viewer,&format);
623:     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
624:       MatInfo info;
625:       MPI_Comm_rank(((PetscObject)mat)->comm,&rank);
626:       MatGetInfo(mat,MAT_LOCAL,&info);
627:       PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Local rows %D nz %D nz alloced %D bs %D mem %D\n",
628:               rank,mat->rmap->N,(PetscInt)info.nz_used*bs,(PetscInt)info.nz_allocated*bs,
629:               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*bs);
632:       MatGetInfo(baij->B,MAT_LOCAL,&info);
633:       PetscViewerASCIISynchronizedPrintf(viewer,"[%d] off-diagonal part: nz %D \n",rank,(PetscInt)info.nz_used*bs);
634:       PetscViewerFlush(viewer);
635:       PetscViewerASCIIPrintf(viewer,"Information on VecScatter used in matrix-vector product: \n");
636:       VecScatterView(baij->Mvctx,viewer);
637:       return(0);
638:     } else if (format == PETSC_VIEWER_ASCII_INFO) {
639:       PetscViewerASCIIPrintf(viewer,"  block size is %D\n",bs);
640:       return(0);
641:     } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
642:       return(0);
643:     }
644:   }

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

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

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

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

680:     for (i=0; i<mbs; i++) {
681:       rvals[0] = bs*(baij->rstartbs + i);
682:       for (j=1; j<bs; j++) { rvals[j] = rvals[j-1] + 1; }
683:       for (j=ai[i]; j<ai[i+1]; j++) {
684:         col = (baij->cstartbs+aj[j])*bs;
685:         for (k=0; k<bs; k++) {
686:           MatSetValues_MPISBAIJ(A,bs,rvals,1,&col,a,INSERT_VALUES);
687:           col++; 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++) {
695: 
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++; a += bs;
703:         }
704:       }
705:     }
706:     PetscFree(rvals);
707:     MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
708:     MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
709:     /* 
710:        Everyone has to call to draw the matrix since the graphics waits are
711:        synchronized across all processors that share the PetscDraw object
712:     */
713:     PetscViewerGetSingleton(viewer,&sviewer);
714:     if (!rank) {
715:       PetscObjectSetName((PetscObject)((Mat_MPISBAIJ*)(A->data))->A,((PetscObject)mat)->name);
716:       MatView(((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:   PetscTruth     iascii,isdraw,issocket,isbinary;

732:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
733:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
734:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);
735:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);
736:   if (iascii || isdraw || issocket || isbinary) {
737:     MatView_MPISBAIJ_ASCIIorDraworSocket(mat,viewer);
738:   } else {
739:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by MPISBAIJ matrices",((PetscObject)viewer)->type_name);
740:   }
741:   return(0);
742: }

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

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

791:   PetscObjectChangeTypeName((PetscObject)mat,0);
792:   PetscObjectComposeFunction((PetscObject)mat,"MatStoreValues_C","",PETSC_NULL);
793:   PetscObjectComposeFunction((PetscObject)mat,"MatRetrieveValues_C","",PETSC_NULL);
794:   PetscObjectComposeFunction((PetscObject)mat,"MatGetDiagonalBlock_C","",PETSC_NULL);
795:   PetscObjectComposeFunction((PetscObject)mat,"MatMPISBAIJSetPreallocation_C","",PETSC_NULL);
796:   return(0);
797: }

801: PetscErrorCode MatMult_MPISBAIJ_Hermitian(Mat A,Vec xx,Vec yy)
802: {
803:   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;
805:   PetscInt       nt,mbs=a->mbs,bs=A->rmap->bs;
806:   PetscScalar    *x,*from;
807: 
809:   VecGetLocalSize(xx,&nt);
810:   if (nt != A->cmap->n) {
811:     SETERRQ(PETSC_ERR_ARG_SIZ,"Incompatible partition of A and xx");
812:   }

814:   /* diagonal part */
815:   (*a->A->ops->mult)(a->A,xx,a->slvec1a);
816:   VecSet(a->slvec1b,0.0);

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

821:   /* copy x into the vec slvec0 */
822:   VecGetArray(a->slvec0,&from);
823:   VecGetArray(xx,&x);

825:   PetscMemcpy(from,x,bs*mbs*sizeof(MatScalar));
826:   VecRestoreArray(a->slvec0,&from);
827:   VecRestoreArray(xx,&x);
828: 
829:   VecScatterBegin(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
830:   VecScatterEnd(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
831:   /* supperdiagonal part */
832:   (*a->B->ops->multadd)(a->B,a->slvec1b,a->slvec1a,yy);
833:   return(0);
834: }

838: PetscErrorCode MatMult_MPISBAIJ(Mat A,Vec xx,Vec yy)
839: {
840:   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;
842:   PetscInt       nt,mbs=a->mbs,bs=A->rmap->bs;
843:   PetscScalar    *x,*from;
844: 
846:   VecGetLocalSize(xx,&nt);
847:   if (nt != A->cmap->n) {
848:     SETERRQ(PETSC_ERR_ARG_SIZ,"Incompatible partition of A and xx");
849:   }

851:   /* diagonal part */
852:   (*a->A->ops->mult)(a->A,xx,a->slvec1a);
853:   VecSet(a->slvec1b,0.0);

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

858:   /* copy x into the vec slvec0 */
859:   VecGetArray(a->slvec0,&from);
860:   VecGetArray(xx,&x);

862:   PetscMemcpy(from,x,bs*mbs*sizeof(MatScalar));
863:   VecRestoreArray(a->slvec0,&from);
864:   VecRestoreArray(xx,&x);
865: 
866:   VecScatterBegin(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
867:   VecScatterEnd(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
868:   /* supperdiagonal part */
869:   (*a->B->ops->multadd)(a->B,a->slvec1b,a->slvec1a,yy);
870:   return(0);
871: }

875: PetscErrorCode MatMult_MPISBAIJ_2comm(Mat A,Vec xx,Vec yy)
876: {
877:   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;
879:   PetscInt       nt;

882:   VecGetLocalSize(xx,&nt);
883:   if (nt != A->cmap->n) {
884:     SETERRQ(PETSC_ERR_ARG_SIZ,"Incompatible partition of A and xx");
885:   }
886:   VecGetLocalSize(yy,&nt);
887:   if (nt != A->rmap->N) {
888:     SETERRQ(PETSC_ERR_ARG_SIZ,"Incompatible parition of A and yy");
889:   }

891:   VecScatterBegin(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);
892:   /* do diagonal part */
893:   (*a->A->ops->mult)(a->A,xx,yy);
894:   /* do supperdiagonal part */
895:   VecScatterEnd(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);
896:   (*a->B->ops->multadd)(a->B,a->lvec,yy,yy);
897:   /* do subdiagonal part */
898:   (*a->B->ops->multtranspose)(a->B,xx,a->lvec);
899:   VecScatterBegin(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);
900:   VecScatterEnd(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);

902:   return(0);
903: }

907: PetscErrorCode MatMultAdd_MPISBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
908: {
909:   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;
911:   PetscInt       mbs=a->mbs,bs=A->rmap->bs;
912:   PetscScalar    *x,*from,zero=0.0;
913: 
915:   /*
916:   PetscSynchronizedPrintf(((PetscObject)A)->comm," MatMultAdd is called ...\n");
917:   PetscSynchronizedFlush(((PetscObject)A)->comm);
918:   */
919:   /* diagonal part */
920:   (*a->A->ops->multadd)(a->A,xx,yy,a->slvec1a);
921:   VecSet(a->slvec1b,zero);

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

926:   /* copy x into the vec slvec0 */
927:   VecGetArray(a->slvec0,&from);
928:   VecGetArray(xx,&x);
929:   PetscMemcpy(from,x,bs*mbs*sizeof(MatScalar));
930:   VecRestoreArray(a->slvec0,&from);
931: 
932:   VecScatterBegin(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
933:   VecRestoreArray(xx,&x);
934:   VecScatterEnd(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);
935: 
936:   /* supperdiagonal part */
937:   (*a->B->ops->multadd)(a->B,a->slvec1b,a->slvec1a,zz);
938: 
939:   return(0);
940: }

944: PetscErrorCode MatMultAdd_MPISBAIJ_2comm(Mat A,Vec xx,Vec yy,Vec zz)
945: {
946:   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;

950:   VecScatterBegin(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);
951:   /* do diagonal part */
952:   (*a->A->ops->multadd)(a->A,xx,yy,zz);
953:   /* do supperdiagonal part */
954:   VecScatterEnd(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);
955:   (*a->B->ops->multadd)(a->B,a->lvec,zz,zz);

957:   /* do subdiagonal part */
958:   (*a->B->ops->multtranspose)(a->B,xx,a->lvec);
959:   VecScatterBegin(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);
960:   VecScatterEnd(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);

962:   return(0);
963: }

965: /*
966:   This only works correctly for square matrices where the subblock A->A is the 
967:    diagonal block
968: */
971: PetscErrorCode MatGetDiagonal_MPISBAIJ(Mat A,Vec v)
972: {
973:   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;

977:   /* if (a->rmap->N != a->cmap->N) SETERRQ(PETSC_ERR_SUP,"Supports only square matrix where A->A is diag block"); */
978:   MatGetDiagonal(a->A,v);
979:   return(0);
980: }

984: PetscErrorCode MatScale_MPISBAIJ(Mat A,PetscScalar aa)
985: {
986:   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;

990:   MatScale(a->A,aa);
991:   MatScale(a->B,aa);
992:   return(0);
993: }

997: PetscErrorCode MatGetRow_MPISBAIJ(Mat matin,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
998: {
999:   Mat_MPISBAIJ   *mat = (Mat_MPISBAIJ*)matin->data;
1000:   PetscScalar    *vworkA,*vworkB,**pvA,**pvB,*v_p;
1002:   PetscInt       bs = matin->rmap->bs,bs2 = mat->bs2,i,*cworkA,*cworkB,**pcA,**pcB;
1003:   PetscInt       nztot,nzA,nzB,lrow,brstart = matin->rmap->rstart,brend = matin->rmap->rend;
1004:   PetscInt       *cmap,*idx_p,cstart = mat->rstartbs;

1007:   if (mat->getrowactive) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Already active");
1008:   mat->getrowactive = PETSC_TRUE;

1010:   if (!mat->rowvalues && (idx || v)) {
1011:     /*
1012:         allocate enough space to hold information from the longest row.
1013:     */
1014:     Mat_SeqSBAIJ *Aa = (Mat_SeqSBAIJ*)mat->A->data;
1015:     Mat_SeqBAIJ  *Ba = (Mat_SeqBAIJ*)mat->B->data;
1016:     PetscInt     max = 1,mbs = mat->mbs,tmp;
1017:     for (i=0; i<mbs; i++) {
1018:       tmp = Aa->i[i+1] - Aa->i[i] + Ba->i[i+1] - Ba->i[i]; /* row length */
1019:       if (max < tmp) { max = tmp; }
1020:     }
1021:     PetscMalloc2(max*bs2,PetscScalar,&mat->rowvalues,max*bs2,PetscInt,&mat->rowindices);
1022:   }
1023: 
1024:   if (row < brstart || row >= brend) SETERRQ(PETSC_ERR_SUP,"Only local rows")
1025:   lrow = row - brstart;  /* local row index */

1027:   pvA = &vworkA; pcA = &cworkA; pvB = &vworkB; pcB = &cworkB;
1028:   if (!v)   {pvA = 0; pvB = 0;}
1029:   if (!idx) {pcA = 0; if (!v) pcB = 0;}
1030:   (*mat->A->ops->getrow)(mat->A,lrow,&nzA,pcA,pvA);
1031:   (*mat->B->ops->getrow)(mat->B,lrow,&nzB,pcB,pvB);
1032:   nztot = nzA + nzB;

1034:   cmap  = mat->garray;
1035:   if (v  || idx) {
1036:     if (nztot) {
1037:       /* Sort by increasing column numbers, assuming A and B already sorted */
1038:       PetscInt imark = -1;
1039:       if (v) {
1040:         *v = v_p = mat->rowvalues;
1041:         for (i=0; i<nzB; i++) {
1042:           if (cmap[cworkB[i]/bs] < cstart)   v_p[i] = vworkB[i];
1043:           else break;
1044:         }
1045:         imark = i;
1046:         for (i=0; i<nzA; i++)     v_p[imark+i] = vworkA[i];
1047:         for (i=imark; i<nzB; i++) v_p[nzA+i]   = vworkB[i];
1048:       }
1049:       if (idx) {
1050:         *idx = idx_p = mat->rowindices;
1051:         if (imark > -1) {
1052:           for (i=0; i<imark; i++) {
1053:             idx_p[i] = cmap[cworkB[i]/bs]*bs + cworkB[i]%bs;
1054:           }
1055:         } else {
1056:           for (i=0; i<nzB; i++) {
1057:             if (cmap[cworkB[i]/bs] < cstart)
1058:               idx_p[i] = cmap[cworkB[i]/bs]*bs + cworkB[i]%bs ;
1059:             else break;
1060:           }
1061:           imark = i;
1062:         }
1063:         for (i=0; i<nzA; i++)     idx_p[imark+i] = cstart*bs + cworkA[i];
1064:         for (i=imark; i<nzB; i++) idx_p[nzA+i]   = cmap[cworkB[i]/bs]*bs + cworkB[i]%bs ;
1065:       }
1066:     } else {
1067:       if (idx) *idx = 0;
1068:       if (v)   *v   = 0;
1069:     }
1070:   }
1071:   *nz = nztot;
1072:   (*mat->A->ops->restorerow)(mat->A,lrow,&nzA,pcA,pvA);
1073:   (*mat->B->ops->restorerow)(mat->B,lrow,&nzB,pcB,pvB);
1074:   return(0);
1075: }

1079: PetscErrorCode MatRestoreRow_MPISBAIJ(Mat mat,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
1080: {
1081:   Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;

1084:   if (!baij->getrowactive) {
1085:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"MatGetRow() must be called first");
1086:   }
1087:   baij->getrowactive = PETSC_FALSE;
1088:   return(0);
1089: }

1093: PetscErrorCode MatGetRowUpperTriangular_MPISBAIJ(Mat A)
1094: {
1095:   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;
1096:   Mat_SeqSBAIJ   *aA = (Mat_SeqSBAIJ*)a->A->data;

1099:   aA->getrow_utriangular = PETSC_TRUE;
1100:   return(0);
1101: }
1104: PetscErrorCode MatRestoreRowUpperTriangular_MPISBAIJ(Mat A)
1105: {
1106:   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;
1107:   Mat_SeqSBAIJ   *aA = (Mat_SeqSBAIJ*)a->A->data;

1110:   aA->getrow_utriangular = PETSC_FALSE;
1111:   return(0);
1112: }

1116: PetscErrorCode MatRealPart_MPISBAIJ(Mat A)
1117: {
1118:   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;

1122:   MatRealPart(a->A);
1123:   MatRealPart(a->B);
1124:   return(0);
1125: }

1129: PetscErrorCode MatImaginaryPart_MPISBAIJ(Mat A)
1130: {
1131:   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;

1135:   MatImaginaryPart(a->A);
1136:   MatImaginaryPart(a->B);
1137:   return(0);
1138: }

1142: PetscErrorCode MatZeroEntries_MPISBAIJ(Mat A)
1143: {
1144:   Mat_MPISBAIJ   *l = (Mat_MPISBAIJ*)A->data;

1148:   MatZeroEntries(l->A);
1149:   MatZeroEntries(l->B);
1150:   return(0);
1151: }

1155: PetscErrorCode MatGetInfo_MPISBAIJ(Mat matin,MatInfoType flag,MatInfo *info)
1156: {
1157:   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)matin->data;
1158:   Mat            A = a->A,B = a->B;
1160:   PetscReal      isend[5],irecv[5];

1163:   info->block_size     = (PetscReal)matin->rmap->bs;
1164:   MatGetInfo(A,MAT_LOCAL,info);
1165:   isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
1166:   isend[3] = info->memory;  isend[4] = info->mallocs;
1167:   MatGetInfo(B,MAT_LOCAL,info);
1168:   isend[0] += info->nz_used; isend[1] += info->nz_allocated; isend[2] += info->nz_unneeded;
1169:   isend[3] += info->memory;  isend[4] += info->mallocs;
1170:   if (flag == MAT_LOCAL) {
1171:     info->nz_used      = isend[0];
1172:     info->nz_allocated = isend[1];
1173:     info->nz_unneeded  = isend[2];
1174:     info->memory       = isend[3];
1175:     info->mallocs      = isend[4];
1176:   } else if (flag == MAT_GLOBAL_MAX) {
1177:     MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_MAX,((PetscObject)matin)->comm);
1178:     info->nz_used      = irecv[0];
1179:     info->nz_allocated = irecv[1];
1180:     info->nz_unneeded  = irecv[2];
1181:     info->memory       = irecv[3];
1182:     info->mallocs      = irecv[4];
1183:   } else if (flag == MAT_GLOBAL_SUM) {
1184:     MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPI_SUM,((PetscObject)matin)->comm);
1185:     info->nz_used      = irecv[0];
1186:     info->nz_allocated = irecv[1];
1187:     info->nz_unneeded  = irecv[2];
1188:     info->memory       = irecv[3];
1189:     info->mallocs      = irecv[4];
1190:   } else {
1191:     SETERRQ1(PETSC_ERR_ARG_WRONG,"Unknown MatInfoType argument %d",(int)flag);
1192:   }
1193:   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
1194:   info->fill_ratio_needed = 0;
1195:   info->factor_mallocs    = 0;
1196:   return(0);
1197: }

1201: PetscErrorCode MatSetOption_MPISBAIJ(Mat A,MatOption op,PetscTruth flg)
1202: {
1203:   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;
1204:   Mat_SeqSBAIJ   *aA = (Mat_SeqSBAIJ*)a->A->data;

1208:   switch (op) {
1209:   case MAT_NEW_NONZERO_LOCATIONS:
1210:   case MAT_NEW_NONZERO_ALLOCATION_ERR:
1211:   case MAT_UNUSED_NONZERO_LOCATION_ERR:
1212:   case MAT_KEEP_NONZERO_PATTERN:
1213:   case MAT_NEW_NONZERO_LOCATION_ERR:
1214:     MatSetOption(a->A,op,flg);
1215:     MatSetOption(a->B,op,flg);
1216:     break;
1217:   case MAT_ROW_ORIENTED:
1218:     a->roworiented = flg;
1219:     MatSetOption(a->A,op,flg);
1220:     MatSetOption(a->B,op,flg);
1221:     break;
1222:   case MAT_NEW_DIAGONALS:
1223:     PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
1224:     break;
1225:   case MAT_IGNORE_OFF_PROC_ENTRIES:
1226:     a->donotstash = flg;
1227:     break;
1228:   case MAT_USE_HASH_TABLE:
1229:     a->ht_flag = flg;
1230:     break;
1231:   case MAT_HERMITIAN:
1232:     if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call MatAssemblyEnd() first");
1233:     MatSetOption(a->A,op,flg);
1234:     A->ops->mult = MatMult_MPISBAIJ_Hermitian;
1235:     break;
1236:   case MAT_SYMMETRIC:
1237:     MatSetOption(a->A,op,flg);
1238:     break;
1239:   case MAT_STRUCTURALLY_SYMMETRIC:
1240:     MatSetOption(a->A,op,flg);
1241:     break;
1242:   case MAT_SYMMETRY_ETERNAL:
1243:     if (!flg) SETERRQ(PETSC_ERR_SUP,"Matrix must be symmetric");
1244:     PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);
1245:     break;
1246:   case MAT_IGNORE_LOWER_TRIANGULAR:
1247:     aA->ignore_ltriangular = flg;
1248:     break;
1249:   case MAT_ERROR_LOWER_TRIANGULAR:
1250:     aA->ignore_ltriangular = flg;
1251:     break;
1252:   case MAT_GETROW_UPPERTRIANGULAR:
1253:     aA->getrow_utriangular = flg;
1254:     break;
1255:   default:
1256:     SETERRQ1(PETSC_ERR_SUP,"unknown option %d",op);
1257:   }
1258:   return(0);
1259: }

1263: PetscErrorCode MatTranspose_MPISBAIJ(Mat A,MatReuse reuse,Mat *B)
1264: {
1267:   if (MAT_INITIAL_MATRIX || *B != A) {
1268:     MatDuplicate(A,MAT_COPY_VALUES,B);
1269:   }
1270:   return(0);
1271: }

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

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

1291:   MatGetLocalSize(mat,&m,&n);
1292:   if (m != n) SETERRQ2(PETSC_ERR_ARG_SIZ,"For symmetric format, local size %d %d must be same",m,n);
1293: 
1294:   VecGetLocalSize(rr,&nv);
1295:   if (nv!=n) SETERRQ(PETSC_ERR_ARG_SIZ,"Left and right vector non-conforming local size");

1297:   VecScatterBegin(baij->Mvctx,rr,baij->lvec,INSERT_VALUES,SCATTER_FORWARD);
1298: 
1299:   /* left diagonalscale the off-diagonal part */
1300:   (*b->ops->diagonalscale)(b,ll,PETSC_NULL);
1301: 
1302:   /* scale the diagonal part */
1303:   (*a->ops->diagonalscale)(a,ll,rr);

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

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

1319:   MatSetUnfactored(a->A);
1320:   return(0);
1321: }

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

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

1335:   a = matA->A; b = matA->B;
1336:   c = matB->A; d = matB->B;

1338:   MatEqual(a,c,&flg);
1339:   if (flg) {
1340:     MatEqual(b,d,&flg);
1341:   }
1342:   MPI_Allreduce(&flg,flag,1,MPI_INT,MPI_LAND,((PetscObject)A)->comm);
1343:   return(0);
1344: }

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

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

1369: PetscErrorCode MatSetUpPreallocation_MPISBAIJ(Mat A)
1370: {

1374:   MatMPISBAIJSetPreallocation(A,-PetscMax(A->rmap->bs,1),PETSC_DEFAULT,0,PETSC_DEFAULT,0);
1375:   return(0);
1376: }

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

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

1409: PetscErrorCode MatSetBlockSize_MPISBAIJ(Mat A,PetscInt bs)
1410: {
1411:   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;
1412:   PetscInt        rbs,cbs;
1413:   PetscErrorCode  ierr;

1416:   MatSetBlockSize(a->A,bs);
1417:   MatSetBlockSize(a->B,bs);
1418:   PetscLayoutGetBlockSize(A->rmap,&rbs);
1419:   PetscLayoutGetBlockSize(A->cmap,&cbs);
1420:   if (rbs != bs) SETERRQ2(PETSC_ERR_ARG_SIZ,"Attempt to set block size %d with SBAIJ %d",bs,rbs);
1421:   if (cbs != bs) SETERRQ2(PETSC_ERR_ARG_SIZ,"Attempt to set block size %d with SBAIJ %d",bs,cbs);
1422:   return(0);
1423: }

1427: PetscErrorCode MatGetSubMatrices_MPISBAIJ(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1428: {
1430:   PetscInt       i;
1431:   PetscTruth     flg;

1434:   for (i=0; i<n; i++) {
1435:     ISEqual(irow[i],icol[i],&flg);
1436:     if (!flg) {
1437:       SETERRQ(PETSC_ERR_SUP,"Can only get symmetric submatrix for MPISBAIJ matrices");
1438:     }
1439:   }
1440:   MatGetSubMatrices_MPIBAIJ(A,n,irow,icol,scall,B);
1441:   return(0);
1442: }
1443: 

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


1562: PetscErrorCode  MatGetDiagonalBlock_MPISBAIJ(Mat A,PetscTruth *iscopy,MatReuse reuse,Mat *a)
1563: {
1565:   *a      = ((Mat_MPISBAIJ *)A->data)->A;
1566:   *iscopy = PETSC_FALSE;
1567:   return(0);
1568: }

1574: PetscErrorCode  MatMPISBAIJSetPreallocation_MPISBAIJ(Mat B,PetscInt bs,PetscInt d_nz,PetscInt *d_nnz,PetscInt o_nz,PetscInt *o_nnz)
1575: {
1576:   Mat_MPISBAIJ   *b;
1578:   PetscInt       i,mbs,Mbs,newbs = PetscAbs(bs);

1581:   if (bs < 0){
1582:     PetscOptionsBegin(((PetscObject)B)->comm,((PetscObject)B)->prefix,"Options for MPISBAIJ matrix","Mat");
1583:       PetscOptionsInt("-mat_block_size","Set the blocksize used to store the matrix","MatMPIBAIJSetPreallocation",newbs,&newbs,PETSC_NULL);
1584:     PetscOptionsEnd();
1585:     bs   = PetscAbs(bs);
1586:   }
1587:   if ((d_nnz || o_nnz) && newbs != bs) {
1588:     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot change blocksize from command line if setting d_nnz or o_nnz");
1589:   }
1590:   bs = newbs;

1592:   if (d_nz == PETSC_DECIDE || d_nz == PETSC_DEFAULT) d_nz = 3;
1593:   if (o_nz == PETSC_DECIDE || o_nz == PETSC_DEFAULT) o_nz = 1;
1594:   if (d_nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"d_nz cannot be less than 0: value %D",d_nz);
1595:   if (o_nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"o_nz cannot be less than 0: value %D",o_nz);

1597:   B->rmap->bs = B->cmap->bs = bs;
1598:   PetscLayoutSetUp(B->rmap);
1599:   PetscLayoutSetUp(B->cmap);

1601:   if (d_nnz) {
1602:     for (i=0; i<B->rmap->n/bs; i++) {
1603:       if (d_nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"d_nnz cannot be less than -1: local row %D value %D",i,d_nnz[i]);
1604:     }
1605:   }
1606:   if (o_nnz) {
1607:     for (i=0; i<B->rmap->n/bs; i++) {
1608:       if (o_nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"o_nnz cannot be less than -1: local row %D value %D",i,o_nnz[i]);
1609:     }
1610:   }

1612:   b   = (Mat_MPISBAIJ*)B->data;
1613:   mbs = B->rmap->n/bs;
1614:   Mbs = B->rmap->N/bs;
1615:   if (mbs*bs != B->rmap->n) {
1616:     SETERRQ2(PETSC_ERR_ARG_SIZ,"No of local rows %D must be divisible by blocksize %D",B->rmap->N,bs);
1617:   }

1619:   B->rmap->bs  = bs;
1620:   b->bs2 = bs*bs;
1621:   b->mbs = mbs;
1622:   b->nbs = mbs;
1623:   b->Mbs = Mbs;
1624:   b->Nbs = Mbs;

1626:   for (i=0; i<=b->size; i++) {
1627:     b->rangebs[i] = B->rmap->range[i]/bs;
1628:   }
1629:   b->rstartbs = B->rmap->rstart/bs;
1630:   b->rendbs   = B->rmap->rend/bs;
1631: 
1632:   b->cstartbs = B->cmap->rstart/bs;
1633:   b->cendbs   = B->cmap->rend/bs;

1635:   if (!B->preallocated) {
1636:     MatCreate(PETSC_COMM_SELF,&b->A);
1637:     MatSetSizes(b->A,B->rmap->n,B->cmap->n,B->rmap->n,B->cmap->n);
1638:     MatSetType(b->A,MATSEQSBAIJ);
1639:     PetscLogObjectParent(B,b->A);
1640:     MatCreate(PETSC_COMM_SELF,&b->B);
1641:     MatSetSizes(b->B,B->rmap->n,B->cmap->N,B->rmap->n,B->cmap->N);
1642:     MatSetType(b->B,MATSEQBAIJ);
1643:     PetscLogObjectParent(B,b->B);
1644:     MatStashCreate_Private(((PetscObject)B)->comm,bs,&B->bstash);
1645:   }

1647:   MatSeqSBAIJSetPreallocation(b->A,bs,d_nz,d_nnz);
1648:   MatSeqBAIJSetPreallocation(b->B,bs,o_nz,o_nnz);
1649:   B->preallocated = PETSC_TRUE;
1650:   return(0);
1651: }

1655: #if defined(PETSC_HAVE_MUMPS)
1657: #endif
1658: #if defined(PETSC_HAVE_SPOOLES)
1660: #endif
1661: #if defined(PETSC_HAVE_PASTIX)
1663: #endif

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

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

1673:   Level: beginner

1675: .seealso: MatCreateMPISBAIJ
1676: M*/

1681: PetscErrorCode  MatCreate_MPISBAIJ(Mat B)
1682: {
1683:   Mat_MPISBAIJ   *b;
1685:   PetscTruth     flg;


1689:   PetscNewLog(B,Mat_MPISBAIJ,&b);
1690:   B->data = (void*)b;
1691:   PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));

1693:   B->ops->destroy    = MatDestroy_MPISBAIJ;
1694:   B->ops->view       = MatView_MPISBAIJ;
1695:   B->mapping         = 0;
1696:   B->assembled       = PETSC_FALSE;

1698:   B->insertmode = NOT_SET_VALUES;
1699:   MPI_Comm_rank(((PetscObject)B)->comm,&b->rank);
1700:   MPI_Comm_size(((PetscObject)B)->comm,&b->size);

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

1705:   /* build cache for off array entries formed */
1706:   MatStashCreate_Private(((PetscObject)B)->comm,1,&B->stash);
1707:   b->donotstash  = PETSC_FALSE;
1708:   b->colmap      = PETSC_NULL;
1709:   b->garray      = PETSC_NULL;
1710:   b->roworiented = PETSC_TRUE;

1712:   /* stuff used in block assembly */
1713:   b->barray       = 0;

1715:   /* stuff used for matrix vector multiply */
1716:   b->lvec         = 0;
1717:   b->Mvctx        = 0;
1718:   b->slvec0       = 0;
1719:   b->slvec0b      = 0;
1720:   b->slvec1       = 0;
1721:   b->slvec1a      = 0;
1722:   b->slvec1b      = 0;
1723:   b->sMvctx       = 0;

1725:   /* stuff for MatGetRow() */
1726:   b->rowindices   = 0;
1727:   b->rowvalues    = 0;
1728:   b->getrowactive = PETSC_FALSE;

1730:   /* hash table stuff */
1731:   b->ht           = 0;
1732:   b->hd           = 0;
1733:   b->ht_size      = 0;
1734:   b->ht_flag      = PETSC_FALSE;
1735:   b->ht_fact      = 0;
1736:   b->ht_total_ct  = 0;
1737:   b->ht_insert_ct = 0;

1739:   b->in_loc       = 0;
1740:   b->v_loc        = 0;
1741:   b->n_loc        = 0;
1742:   PetscOptionsBegin(((PetscObject)B)->comm,PETSC_NULL,"Options for loading MPISBAIJ matrix 1","Mat");
1743:     PetscOptionsTruth("-mat_use_hash_table","Use hash table to save memory in constructing matrix","MatSetOption",PETSC_FALSE,&flg,PETSC_NULL);
1744:     if (flg) {
1745:       PetscReal fact = 1.39;
1746:       MatSetOption(B,MAT_USE_HASH_TABLE,PETSC_TRUE);
1747:       PetscOptionsReal("-mat_use_hash_table","Use hash table factor","MatMPIBAIJSetHashTableFactor",fact,&fact,PETSC_NULL);
1748:       if (fact <= 1.0) fact = 1.39;
1749:       MatMPIBAIJSetHashTableFactor(B,fact);
1750:       PetscInfo1(B,"Hash table Factor used %5.2f\n",fact);
1751:     }
1752:   PetscOptionsEnd();

1754: #if defined(PETSC_HAVE_PASTIX)
1755:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_pastix_C",
1756:                                            "MatGetFactor_mpisbaij_pastix",
1757:                                            MatGetFactor_mpisbaij_pastix);
1758: #endif
1759: #if defined(PETSC_HAVE_MUMPS)
1760:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_mumps_C",
1761:                                      "MatGetFactor_mpisbaij_mumps",
1762:                                      MatGetFactor_mpisbaij_mumps);
1763: #endif
1764: #if defined(PETSC_HAVE_SPOOLES)
1765:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_spooles_C",
1766:                                      "MatGetFactor_mpisbaij_spooles",
1767:                                      MatGetFactor_mpisbaij_spooles);
1768: #endif
1769:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C",
1770:                                      "MatStoreValues_MPISBAIJ",
1771:                                      MatStoreValues_MPISBAIJ);
1772:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C",
1773:                                      "MatRetrieveValues_MPISBAIJ",
1774:                                      MatRetrieveValues_MPISBAIJ);
1775:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetDiagonalBlock_C",
1776:                                      "MatGetDiagonalBlock_MPISBAIJ",
1777:                                      MatGetDiagonalBlock_MPISBAIJ);
1778:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPISBAIJSetPreallocation_C",
1779:                                      "MatMPISBAIJSetPreallocation_MPISBAIJ",
1780:                                      MatMPISBAIJSetPreallocation_MPISBAIJ);
1781:   B->symmetric                  = PETSC_TRUE;
1782:   B->structurally_symmetric     = PETSC_TRUE;
1783:   B->symmetric_set              = PETSC_TRUE;
1784:   B->structurally_symmetric_set = PETSC_TRUE;
1785:   PetscObjectChangeTypeName((PetscObject)B,MATMPISBAIJ);
1786:   return(0);
1787: }

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

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

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

1799:   Level: beginner

1801: .seealso: MatCreateMPISBAIJ,MATSEQSBAIJ,MATMPISBAIJ
1802: M*/

1807: PetscErrorCode  MatCreate_SBAIJ(Mat A)
1808: {
1810:   PetscMPIInt    size;

1813:   MPI_Comm_size(((PetscObject)A)->comm,&size);
1814:   if (size == 1) {
1815:     MatSetType(A,MATSEQSBAIJ);
1816:   } else {
1817:     MatSetType(A,MATMPISBAIJ);
1818:   }
1819:   return(0);
1820: }

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

1831:    Collective on Mat

1833:    Input Parameters:
1834: +  A - the matrix 
1835: .  bs   - size of blockk
1836: .  d_nz  - number of block nonzeros per block row in diagonal portion of local 
1837:            submatrix  (same for all local rows)
1838: .  d_nnz - array containing the number of block nonzeros in the various block rows 
1839:            in the upper triangular and diagonal part of the in diagonal portion of the local
1840:            (possibly different for each block row) or PETSC_NULL.  You must leave room 
1841:            for the diagonal entry even if it is zero.
1842: .  o_nz  - number of block nonzeros per block row in the off-diagonal portion of local
1843:            submatrix (same for all local rows).
1844: -  o_nnz - array containing the number of nonzeros in the various block rows of the
1845:            off-diagonal portion of the local submatrix (possibly different for
1846:            each block row) or PETSC_NULL.


1849:    Options Database Keys:
1850: .   -mat_no_unroll - uses code that does not unroll the loops in the 
1851:                      block calculations (much slower)
1852: .   -mat_block_size - size of the blocks to use

1854:    Notes:

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

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

1861:    Storage Information:
1862:    For a square global matrix we define each processor's diagonal portion 
1863:    to be its local rows and the corresponding columns (a square submatrix);  
1864:    each processor's off-diagonal portion encompasses the remainder of the
1865:    local matrix (a rectangular submatrix). 

1867:    The user can specify preallocated storage for the diagonal part of
1868:    the local submatrix with either d_nz or d_nnz (not both).  Set 
1869:    d_nz=PETSC_DEFAULT and d_nnz=PETSC_NULL for PETSc to control dynamic
1870:    memory allocation.  Likewise, specify preallocated storage for the
1871:    off-diagonal part of the local submatrix with o_nz or o_nnz (not both).

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

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

1881: .vb
1882:            0 1 2 3 4 5 6 7 8 9 10 11
1883:           -------------------
1884:    row 3  |  o o o d d d o o o o o o
1885:    row 4  |  o o o d d d o o o o o o
1886:    row 5  |  o o o d d d o o o o o o
1887:           -------------------
1888: .ve
1889:   
1890:    Thus, any entries in the d locations are stored in the d (diagonal) 
1891:    submatrix, and any entries in the o locations are stored in the
1892:    o (off-diagonal) submatrix.  Note that the d matrix is stored in
1893:    MatSeqSBAIJ format and the o submatrix in MATSEQBAIJ format.

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

1903:    Level: intermediate

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

1907: .seealso: MatCreate(), MatCreateSeqSBAIJ(), MatSetValues(), MatCreateMPIBAIJ()
1908: @*/
1909: PetscErrorCode  MatMPISBAIJSetPreallocation(Mat B,PetscInt bs,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
1910: {
1911:   PetscErrorCode ierr,(*f)(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);

1914:   PetscObjectQueryFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",(void (**)(void))&f);
1915:   if (f) {
1916:     (*f)(B,bs,d_nz,d_nnz,o_nz,o_nnz);
1917:   }
1918:   return(0);
1919: }

1923: /*@C
1924:    MatCreateMPISBAIJ - Creates a sparse parallel matrix in symmetric block AIJ format
1925:    (block compressed row).  For good matrix assembly performance
1926:    the user should preallocate the matrix storage by setting the parameters 
1927:    d_nz (or d_nnz) and o_nz (or o_nnz).  By setting these parameters accurately,
1928:    performance can be increased by more than a factor of 50.

1930:    Collective on MPI_Comm

1932:    Input Parameters:
1933: +  comm - MPI communicator
1934: .  bs   - size of blockk
1935: .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
1936:            This value should be the same as the local size used in creating the 
1937:            y vector for the matrix-vector product y = Ax.
1938: .  n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
1939:            This value should be the same as the local size used in creating the 
1940:            x vector for the matrix-vector product y = Ax.
1941: .  M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
1942: .  N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
1943: .  d_nz  - number of block nonzeros per block row in diagonal portion of local 
1944:            submatrix  (same for all local rows)
1945: .  d_nnz - array containing the number of block nonzeros in the various block rows 
1946:            in the upper triangular portion of the in diagonal portion of the local 
1947:            (possibly different for each block block row) or PETSC_NULL.  
1948:            You must leave room for the diagonal entry even if it is zero.
1949: .  o_nz  - number of block nonzeros per block row in the off-diagonal portion of local
1950:            submatrix (same for all local rows).
1951: -  o_nnz - array containing the number of nonzeros in the various block rows of the
1952:            off-diagonal portion of the local submatrix (possibly different for
1953:            each block row) or PETSC_NULL.

1955:    Output Parameter:
1956: .  A - the matrix 

1958:    Options Database Keys:
1959: .   -mat_no_unroll - uses code that does not unroll the loops in the 
1960:                      block calculations (much slower)
1961: .   -mat_block_size - size of the blocks to use
1962: .   -mat_mpi - use the parallel matrix data structures even on one processor 
1963:                (defaults to using SeqBAIJ format on one processor)

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

1969:    Notes:
1970:    The number of rows and columns must be divisible by blocksize.
1971:    This matrix type does not support complex Hermitian operation.

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

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

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

1981:    Storage Information:
1982:    For a square global matrix we define each processor's diagonal portion 
1983:    to be its local rows and the corresponding columns (a square submatrix);  
1984:    each processor's off-diagonal portion encompasses the remainder of the
1985:    local matrix (a rectangular submatrix). 

1987:    The user can specify preallocated storage for the diagonal part of
1988:    the local submatrix with either d_nz or d_nnz (not both).  Set 
1989:    d_nz=PETSC_DEFAULT and d_nnz=PETSC_NULL for PETSc to control dynamic
1990:    memory allocation.  Likewise, specify preallocated storage for the
1991:    off-diagonal part of the local submatrix with o_nz or o_nnz (not both).

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

1996: .vb
1997:            0 1 2 3 4 5 6 7 8 9 10 11
1998:           -------------------
1999:    row 3  |  o o o d d d o o o o o o
2000:    row 4  |  o o o d d d o o o o o o
2001:    row 5  |  o o o d d d o o o o o o
2002:           -------------------
2003: .ve
2004:   
2005:    Thus, any entries in the d locations are stored in the d (diagonal) 
2006:    submatrix, and any entries in the o locations are stored in the
2007:    o (off-diagonal) submatrix.  Note that the d matrix is stored in
2008:    MatSeqSBAIJ format and the o submatrix in MATSEQBAIJ format.

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

2018:    Level: intermediate

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

2022: .seealso: MatCreate(), MatCreateSeqSBAIJ(), MatSetValues(), MatCreateMPIBAIJ()
2023: @*/

2025: PetscErrorCode  MatCreateMPISBAIJ(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)
2026: {
2028:   PetscMPIInt    size;

2031:   MatCreate(comm,A);
2032:   MatSetSizes(*A,m,n,M,N);
2033:   MPI_Comm_size(comm,&size);
2034:   if (size > 1) {
2035:     MatSetType(*A,MATMPISBAIJ);
2036:     MatMPISBAIJSetPreallocation(*A,bs,d_nz,d_nnz,o_nz,o_nnz);
2037:   } else {
2038:     MatSetType(*A,MATSEQSBAIJ);
2039:     MatSeqSBAIJSetPreallocation(*A,bs,d_nz,d_nnz);
2040:   }
2041:   return(0);
2042: }


2047: static PetscErrorCode MatDuplicate_MPISBAIJ(Mat matin,MatDuplicateOption cpvalues,Mat *newmat)
2048: {
2049:   Mat            mat;
2050:   Mat_MPISBAIJ   *a,*oldmat = (Mat_MPISBAIJ*)matin->data;
2052:   PetscInt       len=0,nt,bs=matin->rmap->bs,mbs=oldmat->mbs;
2053:   PetscScalar    *array;

2056:   *newmat       = 0;
2057:   MatCreate(((PetscObject)matin)->comm,&mat);
2058:   MatSetSizes(mat,matin->rmap->n,matin->cmap->n,matin->rmap->N,matin->cmap->N);
2059:   MatSetType(mat,((PetscObject)matin)->type_name);
2060:   PetscMemcpy(mat->ops,matin->ops,sizeof(struct _MatOps));
2061:   PetscLayoutCopy(matin->rmap,&mat->rmap);
2062:   PetscLayoutCopy(matin->cmap,&mat->cmap);
2063: 
2064:   mat->factor       = matin->factor;
2065:   mat->preallocated = PETSC_TRUE;
2066:   mat->assembled    = PETSC_TRUE;
2067:   mat->insertmode   = NOT_SET_VALUES;

2069:   a = (Mat_MPISBAIJ*)mat->data;
2070:   a->bs2   = oldmat->bs2;
2071:   a->mbs   = oldmat->mbs;
2072:   a->nbs   = oldmat->nbs;
2073:   a->Mbs   = oldmat->Mbs;
2074:   a->Nbs   = oldmat->Nbs;


2077:   a->size         = oldmat->size;
2078:   a->rank         = oldmat->rank;
2079:   a->donotstash   = oldmat->donotstash;
2080:   a->roworiented  = oldmat->roworiented;
2081:   a->rowindices   = 0;
2082:   a->rowvalues    = 0;
2083:   a->getrowactive = PETSC_FALSE;
2084:   a->barray       = 0;
2085:   a->rstartbs    = oldmat->rstartbs;
2086:   a->rendbs      = oldmat->rendbs;
2087:   a->cstartbs    = oldmat->cstartbs;
2088:   a->cendbs      = oldmat->cendbs;

2090:   /* hash table stuff */
2091:   a->ht           = 0;
2092:   a->hd           = 0;
2093:   a->ht_size      = 0;
2094:   a->ht_flag      = oldmat->ht_flag;
2095:   a->ht_fact      = oldmat->ht_fact;
2096:   a->ht_total_ct  = 0;
2097:   a->ht_insert_ct = 0;
2098: 
2099:   PetscMemcpy(a->rangebs,oldmat->rangebs,(a->size+2)*sizeof(PetscInt));
2100:   if (oldmat->colmap) {
2101: #if defined (PETSC_USE_CTABLE)
2102:     PetscTableCreateCopy(oldmat->colmap,&a->colmap);
2103: #else
2104:     PetscMalloc((a->Nbs)*sizeof(PetscInt),&a->colmap);
2105:     PetscLogObjectMemory(mat,(a->Nbs)*sizeof(PetscInt));
2106:     PetscMemcpy(a->colmap,oldmat->colmap,(a->Nbs)*sizeof(PetscInt));
2107: #endif
2108:   } else a->colmap = 0;

2110:   if (oldmat->garray && (len = ((Mat_SeqBAIJ*)(oldmat->B->data))->nbs)) {
2111:     PetscMalloc(len*sizeof(PetscInt),&a->garray);
2112:     PetscLogObjectMemory(mat,len*sizeof(PetscInt));
2113:     PetscMemcpy(a->garray,oldmat->garray,len*sizeof(PetscInt));
2114:   } else a->garray = 0;
2115: 
2116:   MatStashCreate_Private(((PetscObject)matin)->comm,matin->rmap->bs,&mat->bstash);
2117:   VecDuplicate(oldmat->lvec,&a->lvec);
2118:   PetscLogObjectParent(mat,a->lvec);
2119:   VecScatterCopy(oldmat->Mvctx,&a->Mvctx);
2120:   PetscLogObjectParent(mat,a->Mvctx);

2122:    VecDuplicate(oldmat->slvec0,&a->slvec0);
2123:   PetscLogObjectParent(mat,a->slvec0);
2124:    VecDuplicate(oldmat->slvec1,&a->slvec1);
2125:   PetscLogObjectParent(mat,a->slvec1);

2127:   VecGetLocalSize(a->slvec1,&nt);
2128:   VecGetArray(a->slvec1,&array);
2129:   VecCreateSeqWithArray(PETSC_COMM_SELF,bs*mbs,array,&a->slvec1a);
2130:   VecCreateSeqWithArray(PETSC_COMM_SELF,nt-bs*mbs,array+bs*mbs,&a->slvec1b);
2131:   VecRestoreArray(a->slvec1,&array);
2132:   VecGetArray(a->slvec0,&array);
2133:   VecCreateSeqWithArray(PETSC_COMM_SELF,nt-bs*mbs,array+bs*mbs,&a->slvec0b);
2134:   VecRestoreArray(a->slvec0,&array);
2135:   PetscLogObjectParent(mat,a->slvec0);
2136:   PetscLogObjectParent(mat,a->slvec1);
2137:   PetscLogObjectParent(mat,a->slvec0b);
2138:   PetscLogObjectParent(mat,a->slvec1a);
2139:   PetscLogObjectParent(mat,a->slvec1b);

2141:   /*  VecScatterCopy(oldmat->sMvctx,&a->sMvctx); - not written yet, replaced by the lazy trick: */
2142:   PetscObjectReference((PetscObject)oldmat->sMvctx);
2143:   a->sMvctx = oldmat->sMvctx;
2144:   PetscLogObjectParent(mat,a->sMvctx);

2146:    MatDuplicate(oldmat->A,cpvalues,&a->A);
2147:   PetscLogObjectParent(mat,a->A);
2148:    MatDuplicate(oldmat->B,cpvalues,&a->B);
2149:   PetscLogObjectParent(mat,a->B);
2150:   PetscFListDuplicate(((PetscObject)matin)->qlist,&((PetscObject)mat)->qlist);
2151:   *newmat = mat;
2152:   return(0);
2153: }

2157: PetscErrorCode MatLoad_MPISBAIJ(PetscViewer viewer, const MatType type,Mat *newmat)
2158: {
2159:   Mat            A;
2161:   PetscInt       i,nz,j,rstart,rend;
2162:   PetscScalar    *vals,*buf;
2163:   MPI_Comm       comm = ((PetscObject)viewer)->comm;
2164:   MPI_Status     status;
2165:   PetscMPIInt    rank,size,tag = ((PetscObject)viewer)->tag,*sndcounts = 0,*browners,maxnz,*rowners,*locrowlens,mmbs;
2166:   PetscInt       header[4],*rowlengths = 0,M,N,m,*cols;
2167:   PetscInt       *procsnz = 0,jj,*mycols,*ibuf;
2168:   PetscInt       bs=1,Mbs,mbs,extra_rows;
2169:   PetscInt       *dlens,*odlens,*mask,*masked1,*masked2,rowcount,odcount;
2170:   PetscInt       dcount,kmax,k,nzcount,tmp;
2171:   int            fd;
2172: 
2174:   PetscOptionsBegin(comm,PETSC_NULL,"Options for loading MPISBAIJ matrix 2","Mat");
2175:     PetscOptionsInt("-matload_block_size","Set the blocksize used to store the matrix","MatLoad",bs,&bs,PETSC_NULL);
2176:   PetscOptionsEnd();

2178:   MPI_Comm_size(comm,&size);
2179:   MPI_Comm_rank(comm,&rank);
2180:   if (!rank) {
2181:     PetscViewerBinaryGetDescriptor(viewer,&fd);
2182:     PetscBinaryRead(fd,(char *)header,4,PETSC_INT);
2183:     if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object");
2184:     if (header[3] < 0) {
2185:       SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as MPISBAIJ");
2186:     }
2187:   }

2189:   MPI_Bcast(header+1,3,MPIU_INT,0,comm);
2190:   M = header[1]; N = header[2];

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

2194:   /* 
2195:      This code adds extra rows to make sure the number of rows is 
2196:      divisible by the blocksize
2197:   */
2198:   Mbs        = M/bs;
2199:   extra_rows = bs - M + bs*(Mbs);
2200:   if (extra_rows == bs) extra_rows = 0;
2201:   else                  Mbs++;
2202:   if (extra_rows &&!rank) {
2203:     PetscInfo(viewer,"Padding loaded matrix to match blocksize\n");
2204:   }

2206:   /* determine ownership of all rows */
2207:   mbs        = Mbs/size + ((Mbs % size) > rank);
2208:   m          = mbs*bs;
2209:   PetscMalloc2(size+1,PetscMPIInt,&rowners,size+1,PetscMPIInt,&browners);
2210:   mmbs       = PetscMPIIntCast(mbs);
2211:   MPI_Allgather(&mmbs,1,MPI_INT,rowners+1,1,MPI_INT,comm);
2212:   rowners[0] = 0;
2213:   for (i=2; i<=size; i++) rowners[i] += rowners[i-1];
2214:   for (i=0; i<=size;  i++) browners[i] = rowners[i]*bs;
2215:   rstart = rowners[rank];
2216:   rend   = rowners[rank+1];
2217: 
2218:   /* distribute row lengths to all processors */
2219:   PetscMalloc((rend-rstart)*bs*sizeof(PetscMPIInt),&locrowlens);
2220:   if (!rank) {
2221:     PetscMalloc((M+extra_rows)*sizeof(PetscInt),&rowlengths);
2222:     PetscBinaryRead(fd,rowlengths,M,PETSC_INT);
2223:     for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;
2224:     PetscMalloc(size*sizeof(PetscMPIInt),&sndcounts);
2225:     for (i=0; i<size; i++) sndcounts[i] = browners[i+1] - browners[i];
2226:     MPI_Scatterv(rowlengths,sndcounts,browners,MPIU_INT,locrowlens,(rend-rstart)*bs,MPIU_INT,0,comm);
2227:     PetscFree(sndcounts);
2228:   } else {
2229:     MPI_Scatterv(0,0,0,MPIU_INT,locrowlens,(rend-rstart)*bs,MPIU_INT,0,comm);
2230:   }
2231: 
2232:   if (!rank) {   /* procs[0] */
2233:     /* calculate the number of nonzeros on each processor */
2234:     PetscMalloc(size*sizeof(PetscInt),&procsnz);
2235:     PetscMemzero(procsnz,size*sizeof(PetscInt));
2236:     for (i=0; i<size; i++) {
2237:       for (j=rowners[i]*bs; j< rowners[i+1]*bs; j++) {
2238:         procsnz[i] += rowlengths[j];
2239:       }
2240:     }
2241:     PetscFree(rowlengths);
2242: 
2243:     /* determine max buffer needed and allocate it */
2244:     maxnz = 0;
2245:     for (i=0; i<size; i++) {
2246:       maxnz = PetscMax(maxnz,procsnz[i]);
2247:     }
2248:     PetscMalloc(maxnz*sizeof(PetscInt),&cols);

2250:     /* read in my part of the matrix column indices  */
2251:     nz     = procsnz[0];
2252:     PetscMalloc(nz*sizeof(PetscInt),&ibuf);
2253:     mycols = ibuf;
2254:     if (size == 1)  nz -= extra_rows;
2255:     PetscBinaryRead(fd,mycols,nz,PETSC_INT);
2256:     if (size == 1)  for (i=0; i< extra_rows; i++) { mycols[nz+i] = M+i; }

2258:     /* read in every ones (except the last) and ship off */
2259:     for (i=1; i<size-1; i++) {
2260:       nz   = procsnz[i];
2261:       PetscBinaryRead(fd,cols,nz,PETSC_INT);
2262:       MPI_Send(cols,nz,MPIU_INT,i,tag,comm);
2263:     }
2264:     /* read in the stuff for the last proc */
2265:     if (size != 1) {
2266:       nz   = procsnz[size-1] - extra_rows;  /* the extra rows are not on the disk */
2267:       PetscBinaryRead(fd,cols,nz,PETSC_INT);
2268:       for (i=0; i<extra_rows; i++) cols[nz+i] = M+i;
2269:       MPI_Send(cols,nz+extra_rows,MPIU_INT,size-1,tag,comm);
2270:     }
2271:     PetscFree(cols);
2272:   } else {  /* procs[i], i>0 */
2273:     /* determine buffer space needed for message */
2274:     nz = 0;
2275:     for (i=0; i<m; i++) {
2276:       nz += locrowlens[i];
2277:     }
2278:     PetscMalloc(nz*sizeof(PetscInt),&ibuf);
2279:     mycols = ibuf;
2280:     /* receive message of column indices*/
2281:     MPI_Recv(mycols,nz,MPIU_INT,0,tag,comm,&status);
2282:     MPI_Get_count(&status,MPIU_INT,&maxnz);
2283:     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
2284:   }

2286:   /* loop over local rows, determining number of off diagonal entries */
2287:   PetscMalloc2(rend-rstart,PetscInt,&dlens,rend-rstart,PetscInt,&odlens);
2288:   PetscMalloc3(Mbs,PetscInt,&mask,Mbs,PetscInt,&masked1,Mbs,PetscInt,&masked2);
2289:   PetscMemzero(mask,Mbs*sizeof(PetscInt));
2290:   PetscMemzero(masked1,Mbs*sizeof(PetscInt));
2291:   PetscMemzero(masked2,Mbs*sizeof(PetscInt));
2292:   rowcount = 0;
2293:   nzcount  = 0;
2294:   for (i=0; i<mbs; i++) {
2295:     dcount  = 0;
2296:     odcount = 0;
2297:     for (j=0; j<bs; j++) {
2298:       kmax = locrowlens[rowcount];
2299:       for (k=0; k<kmax; k++) {
2300:         tmp = mycols[nzcount++]/bs; /* block col. index */
2301:         if (!mask[tmp]) {
2302:           mask[tmp] = 1;
2303:           if (tmp < rstart || tmp >= rend) masked2[odcount++] = tmp; /* entry in off-diag portion */
2304:           else masked1[dcount++] = tmp; /* entry in diag portion */
2305:         }
2306:       }
2307:       rowcount++;
2308:     }
2309: 
2310:     dlens[i]  = dcount;  /* d_nzz[i] */
2311:     odlens[i] = odcount; /* o_nzz[i] */

2313:     /* zero out the mask elements we set */
2314:     for (j=0; j<dcount; j++) mask[masked1[j]] = 0;
2315:     for (j=0; j<odcount; j++) mask[masked2[j]] = 0;
2316:   }
2317: 
2318:   /* create our matrix */
2319:   MatCreate(comm,&A);
2320:   MatSetSizes(A,m,m,PETSC_DETERMINE,PETSC_DETERMINE);
2321:   MatSetType(A,type);
2322:   MatSetOption(A,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);
2323:   MatMPISBAIJSetPreallocation(A,bs,0,dlens,0,odlens);
2324: 
2325:   if (!rank) {
2326:     PetscMalloc(maxnz*sizeof(PetscScalar),&buf);
2327:     /* read in my part of the matrix numerical values  */
2328:     nz = procsnz[0];
2329:     vals = buf;
2330:     mycols = ibuf;
2331:     if (size == 1)  nz -= extra_rows;
2332:     PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
2333:     if (size == 1)  for (i=0; i< extra_rows; i++) { vals[nz+i] = 1.0; }

2335:     /* insert into matrix */
2336:     jj      = rstart*bs;
2337:     for (i=0; i<m; i++) {
2338:       MatSetValues(A,1,&jj,locrowlens[i],mycols,vals,INSERT_VALUES);
2339:       mycols += locrowlens[i];
2340:       vals   += locrowlens[i];
2341:       jj++;
2342:     }

2344:     /* read in other processors (except the last one) and ship out */
2345:     for (i=1; i<size-1; i++) {
2346:       nz   = procsnz[i];
2347:       vals = buf;
2348:       PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
2349:       MPI_Send(vals,nz,MPIU_SCALAR,i,((PetscObject)A)->tag,comm);
2350:     }
2351:     /* the last proc */
2352:     if (size != 1){
2353:       nz   = procsnz[i] - extra_rows;
2354:       vals = buf;
2355:       PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);
2356:       for (i=0; i<extra_rows; i++) vals[nz+i] = 1.0;
2357:       MPI_Send(vals,nz+extra_rows,MPIU_SCALAR,size-1,((PetscObject)A)->tag,comm);
2358:     }
2359:     PetscFree(procsnz);

2361:   } else {
2362:     /* receive numeric values */
2363:     PetscMalloc(nz*sizeof(PetscScalar),&buf);

2365:     /* receive message of values*/
2366:     vals   = buf;
2367:     mycols = ibuf;
2368:     MPI_Recv(vals,nz,MPIU_SCALAR,0,((PetscObject)A)->tag,comm,&status);
2369:     MPI_Get_count(&status,MPIU_SCALAR,&maxnz);
2370:     if (maxnz != nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");

2372:     /* insert into matrix */
2373:     jj      = rstart*bs;
2374:     for (i=0; i<m; i++) {
2375:       MatSetValues_MPISBAIJ(A,1,&jj,locrowlens[i],mycols,vals,INSERT_VALUES);
2376:       mycols += locrowlens[i];
2377:       vals   += locrowlens[i];
2378:       jj++;
2379:     }
2380:   }

2382:   PetscFree(locrowlens);
2383:   PetscFree(buf);
2384:   PetscFree(ibuf);
2385:   PetscFree2(rowners,browners);
2386:   PetscFree2(dlens,odlens);
2387:   PetscFree3(mask,masked1,masked2);
2388:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
2389:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
2390:   *newmat = A;
2391:   return(0);
2392: }

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

2399:    Input Parameters:
2400: .  mat  - the matrix
2401: .  fact - factor

2403:    Collective on Mat

2405:    Level: advanced

2407:   Notes:
2408:    This can also be set by the command line option: -mat_use_hash_table fact

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

2412: .seealso: MatSetOption()
2413: @XXXXX*/


2418: PetscErrorCode MatGetRowMaxAbs_MPISBAIJ(Mat A,Vec v,PetscInt idx[])
2419: {
2420:   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;
2421:   Mat_SeqBAIJ    *b = (Mat_SeqBAIJ*)(a->B)->data;
2422:   PetscReal      atmp;
2423:   PetscReal      *work,*svalues,*rvalues;
2425:   PetscInt       i,bs,mbs,*bi,*bj,brow,j,ncols,krow,kcol,col,row,Mbs,bcol;
2426:   PetscMPIInt    rank,size;
2427:   PetscInt       *rowners_bs,dest,count,source;
2428:   PetscScalar    *va;
2429:   MatScalar      *ba;
2430:   MPI_Status     stat;

2433:   if (idx) SETERRQ(PETSC_ERR_SUP,"Send email to petsc-maint@mcs.anl.gov");
2434:   MatGetRowMaxAbs(a->A,v,PETSC_NULL);
2435:   VecGetArray(v,&va);

2437:   MPI_Comm_size(((PetscObject)A)->comm,&size);
2438:   MPI_Comm_rank(((PetscObject)A)->comm,&rank);

2440:   bs   = A->rmap->bs;
2441:   mbs  = a->mbs;
2442:   Mbs  = a->Mbs;
2443:   ba   = b->a;
2444:   bi   = b->i;
2445:   bj   = b->j;

2447:   /* find ownerships */
2448:   rowners_bs = A->rmap->range;

2450:   /* each proc creates an array to be distributed */
2451:   PetscMalloc(bs*Mbs*sizeof(PetscReal),&work);
2452:   PetscMemzero(work,bs*Mbs*sizeof(PetscReal));

2454:   /* row_max for B */
2455:   if (rank != size-1){
2456:     for (i=0; i<mbs; i++) {
2457:       ncols = bi[1] - bi[0]; bi++;
2458:       brow  = bs*i;
2459:       for (j=0; j<ncols; j++){
2460:         bcol = bs*(*bj);
2461:         for (kcol=0; kcol<bs; kcol++){
2462:           col = bcol + kcol;                 /* local col index */
2463:           col += rowners_bs[rank+1];      /* global col index */
2464:           for (krow=0; krow<bs; krow++){
2465:             atmp = PetscAbsScalar(*ba); ba++;
2466:             row = brow + krow;    /* local row index */
2467:             if (PetscRealPart(va[row]) < atmp) va[row] = atmp;
2468:             if (work[col] < atmp) work[col] = atmp;
2469:           }
2470:         }
2471:         bj++;
2472:       }
2473:     }

2475:     /* send values to its owners */
2476:     for (dest=rank+1; dest<size; dest++){
2477:       svalues = work + rowners_bs[dest];
2478:       count   = rowners_bs[dest+1]-rowners_bs[dest];
2479:       MPI_Send(svalues,count,MPIU_REAL,dest,rank,((PetscObject)A)->comm);
2480:     }
2481:   }
2482: 
2483:   /* receive values */
2484:   if (rank){
2485:     rvalues = work;
2486:     count   = rowners_bs[rank+1]-rowners_bs[rank];
2487:     for (source=0; source<rank; source++){
2488:       MPI_Recv(rvalues,count,MPIU_REAL,MPI_ANY_SOURCE,MPI_ANY_TAG,((PetscObject)A)->comm,&stat);
2489:       /* process values */
2490:       for (i=0; i<count; i++){
2491:         if (PetscRealPart(va[i]) < rvalues[i]) va[i] = rvalues[i];
2492:       }
2493:     }
2494:   }

2496:   VecRestoreArray(v,&va);
2497:   PetscFree(work);
2498:   return(0);
2499: }

2503: PetscErrorCode MatSOR_MPISBAIJ(Mat matin,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
2504: {
2505:   Mat_MPISBAIJ   *mat = (Mat_MPISBAIJ*)matin->data;
2507:   PetscInt       mbs=mat->mbs,bs=matin->rmap->bs;
2508:   PetscScalar    *x,*b,*ptr,*from;
2509:   Vec            bb1;
2510: 
2512:   if (its <= 0 || lits <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
2513:   if (bs > 1) SETERRQ(PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented");

2515:   if (flag == SOR_APPLY_UPPER) {
2516:     (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx);
2517:     return(0);
2518:   }

2520:   if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP){
2521:     if ( flag & SOR_ZERO_INITIAL_GUESS ) {
2522:       (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,lits,xx);
2523:       its--;
2524:     }

2526:     VecDuplicate(bb,&bb1);
2527:     while (its--){
2528: 
2529:       /* lower triangular part: slvec0b = - B^T*xx */
2530:       (*mat->B->ops->multtranspose)(mat->B,xx,mat->slvec0b);
2531: 
2532:       /* copy xx into slvec0a */
2533:       VecGetArray(mat->slvec0,&ptr);
2534:       VecGetArray(xx,&x);
2535:       PetscMemcpy(ptr,x,bs*mbs*sizeof(MatScalar));
2536:       VecRestoreArray(mat->slvec0,&ptr);

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

2540:       /* copy bb into slvec1a */
2541:       VecGetArray(mat->slvec1,&ptr);
2542:       VecGetArray(bb,&b);
2543:       PetscMemcpy(ptr,b,bs*mbs*sizeof(MatScalar));
2544:       VecRestoreArray(mat->slvec1,&ptr);

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

2549:       VecScatterBegin(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD);
2550:       VecRestoreArray(xx,&x);
2551:       VecRestoreArray(bb,&b);
2552:       VecScatterEnd(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD);

2554:       /* upper triangular part: bb1 = bb1 - B*x */
2555:       (*mat->B->ops->multadd)(mat->B,mat->slvec1b,mat->slvec1a,bb1);
2556: 
2557:       /* local diagonal sweep */
2558:       (*mat->A->ops->sor)(mat->A,bb1,omega,SOR_SYMMETRIC_SWEEP,fshift,lits,lits,xx);
2559:     }
2560:     VecDestroy(bb1);
2561:   } else if ((flag & SOR_LOCAL_FORWARD_SWEEP) && (its == 1) && (flag & SOR_ZERO_INITIAL_GUESS)){
2562:     (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx);
2563:   } else if ((flag & SOR_LOCAL_BACKWARD_SWEEP) && (its == 1) && (flag & SOR_ZERO_INITIAL_GUESS)){
2564:     (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx);
2565:   } else if (flag & SOR_EISENSTAT) {
2566:     Vec               xx1;
2567:     PetscTruth        hasop;
2568:     const PetscScalar *diag;
2569:     PetscScalar       *sl,scale = (omega - 2.0)/omega;
2570:     PetscInt          i,n;

2572:     if (!mat->xx1) {
2573:       VecDuplicate(bb,&mat->xx1);
2574:       VecDuplicate(bb,&mat->bb1);
2575:     }
2576:     xx1 = mat->xx1;
2577:     bb1 = mat->bb1;

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

2581:     if (!mat->diag) {
2582:       /* this is wrong for same matrix with new nonzero values */
2583:       MatGetVecs(matin,&mat->diag,PETSC_NULL);
2584:       MatGetDiagonal(matin,mat->diag);
2585:     }
2586:     MatHasOperation(matin,MATOP_MULT_DIAGONAL_BLOCK,&hasop);

2588:     if (hasop) {
2589:       MatMultDiagonalBlock(matin,xx,bb1);
2590:       VecAYPX(mat->slvec1a,scale,bb);
2591:     } else {
2592:       /*
2593:           These two lines are replaced by code that may be a bit faster for a good compiler
2594:       VecPointwiseMult(mat->slvec1a,mat->diag,xx);
2595:       VecAYPX(mat->slvec1a,scale,bb);
2596:       */
2597:       VecGetArray(mat->slvec1a,&sl);
2598:       VecGetArray(mat->diag,(PetscScalar**)&diag);
2599:       VecGetArray(bb,&b);
2600:       VecGetArray(xx,&x);
2601:       VecGetLocalSize(xx,&n);
2602:       if (omega == 1.0) {
2603:         for (i=0; i<n; i++) {
2604:           sl[i] = b[i] - diag[i]*x[i];
2605:         }
2606:         PetscLogFlops(2.0*n);
2607:       } else {
2608:         for (i=0; i<n; i++) {
2609:           sl[i] = b[i] + scale*diag[i]*x[i];
2610:         }
2611:         PetscLogFlops(3.0*n);
2612:       }
2613:       VecRestoreArray(mat->slvec1a,&sl);
2614:       VecRestoreArray(mat->diag,(PetscScalar**)&diag);
2615:       VecRestoreArray(bb,&b);
2616:       VecRestoreArray(xx,&x);
2617:     }

2619:     /* multiply off-diagonal portion of matrix */
2620:     VecSet(mat->slvec1b,0.0);
2621:     (*mat->B->ops->multtranspose)(mat->B,xx,mat->slvec0b);
2622:     VecGetArray(mat->slvec0,&from);
2623:     VecGetArray(xx,&x);
2624:     PetscMemcpy(from,x,bs*mbs*sizeof(MatScalar));
2625:     VecRestoreArray(mat->slvec0,&from);
2626:     VecRestoreArray(xx,&x);
2627:     VecScatterBegin(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD);
2628:     VecScatterEnd(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD);
2629:     (*mat->B->ops->multadd)(mat->B,mat->slvec1b,mat->slvec1a,mat->slvec1a);

2631:     /* local sweep */
2632:     (*mat->A->ops->sor)(mat->A,mat->slvec1a,omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_FORWARD_SWEEP),fshift,lits,1,xx1);
2633:     VecAXPY(xx,1.0,xx1);
2634:   } else {
2635:     SETERRQ(PETSC_ERR_SUP,"MatSORType is not supported for SBAIJ matrix format");
2636:   }
2637:   return(0);
2638: }

2642: PetscErrorCode MatSOR_MPISBAIJ_2comm(Mat matin,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
2643: {
2644:   Mat_MPISBAIJ   *mat = (Mat_MPISBAIJ*)matin->data;
2646:   Vec            lvec1,bb1;
2647: 
2649:   if (its <= 0 || lits <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
2650:   if (matin->rmap->bs > 1) SETERRQ(PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented");

2652:   if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP){
2653:     if ( flag & SOR_ZERO_INITIAL_GUESS ) {
2654:       (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,lits,xx);
2655:       its--;
2656:     }

2658:     VecDuplicate(mat->lvec,&lvec1);
2659:     VecDuplicate(bb,&bb1);
2660:     while (its--){
2661:       VecScatterBegin(mat->Mvctx,xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD);
2662: 
2663:       /* lower diagonal part: bb1 = bb - B^T*xx */
2664:       (*mat->B->ops->multtranspose)(mat->B,xx,lvec1);
2665:       VecScale(lvec1,-1.0);

2667:       VecScatterEnd(mat->Mvctx,xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD);
2668:       VecCopy(bb,bb1);
2669:       VecScatterBegin(mat->Mvctx,lvec1,bb1,ADD_VALUES,SCATTER_REVERSE);

2671:       /* upper diagonal part: bb1 = bb1 - B*x */
2672:       VecScale(mat->lvec,-1.0);
2673:       (*mat->B->ops->multadd)(mat->B,mat->lvec,bb1,bb1);

2675:       VecScatterEnd(mat->Mvctx,lvec1,bb1,ADD_VALUES,SCATTER_REVERSE);
2676: 
2677:       /* diagonal sweep */
2678:       (*mat->A->ops->sor)(mat->A,bb1,omega,SOR_SYMMETRIC_SWEEP,fshift,lits,lits,xx);
2679:     }
2680:     VecDestroy(lvec1);
2681:     VecDestroy(bb1);
2682:   } else {
2683:     SETERRQ(PETSC_ERR_SUP,"MatSORType is not supported for SBAIJ matrix format");
2684:   }
2685:   return(0);
2686: }