xref: /petsc/src/mat/impls/sbaij/seq/sbaij.c (revision 8cc058d9cd56c1ccb3be12a47760ddfc446aaffc)
1 
2 /*
3     Defines the basic matrix operations for the SBAIJ (compressed row)
4   matrix storage format.
5 */
6 #include <../src/mat/impls/baij/seq/baij.h>         /*I "petscmat.h" I*/
7 #include <../src/mat/impls/sbaij/seq/sbaij.h>
8 #include <petscblaslapack.h>
9 
10 #include <../src/mat/impls/sbaij/seq/relax.h>
11 #define USESHORT
12 #include <../src/mat/impls/sbaij/seq/relax.h>
13 
14 extern PetscErrorCode MatSeqSBAIJSetNumericFactorization_inplace(Mat,PetscBool);
15 
16 /*
17      Checks for missing diagonals
18 */
19 #undef __FUNCT__
20 #define __FUNCT__ "MatMissingDiagonal_SeqSBAIJ"
21 PetscErrorCode MatMissingDiagonal_SeqSBAIJ(Mat A,PetscBool  *missing,PetscInt *dd)
22 {
23   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
24   PetscErrorCode ierr;
25   PetscInt       *diag,*jj = a->j,i;
26 
27   PetscFunctionBegin;
28   ierr     = MatMarkDiagonal_SeqSBAIJ(A);CHKERRQ(ierr);
29   *missing = PETSC_FALSE;
30   if (A->rmap->n > 0 && !jj) {
31     *missing = PETSC_TRUE;
32     if (dd) *dd = 0;
33     PetscInfo(A,"Matrix has no entries therefore is missing diagonal");
34   } else {
35     diag = a->diag;
36     for (i=0; i<a->mbs; i++) {
37       if (jj[diag[i]] != i) {
38         *missing = PETSC_TRUE;
39         if (dd) *dd = i;
40         break;
41       }
42     }
43   }
44   PetscFunctionReturn(0);
45 }
46 
47 #undef __FUNCT__
48 #define __FUNCT__ "MatMarkDiagonal_SeqSBAIJ"
49 PetscErrorCode MatMarkDiagonal_SeqSBAIJ(Mat A)
50 {
51   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
52   PetscErrorCode ierr;
53   PetscInt       i;
54 
55   PetscFunctionBegin;
56   if (!a->diag) {
57     ierr         = PetscMalloc(a->mbs*sizeof(PetscInt),&a->diag);CHKERRQ(ierr);
58     ierr         = PetscLogObjectMemory(A,a->mbs*sizeof(PetscInt));CHKERRQ(ierr);
59     a->free_diag = PETSC_TRUE;
60   }
61   for (i=0; i<a->mbs; i++) a->diag[i] = a->i[i];
62   PetscFunctionReturn(0);
63 }
64 
65 #undef __FUNCT__
66 #define __FUNCT__ "MatGetRowIJ_SeqSBAIJ"
67 static PetscErrorCode MatGetRowIJ_SeqSBAIJ(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *nn,const PetscInt *inia[],const PetscInt *inja[],PetscBool  *done)
68 {
69   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
70   PetscInt       i,j,n = a->mbs,nz = a->i[n],bs = A->rmap->bs;
71   PetscInt       **ia = (PetscInt**)inia,**ja = (PetscInt**)inja;
72   PetscErrorCode ierr;
73 
74   PetscFunctionBegin;
75   *nn = n;
76   if (!ia) PetscFunctionReturn(0);
77   if (!blockcompressed) {
78     /* malloc & create the natural set of indices */
79     ierr = PetscMalloc2((n+1)*bs,PetscInt,ia,nz*bs,PetscInt,ja);CHKERRQ(ierr);
80     for (i=0; i<n+1; i++) {
81       for (j=0; j<bs; j++) {
82         *ia[i*bs+j] = a->i[i]*bs+j+oshift;
83       }
84     }
85     for (i=0; i<nz; i++) {
86       for (j=0; j<bs; j++) {
87         *ja[i*bs+j] = a->j[i]*bs+j+oshift;
88       }
89     }
90   } else { /* blockcompressed */
91     if (oshift == 1) {
92       /* temporarily add 1 to i and j indices */
93       for (i=0; i<nz; i++) a->j[i]++;
94       for (i=0; i<n+1; i++) a->i[i]++;
95     }
96     *ia = a->i; *ja = a->j;
97   }
98   PetscFunctionReturn(0);
99 }
100 
101 #undef __FUNCT__
102 #define __FUNCT__ "MatRestoreRowIJ_SeqSBAIJ"
103 static PetscErrorCode MatRestoreRowIJ_SeqSBAIJ(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *nn,const PetscInt *ia[],const PetscInt *ja[],PetscBool  *done)
104 {
105   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
106   PetscInt       i,n = a->mbs,nz = a->i[n];
107   PetscErrorCode ierr;
108 
109   PetscFunctionBegin;
110   if (!ia) PetscFunctionReturn(0);
111 
112   if (!blockcompressed) {
113     ierr = PetscFree2(*ia,*ja);CHKERRQ(ierr);
114   } else if (oshift == 1) { /* blockcompressed */
115     for (i=0; i<nz; i++) a->j[i]--;
116     for (i=0; i<n+1; i++) a->i[i]--;
117   }
118   PetscFunctionReturn(0);
119 }
120 
121 #undef __FUNCT__
122 #define __FUNCT__ "MatDestroy_SeqSBAIJ"
123 PetscErrorCode MatDestroy_SeqSBAIJ(Mat A)
124 {
125   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
126   PetscErrorCode ierr;
127 
128   PetscFunctionBegin;
129 #if defined(PETSC_USE_LOG)
130   PetscLogObjectState((PetscObject)A,"Rows=%D, NZ=%D",A->rmap->N,a->nz);
131 #endif
132   ierr = MatSeqXAIJFreeAIJ(A,&a->a,&a->j,&a->i);CHKERRQ(ierr);
133   if (a->free_diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);}
134   ierr = ISDestroy(&a->row);CHKERRQ(ierr);
135   ierr = ISDestroy(&a->col);CHKERRQ(ierr);
136   ierr = ISDestroy(&a->icol);CHKERRQ(ierr);
137   ierr = PetscFree(a->idiag);CHKERRQ(ierr);
138   ierr = PetscFree(a->inode.size);CHKERRQ(ierr);
139   if (a->free_imax_ilen) {ierr = PetscFree2(a->imax,a->ilen);CHKERRQ(ierr);}
140   ierr = PetscFree(a->solve_work);CHKERRQ(ierr);
141   ierr = PetscFree(a->sor_work);CHKERRQ(ierr);
142   ierr = PetscFree(a->solves_work);CHKERRQ(ierr);
143   ierr = PetscFree(a->mult_work);CHKERRQ(ierr);
144   ierr = PetscFree(a->saved_values);CHKERRQ(ierr);
145   ierr = PetscFree(a->xtoy);CHKERRQ(ierr);
146   if (a->free_jshort) {ierr = PetscFree(a->jshort);CHKERRQ(ierr);}
147   ierr = PetscFree(a->inew);CHKERRQ(ierr);
148   ierr = MatDestroy(&a->parent);CHKERRQ(ierr);
149   ierr = PetscFree(A->data);CHKERRQ(ierr);
150 
151   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
152   ierr = PetscObjectComposeFunction((PetscObject)A,"MatStoreValues_C","",NULL);CHKERRQ(ierr);
153   ierr = PetscObjectComposeFunction((PetscObject)A,"MatRetrieveValues_C","",NULL);CHKERRQ(ierr);
154   ierr = PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetColumnIndices_C","",NULL);CHKERRQ(ierr);
155   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqaij_C","",NULL);CHKERRQ(ierr);
156   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqbaij_C","",NULL);CHKERRQ(ierr);
157   ierr = PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetPreallocation_C","",NULL);CHKERRQ(ierr);
158   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqsbstrm_C","",NULL);CHKERRQ(ierr);
159   PetscFunctionReturn(0);
160 }
161 
162 #undef __FUNCT__
163 #define __FUNCT__ "MatSetOption_SeqSBAIJ"
164 PetscErrorCode MatSetOption_SeqSBAIJ(Mat A,MatOption op,PetscBool flg)
165 {
166   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
167   PetscErrorCode ierr;
168 
169   PetscFunctionBegin;
170   switch (op) {
171   case MAT_ROW_ORIENTED:
172     a->roworiented = flg;
173     break;
174   case MAT_KEEP_NONZERO_PATTERN:
175     a->keepnonzeropattern = flg;
176     break;
177   case MAT_NEW_NONZERO_LOCATIONS:
178     a->nonew = (flg ? 0 : 1);
179     break;
180   case MAT_NEW_NONZERO_LOCATION_ERR:
181     a->nonew = (flg ? -1 : 0);
182     break;
183   case MAT_NEW_NONZERO_ALLOCATION_ERR:
184     a->nonew = (flg ? -2 : 0);
185     break;
186   case MAT_UNUSED_NONZERO_LOCATION_ERR:
187     a->nounused = (flg ? -1 : 0);
188     break;
189   case MAT_NEW_DIAGONALS:
190   case MAT_IGNORE_OFF_PROC_ENTRIES:
191   case MAT_USE_HASH_TABLE:
192     ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr);
193     break;
194   case MAT_HERMITIAN:
195     if (!A->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call MatAssemblyEnd() first");
196     if (A->cmap->n < 65536 && A->cmap->bs == 1) {
197       A->ops->mult = MatMult_SeqSBAIJ_1_Hermitian_ushort;
198     } else if (A->cmap->bs == 1) {
199       A->ops->mult = MatMult_SeqSBAIJ_1_Hermitian;
200     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for Hermitian with block size greater than 1");
201     break;
202   case MAT_SPD:
203     /* These options are handled directly by MatSetOption() */
204     break;
205   case MAT_SYMMETRIC:
206   case MAT_STRUCTURALLY_SYMMETRIC:
207   case MAT_SYMMETRY_ETERNAL:
208     if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix must be symmetric");
209     ierr = PetscInfo1(A,"Option %s not relevent\n",MatOptions[op]);CHKERRQ(ierr);
210     break;
211   case MAT_IGNORE_LOWER_TRIANGULAR:
212     a->ignore_ltriangular = flg;
213     break;
214   case MAT_ERROR_LOWER_TRIANGULAR:
215     a->ignore_ltriangular = flg;
216     break;
217   case MAT_GETROW_UPPERTRIANGULAR:
218     a->getrow_utriangular = flg;
219     break;
220   default:
221     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %d",op);
222   }
223   PetscFunctionReturn(0);
224 }
225 
226 #undef __FUNCT__
227 #define __FUNCT__ "MatGetRow_SeqSBAIJ"
228 PetscErrorCode MatGetRow_SeqSBAIJ(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **v)
229 {
230   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
231   PetscErrorCode ierr;
232   PetscInt       itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*cols_i,bs2;
233   MatScalar      *aa,*aa_i;
234   PetscScalar    *v_i;
235 
236   PetscFunctionBegin;
237   if (A && !a->getrow_utriangular) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatGetRow is not supported for SBAIJ matrix format. Getting the upper triangular part of row, run with -mat_getrow_uppertriangular, call MatSetOption(mat,MAT_GETROW_UPPERTRIANGULAR,PETSC_TRUE) or MatGetRowUpperTriangular()");
238   /* Get the upper triangular part of the row */
239   bs  = A->rmap->bs;
240   ai  = a->i;
241   aj  = a->j;
242   aa  = a->a;
243   bs2 = a->bs2;
244 
245   if (row < 0 || row >= A->rmap->N) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE, "Row %D out of range", row);
246 
247   bn     = row/bs; /* Block number */
248   bp     = row % bs; /* Block position */
249   M      = ai[bn+1] - ai[bn];
250   *ncols = bs*M;
251 
252   if (v) {
253     *v = 0;
254     if (*ncols) {
255       ierr = PetscMalloc((*ncols+row)*sizeof(PetscScalar),v);CHKERRQ(ierr);
256       for (i=0; i<M; i++) { /* for each block in the block row */
257         v_i  = *v + i*bs;
258         aa_i = aa + bs2*(ai[bn] + i);
259         for (j=bp,k=0; j<bs2; j+=bs,k++) v_i[k] = aa_i[j];
260       }
261     }
262   }
263 
264   if (cols) {
265     *cols = 0;
266     if (*ncols) {
267       ierr = PetscMalloc((*ncols+row)*sizeof(PetscInt),cols);CHKERRQ(ierr);
268       for (i=0; i<M; i++) { /* for each block in the block row */
269         cols_i = *cols + i*bs;
270         itmp   = bs*aj[ai[bn] + i];
271         for (j=0; j<bs; j++) cols_i[j] = itmp++;
272       }
273     }
274   }
275 
276   /*search column A(0:row-1,row) (=A(row,0:row-1)). Could be expensive! */
277   /* this segment is currently removed, so only entries in the upper triangle are obtained */
278 #if defined(column_search)
279   v_i    = *v    + M*bs;
280   cols_i = *cols + M*bs;
281   for (i=0; i<bn; i++) { /* for each block row */
282     M = ai[i+1] - ai[i];
283     for (j=0; j<M; j++) {
284       itmp = aj[ai[i] + j];    /* block column value */
285       if (itmp == bn) {
286         aa_i = aa + bs2*(ai[i] + j) + bs*bp;
287         for (k=0; k<bs; k++) {
288           *cols_i++ = i*bs+k;
289           *v_i++    = aa_i[k];
290         }
291         *ncols += bs;
292         break;
293       }
294     }
295   }
296 #endif
297   PetscFunctionReturn(0);
298 }
299 
300 #undef __FUNCT__
301 #define __FUNCT__ "MatRestoreRow_SeqSBAIJ"
302 PetscErrorCode MatRestoreRow_SeqSBAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
303 {
304   PetscErrorCode ierr;
305 
306   PetscFunctionBegin;
307   if (idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}
308   if (v)   {ierr = PetscFree(*v);CHKERRQ(ierr);}
309   PetscFunctionReturn(0);
310 }
311 
312 #undef __FUNCT__
313 #define __FUNCT__ "MatGetRowUpperTriangular_SeqSBAIJ"
314 PetscErrorCode MatGetRowUpperTriangular_SeqSBAIJ(Mat A)
315 {
316   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
317 
318   PetscFunctionBegin;
319   a->getrow_utriangular = PETSC_TRUE;
320   PetscFunctionReturn(0);
321 }
322 #undef __FUNCT__
323 #define __FUNCT__ "MatRestoreRowUpperTriangular_SeqSBAIJ"
324 PetscErrorCode MatRestoreRowUpperTriangular_SeqSBAIJ(Mat A)
325 {
326   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
327 
328   PetscFunctionBegin;
329   a->getrow_utriangular = PETSC_FALSE;
330   PetscFunctionReturn(0);
331 }
332 
333 #undef __FUNCT__
334 #define __FUNCT__ "MatTranspose_SeqSBAIJ"
335 PetscErrorCode MatTranspose_SeqSBAIJ(Mat A,MatReuse reuse,Mat *B)
336 {
337   PetscErrorCode ierr;
338 
339   PetscFunctionBegin;
340   if (reuse == MAT_INITIAL_MATRIX || *B != A) {
341     ierr = MatDuplicate(A,MAT_COPY_VALUES,B);CHKERRQ(ierr);
342   }
343   PetscFunctionReturn(0);
344 }
345 
346 #undef __FUNCT__
347 #define __FUNCT__ "MatView_SeqSBAIJ_ASCII"
348 static PetscErrorCode MatView_SeqSBAIJ_ASCII(Mat A,PetscViewer viewer)
349 {
350   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
351   PetscErrorCode    ierr;
352   PetscInt          i,j,bs = A->rmap->bs,k,l,bs2=a->bs2;
353   PetscViewerFormat format;
354   PetscInt          *diag;
355 
356   PetscFunctionBegin;
357   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
358   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
359     ierr = PetscViewerASCIIPrintf(viewer,"  block size is %D\n",bs);CHKERRQ(ierr);
360   } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
361     Mat aij;
362     if (A->factortype && bs>1) {
363       ierr = PetscPrintf(PETSC_COMM_SELF,"Warning: matrix is factored with bs>1. MatView() with PETSC_VIEWER_ASCII_MATLAB is not supported and ignored!\n");CHKERRQ(ierr);
364       PetscFunctionReturn(0);
365     }
366     ierr = MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&aij);CHKERRQ(ierr);
367     ierr = MatView(aij,viewer);CHKERRQ(ierr);
368     ierr = MatDestroy(&aij);CHKERRQ(ierr);
369   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
370     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
371     for (i=0; i<a->mbs; i++) {
372       for (j=0; j<bs; j++) {
373         ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);CHKERRQ(ierr);
374         for (k=a->i[i]; k<a->i[i+1]; k++) {
375           for (l=0; l<bs; l++) {
376 #if defined(PETSC_USE_COMPLEX)
377             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
378               ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",bs*a->j[k]+l,
379                                             PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
380             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
381               ierr = PetscViewerASCIIPrintf(viewer," (%D, %G - %G i) ",bs*a->j[k]+l,
382                                             PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
383             } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
384               ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
385             }
386 #else
387             if (a->a[bs2*k + l*bs + j] != 0.0) {
388               ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
389             }
390 #endif
391           }
392         }
393         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
394       }
395     }
396     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
397   } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
398     PetscFunctionReturn(0);
399   } else {
400     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
401     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)A,viewer,"Matrix Object");CHKERRQ(ierr);
402     if (A->factortype) { /* for factored matrix */
403       if (bs>1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"matrix is factored with bs>1. Not implemented yet");
404 
405       diag=a->diag;
406       for (i=0; i<a->mbs; i++) { /* for row block i */
407         ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i);CHKERRQ(ierr);
408         /* diagonal entry */
409 #if defined(PETSC_USE_COMPLEX)
410         if (PetscImaginaryPart(a->a[diag[i]]) > 0.0) {
411           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",a->j[diag[i]],PetscRealPart(1.0/a->a[diag[i]]),PetscImaginaryPart(1.0/a->a[diag[i]]));CHKERRQ(ierr);
412         } else if (PetscImaginaryPart(a->a[diag[i]]) < 0.0) {
413           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G - %G i) ",a->j[diag[i]],PetscRealPart(1.0/a->a[diag[i]]),-PetscImaginaryPart(1.0/a->a[diag[i]]));CHKERRQ(ierr);
414         } else {
415           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",a->j[diag[i]],PetscRealPart(1.0/a->a[diag[i]]));CHKERRQ(ierr);
416         }
417 #else
418         ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",a->j[diag[i]],1.0/a->a[diag[i]]);CHKERRQ(ierr);
419 #endif
420         /* off-diagonal entries */
421         for (k=a->i[i]; k<a->i[i+1]-1; k++) {
422 #if defined(PETSC_USE_COMPLEX)
423           if (PetscImaginaryPart(a->a[k]) > 0.0) {
424             ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",bs*a->j[k],PetscRealPart(a->a[k]),PetscImaginaryPart(a->a[k]));CHKERRQ(ierr);
425           } else if (PetscImaginaryPart(a->a[k]) < 0.0) {
426             ierr = PetscViewerASCIIPrintf(viewer," (%D, %G - %G i) ",bs*a->j[k],PetscRealPart(a->a[k]),-PetscImaginaryPart(a->a[k]));CHKERRQ(ierr);
427           } else {
428             ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k],PetscRealPart(a->a[k]));CHKERRQ(ierr);
429           }
430 #else
431           ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",a->j[k],a->a[k]);CHKERRQ(ierr);
432 #endif
433         }
434         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
435       }
436 
437     } else { /* for non-factored matrix */
438       for (i=0; i<a->mbs; i++) { /* for row block i */
439         for (j=0; j<bs; j++) {   /* for row bs*i + j */
440           ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);CHKERRQ(ierr);
441           for (k=a->i[i]; k<a->i[i+1]; k++) { /* for column block */
442             for (l=0; l<bs; l++) {            /* for column */
443 #if defined(PETSC_USE_COMPLEX)
444               if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) {
445                 ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",bs*a->j[k]+l,
446                                               PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
447               } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) {
448                 ierr = PetscViewerASCIIPrintf(viewer," (%D, %G - %G i) ",bs*a->j[k]+l,
449                                               PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
450               } else {
451                 ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
452               }
453 #else
454               ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
455 #endif
456             }
457           }
458           ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
459         }
460       }
461     }
462     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
463   }
464   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
465   PetscFunctionReturn(0);
466 }
467 
468 #include <petscdraw.h>
469 #undef __FUNCT__
470 #define __FUNCT__ "MatView_SeqSBAIJ_Draw_Zoom"
471 static PetscErrorCode MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw,void *Aa)
472 {
473   Mat            A = (Mat) Aa;
474   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)A->data;
475   PetscErrorCode ierr;
476   PetscInt       row,i,j,k,l,mbs=a->mbs,color,bs=A->rmap->bs,bs2=a->bs2;
477   PetscMPIInt    rank;
478   PetscReal      xl,yl,xr,yr,x_l,x_r,y_l,y_r;
479   MatScalar      *aa;
480   MPI_Comm       comm;
481   PetscViewer    viewer;
482 
483   PetscFunctionBegin;
484   /*
485     This is nasty. If this is called from an originally parallel matrix
486     then all processes call this,but only the first has the matrix so the
487     rest should return immediately.
488   */
489   ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr);
490   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
491   if (rank) PetscFunctionReturn(0);
492 
493   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
494 
495   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
496   PetscDrawString(draw, .3*(xl+xr), .3*(yl+yr), PETSC_DRAW_BLACK, "symmetric");
497 
498   /* loop over matrix elements drawing boxes */
499   color = PETSC_DRAW_BLUE;
500   for (i=0,row=0; i<mbs; i++,row+=bs) {
501     for (j=a->i[i]; j<a->i[i+1]; j++) {
502       y_l = A->rmap->N - row - 1.0; y_r = y_l + 1.0;
503       x_l = a->j[j]*bs; x_r = x_l + 1.0;
504       aa  = a->a + j*bs2;
505       for (k=0; k<bs; k++) {
506         for (l=0; l<bs; l++) {
507           if (PetscRealPart(*aa++) >=  0.) continue;
508           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
509         }
510       }
511     }
512   }
513   color = PETSC_DRAW_CYAN;
514   for (i=0,row=0; i<mbs; i++,row+=bs) {
515     for (j=a->i[i]; j<a->i[i+1]; j++) {
516       y_l = A->rmap->N - row - 1.0; y_r = y_l + 1.0;
517       x_l = a->j[j]*bs; x_r = x_l + 1.0;
518       aa = a->a + j*bs2;
519       for (k=0; k<bs; k++) {
520         for (l=0; l<bs; l++) {
521           if (PetscRealPart(*aa++) != 0.) continue;
522           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
523         }
524       }
525     }
526   }
527 
528   color = PETSC_DRAW_RED;
529   for (i=0,row=0; i<mbs; i++,row+=bs) {
530     for (j=a->i[i]; j<a->i[i+1]; j++) {
531       y_l = A->rmap->N - row - 1.0; y_r = y_l + 1.0;
532       x_l = a->j[j]*bs; x_r = x_l + 1.0;
533       aa = a->a + j*bs2;
534       for (k=0; k<bs; k++) {
535         for (l=0; l<bs; l++) {
536           if (PetscRealPart(*aa++) <= 0.) continue;
537           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
538         }
539       }
540     }
541   }
542   PetscFunctionReturn(0);
543 }
544 
545 #undef __FUNCT__
546 #define __FUNCT__ "MatView_SeqSBAIJ_Draw"
547 static PetscErrorCode MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer)
548 {
549   PetscErrorCode ierr;
550   PetscReal      xl,yl,xr,yr,w,h;
551   PetscDraw      draw;
552   PetscBool      isnull;
553 
554   PetscFunctionBegin;
555   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
556   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
557 
558   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
559   xr   = A->rmap->N; yr = A->rmap->N; h = yr/10.0; w = xr/10.0;
560   xr  += w;    yr += h;  xl = -w;     yl = -h;
561   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
562   ierr = PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);CHKERRQ(ierr);
563   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",NULL);CHKERRQ(ierr);
564   PetscFunctionReturn(0);
565 }
566 
567 #undef __FUNCT__
568 #define __FUNCT__ "MatView_SeqSBAIJ"
569 PetscErrorCode MatView_SeqSBAIJ(Mat A,PetscViewer viewer)
570 {
571   PetscErrorCode ierr;
572   PetscBool      iascii,isdraw;
573   FILE           *file = 0;
574 
575   PetscFunctionBegin;
576   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
577   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
578   if (iascii) {
579     ierr = MatView_SeqSBAIJ_ASCII(A,viewer);CHKERRQ(ierr);
580   } else if (isdraw) {
581     ierr = MatView_SeqSBAIJ_Draw(A,viewer);CHKERRQ(ierr);
582   } else {
583     Mat B;
584     ierr = MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr);
585     ierr = MatView(B,viewer);CHKERRQ(ierr);
586     ierr = MatDestroy(&B);CHKERRQ(ierr);
587     ierr = PetscViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr);
588     if (file) {
589       fprintf(file,"-matload_block_size %d\n",(int)A->rmap->bs);
590     }
591   }
592   PetscFunctionReturn(0);
593 }
594 
595 
596 #undef __FUNCT__
597 #define __FUNCT__ "MatGetValues_SeqSBAIJ"
598 PetscErrorCode MatGetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],PetscScalar v[])
599 {
600   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
601   PetscInt     *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
602   PetscInt     *ai = a->i,*ailen = a->ilen;
603   PetscInt     brow,bcol,ridx,cidx,bs=A->rmap->bs,bs2=a->bs2;
604   MatScalar    *ap,*aa = a->a;
605 
606   PetscFunctionBegin;
607   for (k=0; k<m; k++) { /* loop over rows */
608     row = im[k]; brow = row/bs;
609     if (row < 0) {v += n; continue;} /* SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",row); */
610     if (row >= A->rmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap->N-1);
611     rp   = aj + ai[brow]; ap = aa + bs2*ai[brow];
612     nrow = ailen[brow];
613     for (l=0; l<n; l++) { /* loop over columns */
614       if (in[l] < 0) {v++; continue;} /* SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %D",in[l]); */
615       if (in[l] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->cmap->n-1);
616       col  = in[l];
617       bcol = col/bs;
618       cidx = col%bs;
619       ridx = row%bs;
620       high = nrow;
621       low  = 0; /* assume unsorted */
622       while (high-low > 5) {
623         t = (low+high)/2;
624         if (rp[t] > bcol) high = t;
625         else              low  = t;
626       }
627       for (i=low; i<high; i++) {
628         if (rp[i] > bcol) break;
629         if (rp[i] == bcol) {
630           *v++ = ap[bs2*i+bs*cidx+ridx];
631           goto finished;
632         }
633       }
634       *v++ = 0.0;
635 finished:;
636     }
637   }
638   PetscFunctionReturn(0);
639 }
640 
641 
642 #undef __FUNCT__
643 #define __FUNCT__ "MatSetValuesBlocked_SeqSBAIJ"
644 PetscErrorCode MatSetValuesBlocked_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
645 {
646   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
647   PetscErrorCode    ierr;
648   PetscInt          *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,lastcol = -1;
649   PetscInt          *imax      =a->imax,*ai=a->i,*ailen=a->ilen;
650   PetscInt          *aj        =a->j,nonew=a->nonew,bs2=a->bs2,bs=A->rmap->bs,stepval;
651   PetscBool         roworiented=a->roworiented;
652   const PetscScalar *value     = v;
653   MatScalar         *ap,*aa = a->a,*bap;
654 
655   PetscFunctionBegin;
656   if (roworiented) stepval = (n-1)*bs;
657   else stepval = (m-1)*bs;
658 
659   for (k=0; k<m; k++) { /* loop over added rows */
660     row = im[k];
661     if (row < 0) continue;
662 #if defined(PETSC_USE_DEBUG)
663     if (row >= a->mbs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,a->mbs-1);
664 #endif
665     rp   = aj + ai[row];
666     ap   = aa + bs2*ai[row];
667     rmax = imax[row];
668     nrow = ailen[row];
669     low  = 0;
670     high = nrow;
671     for (l=0; l<n; l++) { /* loop over added columns */
672       if (in[l] < 0) continue;
673       col = in[l];
674 #if defined(PETSC_USE_DEBUG)
675       if (col >= a->nbs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",col,a->nbs-1);
676 #endif
677       if (col < row) {
678         if (a->ignore_ltriangular) continue; /* ignore lower triangular block */
679         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)");
680       }
681       if (roworiented) value = v + k*(stepval+bs)*bs + l*bs;
682       else value = v + l*(stepval+bs)*bs + k*bs;
683 
684       if (col <= lastcol) low = 0;
685       else high = nrow;
686 
687       lastcol = col;
688       while (high-low > 7) {
689         t = (low+high)/2;
690         if (rp[t] > col) high = t;
691         else             low  = t;
692       }
693       for (i=low; i<high; i++) {
694         if (rp[i] > col) break;
695         if (rp[i] == col) {
696           bap = ap +  bs2*i;
697           if (roworiented) {
698             if (is == ADD_VALUES) {
699               for (ii=0; ii<bs; ii++,value+=stepval) {
700                 for (jj=ii; jj<bs2; jj+=bs) {
701                   bap[jj] += *value++;
702                 }
703               }
704             } else {
705               for (ii=0; ii<bs; ii++,value+=stepval) {
706                 for (jj=ii; jj<bs2; jj+=bs) {
707                   bap[jj] = *value++;
708                 }
709                }
710             }
711           } else {
712             if (is == ADD_VALUES) {
713               for (ii=0; ii<bs; ii++,value+=stepval) {
714                 for (jj=0; jj<bs; jj++) {
715                   *bap++ += *value++;
716                 }
717               }
718             } else {
719               for (ii=0; ii<bs; ii++,value+=stepval) {
720                 for (jj=0; jj<bs; jj++) {
721                   *bap++  = *value++;
722                 }
723               }
724             }
725           }
726           goto noinsert2;
727         }
728       }
729       if (nonew == 1) goto noinsert2;
730       if (nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col);
731       MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar);
732       N = nrow++ - 1; high++;
733       /* shift up all the later entries in this row */
734       for (ii=N; ii>=i; ii--) {
735         rp[ii+1] = rp[ii];
736         ierr     = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
737       }
738       if (N >= i) {
739         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
740       }
741       rp[i] = col;
742       bap   = ap +  bs2*i;
743       if (roworiented) {
744         for (ii=0; ii<bs; ii++,value+=stepval) {
745           for (jj=ii; jj<bs2; jj+=bs) {
746             bap[jj] = *value++;
747           }
748         }
749       } else {
750         for (ii=0; ii<bs; ii++,value+=stepval) {
751           for (jj=0; jj<bs; jj++) {
752             *bap++ = *value++;
753           }
754         }
755        }
756     noinsert2:;
757       low = i;
758     }
759     ailen[row] = nrow;
760   }
761   PetscFunctionReturn(0);
762 }
763 
764 /*
765     This is not yet used
766 */
767 #undef __FUNCT__
768 #define __FUNCT__ "MatAssemblyEnd_SeqSBAIJ_SeqAIJ_Inode"
769 PetscErrorCode MatAssemblyEnd_SeqSBAIJ_SeqAIJ_Inode(Mat A)
770 {
771   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
772   PetscErrorCode ierr;
773   const PetscInt *ai = a->i, *aj = a->j,*cols;
774   PetscInt       i   = 0,j,blk_size,m = A->rmap->n,node_count = 0,nzx,nzy,*ns,row,nz,cnt,cnt2,*counts;
775   PetscBool      flag;
776 
777   PetscFunctionBegin;
778   ierr = PetscMalloc(m*sizeof(PetscInt),&ns);CHKERRQ(ierr);
779   while (i < m) {
780     nzx = ai[i+1] - ai[i];       /* Number of nonzeros */
781     /* Limits the number of elements in a node to 'a->inode.limit' */
782     for (j=i+1,blk_size=1; j<m && blk_size <a->inode.limit; ++j,++blk_size) {
783       nzy = ai[j+1] - ai[j];
784       if (nzy != (nzx - j + i)) break;
785       ierr = PetscMemcmp(aj + ai[i] + j - i,aj + ai[j],nzy*sizeof(PetscInt),&flag);CHKERRQ(ierr);
786       if (!flag) break;
787     }
788     ns[node_count++] = blk_size;
789 
790     i = j;
791   }
792   if (!a->inode.size && m && node_count > .9*m) {
793     ierr = PetscFree(ns);CHKERRQ(ierr);
794     ierr = PetscInfo2(A,"Found %D nodes out of %D rows. Not using Inode routines\n",node_count,m);CHKERRQ(ierr);
795   } else {
796     a->inode.node_count = node_count;
797 
798     ierr = PetscMalloc(node_count*sizeof(PetscInt),&a->inode.size);CHKERRQ(ierr);
799     ierr = PetscLogObjectMemory(A,node_count*sizeof(PetscInt));CHKERRQ(ierr);
800     ierr = PetscMemcpy(a->inode.size,ns,node_count*sizeof(PetscInt));CHKERRQ(ierr);
801     ierr = PetscFree(ns);CHKERRQ(ierr);
802     ierr = PetscInfo3(A,"Found %D nodes of %D. Limit used: %D. Using Inode routines\n",node_count,m,a->inode.limit);CHKERRQ(ierr);
803 
804     /* count collections of adjacent columns in each inode */
805     row = 0;
806     cnt = 0;
807     for (i=0; i<node_count; i++) {
808       cols = aj + ai[row] + a->inode.size[i];
809       nz   = ai[row+1] - ai[row] - a->inode.size[i];
810       for (j=1; j<nz; j++) {
811         if (cols[j] != cols[j-1]+1) cnt++;
812       }
813       cnt++;
814       row += a->inode.size[i];
815     }
816     ierr = PetscMalloc(2*cnt*sizeof(PetscInt),&counts);CHKERRQ(ierr);
817     cnt  = 0;
818     row  = 0;
819     for (i=0; i<node_count; i++) {
820       cols = aj + ai[row] + a->inode.size[i];
821       CHKMEMQ;
822       counts[2*cnt] = cols[0];
823       CHKMEMQ;
824       nz   = ai[row+1] - ai[row] - a->inode.size[i];
825       cnt2 = 1;
826       for (j=1; j<nz; j++) {
827         if (cols[j] != cols[j-1]+1) {
828           CHKMEMQ;
829           counts[2*(cnt++)+1] = cnt2;
830           counts[2*cnt]       = cols[j];
831           CHKMEMQ;
832           cnt2 = 1;
833         } else cnt2++;
834       }
835       CHKMEMQ;
836       counts[2*(cnt++)+1] = cnt2;
837       CHKMEMQ;
838       row += a->inode.size[i];
839     }
840     ierr = PetscIntView(2*cnt,counts,0);CHKERRQ(ierr);
841   }
842   PetscFunctionReturn(0);
843 }
844 
845 #undef __FUNCT__
846 #define __FUNCT__ "MatAssemblyEnd_SeqSBAIJ"
847 PetscErrorCode MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode)
848 {
849   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
850   PetscErrorCode ierr;
851   PetscInt       fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
852   PetscInt       m      = A->rmap->N,*ip,N,*ailen = a->ilen;
853   PetscInt       mbs    = a->mbs,bs2 = a->bs2,rmax = 0;
854   MatScalar      *aa    = a->a,*ap;
855 
856   PetscFunctionBegin;
857   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
858 
859   if (m) rmax = ailen[0];
860   for (i=1; i<mbs; i++) {
861     /* move each row back by the amount of empty slots (fshift) before it*/
862     fshift += imax[i-1] - ailen[i-1];
863     rmax    = PetscMax(rmax,ailen[i]);
864     if (fshift) {
865       ip = aj + ai[i]; ap = aa + bs2*ai[i];
866       N  = ailen[i];
867       for (j=0; j<N; j++) {
868         ip[j-fshift] = ip[j];
869         ierr         = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
870       }
871     }
872     ai[i] = ai[i-1] + ailen[i-1];
873   }
874   if (mbs) {
875     fshift += imax[mbs-1] - ailen[mbs-1];
876     ai[mbs] = ai[mbs-1] + ailen[mbs-1];
877   }
878   /* reset ilen and imax for each row */
879   for (i=0; i<mbs; i++) {
880     ailen[i] = imax[i] = ai[i+1] - ai[i];
881   }
882   a->nz = ai[mbs];
883 
884   /* diagonals may have moved, reset it */
885   if (a->diag) {
886     ierr = PetscMemcpy(a->diag,ai,mbs*sizeof(PetscInt));CHKERRQ(ierr);
887   }
888   if (fshift && a->nounused == -1) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_PLIB, "Unused space detected in matrix: %D X %D block size %D, %D unneeded", m, A->cmap->n, A->rmap->bs, fshift*bs2);
889 
890   ierr = PetscInfo5(A,"Matrix size: %D X %D, block size %D; storage space: %D unneeded, %D used\n",m,A->rmap->N,A->rmap->bs,fshift*bs2,a->nz*bs2);CHKERRQ(ierr);
891   ierr = PetscInfo1(A,"Number of mallocs during MatSetValues is %D\n",a->reallocs);CHKERRQ(ierr);
892   ierr = PetscInfo1(A,"Most nonzeros blocks in any row is %D\n",rmax);CHKERRQ(ierr);
893 
894   A->info.mallocs    += a->reallocs;
895   a->reallocs         = 0;
896   A->info.nz_unneeded = (PetscReal)fshift*bs2;
897   a->idiagvalid       = PETSC_FALSE;
898 
899   if (A->cmap->n < 65536 && A->cmap->bs == 1) {
900     if (a->jshort && a->free_jshort) {
901       /* when matrix data structure is changed, previous jshort must be replaced */
902       ierr = PetscFree(a->jshort);CHKERRQ(ierr);
903     }
904     ierr = PetscMalloc(a->i[A->rmap->n]*sizeof(unsigned short),&a->jshort);CHKERRQ(ierr);
905     ierr = PetscLogObjectMemory(A,a->i[A->rmap->n]*sizeof(unsigned short));CHKERRQ(ierr);
906     for (i=0; i<a->i[A->rmap->n]; i++) a->jshort[i] = a->j[i];
907     A->ops->mult   = MatMult_SeqSBAIJ_1_ushort;
908     A->ops->sor    = MatSOR_SeqSBAIJ_ushort;
909     a->free_jshort = PETSC_TRUE;
910   }
911   PetscFunctionReturn(0);
912 }
913 
914 /*
915    This function returns an array of flags which indicate the locations of contiguous
916    blocks that should be zeroed. for eg: if bs = 3  and is = [0,1,2,3,5,6,7,8,9]
917    then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)]
918    Assume: sizes should be long enough to hold all the values.
919 */
920 #undef __FUNCT__
921 #define __FUNCT__ "MatZeroRows_SeqSBAIJ_Check_Blocks"
922 PetscErrorCode MatZeroRows_SeqSBAIJ_Check_Blocks(PetscInt idx[],PetscInt n,PetscInt bs,PetscInt sizes[], PetscInt *bs_max)
923 {
924   PetscInt  i,j,k,row;
925   PetscBool flg;
926 
927   PetscFunctionBegin;
928   for (i=0,j=0; i<n; j++) {
929     row = idx[i];
930     if (row%bs!=0) { /* Not the begining of a block */
931       sizes[j] = 1;
932       i++;
933     } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */
934       sizes[j] = 1;         /* Also makes sure atleast 'bs' values exist for next else */
935       i++;
936     } else { /* Begining of the block, so check if the complete block exists */
937       flg = PETSC_TRUE;
938       for (k=1; k<bs; k++) {
939         if (row+k != idx[i+k]) { /* break in the block */
940           flg = PETSC_FALSE;
941           break;
942         }
943       }
944       if (flg) { /* No break in the bs */
945         sizes[j] = bs;
946         i       += bs;
947       } else {
948         sizes[j] = 1;
949         i++;
950       }
951     }
952   }
953   *bs_max = j;
954   PetscFunctionReturn(0);
955 }
956 
957 
958 /* Only add/insert a(i,j) with i<=j (blocks).
959    Any a(i,j) with i>j input by user is ingored.
960 */
961 
962 #undef __FUNCT__
963 #define __FUNCT__ "MatSetValues_SeqSBAIJ"
964 PetscErrorCode MatSetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
965 {
966   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
967   PetscErrorCode ierr;
968   PetscInt       *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,lastcol = -1;
969   PetscInt       *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented;
970   PetscInt       *aj  =a->j,nonew=a->nonew,bs=A->rmap->bs,brow,bcol;
971   PetscInt       ridx,cidx,bs2=a->bs2;
972   MatScalar      *ap,value,*aa=a->a,*bap;
973 
974   PetscFunctionBegin;
975   if (v) PetscValidScalarPointer(v,6);
976   for (k=0; k<m; k++) { /* loop over added rows */
977     row  = im[k];       /* row number */
978     brow = row/bs;      /* block row number */
979     if (row < 0) continue;
980 #if defined(PETSC_USE_DEBUG)
981     if (row >= A->rmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap->N-1);
982 #endif
983     rp   = aj + ai[brow]; /*ptr to beginning of column value of the row block*/
984     ap   = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/
985     rmax = imax[brow];  /* maximum space allocated for this row */
986     nrow = ailen[brow]; /* actual length of this row */
987     low  = 0;
988 
989     for (l=0; l<n; l++) { /* loop over added columns */
990       if (in[l] < 0) continue;
991 #if defined(PETSC_USE_DEBUG)
992       if (in[l] >= A->rmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->rmap->N-1);
993 #endif
994       col  = in[l];
995       bcol = col/bs;              /* block col number */
996 
997       if (brow > bcol) {
998         if (a->ignore_ltriangular) continue; /* ignore lower triangular values */
999         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)");
1000       }
1001 
1002       ridx = row % bs; cidx = col % bs; /*row and col index inside the block */
1003       if ((brow==bcol && ridx<=cidx) || (brow<bcol)) {
1004         /* element value a(k,l) */
1005         if (roworiented) value = v[l + k*n];
1006         else value = v[k + l*m];
1007 
1008         /* move pointer bap to a(k,l) quickly and add/insert value */
1009         if (col <= lastcol) low = 0;
1010         high = nrow;
1011         lastcol = col;
1012         while (high-low > 7) {
1013           t = (low+high)/2;
1014           if (rp[t] > bcol) high = t;
1015           else              low  = t;
1016         }
1017         for (i=low; i<high; i++) {
1018           if (rp[i] > bcol) break;
1019           if (rp[i] == bcol) {
1020             bap = ap +  bs2*i + bs*cidx + ridx;
1021             if (is == ADD_VALUES) *bap += value;
1022             else                  *bap  = value;
1023             /* for diag block, add/insert its symmetric element a(cidx,ridx) */
1024             if (brow == bcol && ridx < cidx) {
1025               bap = ap +  bs2*i + bs*ridx + cidx;
1026               if (is == ADD_VALUES) *bap += value;
1027               else                  *bap  = value;
1028             }
1029             goto noinsert1;
1030           }
1031         }
1032 
1033         if (nonew == 1) goto noinsert1;
1034         if (nonew == -1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col);
1035         MatSeqXAIJReallocateAIJ(A,a->mbs,bs2,nrow,brow,bcol,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar);
1036 
1037         N = nrow++ - 1; high++;
1038         /* shift up all the later entries in this row */
1039         for (ii=N; ii>=i; ii--) {
1040           rp[ii+1] = rp[ii];
1041           ierr     = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
1042         }
1043         if (N>=i) {
1044           ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1045         }
1046         rp[i]                      = bcol;
1047         ap[bs2*i + bs*cidx + ridx] = value;
1048 noinsert1:;
1049         low = i;
1050       }
1051     }   /* end of loop over added columns */
1052     ailen[brow] = nrow;
1053   }   /* end of loop over added rows */
1054   PetscFunctionReturn(0);
1055 }
1056 
1057 #undef __FUNCT__
1058 #define __FUNCT__ "MatICCFactor_SeqSBAIJ"
1059 PetscErrorCode MatICCFactor_SeqSBAIJ(Mat inA,IS row,const MatFactorInfo *info)
1060 {
1061   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)inA->data;
1062   Mat            outA;
1063   PetscErrorCode ierr;
1064   PetscBool      row_identity;
1065 
1066   PetscFunctionBegin;
1067   if (info->levels != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only levels=0 is supported for in-place icc");
1068   ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr);
1069   if (!row_identity) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix reordering is not supported");
1070   if (inA->rmap->bs != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix block size %D is not supported",inA->rmap->bs); /* Need to replace MatCholeskyFactorSymbolic_SeqSBAIJ_MSR()! */
1071 
1072   outA            = inA;
1073   inA->factortype = MAT_FACTOR_ICC;
1074 
1075   ierr = MatMarkDiagonal_SeqSBAIJ(inA);CHKERRQ(ierr);
1076   ierr = MatSeqSBAIJSetNumericFactorization_inplace(inA,row_identity);CHKERRQ(ierr);
1077 
1078   ierr   = PetscObjectReference((PetscObject)row);CHKERRQ(ierr);
1079   ierr   = ISDestroy(&a->row);CHKERRQ(ierr);
1080   a->row = row;
1081   ierr   = PetscObjectReference((PetscObject)row);CHKERRQ(ierr);
1082   ierr   = ISDestroy(&a->col);CHKERRQ(ierr);
1083   a->col = row;
1084 
1085   /* Create the invert permutation so that it can be used in MatCholeskyFactorNumeric() */
1086   if (a->icol) {ierr = ISInvertPermutation(row,PETSC_DECIDE, &a->icol);CHKERRQ(ierr);}
1087   ierr = PetscLogObjectParent(inA,a->icol);CHKERRQ(ierr);
1088 
1089   if (!a->solve_work) {
1090     ierr = PetscMalloc((inA->rmap->N+inA->rmap->bs)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr);
1091     ierr = PetscLogObjectMemory(inA,(inA->rmap->N+inA->rmap->bs)*sizeof(PetscScalar));CHKERRQ(ierr);
1092   }
1093 
1094   ierr = MatCholeskyFactorNumeric(outA,inA,info);CHKERRQ(ierr);
1095   PetscFunctionReturn(0);
1096 }
1097 
1098 #undef __FUNCT__
1099 #define __FUNCT__ "MatSeqSBAIJSetColumnIndices_SeqSBAIJ"
1100 PetscErrorCode  MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,PetscInt *indices)
1101 {
1102   Mat_SeqSBAIJ   *baij = (Mat_SeqSBAIJ*)mat->data;
1103   PetscInt       i,nz,n;
1104   PetscErrorCode ierr;
1105 
1106   PetscFunctionBegin;
1107   nz = baij->maxnz;
1108   n  = mat->cmap->n;
1109   for (i=0; i<nz; i++) baij->j[i] = indices[i];
1110 
1111   baij->nz = nz;
1112   for (i=0; i<n; i++) baij->ilen[i] = baij->imax[i];
1113 
1114   ierr = MatSetOption(mat,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
1115   PetscFunctionReturn(0);
1116 }
1117 
1118 #undef __FUNCT__
1119 #define __FUNCT__ "MatSeqSBAIJSetColumnIndices"
1120 /*@
1121   MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows
1122   in the matrix.
1123 
1124   Input Parameters:
1125   +  mat     - the SeqSBAIJ matrix
1126   -  indices - the column indices
1127 
1128   Level: advanced
1129 
1130   Notes:
1131   This can be called if you have precomputed the nonzero structure of the
1132   matrix and want to provide it to the matrix object to improve the performance
1133   of the MatSetValues() operation.
1134 
1135   You MUST have set the correct numbers of nonzeros per row in the call to
1136   MatCreateSeqSBAIJ(), and the columns indices MUST be sorted.
1137 
1138   MUST be called before any calls to MatSetValues()
1139 
1140   .seealso: MatCreateSeqSBAIJ
1141 @*/
1142 PetscErrorCode  MatSeqSBAIJSetColumnIndices(Mat mat,PetscInt *indices)
1143 {
1144   PetscErrorCode ierr;
1145 
1146   PetscFunctionBegin;
1147   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
1148   PetscValidPointer(indices,2);
1149   ierr = PetscUseMethod(mat,"MatSeqSBAIJSetColumnIndices_C",(Mat,PetscInt*),(mat,indices));CHKERRQ(ierr);
1150   PetscFunctionReturn(0);
1151 }
1152 
1153 #undef __FUNCT__
1154 #define __FUNCT__ "MatCopy_SeqSBAIJ"
1155 PetscErrorCode MatCopy_SeqSBAIJ(Mat A,Mat B,MatStructure str)
1156 {
1157   PetscErrorCode ierr;
1158 
1159   PetscFunctionBegin;
1160   /* If the two matrices have the same copy implementation, use fast copy. */
1161   if (str == SAME_NONZERO_PATTERN && (A->ops->copy == B->ops->copy)) {
1162     Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1163     Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data;
1164 
1165     if (a->i[A->rmap->N] != b->i[B->rmap->N]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of nonzeros in two matrices are different");
1166     ierr = PetscMemcpy(b->a,a->a,(a->i[A->rmap->N])*sizeof(PetscScalar));CHKERRQ(ierr);
1167   } else {
1168     ierr = MatGetRowUpperTriangular(A);CHKERRQ(ierr);
1169     ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr);
1170     ierr = MatRestoreRowUpperTriangular(A);CHKERRQ(ierr);
1171   }
1172   PetscFunctionReturn(0);
1173 }
1174 
1175 #undef __FUNCT__
1176 #define __FUNCT__ "MatSetUp_SeqSBAIJ"
1177 PetscErrorCode MatSetUp_SeqSBAIJ(Mat A)
1178 {
1179   PetscErrorCode ierr;
1180 
1181   PetscFunctionBegin;
1182   ierr =  MatSeqSBAIJSetPreallocation_SeqSBAIJ(A,A->rmap->bs,PETSC_DEFAULT,0);CHKERRQ(ierr);
1183   PetscFunctionReturn(0);
1184 }
1185 
1186 #undef __FUNCT__
1187 #define __FUNCT__ "MatSeqSBAIJGetArray_SeqSBAIJ"
1188 PetscErrorCode MatSeqSBAIJGetArray_SeqSBAIJ(Mat A,PetscScalar *array[])
1189 {
1190   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1191 
1192   PetscFunctionBegin;
1193   *array = a->a;
1194   PetscFunctionReturn(0);
1195 }
1196 
1197 #undef __FUNCT__
1198 #define __FUNCT__ "MatSeqSBAIJRestoreArray_SeqSBAIJ"
1199 PetscErrorCode MatSeqSBAIJRestoreArray_SeqSBAIJ(Mat A,PetscScalar *array[])
1200 {
1201   PetscFunctionBegin;
1202   PetscFunctionReturn(0);
1203 }
1204 
1205 #undef __FUNCT__
1206 #define __FUNCT__ "MatAXPY_SeqSBAIJ"
1207 PetscErrorCode MatAXPY_SeqSBAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str)
1208 {
1209   Mat_SeqSBAIJ   *x=(Mat_SeqSBAIJ*)X->data, *y=(Mat_SeqSBAIJ*)Y->data;
1210   PetscErrorCode ierr;
1211   PetscInt       i,bs=Y->rmap->bs,bs2=bs*bs,j;
1212   PetscBLASInt   one = 1;
1213 
1214   PetscFunctionBegin;
1215   if (str == SAME_NONZERO_PATTERN) {
1216     PetscScalar  alpha = a;
1217     PetscBLASInt bnz;
1218     ierr = PetscBLASIntCast(x->nz*bs2,&bnz);CHKERRQ(ierr);
1219     PetscStackCall("BLASaxpy",BLASaxpy_(&bnz,&alpha,x->a,&one,y->a,&one));
1220   } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
1221     if (y->xtoy && y->XtoY != X) {
1222       ierr = PetscFree(y->xtoy);CHKERRQ(ierr);
1223       ierr = MatDestroy(&y->XtoY);CHKERRQ(ierr);
1224     }
1225     if (!y->xtoy) { /* get xtoy */
1226       ierr    = MatAXPYGetxtoy_Private(x->mbs,x->i,x->j,NULL, y->i,y->j,NULL, &y->xtoy);CHKERRQ(ierr);
1227       y->XtoY = X;
1228     }
1229     for (i=0; i<x->nz; i++) {
1230       j = 0;
1231       while (j < bs2) {
1232         y->a[bs2*y->xtoy[i]+j] += a*(x->a[bs2*i+j]);
1233         j++;
1234       }
1235     }
1236     ierr = PetscInfo3(Y,"ratio of nnz_s(X)/nnz_s(Y): %D/%D = %G\n",bs2*x->nz,bs2*y->nz,(PetscReal)(bs2*x->nz)/(bs2*y->nz));CHKERRQ(ierr);
1237   } else {
1238     ierr = MatGetRowUpperTriangular(X);CHKERRQ(ierr);
1239     ierr = MatAXPY_Basic(Y,a,X,str);CHKERRQ(ierr);
1240     ierr = MatRestoreRowUpperTriangular(X);CHKERRQ(ierr);
1241   }
1242   PetscFunctionReturn(0);
1243 }
1244 
1245 #undef __FUNCT__
1246 #define __FUNCT__ "MatIsSymmetric_SeqSBAIJ"
1247 PetscErrorCode MatIsSymmetric_SeqSBAIJ(Mat A,PetscReal tol,PetscBool  *flg)
1248 {
1249   PetscFunctionBegin;
1250   *flg = PETSC_TRUE;
1251   PetscFunctionReturn(0);
1252 }
1253 
1254 #undef __FUNCT__
1255 #define __FUNCT__ "MatIsStructurallySymmetric_SeqSBAIJ"
1256 PetscErrorCode MatIsStructurallySymmetric_SeqSBAIJ(Mat A,PetscBool  *flg)
1257 {
1258   PetscFunctionBegin;
1259   *flg = PETSC_TRUE;
1260   PetscFunctionReturn(0);
1261 }
1262 
1263 #undef __FUNCT__
1264 #define __FUNCT__ "MatIsHermitian_SeqSBAIJ"
1265 PetscErrorCode MatIsHermitian_SeqSBAIJ(Mat A,PetscReal tol,PetscBool  *flg)
1266 {
1267   PetscFunctionBegin;
1268   *flg = PETSC_FALSE;
1269   PetscFunctionReturn(0);
1270 }
1271 
1272 #undef __FUNCT__
1273 #define __FUNCT__ "MatRealPart_SeqSBAIJ"
1274 PetscErrorCode MatRealPart_SeqSBAIJ(Mat A)
1275 {
1276   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1277   PetscInt     i,nz = a->bs2*a->i[a->mbs];
1278   MatScalar    *aa = a->a;
1279 
1280   PetscFunctionBegin;
1281   for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]);
1282   PetscFunctionReturn(0);
1283 }
1284 
1285 #undef __FUNCT__
1286 #define __FUNCT__ "MatImaginaryPart_SeqSBAIJ"
1287 PetscErrorCode MatImaginaryPart_SeqSBAIJ(Mat A)
1288 {
1289   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1290   PetscInt     i,nz = a->bs2*a->i[a->mbs];
1291   MatScalar    *aa = a->a;
1292 
1293   PetscFunctionBegin;
1294   for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]);
1295   PetscFunctionReturn(0);
1296 }
1297 
1298 #undef __FUNCT__
1299 #define __FUNCT__ "MatZeroRowsColumns_SeqSBAIJ"
1300 PetscErrorCode MatZeroRowsColumns_SeqSBAIJ(Mat A,PetscInt is_n,const PetscInt is_idx[],PetscScalar diag,Vec x, Vec b)
1301 {
1302   Mat_SeqSBAIJ      *baij=(Mat_SeqSBAIJ*)A->data;
1303   PetscErrorCode    ierr;
1304   PetscInt          i,j,k,count;
1305   PetscInt          bs   =A->rmap->bs,bs2=baij->bs2,row,col;
1306   PetscScalar       zero = 0.0;
1307   MatScalar         *aa;
1308   const PetscScalar *xx;
1309   PetscScalar       *bb;
1310   PetscBool         *zeroed,vecs = PETSC_FALSE;
1311 
1312   PetscFunctionBegin;
1313   /* fix right hand side if needed */
1314   if (x && b) {
1315     ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
1316     ierr = VecGetArray(b,&bb);CHKERRQ(ierr);
1317     vecs = PETSC_TRUE;
1318   }
1319   A->same_nonzero = PETSC_TRUE;
1320 
1321   /* zero the columns */
1322   ierr = PetscMalloc(A->rmap->n*sizeof(PetscBool),&zeroed);CHKERRQ(ierr);
1323   ierr = PetscMemzero(zeroed,A->rmap->n*sizeof(PetscBool));CHKERRQ(ierr);
1324   for (i=0; i<is_n; i++) {
1325     if (is_idx[i] < 0 || is_idx[i] >= A->rmap->N) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range",is_idx[i]);
1326     zeroed[is_idx[i]] = PETSC_TRUE;
1327   }
1328   if (vecs) {
1329     for (i=0; i<A->rmap->N; i++) {
1330       row = i/bs;
1331       for (j=baij->i[row]; j<baij->i[row+1]; j++) {
1332         for (k=0; k<bs; k++) {
1333           col = bs*baij->j[j] + k;
1334           if (col <= i) continue;
1335           aa = ((MatScalar*)(baij->a)) + j*bs2 + (i%bs) + bs*k;
1336           if (!zeroed[i] && zeroed[col]) bb[i]   -= aa[0]*xx[col];
1337           if (zeroed[i] && !zeroed[col]) bb[col] -= aa[0]*xx[i];
1338         }
1339       }
1340     }
1341     for (i=0; i<is_n; i++) bb[is_idx[i]] = diag*xx[is_idx[i]];
1342   }
1343 
1344   for (i=0; i<A->rmap->N; i++) {
1345     if (!zeroed[i]) {
1346       row = i/bs;
1347       for (j=baij->i[row]; j<baij->i[row+1]; j++) {
1348         for (k=0; k<bs; k++) {
1349           col = bs*baij->j[j] + k;
1350           if (zeroed[col]) {
1351             aa = ((MatScalar*)(baij->a)) + j*bs2 + (i%bs) + bs*k;
1352             aa[0] = 0.0;
1353           }
1354         }
1355       }
1356     }
1357   }
1358   ierr = PetscFree(zeroed);CHKERRQ(ierr);
1359   if (vecs) {
1360     ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
1361     ierr = VecRestoreArray(b,&bb);CHKERRQ(ierr);
1362   }
1363 
1364   /* zero the rows */
1365   for (i=0; i<is_n; i++) {
1366     row   = is_idx[i];
1367     count = (baij->i[row/bs +1] - baij->i[row/bs])*bs;
1368     aa    = ((MatScalar*)(baij->a)) + baij->i[row/bs]*bs2 + (row%bs);
1369     for (k=0; k<count; k++) {
1370       aa[0] =  zero;
1371       aa   += bs;
1372     }
1373     if (diag != 0.0) {
1374       ierr = (*A->ops->setvalues)(A,1,&row,1,&row,&diag,INSERT_VALUES);CHKERRQ(ierr);
1375     }
1376   }
1377   ierr = MatAssemblyEnd_SeqSBAIJ(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1378   PetscFunctionReturn(0);
1379 }
1380 
1381 /* -------------------------------------------------------------------*/
1382 static struct _MatOps MatOps_Values = {
1383         MatSetValues_SeqSBAIJ,
1384         MatGetRow_SeqSBAIJ,
1385         MatRestoreRow_SeqSBAIJ,
1386         MatMult_SeqSBAIJ_N,
1387 /*  4*/ MatMultAdd_SeqSBAIJ_N,
1388         MatMult_SeqSBAIJ_N,       /* transpose versions are same as non-transpose versions */
1389         MatMultAdd_SeqSBAIJ_N,
1390         0,
1391         0,
1392         0,
1393 /* 10*/ 0,
1394         0,
1395         MatCholeskyFactor_SeqSBAIJ,
1396         MatSOR_SeqSBAIJ,
1397         MatTranspose_SeqSBAIJ,
1398 /* 15*/ MatGetInfo_SeqSBAIJ,
1399         MatEqual_SeqSBAIJ,
1400         MatGetDiagonal_SeqSBAIJ,
1401         MatDiagonalScale_SeqSBAIJ,
1402         MatNorm_SeqSBAIJ,
1403 /* 20*/ 0,
1404         MatAssemblyEnd_SeqSBAIJ,
1405         MatSetOption_SeqSBAIJ,
1406         MatZeroEntries_SeqSBAIJ,
1407 /* 24*/ 0,
1408         0,
1409         0,
1410         0,
1411         0,
1412 /* 29*/ MatSetUp_SeqSBAIJ,
1413         0,
1414         0,
1415         0,
1416         0,
1417 /* 34*/ MatDuplicate_SeqSBAIJ,
1418         0,
1419         0,
1420         0,
1421         MatICCFactor_SeqSBAIJ,
1422 /* 39*/ MatAXPY_SeqSBAIJ,
1423         MatGetSubMatrices_SeqSBAIJ,
1424         MatIncreaseOverlap_SeqSBAIJ,
1425         MatGetValues_SeqSBAIJ,
1426         MatCopy_SeqSBAIJ,
1427 /* 44*/ 0,
1428         MatScale_SeqSBAIJ,
1429         0,
1430         0,
1431         MatZeroRowsColumns_SeqSBAIJ,
1432 /* 49*/ 0,
1433         MatGetRowIJ_SeqSBAIJ,
1434         MatRestoreRowIJ_SeqSBAIJ,
1435         0,
1436         0,
1437 /* 54*/ 0,
1438         0,
1439         0,
1440         0,
1441         MatSetValuesBlocked_SeqSBAIJ,
1442 /* 59*/ MatGetSubMatrix_SeqSBAIJ,
1443         0,
1444         0,
1445         0,
1446         0,
1447 /* 64*/ 0,
1448         0,
1449         0,
1450         0,
1451         0,
1452 /* 69*/ MatGetRowMaxAbs_SeqSBAIJ,
1453         0,
1454         0,
1455         0,
1456         0,
1457 /* 74*/ 0,
1458         0,
1459         0,
1460         0,
1461         0,
1462 /* 79*/ 0,
1463         0,
1464         0,
1465         MatGetInertia_SeqSBAIJ,
1466         MatLoad_SeqSBAIJ,
1467 /* 84*/ MatIsSymmetric_SeqSBAIJ,
1468         MatIsHermitian_SeqSBAIJ,
1469         MatIsStructurallySymmetric_SeqSBAIJ,
1470         0,
1471         0,
1472 /* 89*/ 0,
1473         0,
1474         0,
1475         0,
1476         0,
1477 /* 94*/ 0,
1478         0,
1479         0,
1480         0,
1481         0,
1482 /* 99*/ 0,
1483         0,
1484         0,
1485         0,
1486         0,
1487 /*104*/ 0,
1488         MatRealPart_SeqSBAIJ,
1489         MatImaginaryPart_SeqSBAIJ,
1490         MatGetRowUpperTriangular_SeqSBAIJ,
1491         MatRestoreRowUpperTriangular_SeqSBAIJ,
1492 /*109*/ 0,
1493         0,
1494         0,
1495         0,
1496         MatMissingDiagonal_SeqSBAIJ,
1497 /*114*/ 0,
1498         0,
1499         0,
1500         0,
1501         0,
1502 /*119*/ 0,
1503         0,
1504         0,
1505         0
1506 };
1507 
1508 #undef __FUNCT__
1509 #define __FUNCT__ "MatStoreValues_SeqSBAIJ"
1510 PetscErrorCode  MatStoreValues_SeqSBAIJ(Mat mat)
1511 {
1512   Mat_SeqSBAIJ   *aij = (Mat_SeqSBAIJ*)mat->data;
1513   PetscInt       nz   = aij->i[mat->rmap->N]*mat->rmap->bs*aij->bs2;
1514   PetscErrorCode ierr;
1515 
1516   PetscFunctionBegin;
1517   if (aij->nonew != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first");
1518 
1519   /* allocate space for values if not already there */
1520   if (!aij->saved_values) {
1521     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr);
1522   }
1523 
1524   /* copy values over */
1525   ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr);
1526   PetscFunctionReturn(0);
1527 }
1528 
1529 #undef __FUNCT__
1530 #define __FUNCT__ "MatRetrieveValues_SeqSBAIJ"
1531 PetscErrorCode  MatRetrieveValues_SeqSBAIJ(Mat mat)
1532 {
1533   Mat_SeqSBAIJ   *aij = (Mat_SeqSBAIJ*)mat->data;
1534   PetscErrorCode ierr;
1535   PetscInt       nz = aij->i[mat->rmap->N]*mat->rmap->bs*aij->bs2;
1536 
1537   PetscFunctionBegin;
1538   if (aij->nonew != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NEW_NONZERO_LOCATIONS,PETSC_FALSE);first");
1539   if (!aij->saved_values) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Must call MatStoreValues(A);first");
1540 
1541   /* copy values over */
1542   ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr);
1543   PetscFunctionReturn(0);
1544 }
1545 
1546 #undef __FUNCT__
1547 #define __FUNCT__ "MatSeqSBAIJSetPreallocation_SeqSBAIJ"
1548 PetscErrorCode  MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat B,PetscInt bs,PetscInt nz,PetscInt *nnz)
1549 {
1550   Mat_SeqSBAIJ   *b = (Mat_SeqSBAIJ*)B->data;
1551   PetscErrorCode ierr;
1552   PetscInt       i,mbs,bs2;
1553   PetscBool      skipallocation = PETSC_FALSE,flg = PETSC_FALSE,realalloc = PETSC_FALSE;
1554 
1555   PetscFunctionBegin;
1556   if (nz >= 0 || nnz) realalloc = PETSC_TRUE;
1557   B->preallocated = PETSC_TRUE;
1558 
1559   ierr = PetscLayoutSetBlockSize(B->rmap,bs);CHKERRQ(ierr);
1560   ierr = PetscLayoutSetBlockSize(B->cmap,bs);CHKERRQ(ierr);
1561   ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr);
1562   ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr);
1563   ierr = PetscLayoutGetBlockSize(B->rmap,&bs);CHKERRQ(ierr);
1564 
1565   mbs = B->rmap->N/bs;
1566   bs2 = bs*bs;
1567 
1568   if (mbs*bs != B->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize");
1569 
1570   if (nz == MAT_SKIP_ALLOCATION) {
1571     skipallocation = PETSC_TRUE;
1572     nz             = 0;
1573   }
1574 
1575   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3;
1576   if (nz < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %D",nz);
1577   if (nnz) {
1578     for (i=0; i<mbs; i++) {
1579       if (nnz[i] < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %D value %D",i,nnz[i]);
1580       if (nnz[i] > mbs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be greater than block row length: local row %D value %D rowlength %D",i,nnz[i],mbs);
1581     }
1582   }
1583 
1584   B->ops->mult             = MatMult_SeqSBAIJ_N;
1585   B->ops->multadd          = MatMultAdd_SeqSBAIJ_N;
1586   B->ops->multtranspose    = MatMult_SeqSBAIJ_N;
1587   B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_N;
1588 
1589   ierr  = PetscOptionsGetBool(((PetscObject)B)->prefix,"-mat_no_unroll",&flg,NULL);CHKERRQ(ierr);
1590   if (!flg) {
1591     switch (bs) {
1592     case 1:
1593       B->ops->mult             = MatMult_SeqSBAIJ_1;
1594       B->ops->multadd          = MatMultAdd_SeqSBAIJ_1;
1595       B->ops->multtranspose    = MatMult_SeqSBAIJ_1;
1596       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_1;
1597       break;
1598     case 2:
1599       B->ops->mult             = MatMult_SeqSBAIJ_2;
1600       B->ops->multadd          = MatMultAdd_SeqSBAIJ_2;
1601       B->ops->multtranspose    = MatMult_SeqSBAIJ_2;
1602       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_2;
1603       break;
1604     case 3:
1605       B->ops->mult             = MatMult_SeqSBAIJ_3;
1606       B->ops->multadd          = MatMultAdd_SeqSBAIJ_3;
1607       B->ops->multtranspose    = MatMult_SeqSBAIJ_3;
1608       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_3;
1609       break;
1610     case 4:
1611       B->ops->mult             = MatMult_SeqSBAIJ_4;
1612       B->ops->multadd          = MatMultAdd_SeqSBAIJ_4;
1613       B->ops->multtranspose    = MatMult_SeqSBAIJ_4;
1614       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_4;
1615       break;
1616     case 5:
1617       B->ops->mult             = MatMult_SeqSBAIJ_5;
1618       B->ops->multadd          = MatMultAdd_SeqSBAIJ_5;
1619       B->ops->multtranspose    = MatMult_SeqSBAIJ_5;
1620       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_5;
1621       break;
1622     case 6:
1623       B->ops->mult             = MatMult_SeqSBAIJ_6;
1624       B->ops->multadd          = MatMultAdd_SeqSBAIJ_6;
1625       B->ops->multtranspose    = MatMult_SeqSBAIJ_6;
1626       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_6;
1627       break;
1628     case 7:
1629       B->ops->mult             = MatMult_SeqSBAIJ_7;
1630       B->ops->multadd          = MatMultAdd_SeqSBAIJ_7;
1631       B->ops->multtranspose    = MatMult_SeqSBAIJ_7;
1632       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_7;
1633       break;
1634     }
1635   }
1636 
1637   b->mbs = mbs;
1638   b->nbs = mbs;
1639   if (!skipallocation) {
1640     if (!b->imax) {
1641       ierr = PetscMalloc2(mbs,PetscInt,&b->imax,mbs,PetscInt,&b->ilen);CHKERRQ(ierr);
1642 
1643       b->free_imax_ilen = PETSC_TRUE;
1644 
1645       ierr = PetscLogObjectMemory(B,2*mbs*sizeof(PetscInt));CHKERRQ(ierr);
1646     }
1647     if (!nnz) {
1648       if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
1649       else if (nz <= 0) nz = 1;
1650       for (i=0; i<mbs; i++) b->imax[i] = nz;
1651       nz = nz*mbs; /* total nz */
1652     } else {
1653       nz = 0;
1654       for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
1655     }
1656     /* b->ilen will count nonzeros in each block row so far. */
1657     for (i=0; i<mbs; i++) b->ilen[i] = 0;
1658     /* nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */
1659 
1660     /* allocate the matrix space */
1661     ierr = MatSeqXAIJFreeAIJ(B,&b->a,&b->j,&b->i);CHKERRQ(ierr);
1662     ierr = PetscMalloc3(bs2*nz,PetscScalar,&b->a,nz,PetscInt,&b->j,B->rmap->N+1,PetscInt,&b->i);CHKERRQ(ierr);
1663     ierr = PetscLogObjectMemory(B,(B->rmap->N+1)*sizeof(PetscInt)+nz*(bs2*sizeof(PetscScalar)+sizeof(PetscInt)));CHKERRQ(ierr);
1664     ierr = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr);
1665     ierr = PetscMemzero(b->j,nz*sizeof(PetscInt));CHKERRQ(ierr);
1666 
1667     b->singlemalloc = PETSC_TRUE;
1668 
1669     /* pointer to beginning of each row */
1670     b->i[0] = 0;
1671     for (i=1; i<mbs+1; i++) b->i[i] = b->i[i-1] + b->imax[i-1];
1672 
1673     b->free_a  = PETSC_TRUE;
1674     b->free_ij = PETSC_TRUE;
1675   } else {
1676     b->free_a  = PETSC_FALSE;
1677     b->free_ij = PETSC_FALSE;
1678   }
1679 
1680   B->rmap->bs = bs;
1681   b->bs2      = bs2;
1682   b->nz       = 0;
1683   b->maxnz    = nz;
1684 
1685   b->inew    = 0;
1686   b->jnew    = 0;
1687   b->anew    = 0;
1688   b->a2anew  = 0;
1689   b->permute = PETSC_FALSE;
1690   if (realalloc) {ierr = MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);}
1691   PetscFunctionReturn(0);
1692 }
1693 
1694 /*
1695    This is used to set the numeric factorization for both Cholesky and ICC symbolic factorization
1696 */
1697 #undef __FUNCT__
1698 #define __FUNCT__ "MatSeqSBAIJSetNumericFactorization_inplace"
1699 PetscErrorCode MatSeqSBAIJSetNumericFactorization_inplace(Mat B,PetscBool natural)
1700 {
1701   PetscErrorCode ierr;
1702   PetscBool      flg = PETSC_FALSE;
1703   PetscInt       bs  = B->rmap->bs;
1704 
1705   PetscFunctionBegin;
1706   ierr = PetscOptionsGetBool(((PetscObject)B)->prefix,"-mat_no_unroll",&flg,NULL);CHKERRQ(ierr);
1707   if (flg) bs = 8;
1708 
1709   if (!natural) {
1710     switch (bs) {
1711     case 1:
1712       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_inplace;
1713       break;
1714     case 2:
1715       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2;
1716       break;
1717     case 3:
1718       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3;
1719       break;
1720     case 4:
1721       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4;
1722       break;
1723     case 5:
1724       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5;
1725       break;
1726     case 6:
1727       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6;
1728       break;
1729     case 7:
1730       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7;
1731       break;
1732     default:
1733       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N;
1734       break;
1735     }
1736   } else {
1737     switch (bs) {
1738     case 1:
1739       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering_inplace;
1740       break;
1741     case 2:
1742       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering;
1743       break;
1744     case 3:
1745       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering;
1746       break;
1747     case 4:
1748       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering;
1749       break;
1750     case 5:
1751       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering;
1752       break;
1753     case 6:
1754       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering;
1755       break;
1756     case 7:
1757       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering;
1758       break;
1759     default:
1760       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering;
1761       break;
1762     }
1763   }
1764   PetscFunctionReturn(0);
1765 }
1766 
1767 PETSC_EXTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqAIJ(Mat, MatType,MatReuse,Mat*);
1768 PETSC_EXTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqBAIJ(Mat, MatType,MatReuse,Mat*);
1769 
1770 #undef __FUNCT__
1771 #define __FUNCT__ "MatGetFactor_seqsbaij_petsc"
1772 PETSC_EXTERN PetscErrorCode MatGetFactor_seqsbaij_petsc(Mat A,MatFactorType ftype,Mat *B)
1773 {
1774   PetscInt       n = A->rmap->n;
1775   PetscErrorCode ierr;
1776 
1777   PetscFunctionBegin;
1778 #if defined(PETSC_USE_COMPLEX)
1779   if (A->hermitian) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Hermitian Factor is not supported");
1780 #endif
1781   ierr = MatCreate(PetscObjectComm((PetscObject)A),B);CHKERRQ(ierr);
1782   ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr);
1783   if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) {
1784     ierr = MatSetType(*B,MATSEQSBAIJ);CHKERRQ(ierr);
1785     ierr = MatSeqSBAIJSetPreallocation(*B,A->rmap->bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr);
1786 
1787     (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqSBAIJ;
1788     (*B)->ops->iccfactorsymbolic      = MatICCFactorSymbolic_SeqSBAIJ;
1789   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
1790   (*B)->factortype = ftype;
1791   PetscFunctionReturn(0);
1792 }
1793 
1794 #undef __FUNCT__
1795 #define __FUNCT__ "MatGetFactorAvailable_seqsbaij_petsc"
1796 PetscErrorCode MatGetFactorAvailable_seqsbaij_petsc(Mat A,MatFactorType ftype,PetscBool  *flg)
1797 {
1798   PetscFunctionBegin;
1799   if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) {
1800     *flg = PETSC_TRUE;
1801   } else {
1802     *flg = PETSC_FALSE;
1803   }
1804   PetscFunctionReturn(0);
1805 }
1806 
1807 #if defined(PETSC_HAVE_MUMPS)
1808 PETSC_EXTERN PetscErrorCode MatGetFactor_sbaij_mumps(Mat,MatFactorType,Mat*);
1809 #endif
1810 #if defined(PETSC_HAVE_PASTIX)
1811 PETSC_EXTERN PetscErrorCode MatGetFactor_seqsbaij_pastix(Mat,MatFactorType,Mat*);
1812 #endif
1813 #if defined(PETSC_HAVE_CHOLMOD)
1814 PETSC_EXTERN PetscErrorCode MatGetFactor_seqsbaij_cholmod(Mat,MatFactorType,Mat*);
1815 #endif
1816 PETSC_EXTERN PetscErrorCode MatGetFactor_seqsbaij_sbstrm(Mat,MatFactorType,Mat*);
1817 
1818 /*MC
1819   MATSEQSBAIJ - MATSEQSBAIJ = "seqsbaij" - A matrix type to be used for sequential symmetric block sparse matrices,
1820   based on block compressed sparse row format.  Only the upper triangular portion of the matrix is stored.
1821 
1822   For complex numbers by default this matrix is symmetric, NOT Hermitian symmetric. To make it Hermitian symmetric you
1823   can call MatSetOption(Mat, MAT_HERMITIAN); after MatAssemblyEnd()
1824 
1825   Options Database Keys:
1826   . -mat_type seqsbaij - sets the matrix type to "seqsbaij" during a call to MatSetFromOptions()
1827 
1828   Notes: By default if you insert values into the lower triangular part of the matrix they are simply ignored (since they are not
1829      stored and it is assumed they symmetric to the upper triangular). If you call MatSetOption(Mat,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_FALSE) or use
1830      the options database -mat_ignore_lower_triangular false it will generate an error if you try to set a value in the lower triangular portion.
1831 
1832 
1833   Level: beginner
1834 
1835   .seealso: MatCreateSeqSBAIJ
1836 M*/
1837 
1838 PETSC_EXTERN PetscErrorCode MatConvert_SeqSBAIJ_SeqSBSTRM(Mat, MatType,MatReuse,Mat*);
1839 
1840 #undef __FUNCT__
1841 #define __FUNCT__ "MatCreate_SeqSBAIJ"
1842 PETSC_EXTERN PetscErrorCode MatCreate_SeqSBAIJ(Mat B)
1843 {
1844   Mat_SeqSBAIJ   *b;
1845   PetscErrorCode ierr;
1846   PetscMPIInt    size;
1847   PetscBool      no_unroll = PETSC_FALSE,no_inode = PETSC_FALSE;
1848 
1849   PetscFunctionBegin;
1850   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)B),&size);CHKERRQ(ierr);
1851   if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
1852 
1853   ierr    = PetscNewLog(B,Mat_SeqSBAIJ,&b);CHKERRQ(ierr);
1854   B->data = (void*)b;
1855   ierr    = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1856 
1857   B->ops->destroy    = MatDestroy_SeqSBAIJ;
1858   B->ops->view       = MatView_SeqSBAIJ;
1859   b->row             = 0;
1860   b->icol            = 0;
1861   b->reallocs        = 0;
1862   b->saved_values    = 0;
1863   b->inode.limit     = 5;
1864   b->inode.max_limit = 5;
1865 
1866   b->roworiented        = PETSC_TRUE;
1867   b->nonew              = 0;
1868   b->diag               = 0;
1869   b->solve_work         = 0;
1870   b->mult_work          = 0;
1871   B->spptr              = 0;
1872   B->info.nz_unneeded   = (PetscReal)b->maxnz*b->bs2;
1873   b->keepnonzeropattern = PETSC_FALSE;
1874   b->xtoy               = 0;
1875   b->XtoY               = 0;
1876 
1877   b->inew    = 0;
1878   b->jnew    = 0;
1879   b->anew    = 0;
1880   b->a2anew  = 0;
1881   b->permute = PETSC_FALSE;
1882 
1883   b->ignore_ltriangular = PETSC_TRUE;
1884 
1885   ierr = PetscOptionsGetBool(((PetscObject)B)->prefix,"-mat_ignore_lower_triangular",&b->ignore_ltriangular,NULL);CHKERRQ(ierr);
1886 
1887   b->getrow_utriangular = PETSC_FALSE;
1888 
1889   ierr = PetscOptionsGetBool(((PetscObject)B)->prefix,"-mat_getrow_uppertriangular",&b->getrow_utriangular,NULL);CHKERRQ(ierr);
1890 
1891 #if defined(PETSC_HAVE_PASTIX)
1892   ierr = PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_pastix_C","MatGetFactor_seqsbaij_pastix",MatGetFactor_seqsbaij_pastix);CHKERRQ(ierr);
1893 #endif
1894 #if defined(PETSC_HAVE_MUMPS)
1895   ierr = PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_mumps_C","MatGetFactor_sbaij_mumps",MatGetFactor_sbaij_mumps);CHKERRQ(ierr);
1896 #endif
1897 #if defined(PETSC_HAVE_CHOLMOD)
1898   ierr = PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_cholmod_C","MatGetFactor_seqsbaij_cholmod",MatGetFactor_seqsbaij_cholmod);CHKERRQ(ierr);
1899 #endif
1900   ierr = PetscObjectComposeFunction((PetscObject)B,"MatGetFactorAvailable_petsc_C","MatGetFactorAvailable_seqsbaij_petsc",MatGetFactorAvailable_seqsbaij_petsc);CHKERRQ(ierr);
1901   ierr = PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_petsc_C","MatGetFactor_seqsbaij_petsc",MatGetFactor_seqsbaij_petsc);CHKERRQ(ierr);
1902   ierr = PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_sbstrm_C","MatGetFactor_seqsbaij_sbstrm",MatGetFactor_seqsbaij_sbstrm);CHKERRQ(ierr);
1903   ierr = PetscObjectComposeFunction((PetscObject)B,"MatStoreValues_C","MatStoreValues_SeqSBAIJ",MatStoreValues_SeqSBAIJ);CHKERRQ(ierr);
1904   ierr = PetscObjectComposeFunction((PetscObject)B,"MatRetrieveValues_C","MatRetrieveValues_SeqSBAIJ",MatRetrieveValues_SeqSBAIJ);CHKERRQ(ierr);
1905   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C","MatSeqSBAIJSetColumnIndices_SeqSBAIJ",MatSeqSBAIJSetColumnIndices_SeqSBAIJ);CHKERRQ(ierr);
1906   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_seqaij_C","MatConvert_SeqSBAIJ_SeqAIJ",MatConvert_SeqSBAIJ_SeqAIJ);CHKERRQ(ierr);
1907   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_seqbaij_C","MatConvert_SeqSBAIJ_SeqBAIJ",MatConvert_SeqSBAIJ_SeqBAIJ);CHKERRQ(ierr);
1908   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSeqSBAIJSetPreallocation_C","MatSeqSBAIJSetPreallocation_SeqSBAIJ",MatSeqSBAIJSetPreallocation_SeqSBAIJ);CHKERRQ(ierr);
1909   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_seqsbstrm_C","MatConvert_SeqSBAIJ_SeqSBSTRM",MatConvert_SeqSBAIJ_SeqSBSTRM);CHKERRQ(ierr);
1910 
1911   B->symmetric                  = PETSC_TRUE;
1912   B->structurally_symmetric     = PETSC_TRUE;
1913   B->symmetric_set              = PETSC_TRUE;
1914   B->structurally_symmetric_set = PETSC_TRUE;
1915 
1916   ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQSBAIJ);CHKERRQ(ierr);
1917 
1918   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)B),((PetscObject)B)->prefix,"Options for SEQSBAIJ matrix","Mat");CHKERRQ(ierr);
1919   ierr = PetscOptionsBool("-mat_no_unroll","Do not optimize for inodes (slower)",NULL,no_unroll,&no_unroll,NULL);CHKERRQ(ierr);
1920   if (no_unroll) {
1921     ierr = PetscInfo(B,"Not using Inode routines due to -mat_no_unroll\n");CHKERRQ(ierr);
1922   }
1923   ierr = PetscOptionsBool("-mat_no_inode","Do not optimize for inodes (slower)",NULL,no_inode,&no_inode,NULL);CHKERRQ(ierr);
1924   if (no_inode) {
1925     ierr = PetscInfo(B,"Not using Inode routines due to -mat_no_inode\n");CHKERRQ(ierr);
1926   }
1927   ierr = PetscOptionsInt("-mat_inode_limit","Do not use inodes larger then this value",NULL,b->inode.limit,&b->inode.limit,NULL);CHKERRQ(ierr);
1928   ierr = PetscOptionsEnd();CHKERRQ(ierr);
1929   b->inode.use = (PetscBool)(!(no_unroll || no_inode));
1930   if (b->inode.limit > b->inode.max_limit) b->inode.limit = b->inode.max_limit;
1931   PetscFunctionReturn(0);
1932 }
1933 
1934 #undef __FUNCT__
1935 #define __FUNCT__ "MatSeqSBAIJSetPreallocation"
1936 /*@C
1937    MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block
1938    compressed row) format.  For good matrix assembly performance the
1939    user should preallocate the matrix storage by setting the parameter nz
1940    (or the array nnz).  By setting these parameters accurately, performance
1941    during matrix assembly can be increased by more than a factor of 50.
1942 
1943    Collective on Mat
1944 
1945    Input Parameters:
1946 +  A - the symmetric matrix
1947 .  bs - size of block
1948 .  nz - number of block nonzeros per block row (same for all rows)
1949 -  nnz - array containing the number of block nonzeros in the upper triangular plus
1950          diagonal portion of each block (possibly different for each block row) or NULL
1951 
1952    Options Database Keys:
1953 .   -mat_no_unroll - uses code that does not unroll the loops in the
1954                      block calculations (much slower)
1955 .    -mat_block_size - size of the blocks to use (only works if a negative bs is passed in
1956 
1957    Level: intermediate
1958 
1959    Notes:
1960    Specify the preallocated storage with either nz or nnz (not both).
1961    Set nz=PETSC_DEFAULT and nnz=NULL for PETSc to control dynamic memory
1962    allocation.  See the <a href="../../docs/manual.pdf#nameddest=ch_mat">Mat chapter of the users manual</a> for details.
1963 
1964    You can call MatGetInfo() to get information on how effective the preallocation was;
1965    for example the fields mallocs,nz_allocated,nz_used,nz_unneeded;
1966    You can also run with the option -info and look for messages with the string
1967    malloc in them to see if additional memory allocation was needed.
1968 
1969    If the nnz parameter is given then the nz parameter is ignored
1970 
1971 
1972 .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateSBAIJ()
1973 @*/
1974 PetscErrorCode  MatSeqSBAIJSetPreallocation(Mat B,PetscInt bs,PetscInt nz,const PetscInt nnz[])
1975 {
1976   PetscErrorCode ierr;
1977 
1978   PetscFunctionBegin;
1979   PetscValidHeaderSpecific(B,MAT_CLASSID,1);
1980   PetscValidType(B,1);
1981   PetscValidLogicalCollectiveInt(B,bs,2);
1982   ierr = PetscTryMethod(B,"MatSeqSBAIJSetPreallocation_C",(Mat,PetscInt,PetscInt,const PetscInt[]),(B,bs,nz,nnz));CHKERRQ(ierr);
1983   PetscFunctionReturn(0);
1984 }
1985 
1986 #undef __FUNCT__
1987 #define __FUNCT__ "MatCreateSeqSBAIJ"
1988 /*@C
1989    MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block
1990    compressed row) format.  For good matrix assembly performance the
1991    user should preallocate the matrix storage by setting the parameter nz
1992    (or the array nnz).  By setting these parameters accurately, performance
1993    during matrix assembly can be increased by more than a factor of 50.
1994 
1995    Collective on MPI_Comm
1996 
1997    Input Parameters:
1998 +  comm - MPI communicator, set to PETSC_COMM_SELF
1999 .  bs - size of block
2000 .  m - number of rows, or number of columns
2001 .  nz - number of block nonzeros per block row (same for all rows)
2002 -  nnz - array containing the number of block nonzeros in the upper triangular plus
2003          diagonal portion of each block (possibly different for each block row) or NULL
2004 
2005    Output Parameter:
2006 .  A - the symmetric matrix
2007 
2008    Options Database Keys:
2009 .   -mat_no_unroll - uses code that does not unroll the loops in the
2010                      block calculations (much slower)
2011 .    -mat_block_size - size of the blocks to use
2012 
2013    Level: intermediate
2014 
2015    It is recommended that one use the MatCreate(), MatSetType() and/or MatSetFromOptions(),
2016    MatXXXXSetPreallocation() paradgm instead of this routine directly.
2017    [MatXXXXSetPreallocation() is, for example, MatSeqAIJSetPreallocation]
2018 
2019    Notes:
2020    The number of rows and columns must be divisible by blocksize.
2021    This matrix type does not support complex Hermitian operation.
2022 
2023    Specify the preallocated storage with either nz or nnz (not both).
2024    Set nz=PETSC_DEFAULT and nnz=NULL for PETSc to control dynamic memory
2025    allocation.  See the <a href="../../docs/manual.pdf#nameddest=ch_mat">Mat chapter of the users manual</a> for details.
2026 
2027    If the nnz parameter is given then the nz parameter is ignored
2028 
2029 .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateSBAIJ()
2030 @*/
2031 PetscErrorCode  MatCreateSeqSBAIJ(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A)
2032 {
2033   PetscErrorCode ierr;
2034 
2035   PetscFunctionBegin;
2036   ierr = MatCreate(comm,A);CHKERRQ(ierr);
2037   ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr);
2038   ierr = MatSetType(*A,MATSEQSBAIJ);CHKERRQ(ierr);
2039   ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(*A,bs,nz,(PetscInt*)nnz);CHKERRQ(ierr);
2040   PetscFunctionReturn(0);
2041 }
2042 
2043 #undef __FUNCT__
2044 #define __FUNCT__ "MatDuplicate_SeqSBAIJ"
2045 PetscErrorCode MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
2046 {
2047   Mat            C;
2048   Mat_SeqSBAIJ   *c,*a = (Mat_SeqSBAIJ*)A->data;
2049   PetscErrorCode ierr;
2050   PetscInt       i,mbs = a->mbs,nz = a->nz,bs2 =a->bs2;
2051 
2052   PetscFunctionBegin;
2053   if (a->i[mbs] != nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Corrupt matrix");
2054 
2055   *B   = 0;
2056   ierr = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr);
2057   ierr = MatSetSizes(C,A->rmap->N,A->cmap->n,A->rmap->N,A->cmap->n);CHKERRQ(ierr);
2058   ierr = MatSetType(C,MATSEQSBAIJ);CHKERRQ(ierr);
2059   ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
2060   c    = (Mat_SeqSBAIJ*)C->data;
2061 
2062   C->preallocated       = PETSC_TRUE;
2063   C->factortype         = A->factortype;
2064   c->row                = 0;
2065   c->icol               = 0;
2066   c->saved_values       = 0;
2067   c->keepnonzeropattern = a->keepnonzeropattern;
2068   C->assembled          = PETSC_TRUE;
2069 
2070   ierr   = PetscLayoutReference(A->rmap,&C->rmap);CHKERRQ(ierr);
2071   ierr   = PetscLayoutReference(A->cmap,&C->cmap);CHKERRQ(ierr);
2072   c->bs2 = a->bs2;
2073   c->mbs = a->mbs;
2074   c->nbs = a->nbs;
2075 
2076   if  (cpvalues == MAT_SHARE_NONZERO_PATTERN) {
2077     c->imax           = a->imax;
2078     c->ilen           = a->ilen;
2079     c->free_imax_ilen = PETSC_FALSE;
2080   } else {
2081     ierr = PetscMalloc2((mbs+1),PetscInt,&c->imax,(mbs+1),PetscInt,&c->ilen);CHKERRQ(ierr);
2082     ierr = PetscLogObjectMemory(C,2*(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr);
2083     for (i=0; i<mbs; i++) {
2084       c->imax[i] = a->imax[i];
2085       c->ilen[i] = a->ilen[i];
2086     }
2087     c->free_imax_ilen = PETSC_TRUE;
2088   }
2089 
2090   /* allocate the matrix space */
2091   if (cpvalues == MAT_SHARE_NONZERO_PATTERN) {
2092     ierr            = PetscMalloc(bs2*nz*sizeof(MatScalar),&c->a);CHKERRQ(ierr);
2093     ierr            = PetscLogObjectMemory(C,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr);
2094     c->i            = a->i;
2095     c->j            = a->j;
2096     c->singlemalloc = PETSC_FALSE;
2097     c->free_a       = PETSC_TRUE;
2098     c->free_ij      = PETSC_FALSE;
2099     c->parent       = A;
2100     ierr            = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
2101     ierr            = MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
2102     ierr            = MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
2103   } else {
2104     ierr            = PetscMalloc3(bs2*nz,MatScalar,&c->a,nz,PetscInt,&c->j,mbs+1,PetscInt,&c->i);CHKERRQ(ierr);
2105     ierr            = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr);
2106     ierr            = PetscLogObjectMemory(C,(mbs+1)*sizeof(PetscInt) + nz*(bs2*sizeof(MatScalar) + sizeof(PetscInt)));CHKERRQ(ierr);
2107     c->singlemalloc = PETSC_TRUE;
2108     c->free_a       = PETSC_TRUE;
2109     c->free_ij      = PETSC_TRUE;
2110   }
2111   if (mbs > 0) {
2112     if (cpvalues != MAT_SHARE_NONZERO_PATTERN) {
2113       ierr = PetscMemcpy(c->j,a->j,nz*sizeof(PetscInt));CHKERRQ(ierr);
2114     }
2115     if (cpvalues == MAT_COPY_VALUES) {
2116       ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
2117     } else {
2118       ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
2119     }
2120     if (a->jshort) {
2121       /* cannot share jshort, it is reallocated in MatAssemblyEnd_SeqSBAIJ() */
2122       /* if the parent matrix is reassembled, this child matrix will never notice */
2123       ierr = PetscMalloc(nz*sizeof(unsigned short),&c->jshort);CHKERRQ(ierr);
2124       ierr = PetscLogObjectMemory(C,nz*sizeof(unsigned short));CHKERRQ(ierr);
2125       ierr = PetscMemcpy(c->jshort,a->jshort,nz*sizeof(unsigned short));CHKERRQ(ierr);
2126 
2127       c->free_jshort = PETSC_TRUE;
2128     }
2129   }
2130 
2131   c->roworiented = a->roworiented;
2132   c->nonew       = a->nonew;
2133 
2134   if (a->diag) {
2135     if (cpvalues == MAT_SHARE_NONZERO_PATTERN) {
2136       c->diag      = a->diag;
2137       c->free_diag = PETSC_FALSE;
2138     } else {
2139       ierr = PetscMalloc(mbs*sizeof(PetscInt),&c->diag);CHKERRQ(ierr);
2140       ierr = PetscLogObjectMemory(C,mbs*sizeof(PetscInt));CHKERRQ(ierr);
2141       for (i=0; i<mbs; i++) c->diag[i] = a->diag[i];
2142       c->free_diag = PETSC_TRUE;
2143     }
2144   }
2145   c->nz         = a->nz;
2146   c->maxnz      = a->nz; /* Since we allocate exactly the right amount */
2147   c->solve_work = 0;
2148   c->mult_work  = 0;
2149 
2150   *B   = C;
2151   ierr = PetscFunctionListDuplicate(((PetscObject)A)->qlist,&((PetscObject)C)->qlist);CHKERRQ(ierr);
2152   PetscFunctionReturn(0);
2153 }
2154 
2155 #undef __FUNCT__
2156 #define __FUNCT__ "MatLoad_SeqSBAIJ"
2157 PetscErrorCode MatLoad_SeqSBAIJ(Mat newmat,PetscViewer viewer)
2158 {
2159   Mat_SeqSBAIJ   *a;
2160   PetscErrorCode ierr;
2161   int            fd;
2162   PetscMPIInt    size;
2163   PetscInt       i,nz,header[4],*rowlengths=0,M,N,bs=1;
2164   PetscInt       *mask,mbs,*jj,j,rowcount,nzcount,k,*s_browlengths,maskcount;
2165   PetscInt       kmax,jcount,block,idx,point,nzcountb,extra_rows,rows,cols;
2166   PetscInt       *masked,nmask,tmp,bs2,ishift;
2167   PetscScalar    *aa;
2168   MPI_Comm       comm;
2169 
2170   PetscFunctionBegin;
2171   ierr = PetscObjectGetComm((PetscObject)viewer,&comm);CHKERRQ(ierr);
2172   ierr = PetscOptionsGetInt(((PetscObject)newmat)->prefix,"-matload_block_size",&bs,NULL);CHKERRQ(ierr);
2173   bs2  = bs*bs;
2174 
2175   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
2176   if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"view must have one processor");
2177   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
2178   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
2179   if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"not Mat object");
2180   M = header[1]; N = header[2]; nz = header[3];
2181 
2182   if (header[3] < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ");
2183 
2184   if (M != N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Can only do square matrices");
2185 
2186   /*
2187      This code adds extra rows to make sure the number of rows is
2188     divisible by the blocksize
2189   */
2190   mbs        = M/bs;
2191   extra_rows = bs - M + bs*(mbs);
2192   if (extra_rows == bs) extra_rows = 0;
2193   else                  mbs++;
2194   if (extra_rows) {
2195     ierr = PetscInfo(viewer,"Padding loaded matrix to match blocksize\n");CHKERRQ(ierr);
2196   }
2197 
2198   /* Set global sizes if not already set */
2199   if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) {
2200     ierr = MatSetSizes(newmat,PETSC_DECIDE,PETSC_DECIDE,M+extra_rows,N+extra_rows);CHKERRQ(ierr);
2201   } else { /* Check if the matrix global sizes are correct */
2202     ierr = MatGetSize(newmat,&rows,&cols);CHKERRQ(ierr);
2203     if (M != rows ||  N != cols) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Matrix in file of different length (%d, %d) than the input matrix (%d, %d)",M,N,rows,cols);
2204   }
2205 
2206   /* read in row lengths */
2207   ierr = PetscMalloc((M+extra_rows)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
2208   ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
2209   for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;
2210 
2211   /* read in column indices */
2212   ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscInt),&jj);CHKERRQ(ierr);
2213   ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr);
2214   for (i=0; i<extra_rows; i++) jj[nz+i] = M+i;
2215 
2216   /* loop over row lengths determining block row lengths */
2217   ierr     = PetscMalloc(mbs*sizeof(PetscInt),&s_browlengths);CHKERRQ(ierr);
2218   ierr     = PetscMemzero(s_browlengths,mbs*sizeof(PetscInt));CHKERRQ(ierr);
2219   ierr     = PetscMalloc2(mbs,PetscInt,&mask,mbs,PetscInt,&masked);CHKERRQ(ierr);
2220   ierr     = PetscMemzero(mask,mbs*sizeof(PetscInt));CHKERRQ(ierr);
2221   rowcount = 0;
2222   nzcount  = 0;
2223   for (i=0; i<mbs; i++) {
2224     nmask = 0;
2225     for (j=0; j<bs; j++) {
2226       kmax = rowlengths[rowcount];
2227       for (k=0; k<kmax; k++) {
2228         tmp = jj[nzcount++]/bs;   /* block col. index */
2229         if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;}
2230       }
2231       rowcount++;
2232     }
2233     s_browlengths[i] += nmask;
2234 
2235     /* zero out the mask elements we set */
2236     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
2237   }
2238 
2239   /* Do preallocation */
2240   ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(newmat,bs,0,s_browlengths);CHKERRQ(ierr);
2241   a    = (Mat_SeqSBAIJ*)newmat->data;
2242 
2243   /* set matrix "i" values */
2244   a->i[0] = 0;
2245   for (i=1; i<= mbs; i++) {
2246     a->i[i]      = a->i[i-1] + s_browlengths[i-1];
2247     a->ilen[i-1] = s_browlengths[i-1];
2248   }
2249   a->nz = a->i[mbs];
2250 
2251   /* read in nonzero values */
2252   ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);CHKERRQ(ierr);
2253   ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr);
2254   for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0;
2255 
2256   /* set "a" and "j" values into matrix */
2257   nzcount = 0; jcount = 0;
2258   for (i=0; i<mbs; i++) {
2259     nzcountb = nzcount;
2260     nmask    = 0;
2261     for (j=0; j<bs; j++) {
2262       kmax = rowlengths[i*bs+j];
2263       for (k=0; k<kmax; k++) {
2264         tmp = jj[nzcount++]/bs; /* block col. index */
2265         if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;}
2266       }
2267     }
2268     /* sort the masked values */
2269     ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr);
2270 
2271     /* set "j" values into matrix */
2272     maskcount = 1;
2273     for (j=0; j<nmask; j++) {
2274       a->j[jcount++]  = masked[j];
2275       mask[masked[j]] = maskcount++;
2276     }
2277 
2278     /* set "a" values into matrix */
2279     ishift = bs2*a->i[i];
2280     for (j=0; j<bs; j++) {
2281       kmax = rowlengths[i*bs+j];
2282       for (k=0; k<kmax; k++) {
2283         tmp = jj[nzcountb]/bs;        /* block col. index */
2284         if (tmp >= i) {
2285           block     = mask[tmp] - 1;
2286           point     = jj[nzcountb] - bs*tmp;
2287           idx       = ishift + bs2*block + j + bs*point;
2288           a->a[idx] = aa[nzcountb];
2289         }
2290         nzcountb++;
2291       }
2292     }
2293     /* zero out the mask elements we set */
2294     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
2295   }
2296   if (jcount != a->nz) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix");
2297 
2298   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
2299   ierr = PetscFree(s_browlengths);CHKERRQ(ierr);
2300   ierr = PetscFree(aa);CHKERRQ(ierr);
2301   ierr = PetscFree(jj);CHKERRQ(ierr);
2302   ierr = PetscFree2(mask,masked);CHKERRQ(ierr);
2303 
2304   ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2305   ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2306   PetscFunctionReturn(0);
2307 }
2308 
2309 #undef __FUNCT__
2310 #define __FUNCT__ "MatCreateSeqSBAIJWithArrays"
2311 /*@
2312      MatCreateSeqSBAIJWithArrays - Creates an sequential SBAIJ matrix using matrix elements
2313               (upper triangular entries in CSR format) provided by the user.
2314 
2315      Collective on MPI_Comm
2316 
2317    Input Parameters:
2318 +  comm - must be an MPI communicator of size 1
2319 .  bs - size of block
2320 .  m - number of rows
2321 .  n - number of columns
2322 .  i - row indices
2323 .  j - column indices
2324 -  a - matrix values
2325 
2326    Output Parameter:
2327 .  mat - the matrix
2328 
2329    Level: advanced
2330 
2331    Notes:
2332        The i, j, and a arrays are not copied by this routine, the user must free these arrays
2333     once the matrix is destroyed
2334 
2335        You cannot set new nonzero locations into this matrix, that will generate an error.
2336 
2337        The i and j indices are 0 based
2338 
2339        When block size is greater than 1 the matrix values must be stored using the SBAIJ storage format (see the SBAIJ code to determine this). For block size of 1
2340        it is the regular CSR format excluding the lower triangular elements.
2341 
2342 .seealso: MatCreate(), MatCreateSBAIJ(), MatCreateSeqSBAIJ()
2343 
2344 @*/
2345 PetscErrorCode  MatCreateSeqSBAIJWithArrays(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt *i,PetscInt *j,PetscScalar *a,Mat *mat)
2346 {
2347   PetscErrorCode ierr;
2348   PetscInt       ii;
2349   Mat_SeqSBAIJ   *sbaij;
2350 
2351   PetscFunctionBegin;
2352   if (bs != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"block size %D > 1 is not supported yet",bs);
2353   if (i[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"i (row indices) must start with 0");
2354 
2355   ierr  = MatCreate(comm,mat);CHKERRQ(ierr);
2356   ierr  = MatSetSizes(*mat,m,n,m,n);CHKERRQ(ierr);
2357   ierr  = MatSetType(*mat,MATSEQSBAIJ);CHKERRQ(ierr);
2358   ierr  = MatSeqSBAIJSetPreallocation_SeqSBAIJ(*mat,bs,MAT_SKIP_ALLOCATION,0);CHKERRQ(ierr);
2359   sbaij = (Mat_SeqSBAIJ*)(*mat)->data;
2360   ierr  = PetscMalloc2(m,PetscInt,&sbaij->imax,m,PetscInt,&sbaij->ilen);CHKERRQ(ierr);
2361   ierr  = PetscLogObjectMemory(*mat,2*m*sizeof(PetscInt));CHKERRQ(ierr);
2362 
2363   sbaij->i = i;
2364   sbaij->j = j;
2365   sbaij->a = a;
2366 
2367   sbaij->singlemalloc = PETSC_FALSE;
2368   sbaij->nonew        = -1;             /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/
2369   sbaij->free_a       = PETSC_FALSE;
2370   sbaij->free_ij      = PETSC_FALSE;
2371 
2372   for (ii=0; ii<m; ii++) {
2373     sbaij->ilen[ii] = sbaij->imax[ii] = i[ii+1] - i[ii];
2374 #if defined(PETSC_USE_DEBUG)
2375     if (i[ii+1] - i[ii] < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row length in i (row indices) row = %d length = %d",ii,i[ii+1] - i[ii]);
2376 #endif
2377   }
2378 #if defined(PETSC_USE_DEBUG)
2379   for (ii=0; ii<sbaij->i[m]; ii++) {
2380     if (j[ii] < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative column index at location = %d index = %d",ii,j[ii]);
2381     if (j[ii] > n - 1) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column index to large at location = %d index = %d",ii,j[ii]);
2382   }
2383 #endif
2384 
2385   ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2386   ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2387   PetscFunctionReturn(0);
2388 }
2389 
2390 
2391 
2392 
2393 
2394