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