xref: /petsc/src/mat/impls/baij/seq/baijfact.c (revision 2e92ee13a8395f820cc1e3fd74a7607ed52efa2a)
1 
2 /*
3     Factorization code for BAIJ format.
4 */
5 #include <../src/mat/impls/baij/seq/baij.h>
6 #include <petsc/private/kernels/blockinvert.h>
7 
8 #undef __FUNCT__
9 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2"
10 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2(Mat B,Mat A,const MatFactorInfo *info)
11 {
12   Mat            C     =B;
13   Mat_SeqBAIJ    *a    =(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data;
14   IS             isrow = b->row,isicol = b->icol;
15   PetscErrorCode ierr;
16   const PetscInt *r,*ic;
17   PetscInt       i,j,k,nz,nzL,row,*pj;
18   const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,bs2=a->bs2;
19   const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag;
20   MatScalar      *rtmp,*pc,*mwork,*pv;
21   MatScalar      *aa=a->a,*v;
22   PetscInt       flg;
23   PetscReal      shift = info->shiftamount;
24   PetscBool      zeropivotdetected;
25 
26   PetscFunctionBegin;
27   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
28   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
29 
30   /* generate work space needed by the factorization */
31   ierr = PetscMalloc2(bs2*n,&rtmp,bs2,&mwork);CHKERRQ(ierr);
32   ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr);
33 
34   for (i=0; i<n; i++) {
35     /* zero rtmp */
36     /* L part */
37     nz    = bi[i+1] - bi[i];
38     bjtmp = bj + bi[i];
39     for  (j=0; j<nz; j++) {
40       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
41     }
42 
43     /* U part */
44     nz    = bdiag[i] - bdiag[i+1];
45     bjtmp = bj + bdiag[i+1]+1;
46     for  (j=0; j<nz; j++) {
47       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
48     }
49 
50     /* load in initial (unfactored row) */
51     nz    = ai[r[i]+1] - ai[r[i]];
52     ajtmp = aj + ai[r[i]];
53     v     = aa + bs2*ai[r[i]];
54     for (j=0; j<nz; j++) {
55       ierr = PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr);
56     }
57 
58     /* elimination */
59     bjtmp = bj + bi[i];
60     nzL   = bi[i+1] - bi[i];
61     for (k=0; k < nzL; k++) {
62       row = bjtmp[k];
63       pc  = rtmp + bs2*row;
64       for (flg=0,j=0; j<bs2; j++) {
65         if (pc[j] != (PetscScalar)0.0) {
66           flg = 1;
67           break;
68         }
69       }
70       if (flg) {
71         pv = b->a + bs2*bdiag[row];
72         /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
73         ierr = PetscKernel_A_gets_A_times_B_2(pc,pv,mwork);CHKERRQ(ierr);
74 
75         pj = b->j + bdiag[row+1]+1; /* begining of U(row,:) */
76         pv = b->a + bs2*(bdiag[row+1]+1);
77         nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries inU(row,:), excluding diag */
78         for (j=0; j<nz; j++) {
79           /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
80           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
81           v    = rtmp + 4*pj[j];
82           ierr = PetscKernel_A_gets_A_minus_B_times_C_2(v,pc,pv);CHKERRQ(ierr);
83           pv  += 4;
84         }
85         ierr = PetscLogFlops(16*nz+12);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
86       }
87     }
88 
89     /* finished row so stick it into b->a */
90     /* L part */
91     pv = b->a + bs2*bi[i];
92     pj = b->j + bi[i];
93     nz = bi[i+1] - bi[i];
94     for (j=0; j<nz; j++) {
95       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
96     }
97 
98     /* Mark diagonal and invert diagonal for simplier triangular solves */
99     pv   = b->a + bs2*bdiag[i];
100     pj   = b->j + bdiag[i];
101     ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr);
102     /* ierr = PetscKernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */
103     ierr = PetscKernel_A_gets_inverse_A_2(pv,shift,!A->erroriffailure,&zeropivotdetected);CHKERRQ(ierr);
104     if (zeropivotdetected) B->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
105 
106     /* U part */
107     pv = b->a + bs2*(bdiag[i+1]+1);
108     pj = b->j + bdiag[i+1]+1;
109     nz = bdiag[i] - bdiag[i+1] - 1;
110     for (j=0; j<nz; j++) {
111       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
112     }
113   }
114 
115   ierr = PetscFree2(rtmp,mwork);CHKERRQ(ierr);
116   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
117   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
118 
119   C->ops->solve          = MatSolve_SeqBAIJ_2;
120   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2;
121   C->assembled           = PETSC_TRUE;
122 
123   ierr = PetscLogFlops(1.333333333333*2*2*2*n);CHKERRQ(ierr); /* from inverting diagonal blocks */
124   PetscFunctionReturn(0);
125 }
126 
127 #undef __FUNCT__
128 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering"
129 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info)
130 {
131   Mat            C =B;
132   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data;
133   PetscErrorCode ierr;
134   PetscInt       i,j,k,nz,nzL,row,*pj;
135   const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,bs2=a->bs2;
136   const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag;
137   MatScalar      *rtmp,*pc,*mwork,*pv;
138   MatScalar      *aa=a->a,*v;
139   PetscInt       flg;
140   PetscReal      shift = info->shiftamount;
141   PetscBool      zeropivotdetected;
142 
143   PetscFunctionBegin;
144   /* generate work space needed by the factorization */
145   ierr = PetscMalloc2(bs2*n,&rtmp,bs2,&mwork);CHKERRQ(ierr);
146   ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr);
147 
148   for (i=0; i<n; i++) {
149     /* zero rtmp */
150     /* L part */
151     nz    = bi[i+1] - bi[i];
152     bjtmp = bj + bi[i];
153     for  (j=0; j<nz; j++) {
154       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
155     }
156 
157     /* U part */
158     nz    = bdiag[i] - bdiag[i+1];
159     bjtmp = bj + bdiag[i+1]+1;
160     for  (j=0; j<nz; j++) {
161       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
162     }
163 
164     /* load in initial (unfactored row) */
165     nz    = ai[i+1] - ai[i];
166     ajtmp = aj + ai[i];
167     v     = aa + bs2*ai[i];
168     for (j=0; j<nz; j++) {
169       ierr = PetscMemcpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr);
170     }
171 
172     /* elimination */
173     bjtmp = bj + bi[i];
174     nzL   = bi[i+1] - bi[i];
175     for (k=0; k < nzL; k++) {
176       row = bjtmp[k];
177       pc  = rtmp + bs2*row;
178       for (flg=0,j=0; j<bs2; j++) {
179         if (pc[j]!=(PetscScalar)0.0) {
180           flg = 1;
181           break;
182         }
183       }
184       if (flg) {
185         pv = b->a + bs2*bdiag[row];
186         /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
187         ierr = PetscKernel_A_gets_A_times_B_2(pc,pv,mwork);CHKERRQ(ierr);
188 
189         pj = b->j + bdiag[row+1]+1; /* beginning of U(row,:) */
190         pv = b->a + bs2*(bdiag[row+1]+1);
191         nz = bdiag[row]-bdiag[row+1] - 1; /* num of entries in U(row,:) excluding diag */
192         for (j=0; j<nz; j++) {
193           /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
194           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
195           v    = rtmp + 4*pj[j];
196           ierr = PetscKernel_A_gets_A_minus_B_times_C_2(v,pc,pv);CHKERRQ(ierr);
197           pv  += 4;
198         }
199         ierr = PetscLogFlops(16*nz+12);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
200       }
201     }
202 
203     /* finished row so stick it into b->a */
204     /* L part */
205     pv = b->a + bs2*bi[i];
206     pj = b->j + bi[i];
207     nz = bi[i+1] - bi[i];
208     for (j=0; j<nz; j++) {
209       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
210     }
211 
212     /* Mark diagonal and invert diagonal for simplier triangular solves */
213     pv   = b->a + bs2*bdiag[i];
214     pj   = b->j + bdiag[i];
215     ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr);
216     /* ierr = PetscKernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */
217     ierr = PetscKernel_A_gets_inverse_A_2(pv,shift,!A->erroriffailure,&zeropivotdetected);CHKERRQ(ierr);
218     if (zeropivotdetected) B->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
219 
220     /* U part */
221     /*
222     pv = b->a + bs2*bi[2*n-i];
223     pj = b->j + bi[2*n-i];
224     nz = bi[2*n-i+1] - bi[2*n-i] - 1;
225     */
226     pv = b->a + bs2*(bdiag[i+1]+1);
227     pj = b->j + bdiag[i+1]+1;
228     nz = bdiag[i] - bdiag[i+1] - 1;
229     for (j=0; j<nz; j++) {
230       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
231     }
232   }
233   ierr = PetscFree2(rtmp,mwork);CHKERRQ(ierr);
234 
235   C->ops->solve          = MatSolve_SeqBAIJ_2_NaturalOrdering;
236   C->ops->forwardsolve   = MatForwardSolve_SeqBAIJ_2_NaturalOrdering;
237   C->ops->backwardsolve  = MatBackwardSolve_SeqBAIJ_2_NaturalOrdering;
238   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_NaturalOrdering;
239   C->assembled           = PETSC_TRUE;
240 
241   ierr = PetscLogFlops(1.333333333333*2*2*2*n);CHKERRQ(ierr); /* from inverting diagonal blocks */
242   PetscFunctionReturn(0);
243 }
244 
245 #undef __FUNCT__
246 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2_inplace"
247 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_inplace(Mat B,Mat A,const MatFactorInfo *info)
248 {
249   Mat            C     = B;
250   Mat_SeqBAIJ    *a    = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)C->data;
251   IS             isrow = b->row,isicol = b->icol;
252   PetscErrorCode ierr;
253   const PetscInt *r,*ic;
254   PetscInt       i,j,n = a->mbs,*bi = b->i,*bj = b->j;
255   PetscInt       *ajtmpold,*ajtmp,nz,row;
256   PetscInt       *diag_offset=b->diag,idx,*ai=a->i,*aj=a->j,*pj;
257   MatScalar      *pv,*v,*rtmp,m1,m2,m3,m4,*pc,*w,*x,x1,x2,x3,x4;
258   MatScalar      p1,p2,p3,p4;
259   MatScalar      *ba   = b->a,*aa = a->a;
260   PetscReal      shift = info->shiftamount;
261   PetscBool      zeropivotdetected;
262 
263   PetscFunctionBegin;
264   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
265   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
266   ierr = PetscMalloc1(4*(n+1),&rtmp);CHKERRQ(ierr);
267 
268   for (i=0; i<n; i++) {
269     nz    = bi[i+1] - bi[i];
270     ajtmp = bj + bi[i];
271     for  (j=0; j<nz; j++) {
272       x = rtmp+4*ajtmp[j]; x[0] = x[1] = x[2] = x[3] = 0.0;
273     }
274     /* load in initial (unfactored row) */
275     idx      = r[i];
276     nz       = ai[idx+1] - ai[idx];
277     ajtmpold = aj + ai[idx];
278     v        = aa + 4*ai[idx];
279     for (j=0; j<nz; j++) {
280       x    = rtmp+4*ic[ajtmpold[j]];
281       x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3];
282       v   += 4;
283     }
284     row = *ajtmp++;
285     while (row < i) {
286       pc = rtmp + 4*row;
287       p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3];
288       if (p1 != (PetscScalar)0.0 || p2 != (PetscScalar)0.0 || p3 != (PetscScalar)0.0 || p4 != (PetscScalar)0.0) {
289         pv    = ba + 4*diag_offset[row];
290         pj    = bj + diag_offset[row] + 1;
291         x1    = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
292         pc[0] = m1 = p1*x1 + p3*x2;
293         pc[1] = m2 = p2*x1 + p4*x2;
294         pc[2] = m3 = p1*x3 + p3*x4;
295         pc[3] = m4 = p2*x3 + p4*x4;
296         nz    = bi[row+1] - diag_offset[row] - 1;
297         pv   += 4;
298         for (j=0; j<nz; j++) {
299           x1    = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
300           x     = rtmp + 4*pj[j];
301           x[0] -= m1*x1 + m3*x2;
302           x[1] -= m2*x1 + m4*x2;
303           x[2] -= m1*x3 + m3*x4;
304           x[3] -= m2*x3 + m4*x4;
305           pv   += 4;
306         }
307         ierr = PetscLogFlops(16.0*nz+12.0);CHKERRQ(ierr);
308       }
309       row = *ajtmp++;
310     }
311     /* finished row so stick it into b->a */
312     pv = ba + 4*bi[i];
313     pj = bj + bi[i];
314     nz = bi[i+1] - bi[i];
315     for (j=0; j<nz; j++) {
316       x     = rtmp+4*pj[j];
317       pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3];
318       pv   += 4;
319     }
320     /* invert diagonal block */
321     w    = ba + 4*diag_offset[i];
322     ierr = PetscKernel_A_gets_inverse_A_2(w,shift,!A->erroriffailure,&zeropivotdetected);CHKERRQ(ierr);
323     if (zeropivotdetected) C->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
324   }
325 
326   ierr = PetscFree(rtmp);CHKERRQ(ierr);
327   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
328   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
329 
330   C->ops->solve          = MatSolve_SeqBAIJ_2_inplace;
331   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_inplace;
332   C->assembled           = PETSC_TRUE;
333 
334   ierr = PetscLogFlops(1.333333333333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
335   PetscFunctionReturn(0);
336 }
337 /*
338       Version for when blocks are 2 by 2 Using natural ordering
339 */
340 #undef __FUNCT__
341 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_inplace"
342 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_inplace(Mat C,Mat A,const MatFactorInfo *info)
343 {
344   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)C->data;
345   PetscErrorCode ierr;
346   PetscInt       i,j,n = a->mbs,*bi = b->i,*bj = b->j;
347   PetscInt       *ajtmpold,*ajtmp,nz,row;
348   PetscInt       *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj;
349   MatScalar      *pv,*v,*rtmp,*pc,*w,*x;
350   MatScalar      p1,p2,p3,p4,m1,m2,m3,m4,x1,x2,x3,x4;
351   MatScalar      *ba   = b->a,*aa = a->a;
352   PetscReal      shift = info->shiftamount;
353   PetscBool      zeropivotdetected;
354 
355   PetscFunctionBegin;
356   ierr = PetscMalloc1(4*(n+1),&rtmp);CHKERRQ(ierr);
357   for (i=0; i<n; i++) {
358     nz    = bi[i+1] - bi[i];
359     ajtmp = bj + bi[i];
360     for  (j=0; j<nz; j++) {
361       x    = rtmp+4*ajtmp[j];
362       x[0] = x[1]  = x[2]  = x[3]  = 0.0;
363     }
364     /* load in initial (unfactored row) */
365     nz       = ai[i+1] - ai[i];
366     ajtmpold = aj + ai[i];
367     v        = aa + 4*ai[i];
368     for (j=0; j<nz; j++) {
369       x    = rtmp+4*ajtmpold[j];
370       x[0] = v[0];  x[1]  = v[1];  x[2]  = v[2];  x[3]  = v[3];
371       v   += 4;
372     }
373     row = *ajtmp++;
374     while (row < i) {
375       pc = rtmp + 4*row;
376       p1 = pc[0];  p2  = pc[1];  p3  = pc[2];  p4  = pc[3];
377       if (p1 != (PetscScalar)0.0 || p2 != (PetscScalar)0.0 || p3 != (PetscScalar)0.0 || p4 != (PetscScalar)0.0) {
378         pv    = ba + 4*diag_offset[row];
379         pj    = bj + diag_offset[row] + 1;
380         x1    = pv[0];  x2  = pv[1];  x3  = pv[2];  x4  = pv[3];
381         pc[0] = m1 = p1*x1 + p3*x2;
382         pc[1] = m2 = p2*x1 + p4*x2;
383         pc[2] = m3 = p1*x3 + p3*x4;
384         pc[3] = m4 = p2*x3 + p4*x4;
385         nz    = bi[row+1] - diag_offset[row] - 1;
386         pv   += 4;
387         for (j=0; j<nz; j++) {
388           x1    = pv[0];  x2  = pv[1];   x3 = pv[2];  x4  = pv[3];
389           x     = rtmp + 4*pj[j];
390           x[0] -= m1*x1 + m3*x2;
391           x[1] -= m2*x1 + m4*x2;
392           x[2] -= m1*x3 + m3*x4;
393           x[3] -= m2*x3 + m4*x4;
394           pv   += 4;
395         }
396         ierr = PetscLogFlops(16.0*nz+12.0);CHKERRQ(ierr);
397       }
398       row = *ajtmp++;
399     }
400     /* finished row so stick it into b->a */
401     pv = ba + 4*bi[i];
402     pj = bj + bi[i];
403     nz = bi[i+1] - bi[i];
404     for (j=0; j<nz; j++) {
405       x     = rtmp+4*pj[j];
406       pv[0] = x[0];  pv[1]  = x[1];  pv[2]  = x[2];  pv[3]  = x[3];
407       /*
408       printf(" col %d:",pj[j]);
409       PetscInt j1;
410       for (j1=0; j1<4; j1++) printf(" %g,",*(pv+j1));
411       printf("\n");
412       */
413       pv += 4;
414     }
415     /* invert diagonal block */
416     w = ba + 4*diag_offset[i];
417     /*
418     printf(" \n%d -th: diag: ",i);
419     for (j=0; j<4; j++) {
420       printf(" %g,",w[j]);
421     }
422     printf("\n----------------------------\n");
423     */
424     ierr = PetscKernel_A_gets_inverse_A_2(w,shift,!A->erroriffailure,&zeropivotdetected);CHKERRQ(ierr);
425     if (zeropivotdetected) C->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
426   }
427 
428   ierr = PetscFree(rtmp);CHKERRQ(ierr);
429 
430   C->ops->solve          = MatSolve_SeqBAIJ_2_NaturalOrdering_inplace;
431   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_NaturalOrdering_inplace;
432   C->assembled           = PETSC_TRUE;
433 
434   ierr = PetscLogFlops(1.333333333333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
435   PetscFunctionReturn(0);
436 }
437 
438 /* ----------------------------------------------------------- */
439 /*
440      Version for when blocks are 1 by 1.
441 */
442 #undef __FUNCT__
443 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_1"
444 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_1(Mat B,Mat A,const MatFactorInfo *info)
445 {
446   Mat             C     =B;
447   Mat_SeqBAIJ     *a    =(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data;
448   IS              isrow = b->row,isicol = b->icol;
449   PetscErrorCode  ierr;
450   const PetscInt  *r,*ic,*ics;
451   const PetscInt  n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bdiag=b->diag;
452   PetscInt        i,j,k,nz,nzL,row,*pj;
453   const PetscInt  *ajtmp,*bjtmp;
454   MatScalar       *rtmp,*pc,multiplier,*pv;
455   const MatScalar *aa=a->a,*v;
456   PetscBool       row_identity,col_identity;
457   FactorShiftCtx  sctx;
458   const PetscInt  *ddiag;
459   PetscReal       rs;
460   MatScalar       d;
461 
462   PetscFunctionBegin;
463   /* MatPivotSetUp(): initialize shift context sctx */
464   ierr = PetscMemzero(&sctx,sizeof(FactorShiftCtx));CHKERRQ(ierr);
465 
466   if (info->shifttype == (PetscReal) MAT_SHIFT_POSITIVE_DEFINITE) { /* set sctx.shift_top=max{rs} */
467     ddiag          = a->diag;
468     sctx.shift_top = info->zeropivot;
469     for (i=0; i<n; i++) {
470       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
471       d  = (aa)[ddiag[i]];
472       rs = -PetscAbsScalar(d) - PetscRealPart(d);
473       v  = aa+ai[i];
474       nz = ai[i+1] - ai[i];
475       for (j=0; j<nz; j++) rs += PetscAbsScalar(v[j]);
476       if (rs>sctx.shift_top) sctx.shift_top = rs;
477     }
478     sctx.shift_top *= 1.1;
479     sctx.nshift_max = 5;
480     sctx.shift_lo   = 0.;
481     sctx.shift_hi   = 1.;
482   }
483 
484   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
485   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
486   ierr = PetscMalloc1(n+1,&rtmp);CHKERRQ(ierr);
487   ics  = ic;
488 
489   do {
490     sctx.newshift = PETSC_FALSE;
491     for (i=0; i<n; i++) {
492       /* zero rtmp */
493       /* L part */
494       nz    = bi[i+1] - bi[i];
495       bjtmp = bj + bi[i];
496       for  (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0;
497 
498       /* U part */
499       nz    = bdiag[i]-bdiag[i+1];
500       bjtmp = bj + bdiag[i+1]+1;
501       for  (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0;
502 
503       /* load in initial (unfactored row) */
504       nz    = ai[r[i]+1] - ai[r[i]];
505       ajtmp = aj + ai[r[i]];
506       v     = aa + ai[r[i]];
507       for (j=0; j<nz; j++) rtmp[ics[ajtmp[j]]] = v[j];
508 
509       /* ZeropivotApply() */
510       rtmp[i] += sctx.shift_amount;  /* shift the diagonal of the matrix */
511 
512       /* elimination */
513       bjtmp = bj + bi[i];
514       row   = *bjtmp++;
515       nzL   = bi[i+1] - bi[i];
516       for (k=0; k < nzL; k++) {
517         pc = rtmp + row;
518         if (*pc != (PetscScalar)0.0) {
519           pv         = b->a + bdiag[row];
520           multiplier = *pc * (*pv);
521           *pc        = multiplier;
522 
523           pj = b->j + bdiag[row+1]+1; /* beginning of U(row,:) */
524           pv = b->a + bdiag[row+1]+1;
525           nz = bdiag[row]-bdiag[row+1]-1; /* num of entries in U(row,:) excluding diag */
526           for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j];
527           ierr = PetscLogFlops(2.0*nz);CHKERRQ(ierr);
528         }
529         row = *bjtmp++;
530       }
531 
532       /* finished row so stick it into b->a */
533       rs = 0.0;
534       /* L part */
535       pv = b->a + bi[i];
536       pj = b->j + bi[i];
537       nz = bi[i+1] - bi[i];
538       for (j=0; j<nz; j++) {
539         pv[j] = rtmp[pj[j]]; rs += PetscAbsScalar(pv[j]);
540       }
541 
542       /* U part */
543       pv = b->a + bdiag[i+1]+1;
544       pj = b->j + bdiag[i+1]+1;
545       nz = bdiag[i] - bdiag[i+1]-1;
546       for (j=0; j<nz; j++) {
547         pv[j] = rtmp[pj[j]]; rs += PetscAbsScalar(pv[j]);
548       }
549 
550       sctx.rs = rs;
551       sctx.pv = rtmp[i];
552       ierr    = MatPivotCheck(B,A,info,&sctx,i);CHKERRQ(ierr);
553       if (sctx.newshift) break; /* break for-loop */
554       rtmp[i] = sctx.pv; /* sctx.pv might be updated in the case of MAT_SHIFT_INBLOCKS */
555 
556       /* Mark diagonal and invert diagonal for simplier triangular solves */
557       pv  = b->a + bdiag[i];
558       *pv = (PetscScalar)1.0/rtmp[i];
559 
560     } /* endof for (i=0; i<n; i++) { */
561 
562     /* MatPivotRefine() */
563     if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE && !sctx.newshift && sctx.shift_fraction>0 && sctx.nshift<sctx.nshift_max) {
564       /*
565        * if no shift in this attempt & shifting & started shifting & can refine,
566        * then try lower shift
567        */
568       sctx.shift_hi       = sctx.shift_fraction;
569       sctx.shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.;
570       sctx.shift_amount   = sctx.shift_fraction * sctx.shift_top;
571       sctx.newshift       = PETSC_TRUE;
572       sctx.nshift++;
573     }
574   } while (sctx.newshift);
575 
576   ierr = PetscFree(rtmp);CHKERRQ(ierr);
577   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
578   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
579 
580   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
581   ierr = ISIdentity(isicol,&col_identity);CHKERRQ(ierr);
582   if (row_identity && col_identity) {
583     C->ops->solve          = MatSolve_SeqBAIJ_1_NaturalOrdering;
584     C->ops->forwardsolve   = MatForwardSolve_SeqBAIJ_1_NaturalOrdering;
585     C->ops->backwardsolve  = MatBackwardSolve_SeqBAIJ_1_NaturalOrdering;
586     C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1_NaturalOrdering;
587   } else {
588     C->ops->solve          = MatSolve_SeqBAIJ_1;
589     C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1;
590   }
591   C->assembled = PETSC_TRUE;
592   ierr         = PetscLogFlops(C->cmap->n);CHKERRQ(ierr);
593 
594   /* MatShiftView(A,info,&sctx) */
595   if (sctx.nshift) {
596     if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
597       ierr = PetscInfo4(A,"number of shift_pd tries %D, shift_amount %g, diagonal shifted up by %e fraction top_value %e\n",sctx.nshift,(double)sctx.shift_amount,(double)sctx.shift_fraction,(double)sctx.shift_top);CHKERRQ(ierr);
598     } else if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
599       ierr = PetscInfo2(A,"number of shift_nz tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr);
600     } else if (info->shifttype == (PetscReal)MAT_SHIFT_INBLOCKS) {
601       ierr = PetscInfo2(A,"number of shift_inblocks applied %D, each shift_amount %g\n",sctx.nshift,(double)info->shiftamount);CHKERRQ(ierr);
602     }
603   }
604   PetscFunctionReturn(0);
605 }
606 
607 #undef __FUNCT__
608 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_1_inplace"
609 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_1_inplace(Mat C,Mat A,const MatFactorInfo *info)
610 {
611   Mat_SeqBAIJ    *a    = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)C->data;
612   IS             isrow = b->row,isicol = b->icol;
613   PetscErrorCode ierr;
614   const PetscInt *r,*ic;
615   PetscInt       i,j,n = a->mbs,*bi = b->i,*bj = b->j;
616   PetscInt       *ajtmpold,*ajtmp,nz,row,*ai = a->i,*aj = a->j;
617   PetscInt       *diag_offset = b->diag,diag,*pj;
618   MatScalar      *pv,*v,*rtmp,multiplier,*pc;
619   MatScalar      *ba = b->a,*aa = a->a;
620   PetscBool      row_identity, col_identity;
621 
622   PetscFunctionBegin;
623   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
624   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
625   ierr = PetscMalloc1(n+1,&rtmp);CHKERRQ(ierr);
626 
627   for (i=0; i<n; i++) {
628     nz    = bi[i+1] - bi[i];
629     ajtmp = bj + bi[i];
630     for  (j=0; j<nz; j++) rtmp[ajtmp[j]] = 0.0;
631 
632     /* load in initial (unfactored row) */
633     nz       = ai[r[i]+1] - ai[r[i]];
634     ajtmpold = aj + ai[r[i]];
635     v        = aa + ai[r[i]];
636     for (j=0; j<nz; j++) rtmp[ic[ajtmpold[j]]] =  v[j];
637 
638     row = *ajtmp++;
639     while (row < i) {
640       pc = rtmp + row;
641       if (*pc != 0.0) {
642         pv         = ba + diag_offset[row];
643         pj         = bj + diag_offset[row] + 1;
644         multiplier = *pc * *pv++;
645         *pc        = multiplier;
646         nz         = bi[row+1] - diag_offset[row] - 1;
647         for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j];
648         ierr = PetscLogFlops(1.0+2.0*nz);CHKERRQ(ierr);
649       }
650       row = *ajtmp++;
651     }
652     /* finished row so stick it into b->a */
653     pv = ba + bi[i];
654     pj = bj + bi[i];
655     nz = bi[i+1] - bi[i];
656     for (j=0; j<nz; j++) pv[j] = rtmp[pj[j]];
657     diag = diag_offset[i] - bi[i];
658     /* check pivot entry for current row */
659     if (pv[diag] == 0.0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot: row in original ordering %D in permuted ordering %D",r[i],i);
660     pv[diag] = 1.0/pv[diag];
661   }
662 
663   ierr = PetscFree(rtmp);CHKERRQ(ierr);
664   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
665   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
666   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
667   ierr = ISIdentity(isicol,&col_identity);CHKERRQ(ierr);
668   if (row_identity && col_identity) {
669     C->ops->solve          = MatSolve_SeqBAIJ_1_NaturalOrdering_inplace;
670     C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1_NaturalOrdering_inplace;
671   } else {
672     C->ops->solve          = MatSolve_SeqBAIJ_1_inplace;
673     C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1_inplace;
674   }
675   C->assembled = PETSC_TRUE;
676   ierr         = PetscLogFlops(C->cmap->n);CHKERRQ(ierr);
677   PetscFunctionReturn(0);
678 }
679 
680 #undef __FUNCT__
681 #define __FUNCT__ "MatGetFactor_seqbaij_petsc"
682 PETSC_EXTERN PetscErrorCode MatGetFactor_seqbaij_petsc(Mat A,const MatFactorType ftype,Mat *B)
683 {
684   PetscInt       n = A->rmap->n;
685   PetscErrorCode ierr;
686 
687   PetscFunctionBegin;
688 #if defined(PETSC_USE_COMPLEX)
689   if (A->hermitian && (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Hermitian Factor is not supported");
690 #endif
691   ierr = MatCreate(PetscObjectComm((PetscObject)A),B);CHKERRQ(ierr);
692   ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr);
693   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU || ftype == MAT_FACTOR_ILUDT) {
694     ierr = MatSetType(*B,MATSEQBAIJ);CHKERRQ(ierr);
695 
696     (*B)->ops->lufactorsymbolic  = MatLUFactorSymbolic_SeqBAIJ;
697     (*B)->ops->ilufactorsymbolic = MatILUFactorSymbolic_SeqBAIJ;
698   } else if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) {
699     ierr = MatSetType(*B,MATSEQSBAIJ);CHKERRQ(ierr);
700     ierr = MatSeqSBAIJSetPreallocation(*B,A->rmap->bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr);
701 
702     (*B)->ops->iccfactorsymbolic      = MatICCFactorSymbolic_SeqBAIJ;
703     (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqBAIJ;
704   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
705   (*B)->factortype = ftype;
706   PetscFunctionReturn(0);
707 }
708 
709 /* ----------------------------------------------------------- */
710 #undef __FUNCT__
711 #define __FUNCT__ "MatLUFactor_SeqBAIJ"
712 PetscErrorCode MatLUFactor_SeqBAIJ(Mat A,IS row,IS col,const MatFactorInfo *info)
713 {
714   PetscErrorCode ierr;
715   Mat            C;
716 
717   PetscFunctionBegin;
718   ierr = MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_LU,&C);CHKERRQ(ierr);
719   ierr = MatLUFactorSymbolic(C,A,row,col,info);CHKERRQ(ierr);
720   ierr = MatLUFactorNumeric(C,A,info);CHKERRQ(ierr);
721 
722   A->ops->solve          = C->ops->solve;
723   A->ops->solvetranspose = C->ops->solvetranspose;
724 
725   ierr = MatHeaderMerge(A,&C);CHKERRQ(ierr);
726   ierr = PetscLogObjectParent((PetscObject)A,(PetscObject)((Mat_SeqBAIJ*)(A->data))->icol);CHKERRQ(ierr);
727   PetscFunctionReturn(0);
728 }
729 
730 #include <../src/mat/impls/sbaij/seq/sbaij.h>
731 #undef __FUNCT__
732 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N"
733 PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N(Mat C,Mat A,const MatFactorInfo *info)
734 {
735   PetscErrorCode ierr;
736   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data;
737   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
738   IS             ip=b->row;
739   const PetscInt *rip;
740   PetscInt       i,j,mbs=a->mbs,bs=A->rmap->bs,*bi=b->i,*bj=b->j,*bcol;
741   PetscInt       *ai=a->i,*aj=a->j;
742   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
743   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
744   PetscReal      rs;
745   FactorShiftCtx sctx;
746 
747   PetscFunctionBegin;
748   if (bs > 1) { /* convert A to a SBAIJ matrix and apply Cholesky factorization from it */
749     if (!a->sbaijMat) {
750       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
751     }
752     ierr = (a->sbaijMat)->ops->choleskyfactornumeric(C,a->sbaijMat,info);CHKERRQ(ierr);
753     ierr = MatDestroy(&a->sbaijMat);CHKERRQ(ierr);
754     PetscFunctionReturn(0);
755   }
756 
757   /* MatPivotSetUp(): initialize shift context sctx */
758   ierr = PetscMemzero(&sctx,sizeof(FactorShiftCtx));CHKERRQ(ierr);
759 
760   ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr);
761   ierr = PetscMalloc3(mbs,&rtmp,mbs,&il,mbs,&jl);CHKERRQ(ierr);
762 
763   sctx.shift_amount = 0.;
764   sctx.nshift       = 0;
765   do {
766     sctx.newshift = PETSC_FALSE;
767     for (i=0; i<mbs; i++) {
768       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
769     }
770 
771     for (k = 0; k<mbs; k++) {
772       bval = ba + bi[k];
773       /* initialize k-th row by the perm[k]-th row of A */
774       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
775       for (j = jmin; j < jmax; j++) {
776         col = rip[aj[j]];
777         if (col >= k) { /* only take upper triangular entry */
778           rtmp[col] = aa[j];
779           *bval++   = 0.0; /* for in-place factorization */
780         }
781       }
782 
783       /* shift the diagonal of the matrix */
784       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
785 
786       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
787       dk = rtmp[k];
788       i  = jl[k]; /* first row to be added to k_th row  */
789 
790       while (i < k) {
791         nexti = jl[i]; /* next row to be added to k_th row */
792 
793         /* compute multiplier, update diag(k) and U(i,k) */
794         ili     = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
795         uikdi   = -ba[ili]*ba[bi[i]]; /* diagonal(k) */
796         dk     += uikdi*ba[ili];
797         ba[ili] = uikdi; /* -U(i,k) */
798 
799         /* add multiple of row i to k-th row */
800         jmin = ili + 1; jmax = bi[i+1];
801         if (jmin < jmax) {
802           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
803           /* update il and jl for row i */
804           il[i] = jmin;
805           j     = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
806         }
807         i = nexti;
808       }
809 
810       /* shift the diagonals when zero pivot is detected */
811       /* compute rs=sum of abs(off-diagonal) */
812       rs   = 0.0;
813       jmin = bi[k]+1;
814       nz   = bi[k+1] - jmin;
815       if (nz) {
816         bcol = bj + jmin;
817         while (nz--) {
818           rs += PetscAbsScalar(rtmp[*bcol]);
819           bcol++;
820         }
821       }
822 
823       sctx.rs = rs;
824       sctx.pv = dk;
825       ierr    = MatPivotCheck(C,A,info,&sctx,k);CHKERRQ(ierr);
826       if (sctx.newshift) break;
827       dk = sctx.pv;
828 
829       /* copy data into U(k,:) */
830       ba[bi[k]] = 1.0/dk; /* U(k,k) */
831       jmin      = bi[k]+1; jmax = bi[k+1];
832       if (jmin < jmax) {
833         for (j=jmin; j<jmax; j++) {
834           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
835         }
836         /* add the k-th row into il and jl */
837         il[k] = jmin;
838         i     = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
839       }
840     }
841   } while (sctx.newshift);
842   ierr = PetscFree3(rtmp,il,jl);CHKERRQ(ierr);
843 
844   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
845 
846   C->assembled    = PETSC_TRUE;
847   C->preallocated = PETSC_TRUE;
848 
849   ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr);
850   if (sctx.nshift) {
851     if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
852       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr);
853     } else if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
854       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr);
855     }
856   }
857   PetscFunctionReturn(0);
858 }
859 
860 #undef __FUNCT__
861 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering"
862 PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info)
863 {
864   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data;
865   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
866   PetscErrorCode ierr;
867   PetscInt       i,j,am=a->mbs;
868   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
869   PetscInt       k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz;
870   MatScalar      *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval;
871   PetscReal      rs;
872   FactorShiftCtx sctx;
873 
874   PetscFunctionBegin;
875   /* MatPivotSetUp(): initialize shift context sctx */
876   ierr = PetscMemzero(&sctx,sizeof(FactorShiftCtx));CHKERRQ(ierr);
877 
878   ierr = PetscMalloc3(am,&rtmp,am,&il,am,&jl);CHKERRQ(ierr);
879 
880   do {
881     sctx.newshift = PETSC_FALSE;
882     for (i=0; i<am; i++) {
883       rtmp[i] = 0.0; jl[i] = am; il[0] = 0;
884     }
885 
886     for (k = 0; k<am; k++) {
887       /* initialize k-th row with elements nonzero in row perm(k) of A */
888       nz   = ai[k+1] - ai[k];
889       acol = aj + ai[k];
890       aval = aa + ai[k];
891       bval = ba + bi[k];
892       while (nz--) {
893         if (*acol < k) { /* skip lower triangular entries */
894           acol++; aval++;
895         } else {
896           rtmp[*acol++] = *aval++;
897           *bval++       = 0.0; /* for in-place factorization */
898         }
899       }
900 
901       /* shift the diagonal of the matrix */
902       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
903 
904       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
905       dk = rtmp[k];
906       i  = jl[k]; /* first row to be added to k_th row  */
907 
908       while (i < k) {
909         nexti = jl[i]; /* next row to be added to k_th row */
910         /* compute multiplier, update D(k) and U(i,k) */
911         ili     = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
912         uikdi   = -ba[ili]*ba[bi[i]];
913         dk     += uikdi*ba[ili];
914         ba[ili] = uikdi; /* -U(i,k) */
915 
916         /* add multiple of row i to k-th row ... */
917         jmin = ili + 1;
918         nz   = bi[i+1] - jmin;
919         if (nz > 0) {
920           bcol = bj + jmin;
921           bval = ba + jmin;
922           while (nz--) rtmp[*bcol++] += uikdi*(*bval++);
923           /* update il and jl for i-th row */
924           il[i] = jmin;
925           j     = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
926         }
927         i = nexti;
928       }
929 
930       /* shift the diagonals when zero pivot is detected */
931       /* compute rs=sum of abs(off-diagonal) */
932       rs   = 0.0;
933       jmin = bi[k]+1;
934       nz   = bi[k+1] - jmin;
935       if (nz) {
936         bcol = bj + jmin;
937         while (nz--) {
938           rs += PetscAbsScalar(rtmp[*bcol]);
939           bcol++;
940         }
941       }
942 
943       sctx.rs = rs;
944       sctx.pv = dk;
945       ierr    = MatPivotCheck(C,A,info,&sctx,k);CHKERRQ(ierr);
946       if (sctx.newshift) break;    /* sctx.shift_amount is updated */
947       dk = sctx.pv;
948 
949       /* copy data into U(k,:) */
950       ba[bi[k]] = 1.0/dk;
951       jmin      = bi[k]+1;
952       nz        = bi[k+1] - jmin;
953       if (nz) {
954         bcol = bj + jmin;
955         bval = ba + jmin;
956         while (nz--) {
957           *bval++       = rtmp[*bcol];
958           rtmp[*bcol++] = 0.0;
959         }
960         /* add k-th row into il and jl */
961         il[k] = jmin;
962         i     = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
963       }
964     }
965   } while (sctx.newshift);
966   ierr = PetscFree3(rtmp,il,jl);CHKERRQ(ierr);
967 
968   C->ops->solve          = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
969   C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
970   C->assembled           = PETSC_TRUE;
971   C->preallocated        = PETSC_TRUE;
972 
973   ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr);
974   if (sctx.nshift) {
975     if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
976       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr);
977     } else if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
978       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr);
979     }
980   }
981   PetscFunctionReturn(0);
982 }
983 
984 #include <petscbt.h>
985 #include <../src/mat/utils/freespace.h>
986 #undef __FUNCT__
987 #define __FUNCT__ "MatICCFactorSymbolic_SeqBAIJ"
988 PetscErrorCode MatICCFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
989 {
990   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data;
991   Mat_SeqSBAIJ       *b;
992   Mat                B;
993   PetscErrorCode     ierr;
994   PetscBool          perm_identity,missing;
995   PetscInt           reallocs=0,i,*ai=a->i,*aj=a->j,am=a->mbs,bs=A->rmap->bs,*ui;
996   const PetscInt     *rip;
997   PetscInt           jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
998   PetscInt           nlnk,*lnk,*lnk_lvl=NULL,ncols,ncols_upper,*cols,*cols_lvl,*uj,**uj_ptr,**uj_lvl_ptr;
999   PetscReal          fill          =info->fill,levels=info->levels;
1000   PetscFreeSpaceList free_space    =NULL,current_space=NULL;
1001   PetscFreeSpaceList free_space_lvl=NULL,current_space_lvl=NULL;
1002   PetscBT            lnkbt;
1003 
1004   PetscFunctionBegin;
1005   ierr = MatMissingDiagonal(A,&missing,&i);CHKERRQ(ierr);
1006   if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i);
1007 
1008   if (bs > 1) {
1009     if (!a->sbaijMat) {
1010       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
1011     }
1012     (fact)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqSBAIJ;  /* undue the change made in MatGetFactor_seqbaij_petsc */
1013 
1014     ierr = MatICCFactorSymbolic(fact,a->sbaijMat,perm,info);CHKERRQ(ierr);
1015     PetscFunctionReturn(0);
1016   }
1017 
1018   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
1019   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
1020 
1021   /* special case that simply copies fill pattern */
1022   if (!levels && perm_identity) {
1023     ierr = PetscMalloc1(am+1,&ui);CHKERRQ(ierr);
1024     for (i=0; i<am; i++) ui[i] = ai[i+1] - a->diag[i]; /* ui: rowlengths - changes when !perm_identity */
1025     B    = fact;
1026     ierr = MatSeqSBAIJSetPreallocation(B,1,0,ui);CHKERRQ(ierr);
1027 
1028 
1029     b  = (Mat_SeqSBAIJ*)B->data;
1030     uj = b->j;
1031     for (i=0; i<am; i++) {
1032       aj = a->j + a->diag[i];
1033       for (j=0; j<ui[i]; j++) *uj++ = *aj++;
1034       b->ilen[i] = ui[i];
1035     }
1036     ierr = PetscFree(ui);CHKERRQ(ierr);
1037 
1038     B->factortype = MAT_FACTOR_NONE;
1039 
1040     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1041     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1042     B->factortype = MAT_FACTOR_ICC;
1043 
1044     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
1045     PetscFunctionReturn(0);
1046   }
1047 
1048   /* initialization */
1049   ierr  = PetscMalloc1(am+1,&ui);CHKERRQ(ierr);
1050   ui[0] = 0;
1051   ierr  = PetscMalloc1(2*am+1,&cols_lvl);CHKERRQ(ierr);
1052 
1053   /* jl: linked list for storing indices of the pivot rows
1054      il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
1055   ierr = PetscMalloc4(am,&uj_ptr,am,&uj_lvl_ptr,am,&il,am,&jl);CHKERRQ(ierr);
1056   for (i=0; i<am; i++) {
1057     jl[i] = am; il[i] = 0;
1058   }
1059 
1060   /* create and initialize a linked list for storing column indices of the active row k */
1061   nlnk = am + 1;
1062   ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
1063 
1064   /* initial FreeSpace size is fill*(ai[am]+am)/2 */
1065   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am]/2,am/2)),&free_space);CHKERRQ(ierr);
1066 
1067   current_space = free_space;
1068 
1069   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am]/2,am/2)),&free_space_lvl);CHKERRQ(ierr);
1070   current_space_lvl = free_space_lvl;
1071 
1072   for (k=0; k<am; k++) {  /* for each active row k */
1073     /* initialize lnk by the column indices of row rip[k] of A */
1074     nzk         = 0;
1075     ncols       = ai[rip[k]+1] - ai[rip[k]];
1076     ncols_upper = 0;
1077     cols        = cols_lvl + am;
1078     for (j=0; j<ncols; j++) {
1079       i = rip[*(aj + ai[rip[k]] + j)];
1080       if (i >= k) { /* only take upper triangular entry */
1081         cols[ncols_upper]     = i;
1082         cols_lvl[ncols_upper] = -1;  /* initialize level for nonzero entries */
1083         ncols_upper++;
1084       }
1085     }
1086     ierr = PetscIncompleteLLAdd(ncols_upper,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
1087     nzk += nlnk;
1088 
1089     /* update lnk by computing fill-in for each pivot row to be merged in */
1090     prow = jl[k]; /* 1st pivot row */
1091 
1092     while (prow < k) {
1093       nextprow = jl[prow];
1094 
1095       /* merge prow into k-th row */
1096       jmin  = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */
1097       jmax  = ui[prow+1];
1098       ncols = jmax-jmin;
1099       i     = jmin - ui[prow];
1100       cols  = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
1101       for (j=0; j<ncols; j++) cols_lvl[j] = *(uj_lvl_ptr[prow] + i + j);
1102       ierr = PetscIncompleteLLAddSorted(ncols,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
1103       nzk += nlnk;
1104 
1105       /* update il and jl for prow */
1106       if (jmin < jmax) {
1107         il[prow] = jmin;
1108 
1109         j = *cols; jl[prow] = jl[j]; jl[j] = prow;
1110       }
1111       prow = nextprow;
1112     }
1113 
1114     /* if free space is not available, make more free space */
1115     if (current_space->local_remaining<nzk) {
1116       i    = am - k + 1; /* num of unfactored rows */
1117       i    = PetscMin(PetscIntMultTruncate(i,nzk), PetscIntMultTruncate(i,i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
1118       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
1119       ierr = PetscFreeSpaceGet(i,&current_space_lvl);CHKERRQ(ierr);
1120       reallocs++;
1121     }
1122 
1123     /* copy data into free_space and free_space_lvl, then initialize lnk */
1124     ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
1125 
1126     /* add the k-th row into il and jl */
1127     if (nzk-1 > 0) {
1128       i     = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
1129       jl[k] = jl[i]; jl[i] = k;
1130       il[k] = ui[k] + 1;
1131     }
1132     uj_ptr[k]     = current_space->array;
1133     uj_lvl_ptr[k] = current_space_lvl->array;
1134 
1135     current_space->array           += nzk;
1136     current_space->local_used      += nzk;
1137     current_space->local_remaining -= nzk;
1138 
1139     current_space_lvl->array           += nzk;
1140     current_space_lvl->local_used      += nzk;
1141     current_space_lvl->local_remaining -= nzk;
1142 
1143     ui[k+1] = ui[k] + nzk;
1144   }
1145 
1146   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
1147   ierr = PetscFree4(uj_ptr,uj_lvl_ptr,il,jl);CHKERRQ(ierr);
1148   ierr = PetscFree(cols_lvl);CHKERRQ(ierr);
1149 
1150   /* copy free_space into uj and free free_space; set uj in new datastructure; */
1151   ierr = PetscMalloc1(ui[am]+1,&uj);CHKERRQ(ierr);
1152   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
1153   ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1154   ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
1155 
1156   /* put together the new matrix in MATSEQSBAIJ format */
1157   B    = fact;
1158   ierr = MatSeqSBAIJSetPreallocation(B,1,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr);
1159 
1160   b                = (Mat_SeqSBAIJ*)B->data;
1161   b->singlemalloc  = PETSC_FALSE;
1162   b->free_a        = PETSC_TRUE;
1163   b->free_ij       = PETSC_TRUE;
1164 
1165   ierr = PetscMalloc1(ui[am]+1,&b->a);CHKERRQ(ierr);
1166 
1167   b->j             = uj;
1168   b->i             = ui;
1169   b->diag          = 0;
1170   b->ilen          = 0;
1171   b->imax          = 0;
1172   b->row           = perm;
1173   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
1174 
1175   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
1176 
1177   b->icol = perm;
1178 
1179   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
1180   ierr    = PetscMalloc1(am+1,&b->solve_work);CHKERRQ(ierr);
1181   ierr    = PetscLogObjectMemory((PetscObject)B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
1182 
1183   b->maxnz = b->nz = ui[am];
1184 
1185   B->info.factor_mallocs   = reallocs;
1186   B->info.fill_ratio_given = fill;
1187   if (ai[am] != 0.) {
1188     /* nonzeros in lower triangular part of A (includign diagonals)= (ai[am]+am)/2 */
1189     B->info.fill_ratio_needed = ((PetscReal)2*ui[am])/(ai[am]+am);
1190   } else {
1191     B->info.fill_ratio_needed = 0.0;
1192   }
1193 #if defined(PETSC_USE_INFO)
1194   if (ai[am] != 0) {
1195     PetscReal af = B->info.fill_ratio_needed;
1196     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)fill,(double)af);CHKERRQ(ierr);
1197     ierr = PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr);
1198     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%g) for best performance.\n",(double)af);CHKERRQ(ierr);
1199   } else {
1200     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
1201   }
1202 #endif
1203   if (perm_identity) {
1204     B->ops->solve                 = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
1205     B->ops->solvetranspose        = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
1206     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
1207   } else {
1208     (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
1209   }
1210   PetscFunctionReturn(0);
1211 }
1212 
1213 #undef __FUNCT__
1214 #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqBAIJ"
1215 PetscErrorCode MatCholeskyFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
1216 {
1217   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data;
1218   Mat_SeqSBAIJ       *b;
1219   Mat                B;
1220   PetscErrorCode     ierr;
1221   PetscBool          perm_identity,missing;
1222   PetscReal          fill = info->fill;
1223   const PetscInt     *rip;
1224   PetscInt           i,mbs=a->mbs,bs=A->rmap->bs,*ai=a->i,*aj=a->j,reallocs=0,prow;
1225   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
1226   PetscInt           nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr;
1227   PetscFreeSpaceList free_space=NULL,current_space=NULL;
1228   PetscBT            lnkbt;
1229 
1230   PetscFunctionBegin;
1231   if (bs > 1) { /* convert to seqsbaij */
1232     if (!a->sbaijMat) {
1233       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
1234     }
1235     (fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqSBAIJ; /* undue the change made in MatGetFactor_seqbaij_petsc */
1236 
1237     ierr = MatCholeskyFactorSymbolic(fact,a->sbaijMat,perm,info);CHKERRQ(ierr);
1238     PetscFunctionReturn(0);
1239   }
1240 
1241   ierr = MatMissingDiagonal(A,&missing,&i);CHKERRQ(ierr);
1242   if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i);
1243 
1244   /* check whether perm is the identity mapping */
1245   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
1246   if (!perm_identity) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix reordering is not supported");
1247   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
1248 
1249   /* initialization */
1250   ierr  = PetscMalloc1(mbs+1,&ui);CHKERRQ(ierr);
1251   ui[0] = 0;
1252 
1253   /* jl: linked list for storing indices of the pivot rows
1254      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
1255   ierr = PetscMalloc4(mbs,&ui_ptr,mbs,&il,mbs,&jl,mbs,&cols);CHKERRQ(ierr);
1256   for (i=0; i<mbs; i++) {
1257     jl[i] = mbs; il[i] = 0;
1258   }
1259 
1260   /* create and initialize a linked list for storing column indices of the active row k */
1261   nlnk = mbs + 1;
1262   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
1263 
1264   /* initial FreeSpace size is fill* (ai[mbs]+mbs)/2 */
1265   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[mbs]/2,mbs/2)),&free_space);CHKERRQ(ierr);
1266 
1267   current_space = free_space;
1268 
1269   for (k=0; k<mbs; k++) {  /* for each active row k */
1270     /* initialize lnk by the column indices of row rip[k] of A */
1271     nzk         = 0;
1272     ncols       = ai[rip[k]+1] - ai[rip[k]];
1273     ncols_upper = 0;
1274     for (j=0; j<ncols; j++) {
1275       i = rip[*(aj + ai[rip[k]] + j)];
1276       if (i >= k) { /* only take upper triangular entry */
1277         cols[ncols_upper] = i;
1278         ncols_upper++;
1279       }
1280     }
1281     ierr = PetscLLAdd(ncols_upper,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
1282     nzk += nlnk;
1283 
1284     /* update lnk by computing fill-in for each pivot row to be merged in */
1285     prow = jl[k]; /* 1st pivot row */
1286 
1287     while (prow < k) {
1288       nextprow = jl[prow];
1289       /* merge prow into k-th row */
1290       jmin   = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
1291       jmax   = ui[prow+1];
1292       ncols  = jmax-jmin;
1293       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
1294       ierr   = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
1295       nzk   += nlnk;
1296 
1297       /* update il and jl for prow */
1298       if (jmin < jmax) {
1299         il[prow] = jmin;
1300         j        = *uj_ptr;
1301         jl[prow] = jl[j];
1302         jl[j]    = prow;
1303       }
1304       prow = nextprow;
1305     }
1306 
1307     /* if free space is not available, make more free space */
1308     if (current_space->local_remaining<nzk) {
1309       i    = mbs - k + 1; /* num of unfactored rows */
1310       i    = PetscMin(PetscIntMultTruncate(i,nzk), PetscIntMultTruncate(i,i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
1311       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
1312       reallocs++;
1313     }
1314 
1315     /* copy data into free space, then initialize lnk */
1316     ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
1317 
1318     /* add the k-th row into il and jl */
1319     if (nzk-1 > 0) {
1320       i     = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
1321       jl[k] = jl[i]; jl[i] = k;
1322       il[k] = ui[k] + 1;
1323     }
1324     ui_ptr[k]                       = current_space->array;
1325     current_space->array           += nzk;
1326     current_space->local_used      += nzk;
1327     current_space->local_remaining -= nzk;
1328 
1329     ui[k+1] = ui[k] + nzk;
1330   }
1331 
1332   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
1333   ierr = PetscFree4(ui_ptr,il,jl,cols);CHKERRQ(ierr);
1334 
1335   /* copy free_space into uj and free free_space; set uj in new datastructure; */
1336   ierr = PetscMalloc1(ui[mbs]+1,&uj);CHKERRQ(ierr);
1337   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
1338   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1339 
1340   /* put together the new matrix in MATSEQSBAIJ format */
1341   B    = fact;
1342   ierr = MatSeqSBAIJSetPreallocation(B,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr);
1343 
1344   b               = (Mat_SeqSBAIJ*)B->data;
1345   b->singlemalloc = PETSC_FALSE;
1346   b->free_a       = PETSC_TRUE;
1347   b->free_ij      = PETSC_TRUE;
1348 
1349   ierr = PetscMalloc1(ui[mbs]+1,&b->a);CHKERRQ(ierr);
1350 
1351   b->j             = uj;
1352   b->i             = ui;
1353   b->diag          = 0;
1354   b->ilen          = 0;
1355   b->imax          = 0;
1356   b->row           = perm;
1357   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
1358 
1359   ierr     = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
1360   b->icol  = perm;
1361   ierr     = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
1362   ierr     = PetscMalloc1(mbs+1,&b->solve_work);CHKERRQ(ierr);
1363   ierr     = PetscLogObjectMemory((PetscObject)B,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
1364   b->maxnz = b->nz = ui[mbs];
1365 
1366   B->info.factor_mallocs   = reallocs;
1367   B->info.fill_ratio_given = fill;
1368   if (ai[mbs] != 0.) {
1369     /* nonzeros in lower triangular part of A = (ai[mbs]+mbs)/2 */
1370     B->info.fill_ratio_needed = ((PetscReal)2*ui[mbs])/(ai[mbs]+mbs);
1371   } else {
1372     B->info.fill_ratio_needed = 0.0;
1373   }
1374 #if defined(PETSC_USE_INFO)
1375   if (ai[mbs] != 0.) {
1376     PetscReal af = B->info.fill_ratio_needed;
1377     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)fill,(double)af);CHKERRQ(ierr);
1378     ierr = PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr);
1379     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%g) for best performance.\n",(double)af);CHKERRQ(ierr);
1380   } else {
1381     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
1382   }
1383 #endif
1384   if (perm_identity) {
1385     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
1386   } else {
1387     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
1388   }
1389   PetscFunctionReturn(0);
1390 }
1391 
1392 #undef __FUNCT__
1393 #define __FUNCT__ "MatSolve_SeqBAIJ_N_NaturalOrdering"
1394 PetscErrorCode MatSolve_SeqBAIJ_N_NaturalOrdering(Mat A,Vec bb,Vec xx)
1395 {
1396   Mat_SeqBAIJ       *a=(Mat_SeqBAIJ*)A->data;
1397   PetscErrorCode    ierr;
1398   const PetscInt    *ai=a->i,*aj=a->j,*adiag=a->diag,*vi;
1399   PetscInt          i,k,n=a->mbs;
1400   PetscInt          nz,bs=A->rmap->bs,bs2=a->bs2;
1401   const MatScalar   *aa=a->a,*v;
1402   PetscScalar       *x,*s,*t,*ls;
1403   const PetscScalar *b;
1404 
1405   PetscFunctionBegin;
1406   ierr = VecGetArrayRead(bb,&b);CHKERRQ(ierr);
1407   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1408   t    = a->solve_work;
1409 
1410   /* forward solve the lower triangular */
1411   ierr = PetscMemcpy(t,b,bs*sizeof(PetscScalar));CHKERRQ(ierr); /* copy 1st block of b to t */
1412 
1413   for (i=1; i<n; i++) {
1414     v    = aa + bs2*ai[i];
1415     vi   = aj + ai[i];
1416     nz   = ai[i+1] - ai[i];
1417     s    = t + bs*i;
1418     ierr = PetscMemcpy(s,b+bs*i,bs*sizeof(PetscScalar));CHKERRQ(ierr); /* copy i_th block of b to t */
1419     for (k=0;k<nz;k++) {
1420       PetscKernel_v_gets_v_minus_A_times_w(bs,s,v,t+bs*vi[k]);
1421       v += bs2;
1422     }
1423   }
1424 
1425   /* backward solve the upper triangular */
1426   ls = a->solve_work + A->cmap->n;
1427   for (i=n-1; i>=0; i--) {
1428     v    = aa + bs2*(adiag[i+1]+1);
1429     vi   = aj + adiag[i+1]+1;
1430     nz   = adiag[i] - adiag[i+1]-1;
1431     ierr = PetscMemcpy(ls,t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr);
1432     for (k=0; k<nz; k++) {
1433       PetscKernel_v_gets_v_minus_A_times_w(bs,ls,v,t+bs*vi[k]);
1434       v += bs2;
1435     }
1436     PetscKernel_w_gets_A_times_v(bs,ls,aa+bs2*adiag[i],t+i*bs); /* *inv(diagonal[i]) */
1437     ierr = PetscMemcpy(x+i*bs,t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr);
1438   }
1439 
1440   ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr);
1441   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1442   ierr = PetscLogFlops(2.0*(a->bs2)*(a->nz) - A->rmap->bs*A->cmap->n);CHKERRQ(ierr);
1443   PetscFunctionReturn(0);
1444 }
1445 
1446 #undef __FUNCT__
1447 #define __FUNCT__ "MatSolve_SeqBAIJ_N"
1448 PetscErrorCode MatSolve_SeqBAIJ_N(Mat A,Vec bb,Vec xx)
1449 {
1450   Mat_SeqBAIJ        *a   =(Mat_SeqBAIJ*)A->data;
1451   IS                 iscol=a->col,isrow=a->row;
1452   PetscErrorCode     ierr;
1453   const PetscInt     *r,*c,*rout,*cout,*ai=a->i,*aj=a->j,*adiag=a->diag,*vi;
1454   PetscInt           i,m,n=a->mbs;
1455   PetscInt           nz,bs=A->rmap->bs,bs2=a->bs2;
1456   const MatScalar    *aa=a->a,*v;
1457   PetscScalar        *x,*s,*t,*ls;
1458   const PetscScalar  *b;
1459 
1460   PetscFunctionBegin;
1461   ierr = VecGetArrayRead(bb,&b);CHKERRQ(ierr);
1462   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1463   t    = a->solve_work;
1464 
1465   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
1466   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
1467 
1468   /* forward solve the lower triangular */
1469   ierr = PetscMemcpy(t,b+bs*r[0],bs*sizeof(PetscScalar));CHKERRQ(ierr);
1470   for (i=1; i<n; i++) {
1471     v    = aa + bs2*ai[i];
1472     vi   = aj + ai[i];
1473     nz   = ai[i+1] - ai[i];
1474     s    = t + bs*i;
1475     ierr = PetscMemcpy(s,b+bs*r[i],bs*sizeof(PetscScalar));CHKERRQ(ierr);
1476     for (m=0; m<nz; m++) {
1477       PetscKernel_v_gets_v_minus_A_times_w(bs,s,v,t+bs*vi[m]);
1478       v += bs2;
1479     }
1480   }
1481 
1482   /* backward solve the upper triangular */
1483   ls = a->solve_work + A->cmap->n;
1484   for (i=n-1; i>=0; i--) {
1485     v    = aa + bs2*(adiag[i+1]+1);
1486     vi   = aj + adiag[i+1]+1;
1487     nz   = adiag[i] - adiag[i+1] - 1;
1488     ierr = PetscMemcpy(ls,t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr);
1489     for (m=0; m<nz; m++) {
1490       PetscKernel_v_gets_v_minus_A_times_w(bs,ls,v,t+bs*vi[m]);
1491       v += bs2;
1492     }
1493     PetscKernel_w_gets_A_times_v(bs,ls,v,t+i*bs); /* *inv(diagonal[i]) */
1494     ierr = PetscMemcpy(x + bs*c[i],t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr);
1495   }
1496   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
1497   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
1498   ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr);
1499   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1500   ierr = PetscLogFlops(2.0*(a->bs2)*(a->nz) - A->rmap->bs*A->cmap->n);CHKERRQ(ierr);
1501   PetscFunctionReturn(0);
1502 }
1503 
1504 #undef __FUNCT__
1505 #define __FUNCT__ "MatBlockAbs_privat"
1506 /*
1507     For each block in an block array saves the largest absolute value in the block into another array
1508 */
1509 static PetscErrorCode MatBlockAbs_private(PetscInt nbs,PetscInt bs2,PetscScalar *blockarray,PetscReal *absarray)
1510 {
1511   PetscErrorCode ierr;
1512   PetscInt       i,j;
1513 
1514   PetscFunctionBegin;
1515   ierr = PetscMemzero(absarray,(nbs+1)*sizeof(PetscReal));CHKERRQ(ierr);
1516   for (i=0; i<nbs; i++) {
1517     for (j=0; j<bs2; j++) {
1518       if (absarray[i] < PetscAbsScalar(blockarray[i*nbs+j])) absarray[i] = PetscAbsScalar(blockarray[i*nbs+j]);
1519     }
1520   }
1521   PetscFunctionReturn(0);
1522 }
1523 
1524 #undef __FUNCT__
1525 #define __FUNCT__ "MatILUDTFactor_SeqBAIJ"
1526 /*
1527      This needs to be renamed and called by the regular MatILUFactor_SeqBAIJ when drop tolerance is used
1528 */
1529 PetscErrorCode MatILUDTFactor_SeqBAIJ(Mat A,IS isrow,IS iscol,const MatFactorInfo *info,Mat *fact)
1530 {
1531   Mat            B = *fact;
1532   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data,*b;
1533   IS             isicol;
1534   PetscErrorCode ierr;
1535   const PetscInt *r,*ic;
1536   PetscInt       i,mbs=a->mbs,bs=A->rmap->bs,bs2=a->bs2,*ai=a->i,*aj=a->j,*ajtmp,*adiag;
1537   PetscInt       *bi,*bj,*bdiag;
1538 
1539   PetscInt  row,nzi,nzi_bl,nzi_bu,*im,dtcount,nzi_al,nzi_au;
1540   PetscInt  nlnk,*lnk;
1541   PetscBT   lnkbt;
1542   PetscBool row_identity,icol_identity;
1543   MatScalar *aatmp,*pv,*batmp,*ba,*rtmp,*pc,*multiplier,*vtmp;
1544   PetscInt  j,nz,*pj,*bjtmp,k,ncut,*jtmp;
1545 
1546   PetscReal dt=info->dt;          /* shift=info->shiftamount; */
1547   PetscInt  nnz_max;
1548   PetscBool missing;
1549   PetscReal *vtmp_abs;
1550   MatScalar *v_work;
1551   PetscInt  *v_pivots;
1552 
1553   PetscFunctionBegin;
1554   /* ------- symbolic factorization, can be reused ---------*/
1555   ierr = MatMissingDiagonal(A,&missing,&i);CHKERRQ(ierr);
1556   if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i);
1557   adiag=a->diag;
1558 
1559   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
1560 
1561   /* bdiag is location of diagonal in factor */
1562   ierr = PetscMalloc1(mbs+1,&bdiag);CHKERRQ(ierr);
1563 
1564   /* allocate row pointers bi */
1565   ierr = PetscMalloc1(2*mbs+2,&bi);CHKERRQ(ierr);
1566 
1567   /* allocate bj and ba; max num of nonzero entries is (ai[n]+2*n*dtcount+2) */
1568   dtcount = (PetscInt)info->dtcount;
1569   if (dtcount > mbs-1) dtcount = mbs-1;
1570   nnz_max = ai[mbs]+2*mbs*dtcount +2;
1571   /* printf("MatILUDTFactor_SeqBAIJ, bs %d, ai[mbs] %d, nnz_max  %d, dtcount %d\n",bs,ai[mbs],nnz_max,dtcount); */
1572   ierr    = PetscMalloc1(nnz_max,&bj);CHKERRQ(ierr);
1573   nnz_max = nnz_max*bs2;
1574   ierr    = PetscMalloc1(nnz_max,&ba);CHKERRQ(ierr);
1575 
1576   /* put together the new matrix */
1577   ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(B,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr);
1578   ierr = PetscLogObjectParent((PetscObject)B,(PetscObject)isicol);CHKERRQ(ierr);
1579 
1580   b               = (Mat_SeqBAIJ*)(B)->data;
1581   b->free_a       = PETSC_TRUE;
1582   b->free_ij      = PETSC_TRUE;
1583   b->singlemalloc = PETSC_FALSE;
1584 
1585   b->a    = ba;
1586   b->j    = bj;
1587   b->i    = bi;
1588   b->diag = bdiag;
1589   b->ilen = 0;
1590   b->imax = 0;
1591   b->row  = isrow;
1592   b->col  = iscol;
1593 
1594   ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1595   ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
1596 
1597   b->icol  = isicol;
1598   ierr     = PetscMalloc1(bs*(mbs+1),&b->solve_work);CHKERRQ(ierr);
1599   ierr     = PetscLogObjectMemory((PetscObject)B,nnz_max*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
1600   b->maxnz = nnz_max/bs2;
1601 
1602   (B)->factortype            = MAT_FACTOR_ILUDT;
1603   (B)->info.factor_mallocs   = 0;
1604   (B)->info.fill_ratio_given = ((PetscReal)nnz_max)/((PetscReal)(ai[mbs]*bs2));
1605   /* ------- end of symbolic factorization ---------*/
1606   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
1607   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
1608 
1609   /* linked list for storing column indices of the active row */
1610   nlnk = mbs + 1;
1611   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
1612 
1613   /* im: used by PetscLLAddSortedLU(); jtmp: working array for column indices of active row */
1614   ierr = PetscMalloc2(mbs,&im,mbs,&jtmp);CHKERRQ(ierr);
1615   /* rtmp, vtmp: working arrays for sparse and contiguous row entries of active row */
1616   ierr = PetscMalloc2(mbs*bs2,&rtmp,mbs*bs2,&vtmp);CHKERRQ(ierr);
1617   ierr = PetscMalloc1(mbs+1,&vtmp_abs);CHKERRQ(ierr);
1618   ierr = PetscMalloc3(bs,&v_work,bs2,&multiplier,bs,&v_pivots);CHKERRQ(ierr);
1619 
1620   bi[0]       = 0;
1621   bdiag[0]    = (nnz_max/bs2)-1; /* location of diagonal in factor B */
1622   bi[2*mbs+1] = bdiag[0]+1; /* endof bj and ba array */
1623   for (i=0; i<mbs; i++) {
1624     /* copy initial fill into linked list */
1625     nzi = 0; /* nonzeros for active row i */
1626     nzi = ai[r[i]+1] - ai[r[i]];
1627     if (!nzi) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i);
1628     nzi_al = adiag[r[i]] - ai[r[i]];
1629     nzi_au = ai[r[i]+1] - adiag[r[i]] -1;
1630     /* printf("row %d, nzi_al/au %d %d\n",i,nzi_al,nzi_au); */
1631 
1632     /* load in initial unfactored row */
1633     ajtmp = aj + ai[r[i]];
1634     ierr  = PetscLLAddPerm(nzi,ajtmp,ic,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
1635     ierr  = PetscMemzero(rtmp,mbs*bs2*sizeof(PetscScalar));CHKERRQ(ierr);
1636     aatmp = a->a + bs2*ai[r[i]];
1637     for (j=0; j<nzi; j++) {
1638       ierr = PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],aatmp+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1639     }
1640 
1641     /* add pivot rows into linked list */
1642     row = lnk[mbs];
1643     while (row < i) {
1644       nzi_bl = bi[row+1] - bi[row] + 1;
1645       bjtmp  = bj + bdiag[row+1]+1; /* points to 1st column next to the diagonal in U */
1646       ierr   = PetscLLAddSortedLU(bjtmp,row,nlnk,lnk,lnkbt,i,nzi_bl,im);CHKERRQ(ierr);
1647       nzi   += nlnk;
1648       row    = lnk[row];
1649     }
1650 
1651     /* copy data from lnk into jtmp, then initialize lnk */
1652     ierr = PetscLLClean(mbs,mbs,nzi,lnk,jtmp,lnkbt);CHKERRQ(ierr);
1653 
1654     /* numerical factorization */
1655     bjtmp = jtmp;
1656     row   = *bjtmp++; /* 1st pivot row */
1657 
1658     while  (row < i) {
1659       pc = rtmp + bs2*row;
1660       pv = ba + bs2*bdiag[row]; /* inv(diag) of the pivot row */
1661       PetscKernel_A_gets_A_times_B(bs,pc,pv,multiplier); /* pc= multiplier = pc*inv(diag[row]) */
1662       ierr = MatBlockAbs_private(1,bs2,pc,vtmp_abs);CHKERRQ(ierr);
1663       if (vtmp_abs[0] > dt) { /* apply tolerance dropping rule */
1664         pj = bj + bdiag[row+1] + 1;         /* point to 1st entry of U(row,:) */
1665         pv = ba + bs2*(bdiag[row+1] + 1);
1666         nz = bdiag[row] - bdiag[row+1] - 1;         /* num of entries in U(row,:), excluding diagonal */
1667         for (j=0; j<nz; j++) {
1668           PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j);
1669         }
1670         /* ierr = PetscLogFlops(bslog*(nz+1.0)-bs);CHKERRQ(ierr); */
1671       }
1672       row = *bjtmp++;
1673     }
1674 
1675     /* copy sparse rtmp into contiguous vtmp; separate L and U part */
1676     nzi_bl = 0; j = 0;
1677     while (jtmp[j] < i) { /* L-part. Note: jtmp is sorted */
1678       ierr = PetscMemcpy(vtmp+bs2*j,rtmp+bs2*jtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
1679       nzi_bl++; j++;
1680     }
1681     nzi_bu = nzi - nzi_bl -1;
1682     /* printf("nzi %d, nzi_bl %d, nzi_bu %d\n",nzi,nzi_bl,nzi_bu); */
1683 
1684     while (j < nzi) { /* U-part */
1685       ierr = PetscMemcpy(vtmp+bs2*j,rtmp+bs2*jtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
1686       /*
1687       printf(" col %d: ",jtmp[j]);
1688       for (j1=0; j1<bs2; j1++) printf(" %g",*(vtmp+bs2*j+j1));
1689       printf(" \n");
1690       */
1691       j++;
1692     }
1693 
1694     ierr = MatBlockAbs_private(nzi,bs2,vtmp,vtmp_abs);CHKERRQ(ierr);
1695     /*
1696     printf(" row %d, nzi %d, vtmp_abs\n",i,nzi);
1697     for (j1=0; j1<nzi; j1++) printf(" (%d %g),",jtmp[j1],vtmp_abs[j1]);
1698     printf(" \n");
1699     */
1700     bjtmp = bj + bi[i];
1701     batmp = ba + bs2*bi[i];
1702     /* apply level dropping rule to L part */
1703     ncut = nzi_al + dtcount;
1704     if (ncut < nzi_bl) {
1705       ierr = PetscSortSplitReal(ncut,nzi_bl,vtmp_abs,jtmp);CHKERRQ(ierr);
1706       ierr = PetscSortIntWithScalarArray(ncut,jtmp,vtmp);CHKERRQ(ierr);
1707     } else {
1708       ncut = nzi_bl;
1709     }
1710     for (j=0; j<ncut; j++) {
1711       bjtmp[j] = jtmp[j];
1712       ierr     = PetscMemcpy(batmp+bs2*j,rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
1713       /*
1714       printf(" col %d: ",bjtmp[j]);
1715       for (j1=0; j1<bs2; j1++) printf(" %g,",*(batmp+bs2*j+j1));
1716       printf("\n");
1717       */
1718     }
1719     bi[i+1] = bi[i] + ncut;
1720     nzi     = ncut + 1;
1721 
1722     /* apply level dropping rule to U part */
1723     ncut = nzi_au + dtcount;
1724     if (ncut < nzi_bu) {
1725       ierr = PetscSortSplitReal(ncut,nzi_bu,vtmp_abs+nzi_bl+1,jtmp+nzi_bl+1);CHKERRQ(ierr);
1726       ierr = PetscSortIntWithScalarArray(ncut,jtmp+nzi_bl+1,vtmp+nzi_bl+1);CHKERRQ(ierr);
1727     } else {
1728       ncut = nzi_bu;
1729     }
1730     nzi += ncut;
1731 
1732     /* mark bdiagonal */
1733     bdiag[i+1]    = bdiag[i] - (ncut + 1);
1734     bi[2*mbs - i] = bi[2*mbs - i +1] - (ncut + 1);
1735 
1736     bjtmp  = bj + bdiag[i];
1737     batmp  = ba + bs2*bdiag[i];
1738     ierr   = PetscMemcpy(batmp,rtmp+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1739     *bjtmp = i;
1740     /*
1741     printf(" diag %d: ",*bjtmp);
1742     for (j=0; j<bs2; j++) {
1743       printf(" %g,",batmp[j]);
1744     }
1745     printf("\n");
1746     */
1747     bjtmp = bj + bdiag[i+1]+1;
1748     batmp = ba + (bdiag[i+1]+1)*bs2;
1749 
1750     for (k=0; k<ncut; k++) {
1751       bjtmp[k] = jtmp[nzi_bl+1+k];
1752       ierr     = PetscMemcpy(batmp+bs2*k,rtmp+bs2*bjtmp[k],bs2*sizeof(MatScalar));CHKERRQ(ierr);
1753       /*
1754       printf(" col %d:",bjtmp[k]);
1755       for (j1=0; j1<bs2; j1++) printf(" %g,",*(batmp+bs2*k+j1));
1756       printf("\n");
1757       */
1758     }
1759 
1760     im[i] = nzi; /* used by PetscLLAddSortedLU() */
1761 
1762     /* invert diagonal block for simplier triangular solves - add shift??? */
1763     batmp = ba + bs2*bdiag[i];
1764     ierr  = PetscKernel_A_gets_inverse_A(bs,batmp,v_pivots,v_work);CHKERRQ(ierr);
1765   } /* for (i=0; i<mbs; i++) */
1766   ierr = PetscFree3(v_work,multiplier,v_pivots);CHKERRQ(ierr);
1767 
1768   /* printf("end of L %d, beginning of U %d\n",bi[mbs],bdiag[mbs]); */
1769   if (bi[mbs] >= bdiag[mbs]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"end of L array %d cannot >= the beginning of U array %d",bi[mbs],bdiag[mbs]);
1770 
1771   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
1772   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
1773 
1774   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1775 
1776   ierr = PetscFree2(im,jtmp);CHKERRQ(ierr);
1777   ierr = PetscFree2(rtmp,vtmp);CHKERRQ(ierr);
1778 
1779   ierr     = PetscLogFlops(bs2*B->cmap->n);CHKERRQ(ierr);
1780   b->maxnz = b->nz = bi[mbs] + bdiag[0] - bdiag[mbs];
1781 
1782   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
1783   ierr = ISIdentity(isicol,&icol_identity);CHKERRQ(ierr);
1784   if (row_identity && icol_identity) {
1785     B->ops->solve = MatSolve_SeqBAIJ_N_NaturalOrdering;
1786   } else {
1787     B->ops->solve = MatSolve_SeqBAIJ_N;
1788   }
1789 
1790   B->ops->solveadd          = 0;
1791   B->ops->solvetranspose    = 0;
1792   B->ops->solvetransposeadd = 0;
1793   B->ops->matsolve          = 0;
1794   B->assembled              = PETSC_TRUE;
1795   B->preallocated           = PETSC_TRUE;
1796   PetscFunctionReturn(0);
1797 }
1798