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