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