xref: /petsc/src/mat/impls/baij/seq/baijfact.c (revision 6024bd2c273c6b4ac03c2bd49d1968cdcd0f9df3)
1 #define PETSCMAT_DLL
2 
3 /*
4     Factorization code for BAIJ format.
5 */
6 #include "src/mat/impls/baij/seq/baij.h"
7 #include "src/inline/ilu.h"
8 
9 /* ------------------------------------------------------------*/
10 /*
11       Version for when blocks are 2 by 2
12 */
13 #undef __FUNCT__
14 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2"
15 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2(Mat A,MatFactorInfo *info,Mat *B)
16 {
17   Mat            C = *B;
18   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data;
19   IS             isrow = b->row,isicol = b->icol;
20   PetscErrorCode ierr;
21   PetscInt       *r,*ic,i,j,n = a->mbs,*bi = b->i,*bj = b->j;
22   PetscInt       *ajtmpold,*ajtmp,nz,row;
23   PetscInt       *diag_offset=b->diag,idx,*ai=a->i,*aj=a->j,*pj;
24   MatScalar      *pv,*v,*rtmp,m1,m2,m3,m4,*pc,*w,*x,x1,x2,x3,x4;
25   MatScalar      p1,p2,p3,p4;
26   MatScalar      *ba = b->a,*aa = a->a;
27 
28   PetscFunctionBegin;
29   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
30   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
31   ierr  = PetscMalloc(4*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
32 
33   for (i=0; i<n; i++) {
34     nz    = bi[i+1] - bi[i];
35     ajtmp = bj + bi[i];
36     for  (j=0; j<nz; j++) {
37       x = rtmp+4*ajtmp[j]; x[0] = x[1] = x[2] = x[3] = 0.0;
38     }
39     /* load in initial (unfactored row) */
40     idx      = r[i];
41     nz       = ai[idx+1] - ai[idx];
42     ajtmpold = aj + ai[idx];
43     v        = aa + 4*ai[idx];
44     for (j=0; j<nz; j++) {
45       x    = rtmp+4*ic[ajtmpold[j]];
46       x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3];
47       v    += 4;
48     }
49     row = *ajtmp++;
50     while (row < i) {
51       pc = rtmp + 4*row;
52       p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3];
53       if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0) {
54         pv = ba + 4*diag_offset[row];
55         pj = bj + diag_offset[row] + 1;
56         x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
57         pc[0] = m1 = p1*x1 + p3*x2;
58         pc[1] = m2 = p2*x1 + p4*x2;
59         pc[2] = m3 = p1*x3 + p3*x4;
60         pc[3] = m4 = p2*x3 + p4*x4;
61         nz = bi[row+1] - diag_offset[row] - 1;
62         pv += 4;
63         for (j=0; j<nz; j++) {
64           x1   = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
65           x    = rtmp + 4*pj[j];
66           x[0] -= m1*x1 + m3*x2;
67           x[1] -= m2*x1 + m4*x2;
68           x[2] -= m1*x3 + m3*x4;
69           x[3] -= m2*x3 + m4*x4;
70           pv   += 4;
71         }
72         ierr = PetscLogFlops(16*nz+12);CHKERRQ(ierr);
73       }
74       row = *ajtmp++;
75     }
76     /* finished row so stick it into b->a */
77     pv = ba + 4*bi[i];
78     pj = bj + bi[i];
79     nz = bi[i+1] - bi[i];
80     for (j=0; j<nz; j++) {
81       x     = rtmp+4*pj[j];
82       pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3];
83       pv   += 4;
84     }
85     /* invert diagonal block */
86     w = ba + 4*diag_offset[i];
87     ierr = Kernel_A_gets_inverse_A_2(w);CHKERRQ(ierr);
88   }
89 
90   ierr = PetscFree(rtmp);CHKERRQ(ierr);
91   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
92   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
93   C->factor = FACTOR_LU;
94   C->assembled = PETSC_TRUE;
95   ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
96   PetscFunctionReturn(0);
97 }
98 /*
99       Version for when blocks are 2 by 2 Using natural ordering
100 */
101 #undef __FUNCT__
102 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering"
103 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering(Mat A,MatFactorInfo *info,Mat *B)
104 {
105   Mat            C = *B;
106   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data;
107   PetscErrorCode ierr;
108   PetscInt       i,j,n = a->mbs,*bi = b->i,*bj = b->j;
109   PetscInt       *ajtmpold,*ajtmp,nz,row;
110   PetscInt       *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj;
111   MatScalar      *pv,*v,*rtmp,*pc,*w,*x;
112   MatScalar      p1,p2,p3,p4,m1,m2,m3,m4,x1,x2,x3,x4;
113   MatScalar      *ba = b->a,*aa = a->a;
114 
115   PetscFunctionBegin;
116   ierr = PetscMalloc(4*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
117 
118   for (i=0; i<n; i++) {
119     nz    = bi[i+1] - bi[i];
120     ajtmp = bj + bi[i];
121     for  (j=0; j<nz; j++) {
122       x = rtmp+4*ajtmp[j];
123       x[0]  = x[1]  = x[2]  = x[3]  = 0.0;
124     }
125     /* load in initial (unfactored row) */
126     nz       = ai[i+1] - ai[i];
127     ajtmpold = aj + ai[i];
128     v        = aa + 4*ai[i];
129     for (j=0; j<nz; j++) {
130       x    = rtmp+4*ajtmpold[j];
131       x[0]  = v[0];  x[1]  = v[1];  x[2]  = v[2];  x[3]  = v[3];
132       v    += 4;
133     }
134     row = *ajtmp++;
135     while (row < i) {
136       pc  = rtmp + 4*row;
137       p1  = pc[0];  p2  = pc[1];  p3  = pc[2];  p4  = pc[3];
138       if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0) {
139         pv = ba + 4*diag_offset[row];
140         pj = bj + diag_offset[row] + 1;
141         x1  = pv[0];  x2  = pv[1];  x3  = pv[2];  x4  = pv[3];
142         pc[0] = m1 = p1*x1 + p3*x2;
143         pc[1] = m2 = p2*x1 + p4*x2;
144         pc[2] = m3 = p1*x3 + p3*x4;
145         pc[3] = m4 = p2*x3 + p4*x4;
146         nz = bi[row+1] - diag_offset[row] - 1;
147         pv += 4;
148         for (j=0; j<nz; j++) {
149           x1   = pv[0];  x2  = pv[1];   x3 = pv[2];  x4  = pv[3];
150           x    = rtmp + 4*pj[j];
151           x[0] -= m1*x1 + m3*x2;
152           x[1] -= m2*x1 + m4*x2;
153           x[2] -= m1*x3 + m3*x4;
154           x[3] -= m2*x3 + m4*x4;
155           pv   += 4;
156         }
157         ierr = PetscLogFlops(16*nz+12);CHKERRQ(ierr);
158       }
159       row = *ajtmp++;
160     }
161     /* finished row so stick it into b->a */
162     pv = ba + 4*bi[i];
163     pj = bj + bi[i];
164     nz = bi[i+1] - bi[i];
165     for (j=0; j<nz; j++) {
166       x      = rtmp+4*pj[j];
167       pv[0]  = x[0];  pv[1]  = x[1];  pv[2]  = x[2];  pv[3]  = x[3];
168       pv   += 4;
169     }
170     /* invert diagonal block */
171     w = ba + 4*diag_offset[i];
172     ierr = Kernel_A_gets_inverse_A_2(w);CHKERRQ(ierr);
173     /*Kernel_A_gets_inverse_A(bs,w,v_pivots,v_work);*/
174   }
175 
176   ierr = PetscFree(rtmp);CHKERRQ(ierr);
177   C->factor    = FACTOR_LU;
178   C->assembled = PETSC_TRUE;
179   ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
180   PetscFunctionReturn(0);
181 }
182 
183 /* ----------------------------------------------------------- */
184 /*
185      Version for when blocks are 1 by 1.
186 */
187 #undef __FUNCT__
188 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_1"
189 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_1(Mat A,MatFactorInfo *info,Mat *B)
190 {
191   Mat            C = *B;
192   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data;
193   IS             isrow = b->row,isicol = b->icol;
194   PetscErrorCode ierr;
195   PetscInt       *r,*ic,i,j,n = a->mbs,*bi = b->i,*bj = b->j;
196   PetscInt       *ajtmpold,*ajtmp,nz,row,*ai = a->i,*aj = a->j;
197   PetscInt       *diag_offset = b->diag,diag,*pj;
198   MatScalar      *pv,*v,*rtmp,multiplier,*pc;
199   MatScalar      *ba = b->a,*aa = a->a;
200 
201   PetscFunctionBegin;
202   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
203   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
204   ierr  = PetscMalloc((n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
205 
206   for (i=0; i<n; i++) {
207     nz    = bi[i+1] - bi[i];
208     ajtmp = bj + bi[i];
209     for  (j=0; j<nz; j++) rtmp[ajtmp[j]] = 0.0;
210 
211     /* load in initial (unfactored row) */
212     nz       = ai[r[i]+1] - ai[r[i]];
213     ajtmpold = aj + ai[r[i]];
214     v        = aa + ai[r[i]];
215     for (j=0; j<nz; j++) rtmp[ic[ajtmpold[j]]] =  v[j];
216 
217     row = *ajtmp++;
218     while (row < i) {
219       pc = rtmp + row;
220       if (*pc != 0.0) {
221         pv         = ba + diag_offset[row];
222         pj         = bj + diag_offset[row] + 1;
223         multiplier = *pc * *pv++;
224         *pc        = multiplier;
225         nz         = bi[row+1] - diag_offset[row] - 1;
226         for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j];
227         ierr = PetscLogFlops(1+2*nz);CHKERRQ(ierr);
228       }
229       row = *ajtmp++;
230     }
231     /* finished row so stick it into b->a */
232     pv = ba + bi[i];
233     pj = bj + bi[i];
234     nz = bi[i+1] - bi[i];
235     for (j=0; j<nz; j++) {pv[j] = rtmp[pj[j]];}
236     diag = diag_offset[i] - bi[i];
237     /* check pivot entry for current row */
238     if (pv[diag] == 0.0) {
239       SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot");
240     }
241     pv[diag] = 1.0/pv[diag];
242   }
243 
244   ierr = PetscFree(rtmp);CHKERRQ(ierr);
245   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
246   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
247   C->factor    = FACTOR_LU;
248   C->assembled = PETSC_TRUE;
249   ierr = PetscLogFlops(C->n);CHKERRQ(ierr);
250   PetscFunctionReturn(0);
251 }
252 
253 
254 /* ----------------------------------------------------------- */
255 #undef __FUNCT__
256 #define __FUNCT__ "MatLUFactor_SeqBAIJ"
257 PetscErrorCode MatLUFactor_SeqBAIJ(Mat A,IS row,IS col,MatFactorInfo *info)
258 {
259   PetscErrorCode ierr;
260   Mat            C;
261 
262   PetscFunctionBegin;
263   ierr = MatLUFactorSymbolic(A,row,col,info,&C);CHKERRQ(ierr);
264   ierr = MatLUFactorNumeric(A,info,&C);CHKERRQ(ierr);
265   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
266   ierr = PetscLogObjectParent(A,((Mat_SeqBAIJ*)(A->data))->icol);CHKERRQ(ierr);
267   PetscFunctionReturn(0);
268 }
269 
270 #include "src/mat/impls/sbaij/seq/sbaij.h"
271 #undef __FUNCT__
272 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N"
273 PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N(Mat A,MatFactorInfo *info,Mat *B)
274 {
275   PetscErrorCode ierr;
276   Mat            C = *B;
277   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data;
278   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
279   IS             ip=b->row;
280   PetscInt       *rip,i,j,mbs=a->mbs,bs=A->bs,*bi=b->i,*bj=b->j,*bcol;
281   PetscInt       *ai=a->i,*aj=a->j;
282   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
283   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
284   PetscReal      zeropivot,rs,shiftnz;
285   PetscTruth     shiftpd;
286   ChShift_Ctx    sctx;
287   PetscInt       newshift;
288 
289   PetscFunctionBegin;
290   if (bs > 1) {
291     if (!a->sbaijMat){
292       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
293     }
294     ierr = (a->sbaijMat)->ops->choleskyfactornumeric(a->sbaijMat,info,B);CHKERRQ(ierr);
295     ierr = MatDestroy(a->sbaijMat);CHKERRQ(ierr);
296     a->sbaijMat = PETSC_NULL;
297     PetscFunctionReturn(0);
298   }
299 
300   /* initialization */
301   shiftnz   = info->shiftnz;
302   shiftpd   = info->shiftpd;
303   zeropivot = info->zeropivot;
304 
305   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
306   nz   = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar);
307   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
308   jl   = il + mbs;
309   rtmp = (MatScalar*)(jl + mbs);
310 
311   sctx.shift_amount = 0;
312   sctx.nshift       = 0;
313   do {
314     sctx.chshift = PETSC_FALSE;
315     for (i=0; i<mbs; i++) {
316       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
317     }
318 
319     for (k = 0; k<mbs; k++){
320       bval = ba + bi[k];
321       /* initialize k-th row by the perm[k]-th row of A */
322       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
323       for (j = jmin; j < jmax; j++){
324         col = rip[aj[j]];
325         if (col >= k){ /* only take upper triangular entry */
326           rtmp[col] = aa[j];
327           *bval++  = 0.0; /* for in-place factorization */
328         }
329       }
330 
331       /* shift the diagonal of the matrix */
332       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
333 
334       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
335       dk = rtmp[k];
336       i = jl[k]; /* first row to be added to k_th row  */
337 
338       while (i < k){
339         nexti = jl[i]; /* next row to be added to k_th row */
340 
341         /* compute multiplier, update diag(k) and U(i,k) */
342         ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
343         uikdi = - ba[ili]*ba[bi[i]];  /* diagonal(k) */
344         dk += uikdi*ba[ili];
345         ba[ili] = uikdi; /* -U(i,k) */
346 
347         /* add multiple of row i to k-th row */
348         jmin = ili + 1; jmax = bi[i+1];
349         if (jmin < jmax){
350           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
351           /* update il and jl for row i */
352           il[i] = jmin;
353           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
354         }
355         i = nexti;
356       }
357 
358       /* shift the diagonals when zero pivot is detected */
359       /* compute rs=sum of abs(off-diagonal) */
360       rs   = 0.0;
361       jmin = bi[k]+1;
362       nz   = bi[k+1] - jmin;
363       if (nz){
364         bcol = bj + jmin;
365         while (nz--){
366           rs += PetscAbsScalar(rtmp[*bcol]);
367           bcol++;
368         }
369       }
370 
371       sctx.rs = rs;
372       sctx.pv = dk;
373       ierr = MatCholeskyCheckShift_inline(info,sctx,newshift);CHKERRQ(ierr);
374       if (newshift == 1){
375         break;    /* sctx.shift_amount is updated */
376       } else if (newshift == -1){
377         SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g * rs %g",k,PetscAbsScalar(dk),zeropivot,rs);
378       }
379 
380       /* copy data into U(k,:) */
381       ba[bi[k]] = 1.0/dk; /* U(k,k) */
382       jmin = bi[k]+1; jmax = bi[k+1];
383       if (jmin < jmax) {
384         for (j=jmin; j<jmax; j++){
385           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
386         }
387         /* add the k-th row into il and jl */
388         il[k] = jmin;
389         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
390       }
391     }
392   } while (sctx.chshift);
393   ierr = PetscFree(il);CHKERRQ(ierr);
394 
395   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
396   C->factor       = FACTOR_CHOLESKY;
397   C->assembled    = PETSC_TRUE;
398   C->preallocated = PETSC_TRUE;
399   ierr = PetscLogFlops(C->m);CHKERRQ(ierr);
400   if (sctx.nshift){
401     if (shiftnz) {
402       ierr = PetscLogInfo((0,"MatCholeskyFactorNumeric_SeqBAIJ_1: number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount));CHKERRQ(ierr);
403     } else if (shiftpd) {
404       ierr = PetscLogInfo((0,"MatCholeskyFactorNumeric_SeqBAIJ_1: number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount));CHKERRQ(ierr);
405     }
406   }
407   PetscFunctionReturn(0);
408 }
409 
410 #undef __FUNCT__
411 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering"
412 PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering(Mat A,MatFactorInfo *info,Mat *fact)
413 {
414   Mat            C = *fact;
415   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data;
416   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
417   PetscErrorCode ierr;
418   PetscInt       i,j,am=a->mbs;
419   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
420   PetscInt       k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz;
421   MatScalar      *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval;
422   PetscReal      zeropivot,rs,shiftnz;
423   PetscTruth     shiftpd;
424   ChShift_Ctx    sctx;
425   PetscInt       newshift;
426 
427   PetscFunctionBegin;
428   /* initialization */
429   shiftnz   = info->shiftnz;
430   shiftpd   = info->shiftpd;
431   zeropivot = info->zeropivot;
432 
433   nz   = (2*am+1)*sizeof(PetscInt)+am*sizeof(MatScalar);
434   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
435   jl   = il + am;
436   rtmp = (MatScalar*)(jl + am);
437 
438   sctx.shift_amount = 0;
439   sctx.nshift       = 0;
440   do {
441     sctx.chshift = PETSC_FALSE;
442     for (i=0; i<am; i++) {
443       rtmp[i] = 0.0; jl[i] = am; il[0] = 0;
444     }
445 
446     for (k = 0; k<am; k++){
447     /* initialize k-th row with elements nonzero in row perm(k) of A */
448       nz   = ai[k+1] - ai[k];
449       acol = aj + ai[k];
450       aval = aa + ai[k];
451       bval = ba + bi[k];
452       while (nz -- ){
453         if (*acol < k) { /* skip lower triangular entries */
454           acol++; aval++;
455         } else {
456           rtmp[*acol++] = *aval++;
457           *bval++       = 0.0; /* for in-place factorization */
458         }
459       }
460 
461       /* shift the diagonal of the matrix */
462       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
463 
464       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
465       dk = rtmp[k];
466       i  = jl[k]; /* first row to be added to k_th row  */
467 
468       while (i < k){
469         nexti = jl[i]; /* next row to be added to k_th row */
470         /* compute multiplier, update D(k) and U(i,k) */
471         ili   = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
472         uikdi = - ba[ili]*ba[bi[i]];
473         dk   += uikdi*ba[ili];
474         ba[ili] = uikdi; /* -U(i,k) */
475 
476         /* add multiple of row i to k-th row ... */
477         jmin = ili + 1;
478         nz   = bi[i+1] - jmin;
479         if (nz > 0){
480           bcol = bj + jmin;
481           bval = ba + jmin;
482           while (nz --) rtmp[*bcol++] += uikdi*(*bval++);
483           /* update il and jl for i-th row */
484           il[i] = jmin;
485           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
486         }
487         i = nexti;
488       }
489 
490       /* shift the diagonals when zero pivot is detected */
491       /* compute rs=sum of abs(off-diagonal) */
492       rs   = 0.0;
493       jmin = bi[k]+1;
494       nz   = bi[k+1] - jmin;
495       if (nz){
496         bcol = bj + jmin;
497         while (nz--){
498           rs += PetscAbsScalar(rtmp[*bcol]);
499           bcol++;
500         }
501       }
502 
503       sctx.rs = rs;
504       sctx.pv = dk;
505       ierr = MatCholeskyCheckShift_inline(info,sctx,newshift);CHKERRQ(ierr);
506       if (newshift == 1){
507         break;    /* sctx.shift_amount is updated */
508       } else if (newshift == -1){
509         SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g * rs %g",k,PetscAbsScalar(dk),zeropivot,rs);
510       }
511 
512       /* copy data into U(k,:) */
513       ba[bi[k]] = 1.0/dk;
514       jmin      = bi[k]+1;
515       nz        = bi[k+1] - jmin;
516       if (nz){
517         bcol = bj + jmin;
518         bval = ba + jmin;
519         while (nz--){
520           *bval++       = rtmp[*bcol];
521           rtmp[*bcol++] = 0.0;
522         }
523         /* add k-th row into il and jl */
524         il[k] = jmin;
525         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
526       }
527     }
528   } while (sctx.chshift);
529   ierr = PetscFree(il);CHKERRQ(ierr);
530 
531   C->factor       = FACTOR_CHOLESKY;
532   C->assembled    = PETSC_TRUE;
533   C->preallocated = PETSC_TRUE;
534   ierr = PetscLogFlops(C->m);CHKERRQ(ierr);
535     if (sctx.nshift){
536     if (shiftnz) {
537       ierr = PetscLogInfo((0,"MatCholeskyFactorNumeric_SeqBAIJ_1_NaturalOrdering: number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount));CHKERRQ(ierr);
538     } else if (shiftpd) {
539       ierr = PetscLogInfo((0,"MatCholeskyFactorNumeric_SeqBAIJ_1_NaturalOrdering: number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount));CHKERRQ(ierr);
540     }
541   }
542   PetscFunctionReturn(0);
543 }
544 
545 #include "petscbt.h"
546 #include "src/mat/utils/freespace.h"
547 #undef __FUNCT__
548 #define __FUNCT__ "MatICCFactorSymbolic_SeqBAIJ"
549 PetscErrorCode MatICCFactorSymbolic_SeqBAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
550 {
551   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
552   Mat_SeqSBAIJ   *b;
553   Mat            B;
554   PetscErrorCode ierr;
555   PetscTruth     perm_identity;
556   PetscInt       reallocs=0,*rip,i,*ai=a->i,*aj=a->j,am=a->mbs,bs=A->bs,*ui;
557   PetscInt       jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
558   PetscInt       nlnk,*lnk,*lnk_lvl=PETSC_NULL,ncols,ncols_upper,*cols,*cols_lvl,*uj,**uj_ptr,**uj_lvl_ptr;
559   PetscReal      fill=info->fill,levels=info->levels;
560   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
561   FreeSpaceList  free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
562   PetscBT        lnkbt;
563 
564   PetscFunctionBegin;
565   if (bs > 1){
566     if (!a->sbaijMat){
567       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
568     }
569     ierr = MatICCFactorSymbolic(a->sbaijMat,perm,info,fact);CHKERRQ(ierr);
570     B = *fact;
571     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
572     PetscFunctionReturn(0);
573   }
574 
575   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
576   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
577 
578   /* special case that simply copies fill pattern */
579   if (!levels && perm_identity) {
580     ierr = MatMarkDiagonal_SeqBAIJ(A);CHKERRQ(ierr);
581     ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
582     for (i=0; i<am; i++) {
583       ui[i] = ai[i+1] - a->diag[i]; /* ui: rowlengths - changes when !perm_identity */
584     }
585     ierr = MatCreate(PETSC_COMM_SELF,am,am,am,am,fact);CHKERRQ(ierr);
586     B = *fact;
587     ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
588     ierr = MatSeqSBAIJSetPreallocation(B,1,0,ui);CHKERRQ(ierr);
589 
590     b  = (Mat_SeqSBAIJ*)B->data;
591     uj = b->j;
592     for (i=0; i<am; i++) {
593       aj = a->j + a->diag[i];
594       for (j=0; j<ui[i]; j++){
595         *uj++ = *aj++;
596       }
597       b->ilen[i] = ui[i];
598     }
599     ierr = PetscFree(ui);CHKERRQ(ierr);
600     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
601     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
602 
603     B->ops->solve                 = MatSolve_SeqSBAIJ_1_NaturalOrdering;
604     B->ops->solvetranspose        = MatSolve_SeqSBAIJ_1_NaturalOrdering;
605     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
606     PetscFunctionReturn(0);
607   }
608 
609   /* initialization */
610   ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
611   ui[0] = 0;
612   ierr  = PetscMalloc((2*am+1)*sizeof(PetscInt),&cols_lvl);CHKERRQ(ierr);
613 
614   /* jl: linked list for storing indices of the pivot rows
615      il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
616   ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt*),&jl);CHKERRQ(ierr);
617   il         = jl + am;
618   uj_ptr     = (PetscInt**)(il + am);
619   uj_lvl_ptr = (PetscInt**)(uj_ptr + am);
620   for (i=0; i<am; i++){
621     jl[i] = am; il[i] = 0;
622   }
623 
624   /* create and initialize a linked list for storing column indices of the active row k */
625   nlnk = am + 1;
626   ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
627 
628   /* initial FreeSpace size is fill*(ai[am]+1) */
629   ierr = GetMoreSpace((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
630   current_space = free_space;
631   ierr = GetMoreSpace((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr);
632   current_space_lvl = free_space_lvl;
633 
634   for (k=0; k<am; k++){  /* for each active row k */
635     /* initialize lnk by the column indices of row rip[k] of A */
636     nzk   = 0;
637     ncols = ai[rip[k]+1] - ai[rip[k]];
638     ncols_upper = 0;
639     cols        = cols_lvl + am;
640     for (j=0; j<ncols; j++){
641       i = rip[*(aj + ai[rip[k]] + j)];
642       if (i >= k){ /* only take upper triangular entry */
643         cols[ncols_upper] = i;
644         cols_lvl[ncols_upper] = -1;  /* initialize level for nonzero entries */
645         ncols_upper++;
646       }
647     }
648     ierr = PetscIncompleteLLAdd(ncols_upper,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
649     nzk += nlnk;
650 
651     /* update lnk by computing fill-in for each pivot row to be merged in */
652     prow = jl[k]; /* 1st pivot row */
653 
654     while (prow < k){
655       nextprow = jl[prow];
656 
657       /* merge prow into k-th row */
658       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
659       jmax = ui[prow+1];
660       ncols = jmax-jmin;
661       i     = jmin - ui[prow];
662       cols = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
663       for (j=0; j<ncols; j++) cols_lvl[j] = *(uj_lvl_ptr[prow] + i + j);
664       ierr = PetscIncompleteLLAddSorted(ncols,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
665       nzk += nlnk;
666 
667       /* update il and jl for prow */
668       if (jmin < jmax){
669         il[prow] = jmin;
670         j = *cols; jl[prow] = jl[j]; jl[j] = prow;
671       }
672       prow = nextprow;
673     }
674 
675     /* if free space is not available, make more free space */
676     if (current_space->local_remaining<nzk) {
677       i = am - k + 1; /* num of unfactored rows */
678       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
679       ierr = GetMoreSpace(i,&current_space);CHKERRQ(ierr);
680       ierr = GetMoreSpace(i,&current_space_lvl);CHKERRQ(ierr);
681       reallocs++;
682     }
683 
684     /* copy data into free_space and free_space_lvl, then initialize lnk */
685     ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
686 
687     /* add the k-th row into il and jl */
688     if (nzk-1 > 0){
689       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
690       jl[k] = jl[i]; jl[i] = k;
691       il[k] = ui[k] + 1;
692     }
693     uj_ptr[k]     = current_space->array;
694     uj_lvl_ptr[k] = current_space_lvl->array;
695 
696     current_space->array           += nzk;
697     current_space->local_used      += nzk;
698     current_space->local_remaining -= nzk;
699 
700     current_space_lvl->array           += nzk;
701     current_space_lvl->local_used      += nzk;
702     current_space_lvl->local_remaining -= nzk;
703 
704     ui[k+1] = ui[k] + nzk;
705   }
706 
707 #if defined(PETSC_USE_DEBUG)
708   if (ai[am] != 0) {
709     PetscReal af = ((PetscReal)(2*ui[am]-am))/((PetscReal)ai[am]);
710     ierr = PetscLogInfo((A,"MatICCFactorSymbolic_SeqBAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,fill,af));CHKERRQ(ierr);
711     ierr = PetscLogInfo((A,"MatICCFactorSymbolic_SeqBAIJ:Run with -pc_cholesky_fill %g or use \n",af));CHKERRQ(ierr);
712     ierr = PetscLogInfo((A,"MatICCFactorSymbolic_SeqBAIJ:PCCholeskySetFill(pc,%g) for best performance.\n",af));CHKERRQ(ierr);
713   } else {
714     ierr = PetscLogInfo((A,"MatICCFactorSymbolic_SeqBAIJ:Empty matrix.\n"));CHKERRQ(ierr);
715   }
716 #endif
717 
718   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
719   ierr = PetscFree(jl);CHKERRQ(ierr);
720   ierr = PetscFree(cols_lvl);CHKERRQ(ierr);
721 
722   /* destroy list of free space and other temporary array(s) */
723   ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
724   ierr = MakeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
725   ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
726   ierr = DestroySpace(free_space_lvl);CHKERRQ(ierr);
727 
728   /* put together the new matrix in MATSEQSBAIJ format */
729   ierr = MatCreate(PETSC_COMM_SELF,am,am,am,am,fact);CHKERRQ(ierr);
730   B = *fact;
731   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
732   ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
733 
734   b = (Mat_SeqSBAIJ*)B->data;
735   ierr = PetscFree(b->imax);CHKERRQ(ierr);
736   b->singlemalloc = PETSC_FALSE;
737   /* the next line frees the default space generated by the Create() */
738   ierr = PetscFree(b->a);CHKERRQ(ierr);
739   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
740   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
741   b->j    = uj;
742   b->i    = ui;
743   b->diag = 0;
744   b->ilen = 0;
745   b->imax = 0;
746   b->row  = perm;
747   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
748   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
749   b->icol = perm;
750   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
751   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
752   ierr    = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
753   b->maxnz = b->nz = ui[am];
754 
755   B->factor                 = FACTOR_CHOLESKY;
756   B->info.factor_mallocs    = reallocs;
757   B->info.fill_ratio_given  = fill;
758   if (ai[am] != 0) {
759     B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
760   } else {
761     B->info.fill_ratio_needed = 0.0;
762   }
763   if (perm_identity){
764     B->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
765     B->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
766     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
767   } else {
768     (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
769   }
770   PetscFunctionReturn(0);
771 }
772 
773 #undef __FUNCT__
774 #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqBAIJ"
775 PetscErrorCode MatCholeskyFactorSymbolic_SeqBAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
776 {
777   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
778   Mat_SeqSBAIJ   *b;
779   Mat            B;
780   PetscErrorCode ierr;
781   PetscTruth     perm_identity;
782   PetscReal      fill = info->fill;
783   PetscInt       *rip,*riip,i,mbs=a->mbs,bs=A->bs,*ai=a->i,*aj=a->j,reallocs=0,prow;
784   PetscInt       *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
785   PetscInt       nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr;
786   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
787   PetscBT        lnkbt;
788   IS             iperm;
789 
790   PetscFunctionBegin;
791   if (bs > 1) { /* convert to seqsbaij */
792     if (!a->sbaijMat){
793       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
794     }
795     ierr = MatCholeskyFactorSymbolic(a->sbaijMat,perm,info,fact);CHKERRQ(ierr);
796     B    = *fact;
797     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
798     PetscFunctionReturn(0);
799   }
800 
801   /* check whether perm is the identity mapping */
802   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
803   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
804 
805   if (!perm_identity){
806     /* check if perm is symmetric! */
807     ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr);
808     ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr);
809     for (i=0; i<mbs; i++) {
810       if (rip[i] != riip[i]) SETERRQ(PETSC_ERR_ARG_INCOMP,"Non-symmetric permutation, must use symmetric permutation");
811     }
812     ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr);
813     ierr = ISDestroy(iperm);CHKERRQ(ierr);
814   }
815 
816   /* initialization */
817   ierr  = PetscMalloc((mbs+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
818   ui[0] = 0;
819 
820   /* jl: linked list for storing indices of the pivot rows
821      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
822   ierr = PetscMalloc((3*mbs+1)*sizeof(PetscInt)+mbs*sizeof(PetscInt*),&jl);CHKERRQ(ierr);
823   il     = jl + mbs;
824   cols   = il + mbs;
825   ui_ptr = (PetscInt**)(cols + mbs);
826   for (i=0; i<mbs; i++){
827     jl[i] = mbs; il[i] = 0;
828   }
829 
830   /* create and initialize a linked list for storing column indices of the active row k */
831   nlnk = mbs + 1;
832   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
833 
834   /* initial FreeSpace size is fill*(ai[mbs]+1) */
835   ierr = GetMoreSpace((PetscInt)(fill*(ai[mbs]+1)),&free_space);CHKERRQ(ierr);
836   current_space = free_space;
837 
838   for (k=0; k<mbs; k++){  /* for each active row k */
839     /* initialize lnk by the column indices of row rip[k] of A */
840     nzk   = 0;
841     ncols = ai[rip[k]+1] - ai[rip[k]];
842     ncols_upper = 0;
843     for (j=0; j<ncols; j++){
844       i = rip[*(aj + ai[rip[k]] + j)];
845       if (i >= k){ /* only take upper triangular entry */
846         cols[ncols_upper] = i;
847         ncols_upper++;
848       }
849     }
850     ierr = PetscLLAdd(ncols_upper,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
851     nzk += nlnk;
852 
853     /* update lnk by computing fill-in for each pivot row to be merged in */
854     prow = jl[k]; /* 1st pivot row */
855 
856     while (prow < k){
857       nextprow = jl[prow];
858       /* merge prow into k-th row */
859       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
860       jmax = ui[prow+1];
861       ncols = jmax-jmin;
862       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
863       ierr = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
864       nzk += nlnk;
865 
866       /* update il and jl for prow */
867       if (jmin < jmax){
868         il[prow] = jmin;
869         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
870       }
871       prow = nextprow;
872     }
873 
874     /* if free space is not available, make more free space */
875     if (current_space->local_remaining<nzk) {
876       i = mbs - k + 1; /* num of unfactored rows */
877       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
878       ierr = GetMoreSpace(i,&current_space);CHKERRQ(ierr);
879       reallocs++;
880     }
881 
882     /* copy data into free space, then initialize lnk */
883     ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
884 
885     /* add the k-th row into il and jl */
886     if (nzk-1 > 0){
887       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
888       jl[k] = jl[i]; jl[i] = k;
889       il[k] = ui[k] + 1;
890     }
891     ui_ptr[k] = current_space->array;
892     current_space->array           += nzk;
893     current_space->local_used      += nzk;
894     current_space->local_remaining -= nzk;
895 
896     ui[k+1] = ui[k] + nzk;
897   }
898 
899 #if defined(PETSC_USE_DEBUG)
900   if (ai[mbs] != 0) {
901     PetscReal af = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
902     ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,fill,af));CHKERRQ(ierr);
903     ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Run with -pc_cholesky_fill %g or use \n",af));CHKERRQ(ierr);
904     ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqSBAIJ:PCCholeskySetFill(pc,%g) for best performance.\n",af));CHKERRQ(ierr);
905   } else {
906     ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Empty matrix.\n"));CHKERRQ(ierr);
907   }
908 #endif
909 
910   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
911   ierr = PetscFree(jl);CHKERRQ(ierr);
912 
913   /* destroy list of free space and other temporary array(s) */
914   ierr = PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
915   ierr = MakeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
916   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
917 
918   /* put together the new matrix in MATSEQSBAIJ format */
919   ierr = MatCreate(PETSC_COMM_SELF,mbs,mbs,mbs,mbs,fact);CHKERRQ(ierr);
920   B    = *fact;
921   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
922   ierr = MatSeqSBAIJSetPreallocation(B,bs,0,PETSC_NULL);CHKERRQ(ierr);
923 
924   b = (Mat_SeqSBAIJ*)B->data;
925   ierr = PetscFree(b->imax);CHKERRQ(ierr);
926   b->singlemalloc = PETSC_FALSE;
927   /* the next line frees the default space generated by the Create() */
928   ierr = PetscFree(b->a);CHKERRQ(ierr);
929   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
930   ierr = PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
931   b->j    = uj;
932   b->i    = ui;
933   b->diag = 0;
934   b->ilen = 0;
935   b->imax = 0;
936   b->row  = perm;
937   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
938   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
939   b->icol = perm;
940   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
941   ierr    = PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
942   ierr    = PetscLogObjectMemory(B,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
943   b->maxnz = b->nz = ui[mbs];
944 
945   B->factor                 = FACTOR_CHOLESKY;
946   B->info.factor_mallocs    = reallocs;
947   B->info.fill_ratio_given  = fill;
948   if (ai[mbs] != 0) {
949     B->info.fill_ratio_needed = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
950   } else {
951     B->info.fill_ratio_needed = 0.0;
952   }
953   if (perm_identity){
954     B->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
955     B->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
956     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
957   } else {
958     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
959   }
960   PetscFunctionReturn(0);
961 }
962