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