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