xref: /petsc/src/mat/impls/baij/seq/baijfact.c (revision 6abefa4f446b2ea91915da4e25ea0defa83db48e)
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->shiftamount;
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->shiftamount;
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->shiftamount;
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->shiftamount;
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;
564   ChShift_Ctx    sctx;
565   PetscInt       newshift;
566 
567   PetscFunctionBegin;
568   if (bs > 1) {
569     if (!a->sbaijMat){
570       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
571     }
572     ierr = (a->sbaijMat)->ops->choleskyfactornumeric(C,a->sbaijMat,info);CHKERRQ(ierr);
573     ierr = MatDestroy(a->sbaijMat);CHKERRQ(ierr);
574     a->sbaijMat = PETSC_NULL;
575     PetscFunctionReturn(0);
576   }
577 
578   /* initialization */
579   zeropivot = info->zeropivot;
580 
581   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
582   ierr = PetscMalloc3(mbs,MatScalar,&rtmp,mbs,PetscInt,&il,mbs,PetscInt,&jl);CHKERRQ(ierr);
583 
584   sctx.shift_amount = 0.;
585   sctx.nshift       = 0;
586   do {
587     sctx.chshift = PETSC_FALSE;
588     for (i=0; i<mbs; i++) {
589       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
590     }
591 
592     for (k = 0; k<mbs; k++){
593       bval = ba + bi[k];
594       /* initialize k-th row by the perm[k]-th row of A */
595       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
596       for (j = jmin; j < jmax; j++){
597         col = rip[aj[j]];
598         if (col >= k){ /* only take upper triangular entry */
599           rtmp[col] = aa[j];
600           *bval++  = 0.0; /* for in-place factorization */
601         }
602       }
603 
604       /* shift the diagonal of the matrix */
605       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
606 
607       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
608       dk = rtmp[k];
609       i = jl[k]; /* first row to be added to k_th row  */
610 
611       while (i < k){
612         nexti = jl[i]; /* next row to be added to k_th row */
613 
614         /* compute multiplier, update diag(k) and U(i,k) */
615         ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
616         uikdi = - ba[ili]*ba[bi[i]];  /* diagonal(k) */
617         dk += uikdi*ba[ili];
618         ba[ili] = uikdi; /* -U(i,k) */
619 
620         /* add multiple of row i to k-th row */
621         jmin = ili + 1; jmax = bi[i+1];
622         if (jmin < jmax){
623           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
624           /* update il and jl for row i */
625           il[i] = jmin;
626           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
627         }
628         i = nexti;
629       }
630 
631       /* shift the diagonals when zero pivot is detected */
632       /* compute rs=sum of abs(off-diagonal) */
633       rs   = 0.0;
634       jmin = bi[k]+1;
635       nz   = bi[k+1] - jmin;
636       if (nz){
637         bcol = bj + jmin;
638         while (nz--){
639           rs += PetscAbsScalar(rtmp[*bcol]);
640           bcol++;
641         }
642       }
643 
644       sctx.rs = rs;
645       sctx.pv = dk;
646       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
647       if (newshift == 1) break;
648 
649       /* copy data into U(k,:) */
650       ba[bi[k]] = 1.0/dk; /* U(k,k) */
651       jmin = bi[k]+1; jmax = bi[k+1];
652       if (jmin < jmax) {
653         for (j=jmin; j<jmax; j++){
654           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
655         }
656         /* add the k-th row into il and jl */
657         il[k] = jmin;
658         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
659       }
660     }
661   } while (sctx.chshift);
662   ierr = PetscFree3(rtmp,il,jl);CHKERRQ(ierr);
663 
664   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
665   C->assembled    = PETSC_TRUE;
666   C->preallocated = PETSC_TRUE;
667   ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr);
668   if (sctx.nshift){
669     if (info->shifttype == MAT_SHIFT_POSITIVE_DEFINITE) {
670       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
671     } else if (info->shifttype == MAT_SHIFT_NONZERO) {
672       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
673     }
674   }
675   PetscFunctionReturn(0);
676 }
677 
678 #undef __FUNCT__
679 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering"
680 PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info)
681 {
682   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data;
683   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
684   PetscErrorCode ierr;
685   PetscInt       i,j,am=a->mbs;
686   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
687   PetscInt       k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz;
688   MatScalar      *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval;
689   PetscReal      zeropivot,rs;
690   ChShift_Ctx    sctx;
691   PetscInt       newshift;
692 
693   PetscFunctionBegin;
694   /* initialization */
695   zeropivot = info->zeropivot;
696 
697   ierr = PetscMalloc3(am,MatScalar,&rtmp,am,PetscInt,&il,am,PetscInt,&jl);CHKERRQ(ierr);
698 
699   sctx.shift_amount = 0.;
700   sctx.nshift       = 0;
701   do {
702     sctx.chshift = PETSC_FALSE;
703     for (i=0; i<am; i++) {
704       rtmp[i] = 0.0; jl[i] = am; il[0] = 0;
705     }
706 
707     for (k = 0; k<am; k++){
708     /* initialize k-th row with elements nonzero in row perm(k) of A */
709       nz   = ai[k+1] - ai[k];
710       acol = aj + ai[k];
711       aval = aa + ai[k];
712       bval = ba + bi[k];
713       while (nz -- ){
714         if (*acol < k) { /* skip lower triangular entries */
715           acol++; aval++;
716         } else {
717           rtmp[*acol++] = *aval++;
718           *bval++       = 0.0; /* for in-place factorization */
719         }
720       }
721 
722       /* shift the diagonal of the matrix */
723       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
724 
725       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
726       dk = rtmp[k];
727       i  = jl[k]; /* first row to be added to k_th row  */
728 
729       while (i < k){
730         nexti = jl[i]; /* next row to be added to k_th row */
731         /* compute multiplier, update D(k) and U(i,k) */
732         ili   = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
733         uikdi = - ba[ili]*ba[bi[i]];
734         dk   += uikdi*ba[ili];
735         ba[ili] = uikdi; /* -U(i,k) */
736 
737         /* add multiple of row i to k-th row ... */
738         jmin = ili + 1;
739         nz   = bi[i+1] - jmin;
740         if (nz > 0){
741           bcol = bj + jmin;
742           bval = ba + jmin;
743           while (nz --) rtmp[*bcol++] += uikdi*(*bval++);
744           /* update il and jl for i-th row */
745           il[i] = jmin;
746           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
747         }
748         i = nexti;
749       }
750 
751       /* shift the diagonals when zero pivot is detected */
752       /* compute rs=sum of abs(off-diagonal) */
753       rs   = 0.0;
754       jmin = bi[k]+1;
755       nz   = bi[k+1] - jmin;
756       if (nz){
757         bcol = bj + jmin;
758         while (nz--){
759           rs += PetscAbsScalar(rtmp[*bcol]);
760           bcol++;
761         }
762       }
763 
764       sctx.rs = rs;
765       sctx.pv = dk;
766       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
767       if (newshift == 1) break;    /* sctx.shift_amount is updated */
768 
769       /* copy data into U(k,:) */
770       ba[bi[k]] = 1.0/dk;
771       jmin      = bi[k]+1;
772       nz        = bi[k+1] - jmin;
773       if (nz){
774         bcol = bj + jmin;
775         bval = ba + jmin;
776         while (nz--){
777           *bval++       = rtmp[*bcol];
778           rtmp[*bcol++] = 0.0;
779         }
780         /* add k-th row into il and jl */
781         il[k] = jmin;
782         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
783       }
784     }
785   } while (sctx.chshift);
786   ierr = PetscFree3(rtmp,il,jl);CHKERRQ(ierr);
787 
788   C->ops->solve                 = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
789   C->ops->solvetranspose        = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
790   C->assembled    = PETSC_TRUE;
791   C->preallocated = PETSC_TRUE;
792   ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr);
793     if (sctx.nshift){
794     if (info->shifttype == MAT_SHIFT_NONZERO) {
795       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
796     } else if (info->shifttype == MAT_SHIFT_POSITIVE_DEFINITE) {
797       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
798     }
799   }
800   PetscFunctionReturn(0);
801 }
802 
803 #include "petscbt.h"
804 #include "../src/mat/utils/freespace.h"
805 #undef __FUNCT__
806 #define __FUNCT__ "MatICCFactorSymbolic_SeqBAIJ"
807 PetscErrorCode MatICCFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
808 {
809   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data;
810   Mat_SeqSBAIJ       *b;
811   Mat                B;
812   PetscErrorCode     ierr;
813   PetscTruth         perm_identity;
814   PetscInt           reallocs=0,i,*ai=a->i,*aj=a->j,am=a->mbs,bs=A->rmap->bs,*ui;
815   const PetscInt     *rip;
816   PetscInt           jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
817   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL,ncols,ncols_upper,*cols,*cols_lvl,*uj,**uj_ptr,**uj_lvl_ptr;
818   PetscReal          fill=info->fill,levels=info->levels;
819   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
820   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
821   PetscBT            lnkbt;
822 
823   PetscFunctionBegin;
824   if (bs > 1){
825     if (!a->sbaijMat){
826       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
827     }
828     (fact)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqSBAIJ;  /* undue the change made in MatGetFactor_seqbaij_petsc */
829     ierr = MatICCFactorSymbolic(fact,a->sbaijMat,perm,info);CHKERRQ(ierr);
830     PetscFunctionReturn(0);
831   }
832 
833   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
834   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
835 
836   /* special case that simply copies fill pattern */
837   if (!levels && perm_identity) {
838     ierr = MatMarkDiagonal_SeqBAIJ(A);CHKERRQ(ierr);
839     ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
840     for (i=0; i<am; i++) {
841       ui[i] = ai[i+1] - a->diag[i]; /* ui: rowlengths - changes when !perm_identity */
842     }
843     B = fact;
844     ierr = MatSeqSBAIJSetPreallocation(B,1,0,ui);CHKERRQ(ierr);
845 
846 
847     b  = (Mat_SeqSBAIJ*)B->data;
848     uj = b->j;
849     for (i=0; i<am; i++) {
850       aj = a->j + a->diag[i];
851       for (j=0; j<ui[i]; j++){
852         *uj++ = *aj++;
853       }
854       b->ilen[i] = ui[i];
855     }
856     ierr = PetscFree(ui);CHKERRQ(ierr);
857     B->factor = MAT_FACTOR_NONE;
858     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
859     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
860     B->factor = MAT_FACTOR_ICC;
861 
862     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
863     PetscFunctionReturn(0);
864   }
865 
866   /* initialization */
867   ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
868   ui[0] = 0;
869   ierr  = PetscMalloc((2*am+1)*sizeof(PetscInt),&cols_lvl);CHKERRQ(ierr);
870 
871   /* jl: linked list for storing indices of the pivot rows
872      il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
873   ierr = PetscMalloc4(am,PetscInt*,&uj_ptr,am,PetscInt*,&uj_lvl_ptr,am,PetscInt,&il,am,PetscInt,&jl);CHKERRQ(ierr);
874   for (i=0; i<am; i++){
875     jl[i] = am; il[i] = 0;
876   }
877 
878   /* create and initialize a linked list for storing column indices of the active row k */
879   nlnk = am + 1;
880   ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
881 
882   /* initial FreeSpace size is fill*(ai[am]+1) */
883   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
884   current_space = free_space;
885   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr);
886   current_space_lvl = free_space_lvl;
887 
888   for (k=0; k<am; k++){  /* for each active row k */
889     /* initialize lnk by the column indices of row rip[k] of A */
890     nzk   = 0;
891     ncols = ai[rip[k]+1] - ai[rip[k]];
892     ncols_upper = 0;
893     cols        = cols_lvl + am;
894     for (j=0; j<ncols; j++){
895       i = rip[*(aj + ai[rip[k]] + j)];
896       if (i >= k){ /* only take upper triangular entry */
897         cols[ncols_upper] = i;
898         cols_lvl[ncols_upper] = -1;  /* initialize level for nonzero entries */
899         ncols_upper++;
900       }
901     }
902     ierr = PetscIncompleteLLAdd(ncols_upper,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
903     nzk += nlnk;
904 
905     /* update lnk by computing fill-in for each pivot row to be merged in */
906     prow = jl[k]; /* 1st pivot row */
907 
908     while (prow < k){
909       nextprow = jl[prow];
910 
911       /* merge prow into k-th row */
912       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
913       jmax = ui[prow+1];
914       ncols = jmax-jmin;
915       i     = jmin - ui[prow];
916       cols = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
917       for (j=0; j<ncols; j++) cols_lvl[j] = *(uj_lvl_ptr[prow] + i + j);
918       ierr = PetscIncompleteLLAddSorted(ncols,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
919       nzk += nlnk;
920 
921       /* update il and jl for prow */
922       if (jmin < jmax){
923         il[prow] = jmin;
924         j = *cols; jl[prow] = jl[j]; jl[j] = prow;
925       }
926       prow = nextprow;
927     }
928 
929     /* if free space is not available, make more free space */
930     if (current_space->local_remaining<nzk) {
931       i = am - k + 1; /* num of unfactored rows */
932       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
933       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
934       ierr = PetscFreeSpaceGet(i,&current_space_lvl);CHKERRQ(ierr);
935       reallocs++;
936     }
937 
938     /* copy data into free_space and free_space_lvl, then initialize lnk */
939     ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
940 
941     /* add the k-th row into il and jl */
942     if (nzk-1 > 0){
943       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
944       jl[k] = jl[i]; jl[i] = k;
945       il[k] = ui[k] + 1;
946     }
947     uj_ptr[k]     = current_space->array;
948     uj_lvl_ptr[k] = current_space_lvl->array;
949 
950     current_space->array           += nzk;
951     current_space->local_used      += nzk;
952     current_space->local_remaining -= nzk;
953 
954     current_space_lvl->array           += nzk;
955     current_space_lvl->local_used      += nzk;
956     current_space_lvl->local_remaining -= nzk;
957 
958     ui[k+1] = ui[k] + nzk;
959   }
960 
961 #if defined(PETSC_USE_INFO)
962   if (ai[am] != 0) {
963     PetscReal af = ((PetscReal)(2*ui[am]-am))/((PetscReal)ai[am]);
964     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
965     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
966     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
967   } else {
968     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
969   }
970 #endif
971 
972   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
973   ierr = PetscFree4(uj_ptr,uj_lvl_ptr,il,jl);CHKERRQ(ierr);
974   ierr = PetscFree(cols_lvl);CHKERRQ(ierr);
975 
976   /* destroy list of free space and other temporary array(s) */
977   ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
978   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
979   ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
980   ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
981 
982   /* put together the new matrix in MATSEQSBAIJ format */
983   B = fact;
984   ierr = MatSeqSBAIJSetPreallocation(B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
985 
986   b = (Mat_SeqSBAIJ*)B->data;
987   b->singlemalloc = PETSC_FALSE;
988   b->free_a       = PETSC_TRUE;
989   b->free_ij       = PETSC_TRUE;
990   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
991   b->j    = uj;
992   b->i    = ui;
993   b->diag = 0;
994   b->ilen = 0;
995   b->imax = 0;
996   b->row  = perm;
997   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
998   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
999   b->icol = perm;
1000   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
1001   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1002   ierr    = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
1003   b->maxnz = b->nz = ui[am];
1004 
1005   B->info.factor_mallocs    = reallocs;
1006   B->info.fill_ratio_given  = fill;
1007   if (ai[am] != 0.) {
1008     B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
1009   } else {
1010     B->info.fill_ratio_needed = 0.0;
1011   }
1012   if (perm_identity){
1013     B->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
1014     B->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
1015     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
1016   } else {
1017     (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
1018   }
1019   PetscFunctionReturn(0);
1020 }
1021 
1022 #undef __FUNCT__
1023 #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqBAIJ"
1024 PetscErrorCode MatCholeskyFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
1025 {
1026   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data;
1027   Mat_SeqSBAIJ       *b;
1028   Mat                B;
1029   PetscErrorCode     ierr;
1030   PetscTruth         perm_identity;
1031   PetscReal          fill = info->fill;
1032   const PetscInt     *rip;
1033   PetscInt           i,mbs=a->mbs,bs=A->rmap->bs,*ai=a->i,*aj=a->j,reallocs=0,prow;
1034   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
1035   PetscInt           nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr;
1036   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
1037   PetscBT            lnkbt;
1038 
1039   PetscFunctionBegin;
1040   if (bs > 1) { /* convert to seqsbaij */
1041     if (!a->sbaijMat){
1042       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
1043     }
1044     (fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqSBAIJ; /* undue the change made in MatGetFactor_seqbaij_petsc */
1045     ierr = MatCholeskyFactorSymbolic(fact,a->sbaijMat,perm,info);CHKERRQ(ierr);
1046     PetscFunctionReturn(0);
1047   }
1048 
1049   /* check whether perm is the identity mapping */
1050   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
1051   if (!perm_identity) SETERRQ(PETSC_ERR_SUP,"Matrix reordering is not supported");
1052   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
1053 
1054   /* initialization */
1055   ierr  = PetscMalloc((mbs+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
1056   ui[0] = 0;
1057 
1058   /* jl: linked list for storing indices of the pivot rows
1059      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
1060   ierr = PetscMalloc4(mbs,PetscInt*,&ui_ptr,mbs,PetscInt,&il,mbs,PetscInt,&jl,mbs,PetscInt,&cols);CHKERRQ(ierr);
1061   for (i=0; i<mbs; i++){
1062     jl[i] = mbs; il[i] = 0;
1063   }
1064 
1065   /* create and initialize a linked list for storing column indices of the active row k */
1066   nlnk = mbs + 1;
1067   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
1068 
1069   /* initial FreeSpace size is fill*(ai[mbs]+1) */
1070   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[mbs]+1)),&free_space);CHKERRQ(ierr);
1071   current_space = free_space;
1072 
1073   for (k=0; k<mbs; k++){  /* for each active row k */
1074     /* initialize lnk by the column indices of row rip[k] of A */
1075     nzk   = 0;
1076     ncols = ai[rip[k]+1] - ai[rip[k]];
1077     ncols_upper = 0;
1078     for (j=0; j<ncols; j++){
1079       i = rip[*(aj + ai[rip[k]] + j)];
1080       if (i >= k){ /* only take upper triangular entry */
1081         cols[ncols_upper] = i;
1082         ncols_upper++;
1083       }
1084     }
1085     ierr = PetscLLAdd(ncols_upper,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
1086     nzk += nlnk;
1087 
1088     /* update lnk by computing fill-in for each pivot row to be merged in */
1089     prow = jl[k]; /* 1st pivot row */
1090 
1091     while (prow < k){
1092       nextprow = jl[prow];
1093       /* merge prow into k-th row */
1094       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
1095       jmax = ui[prow+1];
1096       ncols = jmax-jmin;
1097       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
1098       ierr = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
1099       nzk += nlnk;
1100 
1101       /* update il and jl for prow */
1102       if (jmin < jmax){
1103         il[prow] = jmin;
1104         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
1105       }
1106       prow = nextprow;
1107     }
1108 
1109     /* if free space is not available, make more free space */
1110     if (current_space->local_remaining<nzk) {
1111       i = mbs - k + 1; /* num of unfactored rows */
1112       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
1113       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
1114       reallocs++;
1115     }
1116 
1117     /* copy data into free space, then initialize lnk */
1118     ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
1119 
1120     /* add the k-th row into il and jl */
1121     if (nzk-1 > 0){
1122       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
1123       jl[k] = jl[i]; jl[i] = k;
1124       il[k] = ui[k] + 1;
1125     }
1126     ui_ptr[k] = current_space->array;
1127     current_space->array           += nzk;
1128     current_space->local_used      += nzk;
1129     current_space->local_remaining -= nzk;
1130 
1131     ui[k+1] = ui[k] + nzk;
1132   }
1133 
1134 #if defined(PETSC_USE_INFO)
1135   if (ai[mbs] != 0.) {
1136     PetscReal af = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
1137     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
1138     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
1139     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
1140   } else {
1141     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
1142   }
1143 #endif
1144 
1145   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
1146   ierr = PetscFree4(ui_ptr,il,jl,cols);CHKERRQ(ierr);
1147 
1148   /* destroy list of free space and other temporary array(s) */
1149   ierr = PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
1150   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
1151   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1152 
1153   /* put together the new matrix in MATSEQSBAIJ format */
1154   B    = fact;
1155   ierr = MatSeqSBAIJSetPreallocation(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
1156 
1157   b = (Mat_SeqSBAIJ*)B->data;
1158   b->singlemalloc = PETSC_FALSE;
1159   b->free_a       = PETSC_TRUE;
1160   b->free_ij      = PETSC_TRUE;
1161   ierr = PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
1162   b->j    = uj;
1163   b->i    = ui;
1164   b->diag = 0;
1165   b->ilen = 0;
1166   b->imax = 0;
1167   b->row  = perm;
1168   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
1169   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
1170   b->icol = perm;
1171   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
1172   ierr    = PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1173   ierr    = PetscLogObjectMemory(B,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
1174   b->maxnz = b->nz = ui[mbs];
1175 
1176   B->info.factor_mallocs    = reallocs;
1177   B->info.fill_ratio_given  = fill;
1178   if (ai[mbs] != 0.) {
1179     B->info.fill_ratio_needed = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
1180   } else {
1181     B->info.fill_ratio_needed = 0.0;
1182   }
1183   if (perm_identity){
1184     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
1185   } else {
1186     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
1187   }
1188   PetscFunctionReturn(0);
1189 }
1190 
1191 #undef __FUNCT__
1192 #define __FUNCT__ "MatSolve_SeqBAIJ_N_NaturalOrdering"
1193 PetscErrorCode MatSolve_SeqBAIJ_N_NaturalOrdering(Mat A,Vec bb,Vec xx)
1194 {
1195   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ *)A->data;
1196   PetscErrorCode ierr;
1197   const PetscInt *ai=a->i,*aj=a->j,*adiag=a->diag,*vi;
1198   PetscInt       i,k,n=a->mbs;
1199   PetscInt       nz,bs=A->rmap->bs,bs2=a->bs2;
1200   MatScalar      *aa=a->a,*v;
1201   PetscScalar    *x,*b,*s,*t,*ls;
1202 
1203   PetscFunctionBegin;
1204   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
1205   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1206   t  = a->solve_work;
1207 
1208   /* forward solve the lower triangular */
1209   ierr = PetscMemcpy(t,b,bs*sizeof(PetscScalar));CHKERRQ(ierr); /* copy 1st block of b to t */
1210 
1211   for (i=1; i<n; i++) {
1212     v   = aa + bs2*ai[i];
1213     vi  = aj + ai[i];
1214     nz = ai[i+1] - ai[i];
1215     s = t + bs*i;
1216     ierr = PetscMemcpy(s,b+bs*i,bs*sizeof(PetscScalar));CHKERRQ(ierr); /* copy i_th block of b to t */
1217     for(k=0;k<nz;k++){
1218       Kernel_v_gets_v_minus_A_times_w(bs,s,v,t+bs*vi[k]);
1219       v += bs2;
1220     }
1221   }
1222 
1223   /* backward solve the upper triangular */
1224   ls = a->solve_work + A->cmap->n;
1225   for (i=n-1; i>=0; i--){
1226     v  = aa + bs2*(adiag[i+1]+1);
1227     vi = aj + adiag[i+1]+1;
1228     nz = adiag[i] - adiag[i+1]-1;
1229     ierr = PetscMemcpy(ls,t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr);
1230     for(k=0;k<nz;k++){
1231       Kernel_v_gets_v_minus_A_times_w(bs,ls,v,t+bs*vi[k]);
1232       v += bs2;
1233     }
1234     Kernel_w_gets_A_times_v(bs,ls,aa+bs2*adiag[i],t+i*bs); /* *inv(diagonal[i]) */
1235     ierr = PetscMemcpy(x+i*bs,t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr);
1236   }
1237 
1238   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
1239   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1240   ierr = PetscLogFlops(2.0*(a->bs2)*(a->nz) - A->rmap->bs*A->cmap->n);CHKERRQ(ierr);
1241   PetscFunctionReturn(0);
1242 }
1243 
1244 #undef __FUNCT__
1245 #define __FUNCT__ "MatSolve_SeqBAIJ_N"
1246 PetscErrorCode MatSolve_SeqBAIJ_N(Mat A,Vec bb,Vec xx)
1247 {
1248   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ *)A->data;
1249   IS             iscol=a->col,isrow=a->row;
1250   PetscErrorCode ierr;
1251   const PetscInt *r,*c,*rout,*cout,*ai=a->i,*aj=a->j,*adiag=a->diag,*vi;
1252   PetscInt       i,m,n=a->mbs;
1253   PetscInt       nz,bs=A->rmap->bs,bs2=a->bs2;
1254   MatScalar      *aa=a->a,*v;
1255   PetscScalar    *x,*b,*s,*t,*ls;
1256 
1257   PetscFunctionBegin;
1258   ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
1259   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1260   t  = a->solve_work;
1261 
1262   ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout;
1263   ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout;
1264 
1265   /* forward solve the lower triangular */
1266   ierr = PetscMemcpy(t,b+bs*r[0],bs*sizeof(PetscScalar));CHKERRQ(ierr);
1267   for (i=1; i<n; i++) {
1268     v   = aa + bs2*ai[i];
1269     vi  = aj + ai[i];
1270     nz = ai[i+1] - ai[i];
1271     s = t + bs*i;
1272     ierr = PetscMemcpy(s,b+bs*r[i],bs*sizeof(PetscScalar));CHKERRQ(ierr);
1273     for(m=0;m<nz;m++){
1274       Kernel_v_gets_v_minus_A_times_w(bs,s,v,t+bs*vi[m]);
1275       v += bs2;
1276     }
1277   }
1278 
1279   /* backward solve the upper triangular */
1280   ls = a->solve_work + A->cmap->n;
1281   for (i=n-1; i>=0; i--){
1282     v  = aa + bs2*(adiag[i+1]+1);
1283     vi = aj + adiag[i+1]+1;
1284     nz = adiag[i] - adiag[i+1] - 1;
1285     ierr = PetscMemcpy(ls,t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr);
1286     for(m=0;m<nz;m++){
1287       Kernel_v_gets_v_minus_A_times_w(bs,ls,v,t+bs*vi[m]);
1288       v += bs2;
1289     }
1290     Kernel_w_gets_A_times_v(bs,ls,v,t+i*bs); /* *inv(diagonal[i]) */
1291     ierr = PetscMemcpy(x + bs*c[i],t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr);
1292   }
1293   ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr);
1294   ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr);
1295   ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
1296   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1297   ierr = PetscLogFlops(2.0*(a->bs2)*(a->nz) - A->rmap->bs*A->cmap->n);CHKERRQ(ierr);
1298   PetscFunctionReturn(0);
1299 }
1300 
1301 #undef __FUNCT__
1302 #define __FUNCT__ "BlockAbs_privat"
1303 PetscErrorCode BlockAbs_private(PetscInt nbs,PetscInt bs2,PetscScalar *blockarray,PetscReal *absarray)
1304 {
1305   PetscErrorCode     ierr;
1306   PetscInt           i,j;
1307   PetscFunctionBegin;
1308   ierr = PetscMemzero(absarray,(nbs+1)*sizeof(PetscReal));CHKERRQ(ierr);
1309   for (i=0; i<nbs; i++){
1310     for (j=0; j<bs2; j++){
1311       if (absarray[i] < PetscAbsScalar(blockarray[i*nbs+j])) absarray[i] = PetscAbsScalar(blockarray[i*nbs+j]);
1312     }
1313   }
1314   PetscFunctionReturn(0);
1315 }
1316 
1317 #undef __FUNCT__
1318 #define __FUNCT__ "MatILUDTFactor_SeqBAIJ"
1319 /*
1320      This needs to be renamed and called by the regular MatILUFactor_SeqBAIJ when drop tolerance is used
1321 */
1322 PetscErrorCode MatILUDTFactor_SeqBAIJ(Mat A,IS isrow,IS iscol,const MatFactorInfo *info,Mat *fact)
1323 {
1324   Mat                B = *fact;
1325   Mat_SeqBAIJ        *a=(Mat_SeqBAIJ*)A->data,*b;
1326   IS                 isicol;
1327   PetscErrorCode     ierr;
1328   const PetscInt     *r,*ic;
1329   PetscInt           i,mbs=a->mbs,bs=A->rmap->bs,bs2=a->bs2,*ai=a->i,*aj=a->j,*ajtmp,*adiag;
1330   PetscInt           *bi,*bj,*bdiag;
1331 
1332   PetscInt           row,nzi,nzi_bl,nzi_bu,*im,dtcount,nzi_al,nzi_au;
1333   PetscInt           nlnk,*lnk;
1334   PetscBT            lnkbt;
1335   PetscTruth         row_identity,icol_identity,both_identity;
1336   MatScalar          *aatmp,*pv,*batmp,*ba,*rtmp,*pc,*multiplier,*vtmp;
1337   const PetscInt     *ics;
1338   PetscInt           j,nz,*pj,*bjtmp,k,ncut,*jtmp;
1339 
1340   PetscReal          dt=info->dt; /* shift=info->shiftamount; */
1341   PetscInt           nnz_max;
1342   PetscTruth         missing;
1343   PetscReal          *vtmp_abs;
1344   MatScalar          *v_work;
1345   PetscInt           *v_pivots;
1346 
1347   PetscFunctionBegin;
1348   /* ------- symbolic factorization, can be reused ---------*/
1349   ierr = MatMissingDiagonal(A,&missing,&i);CHKERRQ(ierr);
1350   if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i);
1351   adiag=a->diag;
1352 
1353   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
1354 
1355   /* bdiag is location of diagonal in factor */
1356   ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr);
1357 
1358   /* allocate row pointers bi */
1359   ierr = PetscMalloc((2*mbs+2)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
1360 
1361   /* allocate bj and ba; max num of nonzero entries is (ai[n]+2*n*dtcount+2) */
1362   dtcount = (PetscInt)info->dtcount;
1363   if (dtcount > mbs-1) dtcount = mbs-1;
1364   nnz_max  = ai[mbs]+2*mbs*dtcount +2;
1365   /* printf("MatILUDTFactor_SeqBAIJ, bs %d, ai[mbs] %d, nnz_max  %d, dtcount %d\n",bs,ai[mbs],nnz_max,dtcount); */
1366   ierr = PetscMalloc(nnz_max*sizeof(PetscInt),&bj);CHKERRQ(ierr);
1367   nnz_max = nnz_max*bs2;
1368   ierr = PetscMalloc(nnz_max*sizeof(MatScalar),&ba);CHKERRQ(ierr);
1369 
1370   /* put together the new matrix */
1371   ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
1372   ierr = PetscLogObjectParent(B,isicol);CHKERRQ(ierr);
1373   b    = (Mat_SeqBAIJ*)(B)->data;
1374   b->free_a       = PETSC_TRUE;
1375   b->free_ij      = PETSC_TRUE;
1376   b->singlemalloc = PETSC_FALSE;
1377   b->a          = ba;
1378   b->j          = bj;
1379   b->i          = bi;
1380   b->diag       = bdiag;
1381   b->ilen       = 0;
1382   b->imax       = 0;
1383   b->row        = isrow;
1384   b->col        = iscol;
1385   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
1386   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
1387   b->icol       = isicol;
1388   ierr = PetscMalloc((bs*(mbs+1))*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
1389 
1390   ierr = PetscLogObjectMemory(B,nnz_max*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
1391   b->maxnz = nnz_max/bs2;
1392 
1393   (B)->factor                = MAT_FACTOR_ILUDT;
1394   (B)->info.factor_mallocs   = 0;
1395   (B)->info.fill_ratio_given = ((PetscReal)nnz_max)/((PetscReal)(ai[mbs]*bs2));
1396   CHKMEMQ;
1397   /* ------- end of symbolic factorization ---------*/
1398   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
1399   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
1400   ics  = ic;
1401 
1402   /* linked list for storing column indices of the active row */
1403   nlnk = mbs + 1;
1404   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
1405 
1406   /* im: used by PetscLLAddSortedLU(); jtmp: working array for column indices of active row */
1407   ierr = PetscMalloc2(mbs,PetscInt,&im,mbs,PetscInt,&jtmp);CHKERRQ(ierr);
1408   /* rtmp, vtmp: working arrays for sparse and contiguous row entries of active row */
1409   ierr = PetscMalloc2(mbs*bs2,MatScalar,&rtmp,mbs*bs2,MatScalar,&vtmp);CHKERRQ(ierr);
1410   ierr = PetscMalloc((mbs+1)*sizeof(PetscReal),&vtmp_abs);CHKERRQ(ierr);
1411   ierr = PetscMalloc3(bs,MatScalar,&v_work,bs2,MatScalar,&multiplier,bs,PetscInt,&v_pivots);CHKERRQ(ierr);
1412 
1413   bi[0]    = 0;
1414   bdiag[0] = (nnz_max/bs2)-1; /* location of diagonal in factor B */
1415   bi[2*mbs+1] = bdiag[0]+1; /* endof bj and ba array */
1416   for (i=0; i<mbs; i++) {
1417     /* copy initial fill into linked list */
1418     nzi = 0; /* nonzeros for active row i */
1419     nzi = ai[r[i]+1] - ai[r[i]];
1420     if (!nzi) SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i);
1421     nzi_al = adiag[r[i]] - ai[r[i]];
1422     nzi_au = ai[r[i]+1] - adiag[r[i]] -1;
1423     /* printf("row %d, nzi_al/au %d %d\n",i,nzi_al,nzi_au); */
1424 
1425     /* load in initial unfactored row */
1426     ajtmp = aj + ai[r[i]];
1427     ierr = PetscLLAddPerm(nzi,ajtmp,ic,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
1428     ierr = PetscMemzero(rtmp,mbs*bs2*sizeof(PetscScalar));CHKERRQ(ierr);
1429     aatmp = a->a + bs2*ai[r[i]];
1430     for (j=0; j<nzi; j++) {
1431       ierr = PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],aatmp+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1432     }
1433 
1434     /* add pivot rows into linked list */
1435     row = lnk[mbs];
1436     while (row < i) {
1437       nzi_bl = bi[row+1] - bi[row] + 1;
1438       bjtmp = bj + bdiag[row+1]+1; /* points to 1st column next to the diagonal in U */
1439       ierr  = PetscLLAddSortedLU(bjtmp,row,nlnk,lnk,lnkbt,i,nzi_bl,im);CHKERRQ(ierr);
1440       nzi  += nlnk;
1441       row   = lnk[row];
1442     }
1443 
1444     /* copy data from lnk into jtmp, then initialize lnk */
1445     ierr = PetscLLClean(mbs,mbs,nzi,lnk,jtmp,lnkbt);CHKERRQ(ierr);
1446 
1447     /* numerical factorization */
1448     bjtmp = jtmp;
1449     row   = *bjtmp++; /* 1st pivot row */
1450 
1451     while  (row < i) {
1452       pc = rtmp + bs2*row;
1453       pv = ba + bs2*bdiag[row]; /* inv(diag) of the pivot row */
1454       Kernel_A_gets_A_times_B(bs,pc,pv,multiplier); /* pc= multiplier = pc*inv(diag[row]) */
1455       ierr = BlockAbs_private(1,bs2,pc,vtmp_abs);CHKERRQ(ierr);
1456       if (vtmp_abs[0] > dt){ /* apply tolerance dropping rule */
1457         pj         = bj + bdiag[row+1] + 1; /* point to 1st entry of U(row,:) */
1458         pv         = ba + bs2*(bdiag[row+1] + 1);
1459         nz         = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:), excluding diagonal */
1460         for (j=0; j<nz; j++){
1461           Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j);
1462         }
1463         /* ierr = PetscLogFlops(bslog*(nz+1.0)-bs);CHKERRQ(ierr); */
1464       }
1465       row = *bjtmp++;
1466     }
1467 
1468     /* copy sparse rtmp into contiguous vtmp; separate L and U part */
1469     nzi_bl = 0; j = 0;
1470     while (jtmp[j] < i){ /* L-part. Note: jtmp is sorted */
1471       ierr = PetscMemcpy(vtmp+bs2*j,rtmp+bs2*jtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
1472       nzi_bl++; j++;
1473     }
1474     nzi_bu = nzi - nzi_bl -1;
1475     /* printf("nzi %d, nzi_bl %d, nzi_bu %d\n",nzi,nzi_bl,nzi_bu); */
1476 
1477     while (j < nzi){ /* U-part */
1478       ierr = PetscMemcpy(vtmp+bs2*j,rtmp+bs2*jtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
1479       /*
1480       printf(" col %d: ",jtmp[j]);
1481       for (j1=0; j1<bs2; j1++) printf(" %g",*(vtmp+bs2*j+j1));
1482       printf(" \n");
1483       */
1484       j++;
1485     }
1486 
1487     ierr = BlockAbs_private(nzi,bs2,vtmp,vtmp_abs);CHKERRQ(ierr);
1488     /*
1489     printf(" row %d, nzi %d, vtmp_abs\n",i,nzi);
1490     for (j1=0; j1<nzi; j1++) printf(" (%d %g),",jtmp[j1],vtmp_abs[j1]);
1491     printf(" \n");
1492     */
1493     bjtmp = bj + bi[i];
1494     batmp = ba + bs2*bi[i];
1495     /* apply level dropping rule to L part */
1496     ncut = nzi_al + dtcount;
1497     if (ncut < nzi_bl){
1498       ierr = PetscSortSplitReal(ncut,nzi_bl,vtmp_abs,jtmp);CHKERRQ(ierr);
1499       ierr = PetscSortIntWithScalarArray(ncut,jtmp,vtmp);CHKERRQ(ierr);
1500     } else {
1501       ncut = nzi_bl;
1502     }
1503     for (j=0; j<ncut; j++){
1504       bjtmp[j] = jtmp[j];
1505       ierr = PetscMemcpy(batmp+bs2*j,rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
1506       /*
1507       printf(" col %d: ",bjtmp[j]);
1508       for (j1=0; j1<bs2; j1++) printf(" %g,",*(batmp+bs2*j+j1));
1509       printf("\n");
1510       */
1511     }
1512     bi[i+1] = bi[i] + ncut;
1513     nzi = ncut + 1;
1514 
1515     /* apply level dropping rule to U part */
1516     ncut = nzi_au + dtcount;
1517     if (ncut < nzi_bu){
1518       ierr = PetscSortSplitReal(ncut,nzi_bu,vtmp_abs+nzi_bl+1,jtmp+nzi_bl+1);CHKERRQ(ierr);
1519       ierr = PetscSortIntWithScalarArray(ncut,jtmp+nzi_bl+1,vtmp+nzi_bl+1);CHKERRQ(ierr);
1520     } else {
1521       ncut = nzi_bu;
1522     }
1523     nzi += ncut;
1524 
1525     /* mark bdiagonal */
1526     bdiag[i+1]    = bdiag[i] - (ncut + 1);
1527     bi[2*mbs - i] = bi[2*mbs - i +1] - (ncut + 1);
1528 
1529     bjtmp = bj + bdiag[i];
1530     batmp = ba + bs2*bdiag[i];
1531     ierr = PetscMemcpy(batmp,rtmp+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
1532     *bjtmp = i;
1533     /*
1534     printf(" diag %d: ",*bjtmp);
1535     for (j=0; j<bs2; j++){
1536       printf(" %g,",batmp[j]);
1537     }
1538     printf("\n");
1539     */
1540     bjtmp = bj + bdiag[i+1]+1;
1541     batmp = ba + (bdiag[i+1]+1)*bs2;
1542 
1543     for (k=0; k<ncut; k++){
1544       bjtmp[k] = jtmp[nzi_bl+1+k];
1545       ierr = PetscMemcpy(batmp+bs2*k,rtmp+bs2*bjtmp[k],bs2*sizeof(MatScalar));CHKERRQ(ierr);
1546       /*
1547       printf(" col %d:",bjtmp[k]);
1548       for (j1=0; j1<bs2; j1++) printf(" %g,",*(batmp+bs2*k+j1));
1549       printf("\n");
1550       */
1551     }
1552 
1553     im[i] = nzi; /* used by PetscLLAddSortedLU() */
1554 
1555     /* invert diagonal block for simplier triangular solves - add shift??? */
1556     batmp = ba + bs2*bdiag[i];
1557     ierr = Kernel_A_gets_inverse_A(bs,batmp,v_pivots,v_work);CHKERRQ(ierr);
1558   } /* for (i=0; i<mbs; i++) */
1559   ierr = PetscFree3(v_work,multiplier,v_pivots);CHKERRQ(ierr);
1560 
1561   /* printf("end of L %d, beginning of U %d\n",bi[mbs],bdiag[mbs]); */
1562   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]);
1563 
1564   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
1565   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
1566 
1567   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1568 
1569   ierr = PetscFree2(im,jtmp);CHKERRQ(ierr);
1570   ierr = PetscFree2(rtmp,vtmp);CHKERRQ(ierr);
1571 
1572   ierr = PetscLogFlops(bs2*B->cmap->n);CHKERRQ(ierr);
1573   b->maxnz = b->nz = bi[mbs] + bdiag[0] - bdiag[mbs];
1574 
1575   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
1576   ierr = ISIdentity(isicol,&icol_identity);CHKERRQ(ierr);
1577   both_identity = (PetscTruth) (row_identity && icol_identity);
1578   if (row_identity && icol_identity) {
1579     B->ops->solve = MatSolve_SeqBAIJ_N_NaturalOrdering;
1580   } else {
1581     B->ops->solve = MatSolve_SeqBAIJ_N;
1582   }
1583 
1584   B->ops->solveadd          = 0;
1585   B->ops->solvetranspose    = 0;
1586   B->ops->solvetransposeadd = 0;
1587   B->ops->matsolve          = 0;
1588   B->assembled              = PETSC_TRUE;
1589   B->preallocated           = PETSC_TRUE;
1590   PetscFunctionReturn(0);
1591 }
1592