xref: /petsc/src/mat/impls/sbaij/mpi/mpisbaij.c (revision a57dd577347cf357bc3142e2987c5ffed089475f)
1 
2 #include <../src/mat/impls/baij/mpi/mpibaij.h>    /*I "petscmat.h" I*/
3 #include <../src/mat/impls/sbaij/mpi/mpisbaij.h>
4 #include <../src/mat/impls/sbaij/seq/sbaij.h>
5 #include <petscblaslapack.h>
6 
7 extern PetscErrorCode MatSetUpMultiply_MPISBAIJ(Mat);
8 extern PetscErrorCode MatSetUpMultiply_MPISBAIJ_2comm(Mat);
9 extern PetscErrorCode MatDisAssemble_MPISBAIJ(Mat);
10 extern PetscErrorCode MatIncreaseOverlap_MPISBAIJ(Mat,PetscInt,IS[],PetscInt);
11 extern PetscErrorCode MatGetValues_SeqSBAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar []);
12 extern PetscErrorCode MatGetValues_SeqBAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar []);
13 extern PetscErrorCode MatSetValues_SeqSBAIJ(Mat,PetscInt,const PetscInt [],PetscInt,const PetscInt [],const PetscScalar [],InsertMode);
14 extern PetscErrorCode MatSetValuesBlocked_SeqSBAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
15 extern PetscErrorCode MatSetValuesBlocked_SeqBAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
16 extern PetscErrorCode MatGetRow_SeqSBAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
17 extern PetscErrorCode MatRestoreRow_SeqSBAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
18 extern PetscErrorCode MatZeroRows_SeqSBAIJ(Mat,IS,PetscScalar*,Vec,Vec);
19 extern PetscErrorCode MatZeroRows_SeqBAIJ(Mat,IS,PetscScalar*,Vec,Vec);
20 extern PetscErrorCode MatGetRowMaxAbs_MPISBAIJ(Mat,Vec,PetscInt[]);
21 extern PetscErrorCode MatSOR_MPISBAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
22 
23 EXTERN_C_BEGIN
24 #undef __FUNCT__
25 #define __FUNCT__ "MatStoreValues_MPISBAIJ"
26 PetscErrorCode  MatStoreValues_MPISBAIJ(Mat mat)
27 {
28   Mat_MPISBAIJ   *aij = (Mat_MPISBAIJ*)mat->data;
29   PetscErrorCode ierr;
30 
31   PetscFunctionBegin;
32   ierr = MatStoreValues(aij->A);CHKERRQ(ierr);
33   ierr = MatStoreValues(aij->B);CHKERRQ(ierr);
34   PetscFunctionReturn(0);
35 }
36 EXTERN_C_END
37 
38 EXTERN_C_BEGIN
39 #undef __FUNCT__
40 #define __FUNCT__ "MatRetrieveValues_MPISBAIJ"
41 PetscErrorCode  MatRetrieveValues_MPISBAIJ(Mat mat)
42 {
43   Mat_MPISBAIJ   *aij = (Mat_MPISBAIJ*)mat->data;
44   PetscErrorCode ierr;
45 
46   PetscFunctionBegin;
47   ierr = MatRetrieveValues(aij->A);CHKERRQ(ierr);
48   ierr = MatRetrieveValues(aij->B);CHKERRQ(ierr);
49   PetscFunctionReturn(0);
50 }
51 EXTERN_C_END
52 
53 
54 #define  MatSetValues_SeqSBAIJ_A_Private(row,col,value,addv) \
55   { \
56  \
57     brow = row/bs;  \
58     rp   = aj + ai[brow]; ap = aa + bs2*ai[brow]; \
59     rmax = aimax[brow]; nrow = ailen[brow]; \
60     bcol = col/bs; \
61     ridx = row % bs; cidx = col % bs; \
62     low  = 0; high = nrow; \
63     while (high-low > 3) { \
64       t = (low+high)/2; \
65       if (rp[t] > bcol) high = t; \
66       else              low  = t; \
67     } \
68     for (_i=low; _i<high; _i++) { \
69       if (rp[_i] > bcol) break; \
70       if (rp[_i] == bcol) { \
71         bap = ap + bs2*_i + bs*cidx + ridx; \
72         if (addv == ADD_VALUES) *bap += value;  \
73         else                    *bap  = value;  \
74         goto a_noinsert; \
75       } \
76     } \
77     if (a->nonew == 1) goto a_noinsert; \
78     if (a->nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) into matrix", row, col); \
79     MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,aimax,a->nonew,MatScalar); \
80     N = nrow++ - 1;  \
81     /* shift up all the later entries in this row */ \
82     for (ii=N; ii>=_i; ii--) { \
83       rp[ii+1] = rp[ii]; \
84       ierr     = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); \
85     } \
86     if (N>=_i) { ierr = PetscMemzero(ap+bs2*_i,bs2*sizeof(MatScalar));CHKERRQ(ierr); }  \
87     rp[_i]                      = bcol;  \
88     ap[bs2*_i + bs*cidx + ridx] = value;  \
89 a_noinsert:; \
90     ailen[brow] = nrow; \
91   }
92 
93 #define  MatSetValues_SeqSBAIJ_B_Private(row,col,value,addv) \
94   { \
95     brow = row/bs;  \
96     rp   = bj + bi[brow]; ap = ba + bs2*bi[brow]; \
97     rmax = bimax[brow]; nrow = bilen[brow]; \
98     bcol = col/bs; \
99     ridx = row % bs; cidx = col % bs; \
100     low  = 0; high = nrow; \
101     while (high-low > 3) { \
102       t = (low+high)/2; \
103       if (rp[t] > bcol) high = t; \
104       else              low  = t; \
105     } \
106     for (_i=low; _i<high; _i++) { \
107       if (rp[_i] > bcol) break; \
108       if (rp[_i] == bcol) { \
109         bap = ap + bs2*_i + bs*cidx + ridx; \
110         if (addv == ADD_VALUES) *bap += value;  \
111         else                    *bap  = value;  \
112         goto b_noinsert; \
113       } \
114     } \
115     if (b->nonew == 1) goto b_noinsert; \
116     if (b->nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) into matrix", row, col); \
117     MatSeqXAIJReallocateAIJ(B,b->mbs,bs2,nrow,brow,bcol,rmax,ba,bi,bj,rp,ap,bimax,b->nonew,MatScalar); \
118     N = nrow++ - 1;  \
119     /* shift up all the later entries in this row */ \
120     for (ii=N; ii>=_i; ii--) { \
121       rp[ii+1] = rp[ii]; \
122       ierr     = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr); \
123     } \
124     if (N>=_i) { ierr = PetscMemzero(ap+bs2*_i,bs2*sizeof(MatScalar));CHKERRQ(ierr);}  \
125     rp[_i]                      = bcol;  \
126     ap[bs2*_i + bs*cidx + ridx] = value;  \
127 b_noinsert:; \
128     bilen[brow] = nrow; \
129   }
130 
131 /* Only add/insert a(i,j) with i<=j (blocks).
132    Any a(i,j) with i>j input by user is ingored.
133 */
134 #undef __FUNCT__
135 #define __FUNCT__ "MatSetValues_MPISBAIJ"
136 PetscErrorCode MatSetValues_MPISBAIJ(Mat mat,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode addv)
137 {
138   Mat_MPISBAIJ   *baij = (Mat_MPISBAIJ*)mat->data;
139   MatScalar      value;
140   PetscBool      roworiented = baij->roworiented;
141   PetscErrorCode ierr;
142   PetscInt       i,j,row,col;
143   PetscInt       rstart_orig=mat->rmap->rstart;
144   PetscInt       rend_orig  =mat->rmap->rend,cstart_orig=mat->cmap->rstart;
145   PetscInt       cend_orig  =mat->cmap->rend,bs=mat->rmap->bs;
146 
147   /* Some Variables required in the macro */
148   Mat          A     = baij->A;
149   Mat_SeqSBAIJ *a    = (Mat_SeqSBAIJ*)(A)->data;
150   PetscInt     *aimax=a->imax,*ai=a->i,*ailen=a->ilen,*aj=a->j;
151   MatScalar    *aa   =a->a;
152 
153   Mat         B     = baij->B;
154   Mat_SeqBAIJ *b    = (Mat_SeqBAIJ*)(B)->data;
155   PetscInt    *bimax=b->imax,*bi=b->i,*bilen=b->ilen,*bj=b->j;
156   MatScalar   *ba   =b->a;
157 
158   PetscInt  *rp,ii,nrow,_i,rmax,N,brow,bcol;
159   PetscInt  low,high,t,ridx,cidx,bs2=a->bs2;
160   MatScalar *ap,*bap;
161 
162   /* for stash */
163   PetscInt  n_loc, *in_loc = NULL;
164   MatScalar *v_loc = NULL;
165 
166   PetscFunctionBegin;
167   if (v) PetscValidScalarPointer(v,6);
168   if (!baij->donotstash) {
169     if (n > baij->n_loc) {
170       ierr = PetscFree(baij->in_loc);CHKERRQ(ierr);
171       ierr = PetscFree(baij->v_loc);CHKERRQ(ierr);
172       ierr = PetscMalloc(n*sizeof(PetscInt),&baij->in_loc);CHKERRQ(ierr);
173       ierr = PetscMalloc(n*sizeof(MatScalar),&baij->v_loc);CHKERRQ(ierr);
174 
175       baij->n_loc = n;
176     }
177     in_loc = baij->in_loc;
178     v_loc  = baij->v_loc;
179   }
180 
181   for (i=0; i<m; i++) {
182     if (im[i] < 0) continue;
183 #if defined(PETSC_USE_DEBUG)
184     if (im[i] >= mat->rmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",im[i],mat->rmap->N-1);
185 #endif
186     if (im[i] >= rstart_orig && im[i] < rend_orig) { /* this processor entry */
187       row = im[i] - rstart_orig;              /* local row index */
188       for (j=0; j<n; j++) {
189         if (im[i]/bs > in[j]/bs) {
190           if (a->ignore_ltriangular) {
191             continue;    /* ignore lower triangular blocks */
192           } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Lower triangular value cannot be set for sbaij format. Ignoring these values, run with -mat_ignore_lower_triangular or call MatSetOption(mat,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE)");
193         }
194         if (in[j] >= cstart_orig && in[j] < cend_orig) {  /* diag entry (A) */
195           col  = in[j] - cstart_orig;         /* local col index */
196           brow = row/bs; bcol = col/bs;
197           if (brow > bcol) continue;  /* ignore lower triangular blocks of A */
198           if (roworiented) value = v[i*n+j];
199           else             value = v[i+j*m];
200           MatSetValues_SeqSBAIJ_A_Private(row,col,value,addv);
201           /* ierr = MatSetValues_SeqBAIJ(baij->A,1,&row,1,&col,&value,addv);CHKERRQ(ierr); */
202         } else if (in[j] < 0) continue;
203 #if defined(PETSC_USE_DEBUG)
204         else if (in[j] >= mat->cmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[j],mat->cmap->N-1);
205 #endif
206         else {  /* off-diag entry (B) */
207           if (mat->was_assembled) {
208             if (!baij->colmap) {
209               ierr = MatCreateColmap_MPIBAIJ_Private(mat);CHKERRQ(ierr);
210             }
211 #if defined(PETSC_USE_CTABLE)
212             ierr = PetscTableFind(baij->colmap,in[j]/bs + 1,&col);CHKERRQ(ierr);
213             col  = col - 1;
214 #else
215             col = baij->colmap[in[j]/bs] - 1;
216 #endif
217             if (col < 0 && !((Mat_SeqSBAIJ*)(baij->A->data))->nonew) {
218               ierr = MatDisAssemble_MPISBAIJ(mat);CHKERRQ(ierr);
219               col  =  in[j];
220               /* Reinitialize the variables required by MatSetValues_SeqBAIJ_B_Private() */
221               B    = baij->B;
222               b    = (Mat_SeqBAIJ*)(B)->data;
223               bimax= b->imax;bi=b->i;bilen=b->ilen;bj=b->j;
224               ba   = b->a;
225             } else col += in[j]%bs;
226           } else col = in[j];
227           if (roworiented) value = v[i*n+j];
228           else             value = v[i+j*m];
229           MatSetValues_SeqSBAIJ_B_Private(row,col,value,addv);
230           /* ierr = MatSetValues_SeqBAIJ(baij->B,1,&row,1,&col,&value,addv);CHKERRQ(ierr); */
231         }
232       }
233     } else {  /* off processor entry */
234       if (mat->nooffprocentries) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Setting off process row %D even though MatSetOption(,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE) was set",im[i]);
235       if (!baij->donotstash) {
236         mat->assembled = PETSC_FALSE;
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         ierr = MatStashValuesRow_Private(&mat->stash,im[i],n_loc,in_loc,v_loc,PETSC_FALSE);CHKERRQ(ierr);
249       }
250     }
251   }
252   PetscFunctionReturn(0);
253 }
254 
255 #undef __FUNCT__
256 #define __FUNCT__ "MatSetValuesBlocked_MPISBAIJ"
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   PetscBool       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;
267 
268   PetscFunctionBegin;
269   if (!barray) {
270     ierr         = PetscMalloc(bs2*sizeof(MatScalar),&barray);CHKERRQ(ierr);
271     baij->barray = barray;
272   }
273 
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_COMM_SELF,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_COMM_SELF,PETSC_ERR_USER,"Lower triangular value cannot be set for sbaij format. Ignoring these values, run with -mat_ignore_lower_triangular or call MatSetOption(mat,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE)");
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           ierr = MatSetValuesBlocked_SeqSBAIJ(baij->A,1,&row,1,&col,barray,addv);CHKERRQ(ierr);
313         } else if (in[j] < 0) continue;
314 #if defined(PETSC_USE_DEBUG)
315         else if (in[j] >= baij->Nbs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large, col %D max %D",in[j],baij->Nbs-1);
316 #endif
317         else {
318           if (mat->was_assembled) {
319             if (!baij->colmap) {
320               ierr = MatCreateColmap_MPIBAIJ_Private(mat);CHKERRQ(ierr);
321             }
322 
323 #if defined(PETSC_USE_DEBUG)
324 #if defined(PETSC_USE_CTABLE)
325             { PetscInt data;
326               ierr = PetscTableFind(baij->colmap,in[j]+1,&data);CHKERRQ(ierr);
327               if ((data - 1) % bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect colmap");
328             }
329 #else
330             if ((baij->colmap[in[j]] - 1) % bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Incorrect colmap");
331 #endif
332 #endif
333 #if defined(PETSC_USE_CTABLE)
334             ierr = PetscTableFind(baij->colmap,in[j]+1,&col);CHKERRQ(ierr);
335             col  = (col - 1)/bs;
336 #else
337             col = (baij->colmap[in[j]] - 1)/bs;
338 #endif
339             if (col < 0 && !((Mat_SeqBAIJ*)(baij->A->data))->nonew) {
340               ierr = MatDisAssemble_MPISBAIJ(mat);CHKERRQ(ierr);
341               col  = in[j];
342             }
343           } else col = in[j];
344           ierr = MatSetValuesBlocked_SeqBAIJ(baij->B,1,&row,1,&col,barray,addv);CHKERRQ(ierr);
345         }
346       }
347     } else {
348       if (mat->nooffprocentries) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Setting off process row %D even though MatSetOption(,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE) was set",im[i]);
349       if (!baij->donotstash) {
350         if (roworiented) {
351           ierr = MatStashValuesRowBlocked_Private(&mat->bstash,im[i],n,in,v,m,n,i);CHKERRQ(ierr);
352         } else {
353           ierr = MatStashValuesColBlocked_Private(&mat->bstash,im[i],n,in,v,m,n,i);CHKERRQ(ierr);
354         }
355       }
356     }
357   }
358   PetscFunctionReturn(0);
359 }
360 
361 #undef __FUNCT__
362 #define __FUNCT__ "MatGetValues_MPISBAIJ"
363 PetscErrorCode MatGetValues_MPISBAIJ(Mat mat,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],PetscScalar v[])
364 {
365   Mat_MPISBAIJ   *baij = (Mat_MPISBAIJ*)mat->data;
366   PetscErrorCode ierr;
367   PetscInt       bs       = mat->rmap->bs,i,j,bsrstart = mat->rmap->rstart,bsrend = mat->rmap->rend;
368   PetscInt       bscstart = mat->cmap->rstart,bscend = mat->cmap->rend,row,col,data;
369 
370   PetscFunctionBegin;
371   for (i=0; i<m; i++) {
372     if (idxm[i] < 0) continue; /* SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",idxm[i]); */
373     if (idxm[i] >= mat->rmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm[i],mat->rmap->N-1);
374     if (idxm[i] >= bsrstart && idxm[i] < bsrend) {
375       row = idxm[i] - bsrstart;
376       for (j=0; j<n; j++) {
377         if (idxn[j] < 0) continue; /* SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column %D",idxn[j]); */
378         if (idxn[j] >= mat->cmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",idxn[j],mat->cmap->N-1);
379         if (idxn[j] >= bscstart && idxn[j] < bscend) {
380           col  = idxn[j] - bscstart;
381           ierr = MatGetValues_SeqSBAIJ(baij->A,1,&row,1,&col,v+i*n+j);CHKERRQ(ierr);
382         } else {
383           if (!baij->colmap) {
384             ierr = MatCreateColmap_MPIBAIJ_Private(mat);CHKERRQ(ierr);
385           }
386 #if defined(PETSC_USE_CTABLE)
387           ierr = PetscTableFind(baij->colmap,idxn[j]/bs+1,&data);CHKERRQ(ierr);
388           data--;
389 #else
390           data = baij->colmap[idxn[j]/bs]-1;
391 #endif
392           if ((data < 0) || (baij->garray[data/bs] != idxn[j]/bs)) *(v+i*n+j) = 0.0;
393           else {
394             col  = data + idxn[j]%bs;
395             ierr = MatGetValues_SeqBAIJ(baij->B,1,&row,1,&col,v+i*n+j);CHKERRQ(ierr);
396           }
397         }
398       }
399     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only local values currently supported");
400   }
401   PetscFunctionReturn(0);
402 }
403 
404 #undef __FUNCT__
405 #define __FUNCT__ "MatNorm_MPISBAIJ"
406 PetscErrorCode MatNorm_MPISBAIJ(Mat mat,NormType type,PetscReal *norm)
407 {
408   Mat_MPISBAIJ   *baij = (Mat_MPISBAIJ*)mat->data;
409   PetscErrorCode ierr;
410   PetscReal      sum[2],*lnorm2;
411 
412   PetscFunctionBegin;
413   if (baij->size == 1) {
414     ierr =  MatNorm(baij->A,type,norm);CHKERRQ(ierr);
415   } else {
416     if (type == NORM_FROBENIUS) {
417       ierr    = PetscMalloc(2*sizeof(PetscReal),&lnorm2);CHKERRQ(ierr);
418       ierr    =  MatNorm(baij->A,type,lnorm2);CHKERRQ(ierr);
419       *lnorm2 = (*lnorm2)*(*lnorm2); lnorm2++;            /* squar power of norm(A) */
420       ierr    =  MatNorm(baij->B,type,lnorm2);CHKERRQ(ierr);
421       *lnorm2 = (*lnorm2)*(*lnorm2); lnorm2--;             /* squar power of norm(B) */
422       ierr    = MPI_Allreduce(lnorm2,sum,2,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
423       *norm   = PetscSqrtReal(sum[0] + 2*sum[1]);
424       ierr    = PetscFree(lnorm2);CHKERRQ(ierr);
425     } else if (type == NORM_INFINITY || type == NORM_1) { /* max row/column sum */
426       Mat_SeqSBAIJ *amat=(Mat_SeqSBAIJ*)baij->A->data;
427       Mat_SeqBAIJ  *bmat=(Mat_SeqBAIJ*)baij->B->data;
428       PetscReal    *rsum,*rsum2,vabs;
429       PetscInt     *jj,*garray=baij->garray,rstart=baij->rstartbs,nz;
430       PetscInt     brow,bcol,col,bs=baij->A->rmap->bs,row,grow,gcol,mbs=amat->mbs;
431       MatScalar    *v;
432 
433       ierr = PetscMalloc2(mat->cmap->N,PetscReal,&rsum,mat->cmap->N,PetscReal,&rsum2);CHKERRQ(ierr);
434       ierr = PetscMemzero(rsum,mat->cmap->N*sizeof(PetscReal));CHKERRQ(ierr);
435       /* Amat */
436       v = amat->a; jj = amat->j;
437       for (brow=0; brow<mbs; brow++) {
438         grow = bs*(rstart + brow);
439         nz   = amat->i[brow+1] - amat->i[brow];
440         for (bcol=0; bcol<nz; bcol++) {
441           gcol = bs*(rstart + *jj); jj++;
442           for (col=0; col<bs; col++) {
443             for (row=0; row<bs; row++) {
444               vabs            = PetscAbsScalar(*v); v++;
445               rsum[gcol+col] += vabs;
446               /* non-diagonal block */
447               if (bcol > 0 && vabs > 0.0) rsum[grow+row] += vabs;
448             }
449           }
450         }
451       }
452       /* Bmat */
453       v = bmat->a; jj = bmat->j;
454       for (brow=0; brow<mbs; brow++) {
455         grow = bs*(rstart + brow);
456         nz = bmat->i[brow+1] - bmat->i[brow];
457         for (bcol=0; bcol<nz; bcol++) {
458           gcol = bs*garray[*jj]; jj++;
459           for (col=0; col<bs; col++) {
460             for (row=0; row<bs; row++) {
461               vabs            = PetscAbsScalar(*v); v++;
462               rsum[gcol+col] += vabs;
463               rsum[grow+row] += vabs;
464             }
465           }
466         }
467       }
468       ierr  = MPI_Allreduce(rsum,rsum2,mat->cmap->N,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
469       *norm = 0.0;
470       for (col=0; col<mat->cmap->N; col++) {
471         if (rsum2[col] > *norm) *norm = rsum2[col];
472       }
473       ierr = PetscFree2(rsum,rsum2);CHKERRQ(ierr);
474     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for this norm yet");
475   }
476   PetscFunctionReturn(0);
477 }
478 
479 #undef __FUNCT__
480 #define __FUNCT__ "MatAssemblyBegin_MPISBAIJ"
481 PetscErrorCode MatAssemblyBegin_MPISBAIJ(Mat mat,MatAssemblyType mode)
482 {
483   Mat_MPISBAIJ   *baij = (Mat_MPISBAIJ*)mat->data;
484   PetscErrorCode ierr;
485   PetscInt       nstash,reallocs;
486   InsertMode     addv;
487 
488   PetscFunctionBegin;
489   if (baij->donotstash || mat->nooffprocentries) PetscFunctionReturn(0);
490 
491   /* make sure all processors are either in INSERTMODE or ADDMODE */
492   ierr = MPI_Allreduce((PetscEnum*)&mat->insertmode,(PetscEnum*)&addv,1,MPIU_ENUM,MPI_BOR,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
493   if (addv == (ADD_VALUES|INSERT_VALUES)) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Some processors inserted others added");
494   mat->insertmode = addv; /* in case this processor had no cache */
495 
496   ierr = MatStashScatterBegin_Private(mat,&mat->stash,mat->rmap->range);CHKERRQ(ierr);
497   ierr = MatStashScatterBegin_Private(mat,&mat->bstash,baij->rangebs);CHKERRQ(ierr);
498   ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);CHKERRQ(ierr);
499   ierr = PetscInfo2(mat,"Stash has %D entries,uses %D mallocs.\n",nstash,reallocs);CHKERRQ(ierr);
500   ierr = MatStashGetInfo_Private(&mat->stash,&nstash,&reallocs);CHKERRQ(ierr);
501   ierr = PetscInfo2(mat,"Block-Stash has %D entries, uses %D mallocs.\n",nstash,reallocs);CHKERRQ(ierr);
502   PetscFunctionReturn(0);
503 }
504 
505 #undef __FUNCT__
506 #define __FUNCT__ "MatAssemblyEnd_MPISBAIJ"
507 PetscErrorCode MatAssemblyEnd_MPISBAIJ(Mat mat,MatAssemblyType mode)
508 {
509   Mat_MPISBAIJ   *baij=(Mat_MPISBAIJ*)mat->data;
510   Mat_SeqSBAIJ   *a   =(Mat_SeqSBAIJ*)baij->A->data;
511   PetscErrorCode ierr;
512   PetscInt       i,j,rstart,ncols,flg,bs2=baij->bs2;
513   PetscInt       *row,*col;
514   PetscBool      other_disassembled;
515   PetscMPIInt    n;
516   PetscBool      r1,r2,r3;
517   MatScalar      *val;
518   InsertMode     addv = mat->insertmode;
519 
520   /* do not use 'b=(Mat_SeqBAIJ*)baij->B->data' as B can be reset in disassembly */
521   PetscFunctionBegin;
522   if (!baij->donotstash &&  !mat->nooffprocentries) {
523     while (1) {
524       ierr = MatStashScatterGetMesg_Private(&mat->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr);
525       if (!flg) break;
526 
527       for (i=0; i<n;) {
528         /* Now identify the consecutive vals belonging to the same row */
529         for (j=i,rstart=row[j]; j<n; j++) {
530           if (row[j] != rstart) break;
531         }
532         if (j < n) ncols = j-i;
533         else       ncols = n-i;
534         /* Now assemble all these values with a single function call */
535         ierr = MatSetValues_MPISBAIJ(mat,1,row+i,ncols,col+i,val+i,addv);CHKERRQ(ierr);
536         i    = j;
537       }
538     }
539     ierr = MatStashScatterEnd_Private(&mat->stash);CHKERRQ(ierr);
540     /* Now process the block-stash. Since the values are stashed column-oriented,
541        set the roworiented flag to column oriented, and after MatSetValues()
542        restore the original flags */
543     r1 = baij->roworiented;
544     r2 = a->roworiented;
545     r3 = ((Mat_SeqBAIJ*)baij->B->data)->roworiented;
546 
547     baij->roworiented = PETSC_FALSE;
548     a->roworiented    = PETSC_FALSE;
549 
550     ((Mat_SeqBAIJ*)baij->B->data)->roworiented = PETSC_FALSE; /* b->roworinted */
551     while (1) {
552       ierr = MatStashScatterGetMesg_Private(&mat->bstash,&n,&row,&col,&val,&flg);CHKERRQ(ierr);
553       if (!flg) break;
554 
555       for (i=0; i<n;) {
556         /* Now identify the consecutive vals belonging to the same row */
557         for (j=i,rstart=row[j]; j<n; j++) {
558           if (row[j] != rstart) break;
559         }
560         if (j < n) ncols = j-i;
561         else       ncols = n-i;
562         ierr = MatSetValuesBlocked_MPISBAIJ(mat,1,row+i,ncols,col+i,val+i*bs2,addv);CHKERRQ(ierr);
563         i    = j;
564       }
565     }
566     ierr = MatStashScatterEnd_Private(&mat->bstash);CHKERRQ(ierr);
567 
568     baij->roworiented = r1;
569     a->roworiented    = r2;
570 
571     ((Mat_SeqBAIJ*)baij->B->data)->roworiented = r3; /* b->roworinted */
572   }
573 
574   ierr = MatAssemblyBegin(baij->A,mode);CHKERRQ(ierr);
575   ierr = MatAssemblyEnd(baij->A,mode);CHKERRQ(ierr);
576 
577   /* determine if any processor has disassembled, if so we must
578      also disassemble ourselfs, in order that we may reassemble. */
579   /*
580      if nonzero structure of submatrix B cannot change then we know that
581      no processor disassembled thus we can skip this stuff
582   */
583   if (!((Mat_SeqBAIJ*)baij->B->data)->nonew) {
584     ierr = MPI_Allreduce(&mat->was_assembled,&other_disassembled,1,MPIU_BOOL,MPI_PROD,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
585     if (mat->was_assembled && !other_disassembled) {
586       ierr = MatDisAssemble_MPISBAIJ(mat);CHKERRQ(ierr);
587     }
588   }
589 
590   if (!mat->was_assembled && mode == MAT_FINAL_ASSEMBLY) {
591     ierr = MatSetUpMultiply_MPISBAIJ(mat);CHKERRQ(ierr); /* setup Mvctx and sMvctx */
592   }
593   ierr = MatSetOption(baij->B,MAT_CHECK_COMPRESSED_ROW,PETSC_TRUE);CHKERRQ(ierr);
594   ierr = MatAssemblyBegin(baij->B,mode);CHKERRQ(ierr);
595   ierr = MatAssemblyEnd(baij->B,mode);CHKERRQ(ierr);
596 
597   ierr = PetscFree2(baij->rowvalues,baij->rowindices);CHKERRQ(ierr);
598 
599   baij->rowvalues = 0;
600   PetscFunctionReturn(0);
601 }
602 
603 extern PetscErrorCode MatSetValues_MPIBAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
604 #include <petscdraw.h>
605 #undef __FUNCT__
606 #define __FUNCT__ "MatView_MPISBAIJ_ASCIIorDraworSocket"
607 static PetscErrorCode MatView_MPISBAIJ_ASCIIorDraworSocket(Mat mat,PetscViewer viewer)
608 {
609   Mat_MPISBAIJ      *baij = (Mat_MPISBAIJ*)mat->data;
610   PetscErrorCode    ierr;
611   PetscInt          bs   = mat->rmap->bs;
612   PetscMPIInt       size = baij->size,rank = baij->rank;
613   PetscBool         iascii,isdraw;
614   PetscViewer       sviewer;
615   PetscViewerFormat format;
616 
617   PetscFunctionBegin;
618   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
619   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
620   if (iascii) {
621     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
622     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
623       MatInfo info;
624       ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&rank);CHKERRQ(ierr);
625       ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr);
626       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
627       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Local rows %D nz %D nz alloced %D bs %D mem %D\n",rank,mat->rmap->n,(PetscInt)info.nz_used,(PetscInt)info.nz_allocated,mat->rmap->bs,(PetscInt)info.memory);CHKERRQ(ierr);
628       ierr = MatGetInfo(baij->A,MAT_LOCAL,&info);CHKERRQ(ierr);
629       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] on-diagonal part: nz %D \n",rank,(PetscInt)info.nz_used);CHKERRQ(ierr);
630       ierr = MatGetInfo(baij->B,MAT_LOCAL,&info);CHKERRQ(ierr);
631       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] off-diagonal part: nz %D \n",rank,(PetscInt)info.nz_used);CHKERRQ(ierr);
632       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
633       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
634       ierr = PetscViewerASCIIPrintf(viewer,"Information on VecScatter used in matrix-vector product: \n");CHKERRQ(ierr);
635       ierr = VecScatterView(baij->Mvctx,viewer);CHKERRQ(ierr);
636       PetscFunctionReturn(0);
637     } else if (format == PETSC_VIEWER_ASCII_INFO) {
638       ierr = PetscViewerASCIIPrintf(viewer,"  block size is %D\n",bs);CHKERRQ(ierr);
639       PetscFunctionReturn(0);
640     } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
641       PetscFunctionReturn(0);
642     }
643   }
644 
645   if (isdraw) {
646     PetscDraw draw;
647     PetscBool isnull;
648     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
649     ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
650   }
651 
652   if (size == 1) {
653     ierr = PetscObjectSetName((PetscObject)baij->A,((PetscObject)mat)->name);CHKERRQ(ierr);
654     ierr = MatView(baij->A,viewer);CHKERRQ(ierr);
655   } else {
656     /* assemble the entire matrix onto first processor. */
657     Mat          A;
658     Mat_SeqSBAIJ *Aloc;
659     Mat_SeqBAIJ  *Bloc;
660     PetscInt     M = mat->rmap->N,N = mat->cmap->N,*ai,*aj,col,i,j,k,*rvals,mbs = baij->mbs;
661     MatScalar    *a;
662 
663     /* Should this be the same type as mat? */
664     ierr = MatCreate(PetscObjectComm((PetscObject)mat),&A);CHKERRQ(ierr);
665     if (!rank) {
666       ierr = MatSetSizes(A,M,N,M,N);CHKERRQ(ierr);
667     } else {
668       ierr = MatSetSizes(A,0,0,M,N);CHKERRQ(ierr);
669     }
670     ierr = MatSetType(A,MATMPISBAIJ);CHKERRQ(ierr);
671     ierr = MatMPISBAIJSetPreallocation(A,mat->rmap->bs,0,NULL,0,NULL);CHKERRQ(ierr);
672     ierr = MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
673     ierr = PetscLogObjectParent(mat,A);CHKERRQ(ierr);
674 
675     /* copy over the A part */
676     Aloc = (Mat_SeqSBAIJ*)baij->A->data;
677     ai   = Aloc->i; aj = Aloc->j; a = Aloc->a;
678     ierr = PetscMalloc(bs*sizeof(PetscInt),&rvals);CHKERRQ(ierr);
679 
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           ierr = MatSetValues_MPISBAIJ(A,bs,rvals,1,&col,a,INSERT_VALUES);CHKERRQ(ierr);
687           col++;
688           a += bs;
689         }
690       }
691     }
692     /* copy over the B part */
693     Bloc = (Mat_SeqBAIJ*)baij->B->data;
694     ai   = Bloc->i; aj = Bloc->j; a = Bloc->a;
695     for (i=0; i<mbs; i++) {
696 
697       rvals[0] = bs*(baij->rstartbs + i);
698       for (j=1; j<bs; j++) rvals[j] = rvals[j-1] + 1;
699       for (j=ai[i]; j<ai[i+1]; j++) {
700         col = baij->garray[aj[j]]*bs;
701         for (k=0; k<bs; k++) {
702           ierr = MatSetValues_MPIBAIJ(A,bs,rvals,1,&col,a,INSERT_VALUES);CHKERRQ(ierr);
703           col++;
704           a += bs;
705         }
706       }
707     }
708     ierr = PetscFree(rvals);CHKERRQ(ierr);
709     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
710     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
711     /*
712        Everyone has to call to draw the matrix since the graphics waits are
713        synchronized across all processors that share the PetscDraw object
714     */
715     ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
716     if (!rank) {
717       ierr = PetscObjectSetName((PetscObject)((Mat_MPISBAIJ*)(A->data))->A,((PetscObject)mat)->name);CHKERRQ(ierr);
718       /* Set the type name to MATMPISBAIJ so that the correct type can be printed out by PetscObjectPrintClassNamePrefixType() in MatView_SeqSBAIJ_ASCII()*/
719       PetscStrcpy(((PetscObject)((Mat_MPISBAIJ*)(A->data))->A)->type_name,MATMPISBAIJ);
720       ierr = MatView(((Mat_MPISBAIJ*)(A->data))->A,sviewer);CHKERRQ(ierr);
721     }
722     ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
723     ierr = MatDestroy(&A);CHKERRQ(ierr);
724   }
725   PetscFunctionReturn(0);
726 }
727 
728 #undef __FUNCT__
729 #define __FUNCT__ "MatView_MPISBAIJ"
730 PetscErrorCode MatView_MPISBAIJ(Mat mat,PetscViewer viewer)
731 {
732   PetscErrorCode ierr;
733   PetscBool      iascii,isdraw,issocket,isbinary;
734 
735   PetscFunctionBegin;
736   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
737   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
738   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSOCKET,&issocket);CHKERRQ(ierr);
739   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
740   if (iascii || isdraw || issocket || isbinary) {
741     ierr = MatView_MPISBAIJ_ASCIIorDraworSocket(mat,viewer);CHKERRQ(ierr);
742   }
743   PetscFunctionReturn(0);
744 }
745 
746 #undef __FUNCT__
747 #define __FUNCT__ "MatDestroy_MPISBAIJ"
748 PetscErrorCode MatDestroy_MPISBAIJ(Mat mat)
749 {
750   Mat_MPISBAIJ   *baij = (Mat_MPISBAIJ*)mat->data;
751   PetscErrorCode ierr;
752 
753   PetscFunctionBegin;
754 #if defined(PETSC_USE_LOG)
755   PetscLogObjectState((PetscObject)mat,"Rows=%D,Cols=%D",mat->rmap->N,mat->cmap->N);
756 #endif
757   ierr = MatStashDestroy_Private(&mat->stash);CHKERRQ(ierr);
758   ierr = MatStashDestroy_Private(&mat->bstash);CHKERRQ(ierr);
759   ierr = MatDestroy(&baij->A);CHKERRQ(ierr);
760   ierr = MatDestroy(&baij->B);CHKERRQ(ierr);
761 #if defined(PETSC_USE_CTABLE)
762   ierr = PetscTableDestroy(&baij->colmap);CHKERRQ(ierr);
763 #else
764   ierr = PetscFree(baij->colmap);CHKERRQ(ierr);
765 #endif
766   ierr = PetscFree(baij->garray);CHKERRQ(ierr);
767   ierr = VecDestroy(&baij->lvec);CHKERRQ(ierr);
768   ierr = VecScatterDestroy(&baij->Mvctx);CHKERRQ(ierr);
769   ierr = VecDestroy(&baij->slvec0);CHKERRQ(ierr);
770   ierr = VecDestroy(&baij->slvec0b);CHKERRQ(ierr);
771   ierr = VecDestroy(&baij->slvec1);CHKERRQ(ierr);
772   ierr = VecDestroy(&baij->slvec1a);CHKERRQ(ierr);
773   ierr = VecDestroy(&baij->slvec1b);CHKERRQ(ierr);
774   ierr = VecScatterDestroy(&baij->sMvctx);CHKERRQ(ierr);
775   ierr = PetscFree2(baij->rowvalues,baij->rowindices);CHKERRQ(ierr);
776   ierr = PetscFree(baij->barray);CHKERRQ(ierr);
777   ierr = PetscFree(baij->hd);CHKERRQ(ierr);
778   ierr = VecDestroy(&baij->diag);CHKERRQ(ierr);
779   ierr = VecDestroy(&baij->bb1);CHKERRQ(ierr);
780   ierr = VecDestroy(&baij->xx1);CHKERRQ(ierr);
781 #if defined(PETSC_USE_REAL_MAT_SINGLE)
782   ierr = PetscFree(baij->setvaluescopy);CHKERRQ(ierr);
783 #endif
784   ierr = PetscFree(baij->in_loc);CHKERRQ(ierr);
785   ierr = PetscFree(baij->v_loc);CHKERRQ(ierr);
786   ierr = PetscFree(baij->rangebs);CHKERRQ(ierr);
787   ierr = PetscFree(mat->data);CHKERRQ(ierr);
788 
789   ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr);
790   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatStoreValues_C","",NULL);CHKERRQ(ierr);
791   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatRetrieveValues_C","",NULL);CHKERRQ(ierr);
792   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatGetDiagonalBlock_C","",NULL);CHKERRQ(ierr);
793   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPISBAIJSetPreallocation_C","",NULL);CHKERRQ(ierr);
794   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatConvert_mpisbaij_mpisbstrm_C","",NULL);CHKERRQ(ierr);
795   PetscFunctionReturn(0);
796 }
797 
798 #undef __FUNCT__
799 #define __FUNCT__ "MatMult_MPISBAIJ_Hermitian"
800 PetscErrorCode MatMult_MPISBAIJ_Hermitian(Mat A,Vec xx,Vec yy)
801 {
802   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;
803   PetscErrorCode ierr;
804   PetscInt       nt,mbs=a->mbs,bs=A->rmap->bs;
805   PetscScalar    *x,*from;
806 
807   PetscFunctionBegin;
808   ierr = VecGetLocalSize(xx,&nt);CHKERRQ(ierr);
809   if (nt != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Incompatible partition of A and xx");
810 
811   /* diagonal part */
812   ierr = (*a->A->ops->mult)(a->A,xx,a->slvec1a);CHKERRQ(ierr);
813   ierr = VecSet(a->slvec1b,0.0);CHKERRQ(ierr);
814 
815   /* subdiagonal part */
816   ierr = (*a->B->ops->multhermitiantranspose)(a->B,xx,a->slvec0b);CHKERRQ(ierr);
817 
818   /* copy x into the vec slvec0 */
819   ierr = VecGetArray(a->slvec0,&from);CHKERRQ(ierr);
820   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
821 
822   ierr = PetscMemcpy(from,x,bs*mbs*sizeof(MatScalar));CHKERRQ(ierr);
823   ierr = VecRestoreArray(a->slvec0,&from);CHKERRQ(ierr);
824   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
825 
826   ierr = VecScatterBegin(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
827   ierr = VecScatterEnd(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
828   /* supperdiagonal part */
829   ierr = (*a->B->ops->multadd)(a->B,a->slvec1b,a->slvec1a,yy);CHKERRQ(ierr);
830   PetscFunctionReturn(0);
831 }
832 
833 #undef __FUNCT__
834 #define __FUNCT__ "MatMult_MPISBAIJ"
835 PetscErrorCode MatMult_MPISBAIJ(Mat A,Vec xx,Vec yy)
836 {
837   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;
838   PetscErrorCode ierr;
839   PetscInt       nt,mbs=a->mbs,bs=A->rmap->bs;
840   PetscScalar    *x,*from;
841 
842   PetscFunctionBegin;
843   ierr = VecGetLocalSize(xx,&nt);CHKERRQ(ierr);
844   if (nt != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Incompatible partition of A and xx");
845 
846   /* diagonal part */
847   ierr = (*a->A->ops->mult)(a->A,xx,a->slvec1a);CHKERRQ(ierr);
848   ierr = VecSet(a->slvec1b,0.0);CHKERRQ(ierr);
849 
850   /* subdiagonal part */
851   ierr = (*a->B->ops->multtranspose)(a->B,xx,a->slvec0b);CHKERRQ(ierr);
852 
853   /* copy x into the vec slvec0 */
854   ierr = VecGetArray(a->slvec0,&from);CHKERRQ(ierr);
855   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
856 
857   ierr = PetscMemcpy(from,x,bs*mbs*sizeof(MatScalar));CHKERRQ(ierr);
858   ierr = VecRestoreArray(a->slvec0,&from);CHKERRQ(ierr);
859   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
860 
861   ierr = VecScatterBegin(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
862   ierr = VecScatterEnd(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
863   /* supperdiagonal part */
864   ierr = (*a->B->ops->multadd)(a->B,a->slvec1b,a->slvec1a,yy);CHKERRQ(ierr);
865   PetscFunctionReturn(0);
866 }
867 
868 #undef __FUNCT__
869 #define __FUNCT__ "MatMult_MPISBAIJ_2comm"
870 PetscErrorCode MatMult_MPISBAIJ_2comm(Mat A,Vec xx,Vec yy)
871 {
872   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;
873   PetscErrorCode ierr;
874   PetscInt       nt;
875 
876   PetscFunctionBegin;
877   ierr = VecGetLocalSize(xx,&nt);CHKERRQ(ierr);
878   if (nt != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Incompatible partition of A and xx");
879 
880   ierr = VecGetLocalSize(yy,&nt);CHKERRQ(ierr);
881   if (nt != A->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Incompatible parition of A and yy");
882 
883   ierr = VecScatterBegin(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
884   /* do diagonal part */
885   ierr = (*a->A->ops->mult)(a->A,xx,yy);CHKERRQ(ierr);
886   /* do supperdiagonal part */
887   ierr = VecScatterEnd(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
888   ierr = (*a->B->ops->multadd)(a->B,a->lvec,yy,yy);CHKERRQ(ierr);
889   /* do subdiagonal part */
890   ierr = (*a->B->ops->multtranspose)(a->B,xx,a->lvec);CHKERRQ(ierr);
891   ierr = VecScatterBegin(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
892   ierr = VecScatterEnd(a->Mvctx,a->lvec,yy,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
893   PetscFunctionReturn(0);
894 }
895 
896 #undef __FUNCT__
897 #define __FUNCT__ "MatMultAdd_MPISBAIJ"
898 PetscErrorCode MatMultAdd_MPISBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
899 {
900   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;
901   PetscErrorCode ierr;
902   PetscInt       mbs=a->mbs,bs=A->rmap->bs;
903   PetscScalar    *x,*from,zero=0.0;
904 
905   PetscFunctionBegin;
906   /*
907   PetscSynchronizedPrintf(PetscObjectComm((PetscObject)A)," MatMultAdd is called ...\n");
908   PetscSynchronizedFlush(PetscObjectComm((PetscObject)A));
909   */
910   /* diagonal part */
911   ierr = (*a->A->ops->multadd)(a->A,xx,yy,a->slvec1a);CHKERRQ(ierr);
912   ierr = VecSet(a->slvec1b,zero);CHKERRQ(ierr);
913 
914   /* subdiagonal part */
915   ierr = (*a->B->ops->multtranspose)(a->B,xx,a->slvec0b);CHKERRQ(ierr);
916 
917   /* copy x into the vec slvec0 */
918   ierr = VecGetArray(a->slvec0,&from);CHKERRQ(ierr);
919   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
920   ierr = PetscMemcpy(from,x,bs*mbs*sizeof(MatScalar));CHKERRQ(ierr);
921   ierr = VecRestoreArray(a->slvec0,&from);CHKERRQ(ierr);
922 
923   ierr = VecScatterBegin(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
924   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
925   ierr = VecScatterEnd(a->sMvctx,a->slvec0,a->slvec1,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
926 
927   /* supperdiagonal part */
928   ierr = (*a->B->ops->multadd)(a->B,a->slvec1b,a->slvec1a,zz);CHKERRQ(ierr);
929   PetscFunctionReturn(0);
930 }
931 
932 #undef __FUNCT__
933 #define __FUNCT__ "MatMultAdd_MPISBAIJ_2comm"
934 PetscErrorCode MatMultAdd_MPISBAIJ_2comm(Mat A,Vec xx,Vec yy,Vec zz)
935 {
936   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;
937   PetscErrorCode ierr;
938 
939   PetscFunctionBegin;
940   ierr = VecScatterBegin(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
941   /* do diagonal part */
942   ierr = (*a->A->ops->multadd)(a->A,xx,yy,zz);CHKERRQ(ierr);
943   /* do supperdiagonal part */
944   ierr = VecScatterEnd(a->Mvctx,xx,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
945   ierr = (*a->B->ops->multadd)(a->B,a->lvec,zz,zz);CHKERRQ(ierr);
946 
947   /* do subdiagonal part */
948   ierr = (*a->B->ops->multtranspose)(a->B,xx,a->lvec);CHKERRQ(ierr);
949   ierr = VecScatterBegin(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
950   ierr = VecScatterEnd(a->Mvctx,a->lvec,zz,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
951   PetscFunctionReturn(0);
952 }
953 
954 /*
955   This only works correctly for square matrices where the subblock A->A is the
956    diagonal block
957 */
958 #undef __FUNCT__
959 #define __FUNCT__ "MatGetDiagonal_MPISBAIJ"
960 PetscErrorCode MatGetDiagonal_MPISBAIJ(Mat A,Vec v)
961 {
962   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;
963   PetscErrorCode ierr;
964 
965   PetscFunctionBegin;
966   /* if (a->rmap->N != a->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Supports only square matrix where A->A is diag block"); */
967   ierr = MatGetDiagonal(a->A,v);CHKERRQ(ierr);
968   PetscFunctionReturn(0);
969 }
970 
971 #undef __FUNCT__
972 #define __FUNCT__ "MatScale_MPISBAIJ"
973 PetscErrorCode MatScale_MPISBAIJ(Mat A,PetscScalar aa)
974 {
975   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;
976   PetscErrorCode ierr;
977 
978   PetscFunctionBegin;
979   ierr = MatScale(a->A,aa);CHKERRQ(ierr);
980   ierr = MatScale(a->B,aa);CHKERRQ(ierr);
981   PetscFunctionReturn(0);
982 }
983 
984 #undef __FUNCT__
985 #define __FUNCT__ "MatGetRow_MPISBAIJ"
986 PetscErrorCode MatGetRow_MPISBAIJ(Mat matin,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
987 {
988   Mat_MPISBAIJ   *mat = (Mat_MPISBAIJ*)matin->data;
989   PetscScalar    *vworkA,*vworkB,**pvA,**pvB,*v_p;
990   PetscErrorCode ierr;
991   PetscInt       bs = matin->rmap->bs,bs2 = mat->bs2,i,*cworkA,*cworkB,**pcA,**pcB;
992   PetscInt       nztot,nzA,nzB,lrow,brstart = matin->rmap->rstart,brend = matin->rmap->rend;
993   PetscInt       *cmap,*idx_p,cstart = mat->rstartbs;
994 
995   PetscFunctionBegin;
996   if (mat->getrowactive) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Already active");
997   mat->getrowactive = PETSC_TRUE;
998 
999   if (!mat->rowvalues && (idx || v)) {
1000     /*
1001         allocate enough space to hold information from the longest row.
1002     */
1003     Mat_SeqSBAIJ *Aa = (Mat_SeqSBAIJ*)mat->A->data;
1004     Mat_SeqBAIJ  *Ba = (Mat_SeqBAIJ*)mat->B->data;
1005     PetscInt     max = 1,mbs = mat->mbs,tmp;
1006     for (i=0; i<mbs; i++) {
1007       tmp = Aa->i[i+1] - Aa->i[i] + Ba->i[i+1] - Ba->i[i]; /* row length */
1008       if (max < tmp) max = tmp;
1009     }
1010     ierr = PetscMalloc2(max*bs2,PetscScalar,&mat->rowvalues,max*bs2,PetscInt,&mat->rowindices);CHKERRQ(ierr);
1011   }
1012 
1013   if (row < brstart || row >= brend) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only local rows");
1014   lrow = row - brstart;  /* local row index */
1015 
1016   pvA = &vworkA; pcA = &cworkA; pvB = &vworkB; pcB = &cworkB;
1017   if (!v)   {pvA = 0; pvB = 0;}
1018   if (!idx) {pcA = 0; if (!v) pcB = 0;}
1019   ierr  = (*mat->A->ops->getrow)(mat->A,lrow,&nzA,pcA,pvA);CHKERRQ(ierr);
1020   ierr  = (*mat->B->ops->getrow)(mat->B,lrow,&nzB,pcB,pvB);CHKERRQ(ierr);
1021   nztot = nzA + nzB;
1022 
1023   cmap = mat->garray;
1024   if (v  || idx) {
1025     if (nztot) {
1026       /* Sort by increasing column numbers, assuming A and B already sorted */
1027       PetscInt imark = -1;
1028       if (v) {
1029         *v = v_p = mat->rowvalues;
1030         for (i=0; i<nzB; i++) {
1031           if (cmap[cworkB[i]/bs] < cstart) v_p[i] = vworkB[i];
1032           else break;
1033         }
1034         imark = i;
1035         for (i=0; i<nzA; i++)     v_p[imark+i] = vworkA[i];
1036         for (i=imark; i<nzB; i++) v_p[nzA+i]   = vworkB[i];
1037       }
1038       if (idx) {
1039         *idx = idx_p = mat->rowindices;
1040         if (imark > -1) {
1041           for (i=0; i<imark; i++) {
1042             idx_p[i] = cmap[cworkB[i]/bs]*bs + cworkB[i]%bs;
1043           }
1044         } else {
1045           for (i=0; i<nzB; i++) {
1046             if (cmap[cworkB[i]/bs] < cstart) idx_p[i] = cmap[cworkB[i]/bs]*bs + cworkB[i]%bs;
1047             else break;
1048           }
1049           imark = i;
1050         }
1051         for (i=0; i<nzA; i++)     idx_p[imark+i] = cstart*bs + cworkA[i];
1052         for (i=imark; i<nzB; i++) idx_p[nzA+i]   = cmap[cworkB[i]/bs]*bs + cworkB[i]%bs ;
1053       }
1054     } else {
1055       if (idx) *idx = 0;
1056       if (v)   *v   = 0;
1057     }
1058   }
1059   *nz  = nztot;
1060   ierr = (*mat->A->ops->restorerow)(mat->A,lrow,&nzA,pcA,pvA);CHKERRQ(ierr);
1061   ierr = (*mat->B->ops->restorerow)(mat->B,lrow,&nzB,pcB,pvB);CHKERRQ(ierr);
1062   PetscFunctionReturn(0);
1063 }
1064 
1065 #undef __FUNCT__
1066 #define __FUNCT__ "MatRestoreRow_MPISBAIJ"
1067 PetscErrorCode MatRestoreRow_MPISBAIJ(Mat mat,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
1068 {
1069   Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
1070 
1071   PetscFunctionBegin;
1072   if (!baij->getrowactive) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"MatGetRow() must be called first");
1073   baij->getrowactive = PETSC_FALSE;
1074   PetscFunctionReturn(0);
1075 }
1076 
1077 #undef __FUNCT__
1078 #define __FUNCT__ "MatGetRowUpperTriangular_MPISBAIJ"
1079 PetscErrorCode MatGetRowUpperTriangular_MPISBAIJ(Mat A)
1080 {
1081   Mat_MPISBAIJ *a  = (Mat_MPISBAIJ*)A->data;
1082   Mat_SeqSBAIJ *aA = (Mat_SeqSBAIJ*)a->A->data;
1083 
1084   PetscFunctionBegin;
1085   aA->getrow_utriangular = PETSC_TRUE;
1086   PetscFunctionReturn(0);
1087 }
1088 #undef __FUNCT__
1089 #define __FUNCT__ "MatRestoreRowUpperTriangular_MPISBAIJ"
1090 PetscErrorCode MatRestoreRowUpperTriangular_MPISBAIJ(Mat A)
1091 {
1092   Mat_MPISBAIJ *a  = (Mat_MPISBAIJ*)A->data;
1093   Mat_SeqSBAIJ *aA = (Mat_SeqSBAIJ*)a->A->data;
1094 
1095   PetscFunctionBegin;
1096   aA->getrow_utriangular = PETSC_FALSE;
1097   PetscFunctionReturn(0);
1098 }
1099 
1100 #undef __FUNCT__
1101 #define __FUNCT__ "MatRealPart_MPISBAIJ"
1102 PetscErrorCode MatRealPart_MPISBAIJ(Mat A)
1103 {
1104   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;
1105   PetscErrorCode ierr;
1106 
1107   PetscFunctionBegin;
1108   ierr = MatRealPart(a->A);CHKERRQ(ierr);
1109   ierr = MatRealPart(a->B);CHKERRQ(ierr);
1110   PetscFunctionReturn(0);
1111 }
1112 
1113 #undef __FUNCT__
1114 #define __FUNCT__ "MatImaginaryPart_MPISBAIJ"
1115 PetscErrorCode MatImaginaryPart_MPISBAIJ(Mat A)
1116 {
1117   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;
1118   PetscErrorCode ierr;
1119 
1120   PetscFunctionBegin;
1121   ierr = MatImaginaryPart(a->A);CHKERRQ(ierr);
1122   ierr = MatImaginaryPart(a->B);CHKERRQ(ierr);
1123   PetscFunctionReturn(0);
1124 }
1125 
1126 #undef __FUNCT__
1127 #define __FUNCT__ "MatZeroEntries_MPISBAIJ"
1128 PetscErrorCode MatZeroEntries_MPISBAIJ(Mat A)
1129 {
1130   Mat_MPISBAIJ   *l = (Mat_MPISBAIJ*)A->data;
1131   PetscErrorCode ierr;
1132 
1133   PetscFunctionBegin;
1134   ierr = MatZeroEntries(l->A);CHKERRQ(ierr);
1135   ierr = MatZeroEntries(l->B);CHKERRQ(ierr);
1136   PetscFunctionReturn(0);
1137 }
1138 
1139 #undef __FUNCT__
1140 #define __FUNCT__ "MatGetInfo_MPISBAIJ"
1141 PetscErrorCode MatGetInfo_MPISBAIJ(Mat matin,MatInfoType flag,MatInfo *info)
1142 {
1143   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)matin->data;
1144   Mat            A  = a->A,B = a->B;
1145   PetscErrorCode ierr;
1146   PetscReal      isend[5],irecv[5];
1147 
1148   PetscFunctionBegin;
1149   info->block_size = (PetscReal)matin->rmap->bs;
1150 
1151   ierr = MatGetInfo(A,MAT_LOCAL,info);CHKERRQ(ierr);
1152 
1153   isend[0] = info->nz_used; isend[1] = info->nz_allocated; isend[2] = info->nz_unneeded;
1154   isend[3] = info->memory;  isend[4] = info->mallocs;
1155 
1156   ierr = MatGetInfo(B,MAT_LOCAL,info);CHKERRQ(ierr);
1157 
1158   isend[0] += info->nz_used; isend[1] += info->nz_allocated; isend[2] += info->nz_unneeded;
1159   isend[3] += info->memory;  isend[4] += info->mallocs;
1160   if (flag == MAT_LOCAL) {
1161     info->nz_used      = isend[0];
1162     info->nz_allocated = isend[1];
1163     info->nz_unneeded  = isend[2];
1164     info->memory       = isend[3];
1165     info->mallocs      = isend[4];
1166   } else if (flag == MAT_GLOBAL_MAX) {
1167     ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)matin));CHKERRQ(ierr);
1168 
1169     info->nz_used      = irecv[0];
1170     info->nz_allocated = irecv[1];
1171     info->nz_unneeded  = irecv[2];
1172     info->memory       = irecv[3];
1173     info->mallocs      = irecv[4];
1174   } else if (flag == MAT_GLOBAL_SUM) {
1175     ierr = MPI_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)matin));CHKERRQ(ierr);
1176 
1177     info->nz_used      = irecv[0];
1178     info->nz_allocated = irecv[1];
1179     info->nz_unneeded  = irecv[2];
1180     info->memory       = irecv[3];
1181     info->mallocs      = irecv[4];
1182   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown MatInfoType argument %d",(int)flag);
1183   info->fill_ratio_given  = 0; /* no parallel LU/ILU/Cholesky */
1184   info->fill_ratio_needed = 0;
1185   info->factor_mallocs    = 0;
1186   PetscFunctionReturn(0);
1187 }
1188 
1189 #undef __FUNCT__
1190 #define __FUNCT__ "MatSetOption_MPISBAIJ"
1191 PetscErrorCode MatSetOption_MPISBAIJ(Mat A,MatOption op,PetscBool flg)
1192 {
1193   Mat_MPISBAIJ   *a  = (Mat_MPISBAIJ*)A->data;
1194   Mat_SeqSBAIJ   *aA = (Mat_SeqSBAIJ*)a->A->data;
1195   PetscErrorCode ierr;
1196 
1197   PetscFunctionBegin;
1198   switch (op) {
1199   case MAT_NEW_NONZERO_LOCATIONS:
1200   case MAT_NEW_NONZERO_ALLOCATION_ERR:
1201   case MAT_UNUSED_NONZERO_LOCATION_ERR:
1202   case MAT_KEEP_NONZERO_PATTERN:
1203   case MAT_NEW_NONZERO_LOCATION_ERR:
1204     ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
1205     ierr = MatSetOption(a->B,op,flg);CHKERRQ(ierr);
1206     break;
1207   case MAT_ROW_ORIENTED:
1208     a->roworiented = flg;
1209 
1210     ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
1211     ierr = MatSetOption(a->B,op,flg);CHKERRQ(ierr);
1212     break;
1213   case MAT_NEW_DIAGONALS:
1214     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
1215     break;
1216   case MAT_IGNORE_OFF_PROC_ENTRIES:
1217     a->donotstash = flg;
1218     break;
1219   case MAT_USE_HASH_TABLE:
1220     a->ht_flag = flg;
1221     break;
1222   case MAT_HERMITIAN:
1223     if (!A->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call MatAssemblyEnd() first");
1224     ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
1225 
1226     A->ops->mult = MatMult_MPISBAIJ_Hermitian;
1227     break;
1228   case MAT_SPD:
1229     A->spd_set = PETSC_TRUE;
1230     A->spd     = flg;
1231     if (flg) {
1232       A->symmetric                  = PETSC_TRUE;
1233       A->structurally_symmetric     = PETSC_TRUE;
1234       A->symmetric_set              = PETSC_TRUE;
1235       A->structurally_symmetric_set = PETSC_TRUE;
1236     }
1237     break;
1238   case MAT_SYMMETRIC:
1239     ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
1240     break;
1241   case MAT_STRUCTURALLY_SYMMETRIC:
1242     ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
1243     break;
1244   case MAT_SYMMETRY_ETERNAL:
1245     if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix must be symmetric");
1246     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
1247     break;
1248   case MAT_IGNORE_LOWER_TRIANGULAR:
1249     aA->ignore_ltriangular = flg;
1250     break;
1251   case MAT_ERROR_LOWER_TRIANGULAR:
1252     aA->ignore_ltriangular = flg;
1253     break;
1254   case MAT_GETROW_UPPERTRIANGULAR:
1255     aA->getrow_utriangular = flg;
1256     break;
1257   default:
1258     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %d",op);
1259   }
1260   PetscFunctionReturn(0);
1261 }
1262 
1263 #undef __FUNCT__
1264 #define __FUNCT__ "MatTranspose_MPISBAIJ"
1265 PetscErrorCode MatTranspose_MPISBAIJ(Mat A,MatReuse reuse,Mat *B)
1266 {
1267   PetscErrorCode ierr;
1268 
1269   PetscFunctionBegin;
1270   if (MAT_INITIAL_MATRIX || *B != A) {
1271     ierr = MatDuplicate(A,MAT_COPY_VALUES,B);CHKERRQ(ierr);
1272   }
1273   PetscFunctionReturn(0);
1274 }
1275 
1276 #undef __FUNCT__
1277 #define __FUNCT__ "MatDiagonalScale_MPISBAIJ"
1278 PetscErrorCode MatDiagonalScale_MPISBAIJ(Mat mat,Vec ll,Vec rr)
1279 {
1280   Mat_MPISBAIJ   *baij = (Mat_MPISBAIJ*)mat->data;
1281   Mat            a     = baij->A, b=baij->B;
1282   PetscErrorCode ierr;
1283   PetscInt       nv,m,n;
1284   PetscBool      flg;
1285 
1286   PetscFunctionBegin;
1287   if (ll != rr) {
1288     ierr = VecEqual(ll,rr,&flg);CHKERRQ(ierr);
1289     if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"For symmetric format, left and right scaling vectors must be same\n");
1290   }
1291   if (!ll) PetscFunctionReturn(0);
1292 
1293   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
1294   if (m != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"For symmetric format, local size %d %d must be same",m,n);
1295 
1296   ierr = VecGetLocalSize(rr,&nv);CHKERRQ(ierr);
1297   if (nv!=n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left and right vector non-conforming local size");
1298 
1299   ierr = VecScatterBegin(baij->Mvctx,rr,baij->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1300 
1301   /* left diagonalscale the off-diagonal part */
1302   ierr = (*b->ops->diagonalscale)(b,ll,NULL);CHKERRQ(ierr);
1303 
1304   /* scale the diagonal part */
1305   ierr = (*a->ops->diagonalscale)(a,ll,rr);CHKERRQ(ierr);
1306 
1307   /* right diagonalscale the off-diagonal part */
1308   ierr = VecScatterEnd(baij->Mvctx,rr,baij->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1309   ierr = (*b->ops->diagonalscale)(b,NULL,baij->lvec);CHKERRQ(ierr);
1310   PetscFunctionReturn(0);
1311 }
1312 
1313 #undef __FUNCT__
1314 #define __FUNCT__ "MatSetUnfactored_MPISBAIJ"
1315 PetscErrorCode MatSetUnfactored_MPISBAIJ(Mat A)
1316 {
1317   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;
1318   PetscErrorCode ierr;
1319 
1320   PetscFunctionBegin;
1321   ierr = MatSetUnfactored(a->A);CHKERRQ(ierr);
1322   PetscFunctionReturn(0);
1323 }
1324 
1325 static PetscErrorCode MatDuplicate_MPISBAIJ(Mat,MatDuplicateOption,Mat*);
1326 
1327 #undef __FUNCT__
1328 #define __FUNCT__ "MatEqual_MPISBAIJ"
1329 PetscErrorCode MatEqual_MPISBAIJ(Mat A,Mat B,PetscBool  *flag)
1330 {
1331   Mat_MPISBAIJ   *matB = (Mat_MPISBAIJ*)B->data,*matA = (Mat_MPISBAIJ*)A->data;
1332   Mat            a,b,c,d;
1333   PetscBool      flg;
1334   PetscErrorCode ierr;
1335 
1336   PetscFunctionBegin;
1337   a = matA->A; b = matA->B;
1338   c = matB->A; d = matB->B;
1339 
1340   ierr = MatEqual(a,c,&flg);CHKERRQ(ierr);
1341   if (flg) {
1342     ierr = MatEqual(b,d,&flg);CHKERRQ(ierr);
1343   }
1344   ierr = MPI_Allreduce(&flg,flag,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
1345   PetscFunctionReturn(0);
1346 }
1347 
1348 #undef __FUNCT__
1349 #define __FUNCT__ "MatCopy_MPISBAIJ"
1350 PetscErrorCode MatCopy_MPISBAIJ(Mat A,Mat B,MatStructure str)
1351 {
1352   PetscErrorCode ierr;
1353   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;
1354   Mat_MPISBAIJ   *b = (Mat_MPISBAIJ*)B->data;
1355 
1356   PetscFunctionBegin;
1357   /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */
1358   if ((str != SAME_NONZERO_PATTERN) || (A->ops->copy != B->ops->copy)) {
1359     ierr = MatGetRowUpperTriangular(A);CHKERRQ(ierr);
1360     ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr);
1361     ierr = MatRestoreRowUpperTriangular(A);CHKERRQ(ierr);
1362   } else {
1363     ierr = MatCopy(a->A,b->A,str);CHKERRQ(ierr);
1364     ierr = MatCopy(a->B,b->B,str);CHKERRQ(ierr);
1365   }
1366   PetscFunctionReturn(0);
1367 }
1368 
1369 #undef __FUNCT__
1370 #define __FUNCT__ "MatSetUp_MPISBAIJ"
1371 PetscErrorCode MatSetUp_MPISBAIJ(Mat A)
1372 {
1373   PetscErrorCode ierr;
1374 
1375   PetscFunctionBegin;
1376   ierr = MatMPISBAIJSetPreallocation(A,A->rmap->bs,PETSC_DEFAULT,0,PETSC_DEFAULT,0);CHKERRQ(ierr);
1377   PetscFunctionReturn(0);
1378 }
1379 
1380 #undef __FUNCT__
1381 #define __FUNCT__ "MatAXPY_MPISBAIJ"
1382 PetscErrorCode MatAXPY_MPISBAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str)
1383 {
1384   PetscErrorCode ierr;
1385   Mat_MPISBAIJ   *xx=(Mat_MPISBAIJ*)X->data,*yy=(Mat_MPISBAIJ*)Y->data;
1386   PetscBLASInt   bnz,one=1;
1387   Mat_SeqSBAIJ   *xa,*ya;
1388   Mat_SeqBAIJ    *xb,*yb;
1389 
1390   PetscFunctionBegin;
1391   if (str == SAME_NONZERO_PATTERN) {
1392     PetscScalar alpha = a;
1393     xa   = (Mat_SeqSBAIJ*)xx->A->data;
1394     ya   = (Mat_SeqSBAIJ*)yy->A->data;
1395     ierr = PetscBLASIntCast(xa->nz,&bnz);CHKERRQ(ierr);
1396     PetscStackCall("BLASaxpy",BLASaxpy_(&bnz,&alpha,xa->a,&one,ya->a,&one));
1397     xb   = (Mat_SeqBAIJ*)xx->B->data;
1398     yb   = (Mat_SeqBAIJ*)yy->B->data;
1399     ierr = PetscBLASIntCast(xb->nz,&bnz);CHKERRQ(ierr);
1400     PetscStackCall("BLASaxpy",BLASaxpy_(&bnz,&alpha,xb->a,&one,yb->a,&one));
1401   } else {
1402     ierr = MatGetRowUpperTriangular(X);CHKERRQ(ierr);
1403     ierr = MatAXPY_Basic(Y,a,X,str);CHKERRQ(ierr);
1404     ierr = MatRestoreRowUpperTriangular(X);CHKERRQ(ierr);
1405   }
1406   PetscFunctionReturn(0);
1407 }
1408 
1409 #undef __FUNCT__
1410 #define __FUNCT__ "MatGetSubMatrices_MPISBAIJ"
1411 PetscErrorCode MatGetSubMatrices_MPISBAIJ(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
1412 {
1413   PetscErrorCode ierr;
1414   PetscInt       i;
1415   PetscBool      flg;
1416 
1417   PetscFunctionBegin;
1418   for (i=0; i<n; i++) {
1419     ierr = ISEqual(irow[i],icol[i],&flg);CHKERRQ(ierr);
1420     if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Can only get symmetric submatrix for MPISBAIJ matrices");
1421   }
1422   ierr = MatGetSubMatrices_MPIBAIJ(A,n,irow,icol,scall,B);CHKERRQ(ierr);
1423   PetscFunctionReturn(0);
1424 }
1425 
1426 
1427 /* -------------------------------------------------------------------*/
1428 static struct _MatOps MatOps_Values = {
1429         MatSetValues_MPISBAIJ,
1430         MatGetRow_MPISBAIJ,
1431         MatRestoreRow_MPISBAIJ,
1432         MatMult_MPISBAIJ,
1433 /*  4*/ MatMultAdd_MPISBAIJ,
1434         MatMult_MPISBAIJ,       /* transpose versions are same as non-transpose */
1435         MatMultAdd_MPISBAIJ,
1436         0,
1437         0,
1438         0,
1439 /* 10*/ 0,
1440         0,
1441         0,
1442         MatSOR_MPISBAIJ,
1443         MatTranspose_MPISBAIJ,
1444 /* 15*/ MatGetInfo_MPISBAIJ,
1445         MatEqual_MPISBAIJ,
1446         MatGetDiagonal_MPISBAIJ,
1447         MatDiagonalScale_MPISBAIJ,
1448         MatNorm_MPISBAIJ,
1449 /* 20*/ MatAssemblyBegin_MPISBAIJ,
1450         MatAssemblyEnd_MPISBAIJ,
1451         MatSetOption_MPISBAIJ,
1452         MatZeroEntries_MPISBAIJ,
1453 /* 24*/ 0,
1454         0,
1455         0,
1456         0,
1457         0,
1458 /* 29*/ MatSetUp_MPISBAIJ,
1459         0,
1460         0,
1461         0,
1462         0,
1463 /* 34*/ MatDuplicate_MPISBAIJ,
1464         0,
1465         0,
1466         0,
1467         0,
1468 /* 39*/ MatAXPY_MPISBAIJ,
1469         MatGetSubMatrices_MPISBAIJ,
1470         MatIncreaseOverlap_MPISBAIJ,
1471         MatGetValues_MPISBAIJ,
1472         MatCopy_MPISBAIJ,
1473 /* 44*/ 0,
1474         MatScale_MPISBAIJ,
1475         0,
1476         0,
1477         0,
1478 /* 49*/ 0,
1479         0,
1480         0,
1481         0,
1482         0,
1483 /* 54*/ 0,
1484         0,
1485         MatSetUnfactored_MPISBAIJ,
1486         0,
1487         MatSetValuesBlocked_MPISBAIJ,
1488 /* 59*/ 0,
1489         0,
1490         0,
1491         0,
1492         0,
1493 /* 64*/ 0,
1494         0,
1495         0,
1496         0,
1497         0,
1498 /* 69*/ MatGetRowMaxAbs_MPISBAIJ,
1499         0,
1500         0,
1501         0,
1502         0,
1503 /* 74*/ 0,
1504         0,
1505         0,
1506         0,
1507         0,
1508 /* 79*/ 0,
1509         0,
1510         0,
1511         0,
1512         MatLoad_MPISBAIJ,
1513 /* 84*/ 0,
1514         0,
1515         0,
1516         0,
1517         0,
1518 /* 89*/ 0,
1519         0,
1520         0,
1521         0,
1522         0,
1523 /* 94*/ 0,
1524         0,
1525         0,
1526         0,
1527         0,
1528 /* 99*/ 0,
1529         0,
1530         0,
1531         0,
1532         0,
1533 /*104*/ 0,
1534         MatRealPart_MPISBAIJ,
1535         MatImaginaryPart_MPISBAIJ,
1536         MatGetRowUpperTriangular_MPISBAIJ,
1537         MatRestoreRowUpperTriangular_MPISBAIJ,
1538 /*109*/ 0,
1539         0,
1540         0,
1541         0,
1542         0,
1543 /*114*/ 0,
1544         0,
1545         0,
1546         0,
1547         0,
1548 /*119*/ 0,
1549         0,
1550         0,
1551         0
1552 };
1553 
1554 
1555 EXTERN_C_BEGIN
1556 #undef __FUNCT__
1557 #define __FUNCT__ "MatGetDiagonalBlock_MPISBAIJ"
1558 PetscErrorCode  MatGetDiagonalBlock_MPISBAIJ(Mat A,Mat *a)
1559 {
1560   PetscFunctionBegin;
1561   *a = ((Mat_MPISBAIJ*)A->data)->A;
1562   PetscFunctionReturn(0);
1563 }
1564 EXTERN_C_END
1565 
1566 EXTERN_C_BEGIN
1567 #undef __FUNCT__
1568 #define __FUNCT__ "MatMPISBAIJSetPreallocation_MPISBAIJ"
1569 PetscErrorCode  MatMPISBAIJSetPreallocation_MPISBAIJ(Mat B,PetscInt bs,PetscInt d_nz,PetscInt *d_nnz,PetscInt o_nz,PetscInt *o_nnz)
1570 {
1571   Mat_MPISBAIJ   *b;
1572   PetscErrorCode ierr;
1573   PetscInt       i,mbs,Mbs;
1574 
1575   PetscFunctionBegin;
1576   if (d_nz == PETSC_DECIDE || d_nz == PETSC_DEFAULT) d_nz = 3;
1577   if (o_nz == PETSC_DECIDE || o_nz == PETSC_DEFAULT) o_nz = 1;
1578   if (d_nz < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"d_nz cannot be less than 0: value %D",d_nz);
1579   if (o_nz < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"o_nz cannot be less than 0: value %D",o_nz);
1580 
1581   ierr = PetscLayoutSetBlockSize(B->rmap,bs);CHKERRQ(ierr);
1582   ierr = PetscLayoutSetBlockSize(B->cmap,bs);CHKERRQ(ierr);
1583   ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr);
1584   ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr);
1585   ierr = PetscLayoutGetBlockSize(B->rmap,&bs);CHKERRQ(ierr);
1586 
1587   if (d_nnz) {
1588     for (i=0; i<B->rmap->n/bs; i++) {
1589       if (d_nnz[i] < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"d_nnz cannot be less than -1: local row %D value %D",i,d_nnz[i]);
1590     }
1591   }
1592   if (o_nnz) {
1593     for (i=0; i<B->rmap->n/bs; i++) {
1594       if (o_nnz[i] < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"o_nnz cannot be less than -1: local row %D value %D",i,o_nnz[i]);
1595     }
1596   }
1597 
1598   b   = (Mat_MPISBAIJ*)B->data;
1599   mbs = B->rmap->n/bs;
1600   Mbs = B->rmap->N/bs;
1601   if (mbs*bs != B->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"No of local rows %D must be divisible by blocksize %D",B->rmap->N,bs);
1602 
1603   B->rmap->bs = bs;
1604   b->bs2      = bs*bs;
1605   b->mbs      = mbs;
1606   b->nbs      = mbs;
1607   b->Mbs      = Mbs;
1608   b->Nbs      = Mbs;
1609 
1610   for (i=0; i<=b->size; i++) {
1611     b->rangebs[i] = B->rmap->range[i]/bs;
1612   }
1613   b->rstartbs = B->rmap->rstart/bs;
1614   b->rendbs   = B->rmap->rend/bs;
1615 
1616   b->cstartbs = B->cmap->rstart/bs;
1617   b->cendbs   = B->cmap->rend/bs;
1618 
1619   if (!B->preallocated) {
1620     ierr = MatCreate(PETSC_COMM_SELF,&b->A);CHKERRQ(ierr);
1621     ierr = MatSetSizes(b->A,B->rmap->n,B->cmap->n,B->rmap->n,B->cmap->n);CHKERRQ(ierr);
1622     ierr = MatSetType(b->A,MATSEQSBAIJ);CHKERRQ(ierr);
1623     ierr = PetscLogObjectParent(B,b->A);CHKERRQ(ierr);
1624     ierr = MatCreate(PETSC_COMM_SELF,&b->B);CHKERRQ(ierr);
1625     ierr = MatSetSizes(b->B,B->rmap->n,B->cmap->N,B->rmap->n,B->cmap->N);CHKERRQ(ierr);
1626     ierr = MatSetType(b->B,MATSEQBAIJ);CHKERRQ(ierr);
1627     ierr = PetscLogObjectParent(B,b->B);CHKERRQ(ierr);
1628     ierr = MatStashCreate_Private(PetscObjectComm((PetscObject)B),bs,&B->bstash);CHKERRQ(ierr);
1629   }
1630 
1631   ierr = MatSeqSBAIJSetPreallocation(b->A,bs,d_nz,d_nnz);CHKERRQ(ierr);
1632   ierr = MatSeqBAIJSetPreallocation(b->B,bs,o_nz,o_nnz);CHKERRQ(ierr);
1633 
1634   B->preallocated = PETSC_TRUE;
1635   PetscFunctionReturn(0);
1636 }
1637 EXTERN_C_END
1638 
1639 EXTERN_C_BEGIN
1640 #undef __FUNCT__
1641 #define __FUNCT__ "MatMPISBAIJSetPreallocationCSR_MPISBAIJ"
1642 PetscErrorCode MatMPISBAIJSetPreallocationCSR_MPISBAIJ(Mat B,PetscInt bs,const PetscInt ii[],const PetscInt jj[],const PetscScalar V[])
1643 {
1644   PetscInt       m,rstart,cstart,cend;
1645   PetscInt       i,j,d,nz,nz_max=0,*d_nnz=0,*o_nnz=0;
1646   const PetscInt *JJ    =0;
1647   PetscScalar    *values=0;
1648   PetscErrorCode ierr;
1649 
1650   PetscFunctionBegin;
1651   if (bs < 1) SETERRQ1(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_OUTOFRANGE,"Invalid block size specified, must be positive but it is %D",bs);
1652   ierr   = PetscLayoutSetBlockSize(B->rmap,bs);CHKERRQ(ierr);
1653   ierr   = PetscLayoutSetBlockSize(B->cmap,bs);CHKERRQ(ierr);
1654   ierr   = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr);
1655   ierr   = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr);
1656   ierr   = PetscLayoutGetBlockSize(B->rmap,&bs);CHKERRQ(ierr);
1657   m      = B->rmap->n/bs;
1658   rstart = B->rmap->rstart/bs;
1659   cstart = B->cmap->rstart/bs;
1660   cend   = B->cmap->rend/bs;
1661 
1662   if (ii[0]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"ii[0] must be 0 but it is %D",ii[0]);
1663   ierr = PetscMalloc2(m,PetscInt,&d_nnz,m,PetscInt,&o_nnz);CHKERRQ(ierr);
1664   for (i=0; i<m; i++) {
1665     nz = ii[i+1] - ii[i];
1666     if (nz < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local row %D has a negative number of columns %D",i,nz);
1667     nz_max = PetscMax(nz_max,nz);
1668     JJ     = jj + ii[i];
1669     for (j=0; j<nz; j++) {
1670       if (*JJ >= cstart) break;
1671       JJ++;
1672     }
1673     d = 0;
1674     for (; j<nz; j++) {
1675       if (*JJ++ >= cend) break;
1676       d++;
1677     }
1678     d_nnz[i] = d;
1679     o_nnz[i] = nz - d;
1680   }
1681   ierr = MatMPISBAIJSetPreallocation(B,bs,0,d_nnz,0,o_nnz);CHKERRQ(ierr);
1682   ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr);
1683 
1684   values = (PetscScalar*)V;
1685   if (!values) {
1686     ierr = PetscMalloc(bs*bs*nz_max*sizeof(PetscScalar),&values);CHKERRQ(ierr);
1687     ierr = PetscMemzero(values,bs*bs*nz_max*sizeof(PetscScalar));CHKERRQ(ierr);
1688   }
1689   for (i=0; i<m; i++) {
1690     PetscInt          row    = i + rstart;
1691     PetscInt          ncols  = ii[i+1] - ii[i];
1692     const PetscInt    *icols = jj + ii[i];
1693     const PetscScalar *svals = values + (V ? (bs*bs*ii[i]) : 0);
1694     ierr = MatSetValuesBlocked_MPISBAIJ(B,1,&row,ncols,icols,svals,INSERT_VALUES);CHKERRQ(ierr);
1695   }
1696 
1697   if (!V) { ierr = PetscFree(values);CHKERRQ(ierr); }
1698   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1699   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1700   ierr = MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
1701   PetscFunctionReturn(0);
1702 }
1703 EXTERN_C_END
1704 
1705 EXTERN_C_BEGIN
1706 #if defined(PETSC_HAVE_MUMPS)
1707 extern PetscErrorCode  MatGetFactor_sbaij_mumps(Mat,MatFactorType,Mat*);
1708 #endif
1709 #if defined(PETSC_HAVE_PASTIX)
1710 extern PetscErrorCode MatGetFactor_mpisbaij_pastix(Mat,MatFactorType,Mat*);
1711 #endif
1712 EXTERN_C_END
1713 
1714 /*MC
1715    MATMPISBAIJ - MATMPISBAIJ = "mpisbaij" - A matrix type to be used for distributed symmetric sparse block matrices,
1716    based on block compressed sparse row format.  Only the upper triangular portion of the "diagonal" portion of
1717    the matrix is stored.
1718 
1719   For complex numbers by default this matrix is symmetric, NOT Hermitian symmetric. To make it Hermitian symmetric you
1720   can call MatSetOption(Mat, MAT_HERMITIAN);
1721 
1722    Options Database Keys:
1723 . -mat_type mpisbaij - sets the matrix type to "mpisbaij" during a call to MatSetFromOptions()
1724 
1725   Level: beginner
1726 
1727 .seealso: MatCreateMPISBAIJ
1728 M*/
1729 
1730 EXTERN_C_BEGIN
1731 extern PetscErrorCode MatConvert_MPISBAIJ_MPISBSTRM(Mat,MatType,MatReuse,Mat*);
1732 EXTERN_C_END
1733 
1734 EXTERN_C_BEGIN
1735 #undef __FUNCT__
1736 #define __FUNCT__ "MatCreate_MPISBAIJ"
1737 PetscErrorCode  MatCreate_MPISBAIJ(Mat B)
1738 {
1739   Mat_MPISBAIJ   *b;
1740   PetscErrorCode ierr;
1741   PetscBool      flg;
1742 
1743   PetscFunctionBegin;
1744   ierr    = PetscNewLog(B,Mat_MPISBAIJ,&b);CHKERRQ(ierr);
1745   B->data = (void*)b;
1746   ierr    = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1747 
1748   B->ops->destroy = MatDestroy_MPISBAIJ;
1749   B->ops->view    = MatView_MPISBAIJ;
1750   B->assembled    = PETSC_FALSE;
1751   B->insertmode   = NOT_SET_VALUES;
1752 
1753   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)B),&b->rank);CHKERRQ(ierr);
1754   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)B),&b->size);CHKERRQ(ierr);
1755 
1756   /* build local table of row and column ownerships */
1757   ierr = PetscMalloc((b->size+2)*sizeof(PetscInt),&b->rangebs);CHKERRQ(ierr);
1758 
1759   /* build cache for off array entries formed */
1760   ierr = MatStashCreate_Private(PetscObjectComm((PetscObject)B),1,&B->stash);CHKERRQ(ierr);
1761 
1762   b->donotstash  = PETSC_FALSE;
1763   b->colmap      = NULL;
1764   b->garray      = NULL;
1765   b->roworiented = PETSC_TRUE;
1766 
1767   /* stuff used in block assembly */
1768   b->barray = 0;
1769 
1770   /* stuff used for matrix vector multiply */
1771   b->lvec    = 0;
1772   b->Mvctx   = 0;
1773   b->slvec0  = 0;
1774   b->slvec0b = 0;
1775   b->slvec1  = 0;
1776   b->slvec1a = 0;
1777   b->slvec1b = 0;
1778   b->sMvctx  = 0;
1779 
1780   /* stuff for MatGetRow() */
1781   b->rowindices   = 0;
1782   b->rowvalues    = 0;
1783   b->getrowactive = PETSC_FALSE;
1784 
1785   /* hash table stuff */
1786   b->ht           = 0;
1787   b->hd           = 0;
1788   b->ht_size      = 0;
1789   b->ht_flag      = PETSC_FALSE;
1790   b->ht_fact      = 0;
1791   b->ht_total_ct  = 0;
1792   b->ht_insert_ct = 0;
1793 
1794   /* stuff for MatGetSubMatrices_MPIBAIJ_local() */
1795   b->ijonly = PETSC_FALSE;
1796 
1797   b->in_loc = 0;
1798   b->v_loc  = 0;
1799   b->n_loc  = 0;
1800   ierr      = PetscOptionsBegin(PetscObjectComm((PetscObject)B),NULL,"Options for loading MPISBAIJ matrix 1","Mat");CHKERRQ(ierr);
1801   ierr      = PetscOptionsBool("-mat_use_hash_table","Use hash table to save memory in constructing matrix","MatSetOption",PETSC_FALSE,&flg,NULL);CHKERRQ(ierr);
1802   if (flg) {
1803     PetscReal fact = 1.39;
1804     ierr = MatSetOption(B,MAT_USE_HASH_TABLE,PETSC_TRUE);CHKERRQ(ierr);
1805     ierr = PetscOptionsReal("-mat_use_hash_table","Use hash table factor","MatMPIBAIJSetHashTableFactor",fact,&fact,NULL);CHKERRQ(ierr);
1806     if (fact <= 1.0) fact = 1.39;
1807     ierr = MatMPIBAIJSetHashTableFactor(B,fact);CHKERRQ(ierr);
1808     ierr = PetscInfo1(B,"Hash table Factor used %5.2f\n",fact);CHKERRQ(ierr);
1809   }
1810   ierr = PetscOptionsEnd();CHKERRQ(ierr);
1811 
1812 #if defined(PETSC_HAVE_PASTIX)
1813   ierr = PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_pastix_C","MatGetFactor_mpisbaij_pastix",MatGetFactor_mpisbaij_pastix);CHKERRQ(ierr);
1814 #endif
1815 #if defined(PETSC_HAVE_MUMPS)
1816   ierr = PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_mumps_C","MatGetFactor_sbaij_mumps",MatGetFactor_sbaij_mumps);CHKERRQ(ierr);
1817 #endif
1818   ierr = PetscObjectComposeFunction((PetscObject)B,"MatStoreValues_C","MatStoreValues_MPISBAIJ",MatStoreValues_MPISBAIJ);CHKERRQ(ierr);
1819   ierr = PetscObjectComposeFunction((PetscObject)B,"MatRetrieveValues_C","MatRetrieveValues_MPISBAIJ",MatRetrieveValues_MPISBAIJ);CHKERRQ(ierr);
1820   ierr = PetscObjectComposeFunction((PetscObject)B,"MatGetDiagonalBlock_C","MatGetDiagonalBlock_MPISBAIJ",MatGetDiagonalBlock_MPISBAIJ);CHKERRQ(ierr);
1821   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C","MatMPISBAIJSetPreallocation_MPISBAIJ",MatMPISBAIJSetPreallocation_MPISBAIJ);CHKERRQ(ierr);
1822   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPISBAIJSetPreallocationCSR_C","MatMPISBAIJSetPreallocationCSR_MPISBAIJ",MatMPISBAIJSetPreallocationCSR_MPISBAIJ);CHKERRQ(ierr);
1823   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpisbaij_mpisbstrm_C","MatConvert_MPISBAIJ_MPISBSTRM",MatConvert_MPISBAIJ_MPISBSTRM);CHKERRQ(ierr);
1824 
1825   B->symmetric                  = PETSC_TRUE;
1826   B->structurally_symmetric     = PETSC_TRUE;
1827   B->symmetric_set              = PETSC_TRUE;
1828   B->structurally_symmetric_set = PETSC_TRUE;
1829 
1830   ierr = PetscObjectChangeTypeName((PetscObject)B,MATMPISBAIJ);CHKERRQ(ierr);
1831   PetscFunctionReturn(0);
1832 }
1833 EXTERN_C_END
1834 
1835 /*MC
1836    MATSBAIJ - MATSBAIJ = "sbaij" - A matrix type to be used for symmetric block sparse matrices.
1837 
1838    This matrix type is identical to MATSEQSBAIJ when constructed with a single process communicator,
1839    and MATMPISBAIJ otherwise.
1840 
1841    Options Database Keys:
1842 . -mat_type sbaij - sets the matrix type to "sbaij" during a call to MatSetFromOptions()
1843 
1844   Level: beginner
1845 
1846 .seealso: MatCreateMPISBAIJ,MATSEQSBAIJ,MATMPISBAIJ
1847 M*/
1848 
1849 #undef __FUNCT__
1850 #define __FUNCT__ "MatMPISBAIJSetPreallocation"
1851 /*@C
1852    MatMPISBAIJSetPreallocation - For good matrix assembly performance
1853    the user should preallocate the matrix storage by setting the parameters
1854    d_nz (or d_nnz) and o_nz (or o_nnz).  By setting these parameters accurately,
1855    performance can be increased by more than a factor of 50.
1856 
1857    Collective on Mat
1858 
1859    Input Parameters:
1860 +  A - the matrix
1861 .  bs   - size of blockk
1862 .  d_nz  - number of block nonzeros per block row in diagonal portion of local
1863            submatrix  (same for all local rows)
1864 .  d_nnz - array containing the number of block nonzeros in the various block rows
1865            in the upper triangular and diagonal part of the in diagonal portion of the local
1866            (possibly different for each block row) or NULL.  If you plan to factor the matrix you must leave room
1867            for the diagonal entry and set a value even if it is zero.
1868 .  o_nz  - number of block nonzeros per block row in the off-diagonal portion of local
1869            submatrix (same for all local rows).
1870 -  o_nnz - array containing the number of nonzeros in the various block rows of the
1871            off-diagonal portion of the local submatrix that is right of the diagonal
1872            (possibly different for each block row) or NULL.
1873 
1874 
1875    Options Database Keys:
1876 .   -mat_no_unroll - uses code that does not unroll the loops in the
1877                      block calculations (much slower)
1878 .   -mat_block_size - size of the blocks to use
1879 
1880    Notes:
1881 
1882    If PETSC_DECIDE or  PETSC_DETERMINE is used for a particular argument on one processor
1883    than it must be used on all processors that share the object for that argument.
1884 
1885    If the *_nnz parameter is given then the *_nz parameter is ignored
1886 
1887    Storage Information:
1888    For a square global matrix we define each processor's diagonal portion
1889    to be its local rows and the corresponding columns (a square submatrix);
1890    each processor's off-diagonal portion encompasses the remainder of the
1891    local matrix (a rectangular submatrix).
1892 
1893    The user can specify preallocated storage for the diagonal part of
1894    the local submatrix with either d_nz or d_nnz (not both).  Set
1895    d_nz=PETSC_DEFAULT and d_nnz=NULL for PETSc to control dynamic
1896    memory allocation.  Likewise, specify preallocated storage for the
1897    off-diagonal part of the local submatrix with o_nz or o_nnz (not both).
1898 
1899    You can call MatGetInfo() to get information on how effective the preallocation was;
1900    for example the fields mallocs,nz_allocated,nz_used,nz_unneeded;
1901    You can also run with the option -info and look for messages with the string
1902    malloc in them to see if additional memory allocation was needed.
1903 
1904    Consider a processor that owns rows 3, 4 and 5 of a parallel matrix. In
1905    the figure below we depict these three local rows and all columns (0-11).
1906 
1907 .vb
1908            0 1 2 3 4 5 6 7 8 9 10 11
1909           -------------------
1910    row 3  |  . . . d d d o o o o o o
1911    row 4  |  . . . d d d o o o o o o
1912    row 5  |  . . . d d d o o o o o o
1913           -------------------
1914 .ve
1915 
1916    Thus, any entries in the d locations are stored in the d (diagonal)
1917    submatrix, and any entries in the o locations are stored in the
1918    o (off-diagonal) submatrix.  Note that the d matrix is stored in
1919    MatSeqSBAIJ format and the o submatrix in MATSEQBAIJ format.
1920 
1921    Now d_nz should indicate the number of block nonzeros per row in the upper triangular
1922    plus the diagonal part of the d matrix,
1923    and o_nz should indicate the number of block nonzeros per row in the o matrix
1924 
1925    In general, for PDE problems in which most nonzeros are near the diagonal,
1926    one expects d_nz >> o_nz.   For large problems you MUST preallocate memory
1927    or you will get TERRIBLE performance; see the users' manual chapter on
1928    matrices.
1929 
1930    Level: intermediate
1931 
1932 .keywords: matrix, block, aij, compressed row, sparse, parallel
1933 
1934 .seealso: MatCreate(), MatCreateSeqSBAIJ(), MatSetValues(), MatCreateBAIJ(), PetscSplitOwnership()
1935 @*/
1936 PetscErrorCode  MatMPISBAIJSetPreallocation(Mat B,PetscInt bs,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
1937 {
1938   PetscErrorCode ierr;
1939 
1940   PetscFunctionBegin;
1941   PetscValidHeaderSpecific(B,MAT_CLASSID,1);
1942   PetscValidType(B,1);
1943   PetscValidLogicalCollectiveInt(B,bs,2);
1944   ierr = PetscTryMethod(B,"MatMPISBAIJSetPreallocation_C",(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,bs,d_nz,d_nnz,o_nz,o_nnz));CHKERRQ(ierr);
1945   PetscFunctionReturn(0);
1946 }
1947 
1948 #undef __FUNCT__
1949 #define __FUNCT__ "MatCreateSBAIJ"
1950 /*@C
1951    MatCreateSBAIJ - Creates a sparse parallel matrix in symmetric block AIJ format
1952    (block compressed row).  For good matrix assembly performance
1953    the user should preallocate the matrix storage by setting the parameters
1954    d_nz (or d_nnz) and o_nz (or o_nnz).  By setting these parameters accurately,
1955    performance can be increased by more than a factor of 50.
1956 
1957    Collective on MPI_Comm
1958 
1959    Input Parameters:
1960 +  comm - MPI communicator
1961 .  bs   - size of blockk
1962 .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
1963            This value should be the same as the local size used in creating the
1964            y vector for the matrix-vector product y = Ax.
1965 .  n - number of local columns (or PETSC_DECIDE to have calculated if N is given)
1966            This value should be the same as the local size used in creating the
1967            x vector for the matrix-vector product y = Ax.
1968 .  M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
1969 .  N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
1970 .  d_nz  - number of block nonzeros per block row in diagonal portion of local
1971            submatrix  (same for all local rows)
1972 .  d_nnz - array containing the number of block nonzeros in the various block rows
1973            in the upper triangular portion of the in diagonal portion of the local
1974            (possibly different for each block block row) or NULL.
1975            If you plan to factor the matrix you must leave room for the diagonal entry and
1976            set its value even if it is zero.
1977 .  o_nz  - number of block nonzeros per block row in the off-diagonal portion of local
1978            submatrix (same for all local rows).
1979 -  o_nnz - array containing the number of nonzeros in the various block rows of the
1980            off-diagonal portion of the local submatrix (possibly different for
1981            each block row) or NULL.
1982 
1983    Output Parameter:
1984 .  A - the matrix
1985 
1986    Options Database Keys:
1987 .   -mat_no_unroll - uses code that does not unroll the loops in the
1988                      block calculations (much slower)
1989 .   -mat_block_size - size of the blocks to use
1990 .   -mat_mpi - use the parallel matrix data structures even on one processor
1991                (defaults to using SeqBAIJ format on one processor)
1992 
1993    It is recommended that one use the MatCreate(), MatSetType() and/or MatSetFromOptions(),
1994    MatXXXXSetPreallocation() paradgm instead of this routine directly.
1995    [MatXXXXSetPreallocation() is, for example, MatSeqAIJSetPreallocation]
1996 
1997    Notes:
1998    The number of rows and columns must be divisible by blocksize.
1999    This matrix type does not support complex Hermitian operation.
2000 
2001    The user MUST specify either the local or global matrix dimensions
2002    (possibly both).
2003 
2004    If PETSC_DECIDE or  PETSC_DETERMINE is used for a particular argument on one processor
2005    than it must be used on all processors that share the object for that argument.
2006 
2007    If the *_nnz parameter is given then the *_nz parameter is ignored
2008 
2009    Storage Information:
2010    For a square global matrix we define each processor's diagonal portion
2011    to be its local rows and the corresponding columns (a square submatrix);
2012    each processor's off-diagonal portion encompasses the remainder of the
2013    local matrix (a rectangular submatrix).
2014 
2015    The user can specify preallocated storage for the diagonal part of
2016    the local submatrix with either d_nz or d_nnz (not both).  Set
2017    d_nz=PETSC_DEFAULT and d_nnz=NULL for PETSc to control dynamic
2018    memory allocation.  Likewise, specify preallocated storage for the
2019    off-diagonal part of the local submatrix with o_nz or o_nnz (not both).
2020 
2021    Consider a processor that owns rows 3, 4 and 5 of a parallel matrix. In
2022    the figure below we depict these three local rows and all columns (0-11).
2023 
2024 .vb
2025            0 1 2 3 4 5 6 7 8 9 10 11
2026           -------------------
2027    row 3  |  . . . d d d o o o o o o
2028    row 4  |  . . . d d d o o o o o o
2029    row 5  |  . . . d d d o o o o o o
2030           -------------------
2031 .ve
2032 
2033    Thus, any entries in the d locations are stored in the d (diagonal)
2034    submatrix, and any entries in the o locations are stored in the
2035    o (off-diagonal) submatrix.  Note that the d matrix is stored in
2036    MatSeqSBAIJ format and the o submatrix in MATSEQBAIJ format.
2037 
2038    Now d_nz should indicate the number of block nonzeros per row in the upper triangular
2039    plus the diagonal part of the d matrix,
2040    and o_nz should indicate the number of block nonzeros per row in the o matrix.
2041    In general, for PDE problems in which most nonzeros are near the diagonal,
2042    one expects d_nz >> o_nz.   For large problems you MUST preallocate memory
2043    or you will get TERRIBLE performance; see the users' manual chapter on
2044    matrices.
2045 
2046    Level: intermediate
2047 
2048 .keywords: matrix, block, aij, compressed row, sparse, parallel
2049 
2050 .seealso: MatCreate(), MatCreateSeqSBAIJ(), MatSetValues(), MatCreateBAIJ()
2051 @*/
2052 
2053 PetscErrorCode  MatCreateSBAIJ(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[],Mat *A)
2054 {
2055   PetscErrorCode ierr;
2056   PetscMPIInt    size;
2057 
2058   PetscFunctionBegin;
2059   ierr = MatCreate(comm,A);CHKERRQ(ierr);
2060   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
2061   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
2062   if (size > 1) {
2063     ierr = MatSetType(*A,MATMPISBAIJ);CHKERRQ(ierr);
2064     ierr = MatMPISBAIJSetPreallocation(*A,bs,d_nz,d_nnz,o_nz,o_nnz);CHKERRQ(ierr);
2065   } else {
2066     ierr = MatSetType(*A,MATSEQSBAIJ);CHKERRQ(ierr);
2067     ierr = MatSeqSBAIJSetPreallocation(*A,bs,d_nz,d_nnz);CHKERRQ(ierr);
2068   }
2069   PetscFunctionReturn(0);
2070 }
2071 
2072 
2073 #undef __FUNCT__
2074 #define __FUNCT__ "MatDuplicate_MPISBAIJ"
2075 static PetscErrorCode MatDuplicate_MPISBAIJ(Mat matin,MatDuplicateOption cpvalues,Mat *newmat)
2076 {
2077   Mat            mat;
2078   Mat_MPISBAIJ   *a,*oldmat = (Mat_MPISBAIJ*)matin->data;
2079   PetscErrorCode ierr;
2080   PetscInt       len=0,nt,bs=matin->rmap->bs,mbs=oldmat->mbs;
2081   PetscScalar    *array;
2082 
2083   PetscFunctionBegin;
2084   *newmat = 0;
2085 
2086   ierr = MatCreate(PetscObjectComm((PetscObject)matin),&mat);CHKERRQ(ierr);
2087   ierr = MatSetSizes(mat,matin->rmap->n,matin->cmap->n,matin->rmap->N,matin->cmap->N);CHKERRQ(ierr);
2088   ierr = MatSetType(mat,((PetscObject)matin)->type_name);CHKERRQ(ierr);
2089   ierr = PetscMemcpy(mat->ops,matin->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
2090   ierr = PetscLayoutReference(matin->rmap,&mat->rmap);CHKERRQ(ierr);
2091   ierr = PetscLayoutReference(matin->cmap,&mat->cmap);CHKERRQ(ierr);
2092 
2093   mat->factortype   = matin->factortype;
2094   mat->preallocated = PETSC_TRUE;
2095   mat->assembled    = PETSC_TRUE;
2096   mat->insertmode   = NOT_SET_VALUES;
2097 
2098   a      = (Mat_MPISBAIJ*)mat->data;
2099   a->bs2 = oldmat->bs2;
2100   a->mbs = oldmat->mbs;
2101   a->nbs = oldmat->nbs;
2102   a->Mbs = oldmat->Mbs;
2103   a->Nbs = oldmat->Nbs;
2104 
2105 
2106   a->size         = oldmat->size;
2107   a->rank         = oldmat->rank;
2108   a->donotstash   = oldmat->donotstash;
2109   a->roworiented  = oldmat->roworiented;
2110   a->rowindices   = 0;
2111   a->rowvalues    = 0;
2112   a->getrowactive = PETSC_FALSE;
2113   a->barray       = 0;
2114   a->rstartbs     = oldmat->rstartbs;
2115   a->rendbs       = oldmat->rendbs;
2116   a->cstartbs     = oldmat->cstartbs;
2117   a->cendbs       = oldmat->cendbs;
2118 
2119   /* hash table stuff */
2120   a->ht           = 0;
2121   a->hd           = 0;
2122   a->ht_size      = 0;
2123   a->ht_flag      = oldmat->ht_flag;
2124   a->ht_fact      = oldmat->ht_fact;
2125   a->ht_total_ct  = 0;
2126   a->ht_insert_ct = 0;
2127 
2128   ierr = PetscMemcpy(a->rangebs,oldmat->rangebs,(a->size+2)*sizeof(PetscInt));CHKERRQ(ierr);
2129   if (oldmat->colmap) {
2130 #if defined(PETSC_USE_CTABLE)
2131     ierr = PetscTableCreateCopy(oldmat->colmap,&a->colmap);CHKERRQ(ierr);
2132 #else
2133     ierr = PetscMalloc((a->Nbs)*sizeof(PetscInt),&a->colmap);CHKERRQ(ierr);
2134     ierr = PetscLogObjectMemory(mat,(a->Nbs)*sizeof(PetscInt));CHKERRQ(ierr);
2135     ierr = PetscMemcpy(a->colmap,oldmat->colmap,(a->Nbs)*sizeof(PetscInt));CHKERRQ(ierr);
2136 #endif
2137   } else a->colmap = 0;
2138 
2139   if (oldmat->garray && (len = ((Mat_SeqBAIJ*)(oldmat->B->data))->nbs)) {
2140     ierr = PetscMalloc(len*sizeof(PetscInt),&a->garray);CHKERRQ(ierr);
2141     ierr = PetscLogObjectMemory(mat,len*sizeof(PetscInt));CHKERRQ(ierr);
2142     ierr = PetscMemcpy(a->garray,oldmat->garray,len*sizeof(PetscInt));CHKERRQ(ierr);
2143   } else a->garray = 0;
2144 
2145   ierr = MatStashCreate_Private(PetscObjectComm((PetscObject)matin),matin->rmap->bs,&mat->bstash);CHKERRQ(ierr);
2146   ierr = VecDuplicate(oldmat->lvec,&a->lvec);CHKERRQ(ierr);
2147   ierr = PetscLogObjectParent(mat,a->lvec);CHKERRQ(ierr);
2148   ierr = VecScatterCopy(oldmat->Mvctx,&a->Mvctx);CHKERRQ(ierr);
2149   ierr = PetscLogObjectParent(mat,a->Mvctx);CHKERRQ(ierr);
2150 
2151   ierr =  VecDuplicate(oldmat->slvec0,&a->slvec0);CHKERRQ(ierr);
2152   ierr = PetscLogObjectParent(mat,a->slvec0);CHKERRQ(ierr);
2153   ierr =  VecDuplicate(oldmat->slvec1,&a->slvec1);CHKERRQ(ierr);
2154   ierr = PetscLogObjectParent(mat,a->slvec1);CHKERRQ(ierr);
2155 
2156   ierr = VecGetLocalSize(a->slvec1,&nt);CHKERRQ(ierr);
2157   ierr = VecGetArray(a->slvec1,&array);CHKERRQ(ierr);
2158   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,bs*mbs,array,&a->slvec1a);CHKERRQ(ierr);
2159   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,nt-bs*mbs,array+bs*mbs,&a->slvec1b);CHKERRQ(ierr);
2160   ierr = VecRestoreArray(a->slvec1,&array);CHKERRQ(ierr);
2161   ierr = VecGetArray(a->slvec0,&array);CHKERRQ(ierr);
2162   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,nt-bs*mbs,array+bs*mbs,&a->slvec0b);CHKERRQ(ierr);
2163   ierr = VecRestoreArray(a->slvec0,&array);CHKERRQ(ierr);
2164   ierr = PetscLogObjectParent(mat,a->slvec0);CHKERRQ(ierr);
2165   ierr = PetscLogObjectParent(mat,a->slvec1);CHKERRQ(ierr);
2166   ierr = PetscLogObjectParent(mat,a->slvec0b);CHKERRQ(ierr);
2167   ierr = PetscLogObjectParent(mat,a->slvec1a);CHKERRQ(ierr);
2168   ierr = PetscLogObjectParent(mat,a->slvec1b);CHKERRQ(ierr);
2169 
2170   /* ierr =  VecScatterCopy(oldmat->sMvctx,&a->sMvctx); - not written yet, replaced by the lazy trick: */
2171   ierr      = PetscObjectReference((PetscObject)oldmat->sMvctx);CHKERRQ(ierr);
2172   a->sMvctx = oldmat->sMvctx;
2173   ierr      = PetscLogObjectParent(mat,a->sMvctx);CHKERRQ(ierr);
2174 
2175   ierr    =  MatDuplicate(oldmat->A,cpvalues,&a->A);CHKERRQ(ierr);
2176   ierr    = PetscLogObjectParent(mat,a->A);CHKERRQ(ierr);
2177   ierr    =  MatDuplicate(oldmat->B,cpvalues,&a->B);CHKERRQ(ierr);
2178   ierr    = PetscLogObjectParent(mat,a->B);CHKERRQ(ierr);
2179   ierr    = PetscFunctionListDuplicate(((PetscObject)matin)->qlist,&((PetscObject)mat)->qlist);CHKERRQ(ierr);
2180   *newmat = mat;
2181   PetscFunctionReturn(0);
2182 }
2183 
2184 #undef __FUNCT__
2185 #define __FUNCT__ "MatLoad_MPISBAIJ"
2186 PetscErrorCode MatLoad_MPISBAIJ(Mat newmat,PetscViewer viewer)
2187 {
2188   PetscErrorCode ierr;
2189   PetscInt       i,nz,j,rstart,rend;
2190   PetscScalar    *vals,*buf;
2191   MPI_Comm       comm;
2192   MPI_Status     status;
2193   PetscMPIInt    rank,size,tag = ((PetscObject)viewer)->tag,*sndcounts = 0,*browners,maxnz,*rowners,mmbs;
2194   PetscInt       header[4],*rowlengths = 0,M,N,m,*cols,*locrowlens;
2195   PetscInt       *procsnz = 0,jj,*mycols,*ibuf;
2196   PetscInt       bs       =1,Mbs,mbs,extra_rows;
2197   PetscInt       *dlens,*odlens,*mask,*masked1,*masked2,rowcount,odcount;
2198   PetscInt       dcount,kmax,k,nzcount,tmp,sizesset=1,grows,gcols;
2199   int            fd;
2200 
2201   PetscFunctionBegin;
2202   ierr = PetscObjectGetComm((PetscObject)viewer,&comm);CHKERRQ(ierr);
2203   ierr = PetscOptionsBegin(comm,NULL,"Options for loading MPISBAIJ matrix 2","Mat");CHKERRQ(ierr);
2204   ierr = PetscOptionsInt("-matload_block_size","Set the blocksize used to store the matrix","MatLoad",bs,&bs,NULL);CHKERRQ(ierr);
2205   ierr = PetscOptionsEnd();CHKERRQ(ierr);
2206 
2207   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
2208   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
2209   if (!rank) {
2210     ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
2211     ierr = PetscBinaryRead(fd,(char*)header,4,PETSC_INT);CHKERRQ(ierr);
2212     if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"not matrix object");
2213     if (header[3] < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as MPISBAIJ");
2214   }
2215 
2216   if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) sizesset = 0;
2217 
2218   ierr = MPI_Bcast(header+1,3,MPIU_INT,0,comm);CHKERRQ(ierr);
2219   M    = header[1];
2220   N    = header[2];
2221 
2222   /* If global rows/cols are set to PETSC_DECIDE, set it to the sizes given in the file */
2223   if (sizesset && newmat->rmap->N < 0) newmat->rmap->N = M;
2224   if (sizesset && newmat->cmap->N < 0) newmat->cmap->N = N;
2225 
2226   /* If global sizes are set, check if they are consistent with that given in the file */
2227   if (sizesset) {
2228     ierr = MatGetSize(newmat,&grows,&gcols);CHKERRQ(ierr);
2229   }
2230   if (sizesset && newmat->rmap->N != grows) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Inconsistent # of rows:Matrix in file has (%d) and input matrix has (%d)",M,grows);
2231   if (sizesset && newmat->cmap->N != gcols) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Inconsistent # of cols:Matrix in file has (%d) and input matrix has (%d)",N,gcols);
2232 
2233   if (M != N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Can only do square matrices");
2234 
2235   /*
2236      This code adds extra rows to make sure the number of rows is
2237      divisible by the blocksize
2238   */
2239   Mbs        = M/bs;
2240   extra_rows = bs - M + bs*(Mbs);
2241   if (extra_rows == bs) extra_rows = 0;
2242   else                  Mbs++;
2243   if (extra_rows &&!rank) {
2244     ierr = PetscInfo(viewer,"Padding loaded matrix to match blocksize\n");CHKERRQ(ierr);
2245   }
2246 
2247   /* determine ownership of all rows */
2248   if (newmat->rmap->n < 0) { /* PETSC_DECIDE */
2249     mbs = Mbs/size + ((Mbs % size) > rank);
2250     m   = mbs*bs;
2251   } else { /* User Set */
2252     m   = newmat->rmap->n;
2253     mbs = m/bs;
2254   }
2255   ierr       = PetscMalloc2(size+1,PetscMPIInt,&rowners,size+1,PetscMPIInt,&browners);CHKERRQ(ierr);
2256   ierr       = PetscMPIIntCast(mbs,&mmbs);CHKERRQ(ierr);
2257   ierr       = MPI_Allgather(&mmbs,1,MPI_INT,rowners+1,1,MPI_INT,comm);CHKERRQ(ierr);
2258   rowners[0] = 0;
2259   for (i=2; i<=size; i++) rowners[i] += rowners[i-1];
2260   for (i=0; i<=size; i++) browners[i] = rowners[i]*bs;
2261   rstart = rowners[rank];
2262   rend   = rowners[rank+1];
2263 
2264   /* distribute row lengths to all processors */
2265   ierr = PetscMalloc((rend-rstart)*bs*sizeof(PetscInt),&locrowlens);CHKERRQ(ierr);
2266   if (!rank) {
2267     ierr = PetscMalloc((M+extra_rows)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
2268     ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
2269     for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;
2270     ierr = PetscMalloc(size*sizeof(PetscMPIInt),&sndcounts);CHKERRQ(ierr);
2271     for (i=0; i<size; i++) sndcounts[i] = browners[i+1] - browners[i];
2272     ierr = MPI_Scatterv(rowlengths,sndcounts,browners,MPIU_INT,locrowlens,(rend-rstart)*bs,MPIU_INT,0,comm);CHKERRQ(ierr);
2273     ierr = PetscFree(sndcounts);CHKERRQ(ierr);
2274   } else {
2275     ierr = MPI_Scatterv(0,0,0,MPIU_INT,locrowlens,(rend-rstart)*bs,MPIU_INT,0,comm);CHKERRQ(ierr);
2276   }
2277 
2278   if (!rank) {   /* procs[0] */
2279     /* calculate the number of nonzeros on each processor */
2280     ierr = PetscMalloc(size*sizeof(PetscInt),&procsnz);CHKERRQ(ierr);
2281     ierr = PetscMemzero(procsnz,size*sizeof(PetscInt));CHKERRQ(ierr);
2282     for (i=0; i<size; i++) {
2283       for (j=rowners[i]*bs; j< rowners[i+1]*bs; j++) {
2284         procsnz[i] += rowlengths[j];
2285       }
2286     }
2287     ierr = PetscFree(rowlengths);CHKERRQ(ierr);
2288 
2289     /* determine max buffer needed and allocate it */
2290     maxnz = 0;
2291     for (i=0; i<size; i++) {
2292       maxnz = PetscMax(maxnz,procsnz[i]);
2293     }
2294     ierr = PetscMalloc(maxnz*sizeof(PetscInt),&cols);CHKERRQ(ierr);
2295 
2296     /* read in my part of the matrix column indices  */
2297     nz     = procsnz[0];
2298     ierr   = PetscMalloc(nz*sizeof(PetscInt),&ibuf);CHKERRQ(ierr);
2299     mycols = ibuf;
2300     if (size == 1) nz -= extra_rows;
2301     ierr = PetscBinaryRead(fd,mycols,nz,PETSC_INT);CHKERRQ(ierr);
2302     if (size == 1) {
2303       for (i=0; i< extra_rows; i++) mycols[nz+i] = M+i;
2304     }
2305 
2306     /* read in every ones (except the last) and ship off */
2307     for (i=1; i<size-1; i++) {
2308       nz   = procsnz[i];
2309       ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
2310       ierr = MPI_Send(cols,nz,MPIU_INT,i,tag,comm);CHKERRQ(ierr);
2311     }
2312     /* read in the stuff for the last proc */
2313     if (size != 1) {
2314       nz   = procsnz[size-1] - extra_rows;  /* the extra rows are not on the disk */
2315       ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr);
2316       for (i=0; i<extra_rows; i++) cols[nz+i] = M+i;
2317       ierr = MPI_Send(cols,nz+extra_rows,MPIU_INT,size-1,tag,comm);CHKERRQ(ierr);
2318     }
2319     ierr = PetscFree(cols);CHKERRQ(ierr);
2320   } else {  /* procs[i], i>0 */
2321     /* determine buffer space needed for message */
2322     nz = 0;
2323     for (i=0; i<m; i++) nz += locrowlens[i];
2324     ierr   = PetscMalloc(nz*sizeof(PetscInt),&ibuf);CHKERRQ(ierr);
2325     mycols = ibuf;
2326     /* receive message of column indices*/
2327     ierr = MPI_Recv(mycols,nz,MPIU_INT,0,tag,comm,&status);CHKERRQ(ierr);
2328     ierr = MPI_Get_count(&status,MPIU_INT,&maxnz);CHKERRQ(ierr);
2329     if (maxnz != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
2330   }
2331 
2332   /* loop over local rows, determining number of off diagonal entries */
2333   ierr     = PetscMalloc2(rend-rstart,PetscInt,&dlens,rend-rstart,PetscInt,&odlens);CHKERRQ(ierr);
2334   ierr     = PetscMalloc3(Mbs,PetscInt,&mask,Mbs,PetscInt,&masked1,Mbs,PetscInt,&masked2);CHKERRQ(ierr);
2335   ierr     = PetscMemzero(mask,Mbs*sizeof(PetscInt));CHKERRQ(ierr);
2336   ierr     = PetscMemzero(masked1,Mbs*sizeof(PetscInt));CHKERRQ(ierr);
2337   ierr     = PetscMemzero(masked2,Mbs*sizeof(PetscInt));CHKERRQ(ierr);
2338   rowcount = 0;
2339   nzcount  = 0;
2340   for (i=0; i<mbs; i++) {
2341     dcount  = 0;
2342     odcount = 0;
2343     for (j=0; j<bs; j++) {
2344       kmax = locrowlens[rowcount];
2345       for (k=0; k<kmax; k++) {
2346         tmp = mycols[nzcount++]/bs; /* block col. index */
2347         if (!mask[tmp]) {
2348           mask[tmp] = 1;
2349           if (tmp < rstart || tmp >= rend) masked2[odcount++] = tmp; /* entry in off-diag portion */
2350           else masked1[dcount++] = tmp; /* entry in diag portion */
2351         }
2352       }
2353       rowcount++;
2354     }
2355 
2356     dlens[i]  = dcount;  /* d_nzz[i] */
2357     odlens[i] = odcount; /* o_nzz[i] */
2358 
2359     /* zero out the mask elements we set */
2360     for (j=0; j<dcount; j++) mask[masked1[j]] = 0;
2361     for (j=0; j<odcount; j++) mask[masked2[j]] = 0;
2362   }
2363   if (!sizesset) {
2364     ierr = MatSetSizes(newmat,m,m,M+extra_rows,N+extra_rows);CHKERRQ(ierr);
2365   }
2366   ierr = MatMPISBAIJSetPreallocation(newmat,bs,0,dlens,0,odlens);CHKERRQ(ierr);
2367   ierr = MatSetOption(newmat,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE);CHKERRQ(ierr);
2368 
2369   if (!rank) {
2370     ierr = PetscMalloc(maxnz*sizeof(PetscScalar),&buf);CHKERRQ(ierr);
2371     /* read in my part of the matrix numerical values  */
2372     nz     = procsnz[0];
2373     vals   = buf;
2374     mycols = ibuf;
2375     if (size == 1) nz -= extra_rows;
2376     ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
2377     if (size == 1) {
2378       for (i=0; i< extra_rows; i++) vals[nz+i] = 1.0;
2379     }
2380 
2381     /* insert into matrix */
2382     jj = rstart*bs;
2383     for (i=0; i<m; i++) {
2384       ierr    = MatSetValues(newmat,1,&jj,locrowlens[i],mycols,vals,INSERT_VALUES);CHKERRQ(ierr);
2385       mycols += locrowlens[i];
2386       vals   += locrowlens[i];
2387       jj++;
2388     }
2389 
2390     /* read in other processors (except the last one) and ship out */
2391     for (i=1; i<size-1; i++) {
2392       nz   = procsnz[i];
2393       vals = buf;
2394       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
2395       ierr = MPI_Send(vals,nz,MPIU_SCALAR,i,((PetscObject)newmat)->tag,comm);CHKERRQ(ierr);
2396     }
2397     /* the last proc */
2398     if (size != 1) {
2399       nz   = procsnz[i] - extra_rows;
2400       vals = buf;
2401       ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr);
2402       for (i=0; i<extra_rows; i++) vals[nz+i] = 1.0;
2403       ierr = MPI_Send(vals,nz+extra_rows,MPIU_SCALAR,size-1,((PetscObject)newmat)->tag,comm);CHKERRQ(ierr);
2404     }
2405     ierr = PetscFree(procsnz);CHKERRQ(ierr);
2406 
2407   } else {
2408     /* receive numeric values */
2409     ierr = PetscMalloc(nz*sizeof(PetscScalar),&buf);CHKERRQ(ierr);
2410 
2411     /* receive message of values*/
2412     vals   = buf;
2413     mycols = ibuf;
2414     ierr   = MPI_Recv(vals,nz,MPIU_SCALAR,0,((PetscObject)newmat)->tag,comm,&status);CHKERRQ(ierr);
2415     ierr   = MPI_Get_count(&status,MPIU_SCALAR,&maxnz);CHKERRQ(ierr);
2416     if (maxnz != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"something is wrong with file");
2417 
2418     /* insert into matrix */
2419     jj = rstart*bs;
2420     for (i=0; i<m; i++) {
2421       ierr    = MatSetValues_MPISBAIJ(newmat,1,&jj,locrowlens[i],mycols,vals,INSERT_VALUES);CHKERRQ(ierr);
2422       mycols += locrowlens[i];
2423       vals   += locrowlens[i];
2424       jj++;
2425     }
2426   }
2427 
2428   ierr = PetscFree(locrowlens);CHKERRQ(ierr);
2429   ierr = PetscFree(buf);CHKERRQ(ierr);
2430   ierr = PetscFree(ibuf);CHKERRQ(ierr);
2431   ierr = PetscFree2(rowners,browners);CHKERRQ(ierr);
2432   ierr = PetscFree2(dlens,odlens);CHKERRQ(ierr);
2433   ierr = PetscFree3(mask,masked1,masked2);CHKERRQ(ierr);
2434   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2435   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2436   PetscFunctionReturn(0);
2437 }
2438 
2439 #undef __FUNCT__
2440 #define __FUNCT__ "MatMPISBAIJSetHashTableFactor"
2441 /*XXXXX@
2442    MatMPISBAIJSetHashTableFactor - Sets the factor required to compute the size of the HashTable.
2443 
2444    Input Parameters:
2445 .  mat  - the matrix
2446 .  fact - factor
2447 
2448    Not Collective on Mat, each process can have a different hash factor
2449 
2450    Level: advanced
2451 
2452   Notes:
2453    This can also be set by the command line option: -mat_use_hash_table fact
2454 
2455 .keywords: matrix, hashtable, factor, HT
2456 
2457 .seealso: MatSetOption()
2458 @XXXXX*/
2459 
2460 
2461 #undef __FUNCT__
2462 #define __FUNCT__ "MatGetRowMaxAbs_MPISBAIJ"
2463 PetscErrorCode MatGetRowMaxAbs_MPISBAIJ(Mat A,Vec v,PetscInt idx[])
2464 {
2465   Mat_MPISBAIJ   *a = (Mat_MPISBAIJ*)A->data;
2466   Mat_SeqBAIJ    *b = (Mat_SeqBAIJ*)(a->B)->data;
2467   PetscReal      atmp;
2468   PetscReal      *work,*svalues,*rvalues;
2469   PetscErrorCode ierr;
2470   PetscInt       i,bs,mbs,*bi,*bj,brow,j,ncols,krow,kcol,col,row,Mbs,bcol;
2471   PetscMPIInt    rank,size;
2472   PetscInt       *rowners_bs,dest,count,source;
2473   PetscScalar    *va;
2474   MatScalar      *ba;
2475   MPI_Status     stat;
2476 
2477   PetscFunctionBegin;
2478   if (idx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Send email to petsc-maint@mcs.anl.gov");
2479   ierr = MatGetRowMaxAbs(a->A,v,NULL);CHKERRQ(ierr);
2480   ierr = VecGetArray(v,&va);CHKERRQ(ierr);
2481 
2482   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
2483   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)A),&rank);CHKERRQ(ierr);
2484 
2485   bs  = A->rmap->bs;
2486   mbs = a->mbs;
2487   Mbs = a->Mbs;
2488   ba  = b->a;
2489   bi  = b->i;
2490   bj  = b->j;
2491 
2492   /* find ownerships */
2493   rowners_bs = A->rmap->range;
2494 
2495   /* each proc creates an array to be distributed */
2496   ierr = PetscMalloc(bs*Mbs*sizeof(PetscReal),&work);CHKERRQ(ierr);
2497   ierr = PetscMemzero(work,bs*Mbs*sizeof(PetscReal));CHKERRQ(ierr);
2498 
2499   /* row_max for B */
2500   if (rank != size-1) {
2501     for (i=0; i<mbs; i++) {
2502       ncols = bi[1] - bi[0]; bi++;
2503       brow  = bs*i;
2504       for (j=0; j<ncols; j++) {
2505         bcol = bs*(*bj);
2506         for (kcol=0; kcol<bs; kcol++) {
2507           col  = bcol + kcol;                /* local col index */
2508           col += rowners_bs[rank+1];      /* global col index */
2509           for (krow=0; krow<bs; krow++) {
2510             atmp = PetscAbsScalar(*ba); ba++;
2511             row  = brow + krow;   /* local row index */
2512             if (PetscRealPart(va[row]) < atmp) va[row] = atmp;
2513             if (work[col] < atmp) work[col] = atmp;
2514           }
2515         }
2516         bj++;
2517       }
2518     }
2519 
2520     /* send values to its owners */
2521     for (dest=rank+1; dest<size; dest++) {
2522       svalues = work + rowners_bs[dest];
2523       count   = rowners_bs[dest+1]-rowners_bs[dest];
2524       ierr    = MPI_Send(svalues,count,MPIU_REAL,dest,rank,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
2525     }
2526   }
2527 
2528   /* receive values */
2529   if (rank) {
2530     rvalues = work;
2531     count   = rowners_bs[rank+1]-rowners_bs[rank];
2532     for (source=0; source<rank; source++) {
2533       ierr = MPI_Recv(rvalues,count,MPIU_REAL,MPI_ANY_SOURCE,MPI_ANY_TAG,PetscObjectComm((PetscObject)A),&stat);CHKERRQ(ierr);
2534       /* process values */
2535       for (i=0; i<count; i++) {
2536         if (PetscRealPart(va[i]) < rvalues[i]) va[i] = rvalues[i];
2537       }
2538     }
2539   }
2540 
2541   ierr = VecRestoreArray(v,&va);CHKERRQ(ierr);
2542   ierr = PetscFree(work);CHKERRQ(ierr);
2543   PetscFunctionReturn(0);
2544 }
2545 
2546 #undef __FUNCT__
2547 #define __FUNCT__ "MatSOR_MPISBAIJ"
2548 PetscErrorCode MatSOR_MPISBAIJ(Mat matin,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
2549 {
2550   Mat_MPISBAIJ      *mat = (Mat_MPISBAIJ*)matin->data;
2551   PetscErrorCode    ierr;
2552   PetscInt          mbs=mat->mbs,bs=matin->rmap->bs;
2553   PetscScalar       *x,*ptr,*from;
2554   Vec               bb1;
2555   const PetscScalar *b;
2556 
2557   PetscFunctionBegin;
2558   if (its <= 0 || lits <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
2559   if (bs > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented");
2560 
2561   if (flag == SOR_APPLY_UPPER) {
2562     ierr = (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx);CHKERRQ(ierr);
2563     PetscFunctionReturn(0);
2564   }
2565 
2566   if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP) {
2567     if (flag & SOR_ZERO_INITIAL_GUESS) {
2568       ierr = (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,lits,xx);CHKERRQ(ierr);
2569       its--;
2570     }
2571 
2572     ierr = VecDuplicate(bb,&bb1);CHKERRQ(ierr);
2573     while (its--) {
2574 
2575       /* lower triangular part: slvec0b = - B^T*xx */
2576       ierr = (*mat->B->ops->multtranspose)(mat->B,xx,mat->slvec0b);CHKERRQ(ierr);
2577 
2578       /* copy xx into slvec0a */
2579       ierr = VecGetArray(mat->slvec0,&ptr);CHKERRQ(ierr);
2580       ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2581       ierr = PetscMemcpy(ptr,x,bs*mbs*sizeof(MatScalar));CHKERRQ(ierr);
2582       ierr = VecRestoreArray(mat->slvec0,&ptr);CHKERRQ(ierr);
2583 
2584       ierr = VecScale(mat->slvec0,-1.0);CHKERRQ(ierr);
2585 
2586       /* copy bb into slvec1a */
2587       ierr = VecGetArray(mat->slvec1,&ptr);CHKERRQ(ierr);
2588       ierr = VecGetArrayRead(bb,&b);CHKERRQ(ierr);
2589       ierr = PetscMemcpy(ptr,b,bs*mbs*sizeof(MatScalar));CHKERRQ(ierr);
2590       ierr = VecRestoreArray(mat->slvec1,&ptr);CHKERRQ(ierr);
2591 
2592       /* set slvec1b = 0 */
2593       ierr = VecSet(mat->slvec1b,0.0);CHKERRQ(ierr);
2594 
2595       ierr = VecScatterBegin(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2596       ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2597       ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr);
2598       ierr = VecScatterEnd(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2599 
2600       /* upper triangular part: bb1 = bb1 - B*x */
2601       ierr = (*mat->B->ops->multadd)(mat->B,mat->slvec1b,mat->slvec1a,bb1);CHKERRQ(ierr);
2602 
2603       /* local diagonal sweep */
2604       ierr = (*mat->A->ops->sor)(mat->A,bb1,omega,SOR_SYMMETRIC_SWEEP,fshift,lits,lits,xx);CHKERRQ(ierr);
2605     }
2606     ierr = VecDestroy(&bb1);CHKERRQ(ierr);
2607   } else if ((flag & SOR_LOCAL_FORWARD_SWEEP) && (its == 1) && (flag & SOR_ZERO_INITIAL_GUESS)) {
2608     ierr = (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx);CHKERRQ(ierr);
2609   } else if ((flag & SOR_LOCAL_BACKWARD_SWEEP) && (its == 1) && (flag & SOR_ZERO_INITIAL_GUESS)) {
2610     ierr = (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,1,xx);CHKERRQ(ierr);
2611   } else if (flag & SOR_EISENSTAT) {
2612     Vec               xx1;
2613     PetscBool         hasop;
2614     const PetscScalar *diag;
2615     PetscScalar       *sl,scale = (omega - 2.0)/omega;
2616     PetscInt          i,n;
2617 
2618     if (!mat->xx1) {
2619       ierr = VecDuplicate(bb,&mat->xx1);CHKERRQ(ierr);
2620       ierr = VecDuplicate(bb,&mat->bb1);CHKERRQ(ierr);
2621     }
2622     xx1 = mat->xx1;
2623     bb1 = mat->bb1;
2624 
2625     ierr = (*mat->A->ops->sor)(mat->A,bb,omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_BACKWARD_SWEEP),fshift,lits,1,xx);CHKERRQ(ierr);
2626 
2627     if (!mat->diag) {
2628       /* this is wrong for same matrix with new nonzero values */
2629       ierr = MatGetVecs(matin,&mat->diag,NULL);CHKERRQ(ierr);
2630       ierr = MatGetDiagonal(matin,mat->diag);CHKERRQ(ierr);
2631     }
2632     ierr = MatHasOperation(matin,MATOP_MULT_DIAGONAL_BLOCK,&hasop);CHKERRQ(ierr);
2633 
2634     if (hasop) {
2635       ierr = MatMultDiagonalBlock(matin,xx,bb1);CHKERRQ(ierr);
2636       ierr = VecAYPX(mat->slvec1a,scale,bb);CHKERRQ(ierr);
2637     } else {
2638       /*
2639           These two lines are replaced by code that may be a bit faster for a good compiler
2640       ierr = VecPointwiseMult(mat->slvec1a,mat->diag,xx);CHKERRQ(ierr);
2641       ierr = VecAYPX(mat->slvec1a,scale,bb);CHKERRQ(ierr);
2642       */
2643       ierr = VecGetArray(mat->slvec1a,&sl);CHKERRQ(ierr);
2644       ierr = VecGetArrayRead(mat->diag,&diag);CHKERRQ(ierr);
2645       ierr = VecGetArrayRead(bb,&b);CHKERRQ(ierr);
2646       ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2647       ierr = VecGetLocalSize(xx,&n);CHKERRQ(ierr);
2648       if (omega == 1.0) {
2649         for (i=0; i<n; i++) sl[i] = b[i] - diag[i]*x[i];
2650         ierr = PetscLogFlops(2.0*n);CHKERRQ(ierr);
2651       } else {
2652         for (i=0; i<n; i++) sl[i] = b[i] + scale*diag[i]*x[i];
2653         ierr = PetscLogFlops(3.0*n);CHKERRQ(ierr);
2654       }
2655       ierr = VecRestoreArray(mat->slvec1a,&sl);CHKERRQ(ierr);
2656       ierr = VecRestoreArrayRead(mat->diag,&diag);CHKERRQ(ierr);
2657       ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr);
2658       ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2659     }
2660 
2661     /* multiply off-diagonal portion of matrix */
2662     ierr = VecSet(mat->slvec1b,0.0);CHKERRQ(ierr);
2663     ierr = (*mat->B->ops->multtranspose)(mat->B,xx,mat->slvec0b);CHKERRQ(ierr);
2664     ierr = VecGetArray(mat->slvec0,&from);CHKERRQ(ierr);
2665     ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2666     ierr = PetscMemcpy(from,x,bs*mbs*sizeof(MatScalar));CHKERRQ(ierr);
2667     ierr = VecRestoreArray(mat->slvec0,&from);CHKERRQ(ierr);
2668     ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2669     ierr = VecScatterBegin(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2670     ierr = VecScatterEnd(mat->sMvctx,mat->slvec0,mat->slvec1,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2671     ierr = (*mat->B->ops->multadd)(mat->B,mat->slvec1b,mat->slvec1a,mat->slvec1a);CHKERRQ(ierr);
2672 
2673     /* local sweep */
2674     ierr = (*mat->A->ops->sor)(mat->A,mat->slvec1a,omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_FORWARD_SWEEP),fshift,lits,1,xx1);CHKERRQ(ierr);
2675     ierr = VecAXPY(xx,1.0,xx1);CHKERRQ(ierr);
2676   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatSORType is not supported for SBAIJ matrix format");
2677   PetscFunctionReturn(0);
2678 }
2679 
2680 #undef __FUNCT__
2681 #define __FUNCT__ "MatSOR_MPISBAIJ_2comm"
2682 PetscErrorCode MatSOR_MPISBAIJ_2comm(Mat matin,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
2683 {
2684   Mat_MPISBAIJ   *mat = (Mat_MPISBAIJ*)matin->data;
2685   PetscErrorCode ierr;
2686   Vec            lvec1,bb1;
2687 
2688   PetscFunctionBegin;
2689   if (its <= 0 || lits <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
2690   if (matin->rmap->bs > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented");
2691 
2692   if ((flag & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP) {
2693     if (flag & SOR_ZERO_INITIAL_GUESS) {
2694       ierr = (*mat->A->ops->sor)(mat->A,bb,omega,flag,fshift,lits,lits,xx);CHKERRQ(ierr);
2695       its--;
2696     }
2697 
2698     ierr = VecDuplicate(mat->lvec,&lvec1);CHKERRQ(ierr);
2699     ierr = VecDuplicate(bb,&bb1);CHKERRQ(ierr);
2700     while (its--) {
2701       ierr = VecScatterBegin(mat->Mvctx,xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2702 
2703       /* lower diagonal part: bb1 = bb - B^T*xx */
2704       ierr = (*mat->B->ops->multtranspose)(mat->B,xx,lvec1);CHKERRQ(ierr);
2705       ierr = VecScale(lvec1,-1.0);CHKERRQ(ierr);
2706 
2707       ierr = VecScatterEnd(mat->Mvctx,xx,mat->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2708       ierr = VecCopy(bb,bb1);CHKERRQ(ierr);
2709       ierr = VecScatterBegin(mat->Mvctx,lvec1,bb1,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2710 
2711       /* upper diagonal part: bb1 = bb1 - B*x */
2712       ierr = VecScale(mat->lvec,-1.0);CHKERRQ(ierr);
2713       ierr = (*mat->B->ops->multadd)(mat->B,mat->lvec,bb1,bb1);CHKERRQ(ierr);
2714 
2715       ierr = VecScatterEnd(mat->Mvctx,lvec1,bb1,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
2716 
2717       /* diagonal sweep */
2718       ierr = (*mat->A->ops->sor)(mat->A,bb1,omega,SOR_SYMMETRIC_SWEEP,fshift,lits,lits,xx);CHKERRQ(ierr);
2719     }
2720     ierr = VecDestroy(&lvec1);CHKERRQ(ierr);
2721     ierr = VecDestroy(&bb1);CHKERRQ(ierr);
2722   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatSORType is not supported for SBAIJ matrix format");
2723   PetscFunctionReturn(0);
2724 }
2725 
2726 #undef __FUNCT__
2727 #define __FUNCT__ "MatCreateMPISBAIJWithArrays"
2728 /*@
2729      MatCreateMPISBAIJWithArrays - creates a MPI SBAIJ matrix using arrays that contain in standard
2730          CSR format the local rows.
2731 
2732    Collective on MPI_Comm
2733 
2734    Input Parameters:
2735 +  comm - MPI communicator
2736 .  bs - the block size, only a block size of 1 is supported
2737 .  m - number of local rows (Cannot be PETSC_DECIDE)
2738 .  n - This value should be the same as the local size used in creating the
2739        x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have
2740        calculated if N is given) For square matrices n is almost always m.
2741 .  M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
2742 .  N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
2743 .   i - row indices
2744 .   j - column indices
2745 -   a - matrix values
2746 
2747    Output Parameter:
2748 .   mat - the matrix
2749 
2750    Level: intermediate
2751 
2752    Notes:
2753        The i, j, and a arrays ARE copied by this routine into the internal format used by PETSc;
2754      thus you CANNOT change the matrix entries by changing the values of a[] after you have
2755      called this routine. Use MatCreateMPIAIJWithSplitArrays() to avoid needing to copy the arrays.
2756 
2757        The i and j indices are 0 based, and i indices are indices corresponding to the local j array.
2758 
2759 .keywords: matrix, aij, compressed row, sparse, parallel
2760 
2761 .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatMPIAIJSetPreallocation(), MatMPIAIJSetPreallocationCSR(),
2762           MPIAIJ, MatCreateAIJ(), MatCreateMPIAIJWithSplitArrays()
2763 @*/
2764 PetscErrorCode  MatCreateMPISBAIJWithArrays(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt i[],const PetscInt j[],const PetscScalar a[],Mat *mat)
2765 {
2766   PetscErrorCode ierr;
2767 
2768 
2769   PetscFunctionBegin;
2770   if (i[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"i (row indices) must start with 0");
2771   if (m < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"local number of rows (m) cannot be PETSC_DECIDE, or negative");
2772   ierr = MatCreate(comm,mat);CHKERRQ(ierr);
2773   ierr = MatSetSizes(*mat,m,n,M,N);CHKERRQ(ierr);
2774   ierr = MatSetType(*mat,MATMPISBAIJ);CHKERRQ(ierr);
2775   ierr = MatMPISBAIJSetPreallocationCSR(*mat,bs,i,j,a);CHKERRQ(ierr);
2776   PetscFunctionReturn(0);
2777 }
2778 
2779 
2780 #undef __FUNCT__
2781 #define __FUNCT__ "MatMPISBAIJSetPreallocationCSR"
2782 /*@C
2783    MatMPISBAIJSetPreallocationCSR - Allocates memory for a sparse parallel matrix in BAIJ format
2784    (the default parallel PETSc format).
2785 
2786    Collective on MPI_Comm
2787 
2788    Input Parameters:
2789 +  A - the matrix
2790 .  bs - the block size
2791 .  i - the indices into j for the start of each local row (starts with zero)
2792 .  j - the column indices for each local row (starts with zero) these must be sorted for each row
2793 -  v - optional values in the matrix
2794 
2795    Level: developer
2796 
2797 .keywords: matrix, aij, compressed row, sparse, parallel
2798 
2799 .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatMPIBAIJSetPreallocation(), MatCreateAIJ(), MPIAIJ
2800 @*/
2801 PetscErrorCode  MatMPISBAIJSetPreallocationCSR(Mat B,PetscInt bs,const PetscInt i[],const PetscInt j[], const PetscScalar v[])
2802 {
2803   PetscErrorCode ierr;
2804 
2805   PetscFunctionBegin;
2806   ierr = PetscTryMethod(B,"MatMPISBAIJSetPreallocationCSR_C",(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]),(B,bs,i,j,v));CHKERRQ(ierr);
2807   PetscFunctionReturn(0);
2808 }
2809 
2810 
2811