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