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