xref: /petsc/src/mat/impls/baij/seq/baijfact2.c (revision 4f02bc6a7df4e4b03c783003a8e0899ee35fcb05)
1 
2 /*
3     Factorization code for BAIJ format.
4 */
5 
6 #include <../src/mat/impls/baij/seq/baij.h>
7 #include <petsc/private/kernels/blockinvert.h>
8 #include <petscbt.h>
9 #include <../src/mat/utils/freespace.h>
10 
11 /* ----------------------------------------------------------------*/
12 extern PetscErrorCode MatDuplicateNoCreate_SeqBAIJ(Mat,Mat,MatDuplicateOption,PetscBool);
13 
14 #undef __FUNCT__
15 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_15_NaturalOrdering"
16 /*
17    This is not much faster than MatLUFactorNumeric_SeqBAIJ_N() but the solve is faster at least sometimes
18 */
19 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_15_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info)
20 {
21   Mat             C =B;
22   Mat_SeqBAIJ     *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data;
23   PetscErrorCode  ierr;
24   PetscInt        i,j,k,ipvt[15];
25   const PetscInt  n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ajtmp,*bjtmp,*bdiag=b->diag,*pj;
26   PetscInt        nz,nzL,row;
27   MatScalar       *rtmp,*pc,*mwork,*pv,*vv,work[225];
28   const MatScalar *v,*aa=a->a;
29   PetscInt        bs2 = a->bs2,bs=A->rmap->bs,flg;
30   PetscInt        sol_ver;
31 
32   PetscFunctionBegin;
33   ierr = PetscOptionsGetInt(NULL,((PetscObject)A)->prefix,"-sol_ver",&sol_ver,NULL);CHKERRQ(ierr);
34 
35   /* generate work space needed by the factorization */
36   ierr = PetscMalloc2(bs2*n,&rtmp,bs2,&mwork);CHKERRQ(ierr);
37   ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr);
38 
39   for (i=0; i<n; i++) {
40     /* zero rtmp */
41     /* L part */
42     nz    = bi[i+1] - bi[i];
43     bjtmp = bj + bi[i];
44     for  (j=0; j<nz; j++) {
45       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
46     }
47 
48     /* U part */
49     nz    = bdiag[i] - bdiag[i+1];
50     bjtmp = bj + bdiag[i+1]+1;
51     for  (j=0; j<nz; j++) {
52       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
53     }
54 
55     /* load in initial (unfactored row) */
56     nz    = ai[i+1] - ai[i];
57     ajtmp = aj + ai[i];
58     v     = aa + bs2*ai[i];
59     for (j=0; j<nz; j++) {
60       ierr = PetscMemcpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr);
61     }
62 
63     /* elimination */
64     bjtmp = bj + bi[i];
65     nzL   = bi[i+1] - bi[i];
66     for (k=0; k < nzL; k++) {
67       row = bjtmp[k];
68       pc  = rtmp + bs2*row;
69       for (flg=0,j=0; j<bs2; j++) {
70         if (pc[j]!=0.0) {
71           flg = 1;
72           break;
73         }
74       }
75       if (flg) {
76         pv = b->a + bs2*bdiag[row];
77         PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork);
78         /*ierr = PetscKernel_A_gets_A_times_B_15(pc,pv,mwork);CHKERRQ(ierr);*/
79         pj = b->j + bdiag[row+1]+1; /* begining of U(row,:) */
80         pv = b->a + bs2*(bdiag[row+1]+1);
81         nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries inU(row,:), excluding diag */
82         for (j=0; j<nz; j++) {
83           vv = rtmp + bs2*pj[j];
84           PetscKernel_A_gets_A_minus_B_times_C(bs,vv,pc,pv);
85           /* ierr = PetscKernel_A_gets_A_minus_B_times_C_15(vv,pc,pv);CHKERRQ(ierr); */
86           pv += bs2;
87         }
88         ierr = PetscLogFlops(2*bs2*bs*(nz+1)-bs2);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
89       }
90     }
91 
92     /* finished row so stick it into b->a */
93     /* L part */
94     pv = b->a + bs2*bi[i];
95     pj = b->j + bi[i];
96     nz = bi[i+1] - bi[i];
97     for (j=0; j<nz; j++) {
98       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
99     }
100 
101     /* Mark diagonal and invert diagonal for simplier triangular solves */
102     pv   = b->a + bs2*bdiag[i];
103     pj   = b->j + bdiag[i];
104     ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr);
105     /* PetscKernel_A_gets_inverse_A(bs,pv,pivots,work); */
106     ierr = PetscKernel_A_gets_inverse_A_15(pv,ipvt,work,info->shiftamount);CHKERRQ(ierr);
107 
108     /* U part */
109     pv = b->a + bs2*(bdiag[i+1]+1);
110     pj = b->j + bdiag[i+1]+1;
111     nz = bdiag[i] - bdiag[i+1] - 1;
112     for (j=0; j<nz; j++) {
113       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
114     }
115   }
116 
117   ierr = PetscFree2(rtmp,mwork);CHKERRQ(ierr);
118 
119   C->ops->solve          = MatSolve_SeqBAIJ_15_NaturalOrdering_ver1;
120   C->ops->solvetranspose = MatSolve_SeqBAIJ_N_NaturalOrdering;
121   C->assembled           = PETSC_TRUE;
122 
123   ierr = PetscLogFlops(1.333333333333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
124   PetscFunctionReturn(0);
125 }
126 
127 #undef __FUNCT__
128 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_N"
129 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_N(Mat B,Mat A,const MatFactorInfo *info)
130 {
131   Mat            C     =B;
132   Mat_SeqBAIJ    *a    =(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data;
133   IS             isrow = b->row,isicol = b->icol;
134   PetscErrorCode ierr;
135   const PetscInt *r,*ic;
136   PetscInt       i,j,k,n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
137   PetscInt       *ajtmp,*bjtmp,nz,nzL,row,*bdiag=b->diag,*pj;
138   MatScalar      *rtmp,*pc,*mwork,*v,*pv,*aa=a->a;
139   PetscInt       bs=A->rmap->bs,bs2 = a->bs2,*v_pivots,flg;
140   MatScalar      *v_work;
141   PetscBool      col_identity,row_identity,both_identity;
142 
143   PetscFunctionBegin;
144   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
145   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
146 
147   ierr = PetscMalloc1(bs2*n,&rtmp);CHKERRQ(ierr);
148   ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr);
149 
150   /* generate work space needed by dense LU factorization */
151   ierr = PetscMalloc3(bs,&v_work,bs2,&mwork,bs,&v_pivots);CHKERRQ(ierr);
152 
153   for (i=0; i<n; i++) {
154     /* zero rtmp */
155     /* L part */
156     nz    = bi[i+1] - bi[i];
157     bjtmp = bj + bi[i];
158     for  (j=0; j<nz; j++) {
159       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
160     }
161 
162     /* U part */
163     nz    = bdiag[i] - bdiag[i+1];
164     bjtmp = bj + bdiag[i+1]+1;
165     for  (j=0; j<nz; j++) {
166       ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
167     }
168 
169     /* load in initial (unfactored row) */
170     nz    = ai[r[i]+1] - ai[r[i]];
171     ajtmp = aj + ai[r[i]];
172     v     = aa + bs2*ai[r[i]];
173     for (j=0; j<nz; j++) {
174       ierr = PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr);
175     }
176 
177     /* elimination */
178     bjtmp = bj + bi[i];
179     nzL   = bi[i+1] - bi[i];
180     for (k=0; k < nzL; k++) {
181       row = bjtmp[k];
182       pc  = rtmp + bs2*row;
183       for (flg=0,j=0; j<bs2; j++) {
184         if (pc[j]!=0.0) {
185           flg = 1;
186           break;
187         }
188       }
189       if (flg) {
190         pv = b->a + bs2*bdiag[row];
191         PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); /* *pc = *pc * (*pv); */
192         pj = b->j + bdiag[row+1]+1;         /* begining of U(row,:) */
193         pv = b->a + bs2*(bdiag[row+1]+1);
194         nz = bdiag[row] - bdiag[row+1] - 1;         /* num of entries inU(row,:), excluding diag */
195         for (j=0; j<nz; j++) {
196           PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j);
197         }
198         ierr = PetscLogFlops(2*bs2*bs*(nz+1)-bs2);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
199       }
200     }
201 
202     /* finished row so stick it into b->a */
203     /* L part */
204     pv = b->a + bs2*bi[i];
205     pj = b->j + bi[i];
206     nz = bi[i+1] - bi[i];
207     for (j=0; j<nz; j++) {
208       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
209     }
210 
211     /* Mark diagonal and invert diagonal for simplier triangular solves */
212     pv = b->a + bs2*bdiag[i];
213     pj = b->j + bdiag[i];
214     /* if (*pj != i)SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"row %d != *pj %d",i,*pj); */
215     ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr);
216     ierr = PetscKernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr);
217 
218     /* U part */
219     pv = b->a + bs2*(bdiag[i+1]+1);
220     pj = b->j + bdiag[i+1]+1;
221     nz = bdiag[i] - bdiag[i+1] - 1;
222     for (j=0; j<nz; j++) {
223       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
224     }
225   }
226 
227   ierr = PetscFree(rtmp);CHKERRQ(ierr);
228   ierr = PetscFree3(v_work,mwork,v_pivots);CHKERRQ(ierr);
229   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
230   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
231 
232   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
233   ierr = ISIdentity(isicol,&col_identity);CHKERRQ(ierr);
234 
235   both_identity = (PetscBool) (row_identity && col_identity);
236   if (both_identity) {
237     C->ops->solve = MatSolve_SeqBAIJ_N_NaturalOrdering;
238   } else {
239     C->ops->solve = MatSolve_SeqBAIJ_N;
240   }
241   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_N;
242 
243   C->assembled = PETSC_TRUE;
244 
245   ierr = PetscLogFlops(1.333333333333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
246   PetscFunctionReturn(0);
247 }
248 
249 /*
250    ilu(0) with natural ordering under new data structure.
251    See MatILUFactorSymbolic_SeqAIJ_ilu0() for detailed description
252    because this code is almost identical to MatILUFactorSymbolic_SeqAIJ_ilu0_inplace().
253 */
254 
255 #undef __FUNCT__
256 #define __FUNCT__ "MatILUFactorSymbolic_SeqBAIJ_ilu0"
257 PetscErrorCode MatILUFactorSymbolic_SeqBAIJ_ilu0(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
258 {
259 
260   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b;
261   PetscErrorCode ierr;
262   PetscInt       n=a->mbs,*ai=a->i,*aj,*adiag=a->diag,bs2 = a->bs2;
263   PetscInt       i,j,nz,*bi,*bj,*bdiag,bi_temp;
264 
265   PetscFunctionBegin;
266   ierr = MatDuplicateNoCreate_SeqBAIJ(fact,A,MAT_DO_NOT_COPY_VALUES,PETSC_FALSE);CHKERRQ(ierr);
267   b    = (Mat_SeqBAIJ*)(fact)->data;
268 
269   /* allocate matrix arrays for new data structure */
270   ierr = PetscMalloc3(bs2*ai[n]+1,&b->a,ai[n]+1,&b->j,n+1,&b->i);CHKERRQ(ierr);
271   ierr = PetscLogObjectMemory((PetscObject)fact,ai[n]*(bs2*sizeof(PetscScalar)+sizeof(PetscInt))+(n+1)*sizeof(PetscInt));CHKERRQ(ierr);
272 
273   b->singlemalloc    = PETSC_TRUE;
274   b->free_a          = PETSC_TRUE;
275   b->free_ij         = PETSC_TRUE;
276   fact->preallocated = PETSC_TRUE;
277   fact->assembled    = PETSC_TRUE;
278   if (!b->diag) {
279     ierr = PetscMalloc1(n+1,&b->diag);CHKERRQ(ierr);
280     ierr = PetscLogObjectMemory((PetscObject)fact,(n+1)*sizeof(PetscInt));CHKERRQ(ierr);
281   }
282   bdiag = b->diag;
283 
284   if (n > 0) {
285     ierr = PetscMemzero(b->a,bs2*ai[n]*sizeof(MatScalar));CHKERRQ(ierr);
286   }
287 
288   /* set bi and bj with new data structure */
289   bi = b->i;
290   bj = b->j;
291 
292   /* L part */
293   bi[0] = 0;
294   for (i=0; i<n; i++) {
295     nz      = adiag[i] - ai[i];
296     bi[i+1] = bi[i] + nz;
297     aj      = a->j + ai[i];
298     for (j=0; j<nz; j++) {
299       *bj = aj[j]; bj++;
300     }
301   }
302 
303   /* U part */
304   bi_temp  = bi[n];
305   bdiag[n] = bi[n]-1;
306   for (i=n-1; i>=0; i--) {
307     nz      = ai[i+1] - adiag[i] - 1;
308     bi_temp = bi_temp + nz + 1;
309     aj      = a->j + adiag[i] + 1;
310     for (j=0; j<nz; j++) {
311       *bj = aj[j]; bj++;
312     }
313     /* diag[i] */
314     *bj      = i; bj++;
315     bdiag[i] = bi_temp - 1;
316   }
317   PetscFunctionReturn(0);
318 }
319 
320 #undef __FUNCT__
321 #define __FUNCT__ "MatILUFactorSymbolic_SeqBAIJ"
322 PetscErrorCode MatILUFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
323 {
324   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data,*b;
325   IS                 isicol;
326   PetscErrorCode     ierr;
327   const PetscInt     *r,*ic;
328   PetscInt           n=a->mbs,*ai=a->i,*aj=a->j,d;
329   PetscInt           *bi,*cols,nnz,*cols_lvl;
330   PetscInt           *bdiag,prow,fm,nzbd,reallocs=0,dcount=0;
331   PetscInt           i,levels,diagonal_fill;
332   PetscBool          col_identity,row_identity,both_identity;
333   PetscReal          f;
334   PetscInt           nlnk,*lnk,*lnk_lvl=NULL;
335   PetscBT            lnkbt;
336   PetscInt           nzi,*bj,**bj_ptr,**bjlvl_ptr;
337   PetscFreeSpaceList free_space    =NULL,current_space=NULL;
338   PetscFreeSpaceList free_space_lvl=NULL,current_space_lvl=NULL;
339   PetscBool          missing;
340   PetscInt           bs=A->rmap->bs,bs2=a->bs2;
341 
342   PetscFunctionBegin;
343   if (A->rmap->n != A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Must be square matrix, rows %D columns %D",A->rmap->n,A->cmap->n);
344   if (bs>1) {  /* check shifttype */
345     if (info->shifttype == MAT_SHIFT_NONZERO || info->shifttype == MAT_SHIFT_POSITIVE_DEFINITE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only MAT_SHIFT_NONE and MAT_SHIFT_INBLOCKS are supported for BAIJ matrix");
346   }
347 
348   ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr);
349   if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
350 
351   f             = info->fill;
352   levels        = (PetscInt)info->levels;
353   diagonal_fill = (PetscInt)info->diagonal_fill;
354 
355   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
356 
357   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
358   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
359 
360   both_identity = (PetscBool) (row_identity && col_identity);
361 
362   if (!levels && both_identity) {
363     /* special case: ilu(0) with natural ordering */
364     ierr = MatILUFactorSymbolic_SeqBAIJ_ilu0(fact,A,isrow,iscol,info);CHKERRQ(ierr);
365     ierr = MatSeqBAIJSetNumericFactorization(fact,both_identity);CHKERRQ(ierr);
366 
367     fact->factortype               = MAT_FACTOR_ILU;
368     (fact)->info.factor_mallocs    = 0;
369     (fact)->info.fill_ratio_given  = info->fill;
370     (fact)->info.fill_ratio_needed = 1.0;
371 
372     b                = (Mat_SeqBAIJ*)(fact)->data;
373     b->row           = isrow;
374     b->col           = iscol;
375     b->icol          = isicol;
376     ierr             = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
377     ierr             = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
378     b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE;
379 
380     ierr = PetscMalloc1((n+1)*bs,&b->solve_work);CHKERRQ(ierr);
381     PetscFunctionReturn(0);
382   }
383 
384   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
385   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
386 
387   /* get new row pointers */
388   ierr  = PetscMalloc1(n+1,&bi);CHKERRQ(ierr);
389   bi[0] = 0;
390   /* bdiag is location of diagonal in factor */
391   ierr     = PetscMalloc1(n+1,&bdiag);CHKERRQ(ierr);
392   bdiag[0] = 0;
393 
394   ierr = PetscMalloc2(n,&bj_ptr,n,&bjlvl_ptr);CHKERRQ(ierr);
395 
396   /* create a linked list for storing column indices of the active row */
397   nlnk = n + 1;
398   ierr = PetscIncompleteLLCreate(n,n,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
399 
400   /* initial FreeSpace size is f*(ai[n]+1) */
401   ierr              = PetscFreeSpaceGet(PetscRealIntMultTruncate(f,ai[n]+1),&free_space);CHKERRQ(ierr);
402   current_space     = free_space;
403   ierr              = PetscFreeSpaceGet(PetscRealIntMultTruncate(f,ai[n]+1),&free_space_lvl);CHKERRQ(ierr);
404   current_space_lvl = free_space_lvl;
405 
406   for (i=0; i<n; i++) {
407     nzi = 0;
408     /* copy current row into linked list */
409     nnz = ai[r[i]+1] - ai[r[i]];
410     if (!nnz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i);
411     cols   = aj + ai[r[i]];
412     lnk[i] = -1; /* marker to indicate if diagonal exists */
413     ierr   = PetscIncompleteLLInit(nnz,cols,n,ic,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
414     nzi   += nlnk;
415 
416     /* make sure diagonal entry is included */
417     if (diagonal_fill && lnk[i] == -1) {
418       fm = n;
419       while (lnk[fm] < i) fm = lnk[fm];
420       lnk[i]     = lnk[fm]; /* insert diagonal into linked list */
421       lnk[fm]    = i;
422       lnk_lvl[i] = 0;
423       nzi++; dcount++;
424     }
425 
426     /* add pivot rows into the active row */
427     nzbd = 0;
428     prow = lnk[n];
429     while (prow < i) {
430       nnz      = bdiag[prow];
431       cols     = bj_ptr[prow] + nnz + 1;
432       cols_lvl = bjlvl_ptr[prow] + nnz + 1;
433       nnz      = bi[prow+1] - bi[prow] - nnz - 1;
434 
435       ierr = PetscILULLAddSorted(nnz,cols,levels,cols_lvl,prow,nlnk,lnk,lnk_lvl,lnkbt,prow);CHKERRQ(ierr);
436       nzi += nlnk;
437       prow = lnk[prow];
438       nzbd++;
439     }
440     bdiag[i] = nzbd;
441     bi[i+1]  = bi[i] + nzi;
442 
443     /* if free space is not available, make more free space */
444     if (current_space->local_remaining<nzi) {
445       nnz  = PetscIntMultTruncate(2,PetscIntMultTruncate(nzi,(n - i))); /* estimated and max additional space needed */
446       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
447       ierr = PetscFreeSpaceGet(nnz,&current_space_lvl);CHKERRQ(ierr);
448       reallocs++;
449     }
450 
451     /* copy data into free_space and free_space_lvl, then initialize lnk */
452     ierr = PetscIncompleteLLClean(n,n,nzi,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
453 
454     bj_ptr[i]    = current_space->array;
455     bjlvl_ptr[i] = current_space_lvl->array;
456 
457     /* make sure the active row i has diagonal entry */
458     if (*(bj_ptr[i]+bdiag[i]) != i) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Row %D has missing diagonal in factored matrix\ntry running with -pc_factor_nonzeros_along_diagonal or -pc_factor_diagonal_fill",i);
459 
460     current_space->array           += nzi;
461     current_space->local_used      += nzi;
462     current_space->local_remaining -= nzi;
463 
464     current_space_lvl->array           += nzi;
465     current_space_lvl->local_used      += nzi;
466     current_space_lvl->local_remaining -= nzi;
467   }
468 
469   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
470   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
471 
472   /* copy free_space into bj and free free_space; set bi, bj, bdiag in new datastructure; */
473   ierr = PetscMalloc1(bi[n]+1,&bj);CHKERRQ(ierr);
474   ierr = PetscFreeSpaceContiguous_LU(&free_space,bj,n,bi,bdiag);CHKERRQ(ierr);
475 
476   ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
477   ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
478   ierr = PetscFree2(bj_ptr,bjlvl_ptr);CHKERRQ(ierr);
479 
480 #if defined(PETSC_USE_INFO)
481   {
482     PetscReal af = ((PetscReal)(bdiag[0]+1))/((PetscReal)ai[n]);
483     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)f,(double)af);CHKERRQ(ierr);
484     ierr = PetscInfo1(A,"Run with -[sub_]pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr);
485     ierr = PetscInfo1(A,"PCFactorSetFill([sub]pc,%g);\n",(double)af);CHKERRQ(ierr);
486     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
487     if (diagonal_fill) {
488       ierr = PetscInfo1(A,"Detected and replaced %D missing diagonals\n",dcount);CHKERRQ(ierr);
489     }
490   }
491 #endif
492 
493   /* put together the new matrix */
494   ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(fact,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr);
495   ierr = PetscLogObjectParent((PetscObject)fact,(PetscObject)isicol);CHKERRQ(ierr);
496 
497   b               = (Mat_SeqBAIJ*)(fact)->data;
498   b->free_a       = PETSC_TRUE;
499   b->free_ij      = PETSC_TRUE;
500   b->singlemalloc = PETSC_FALSE;
501 
502   ierr = PetscMalloc1(bs2*(bdiag[0]+1),&b->a);CHKERRQ(ierr);
503 
504   b->j          = bj;
505   b->i          = bi;
506   b->diag       = bdiag;
507   b->free_diag  = PETSC_TRUE;
508   b->ilen       = 0;
509   b->imax       = 0;
510   b->row        = isrow;
511   b->col        = iscol;
512   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
513   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
514   b->icol       = isicol;
515 
516   ierr = PetscMalloc1(bs*n+bs,&b->solve_work);CHKERRQ(ierr);
517   /* In b structure:  Free imax, ilen, old a, old j.
518      Allocate bdiag, solve_work, new a, new j */
519   ierr     = PetscLogObjectMemory((PetscObject)fact,(bdiag[0]+1) * (sizeof(PetscInt)+bs2*sizeof(PetscScalar)));CHKERRQ(ierr);
520   b->maxnz = b->nz = bdiag[0]+1;
521 
522   fact->info.factor_mallocs    = reallocs;
523   fact->info.fill_ratio_given  = f;
524   fact->info.fill_ratio_needed = ((PetscReal)(bdiag[0]+1))/((PetscReal)ai[n]);
525 
526   ierr = MatSeqBAIJSetNumericFactorization(fact,both_identity);CHKERRQ(ierr);
527   PetscFunctionReturn(0);
528 }
529 
530 /*
531      This code is virtually identical to MatILUFactorSymbolic_SeqAIJ
532    except that the data structure of Mat_SeqAIJ is slightly different.
533    Not a good example of code reuse.
534 */
535 #undef __FUNCT__
536 #define __FUNCT__ "MatILUFactorSymbolic_SeqBAIJ_inplace"
537 PetscErrorCode MatILUFactorSymbolic_SeqBAIJ_inplace(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
538 {
539   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b;
540   IS             isicol;
541   PetscErrorCode ierr;
542   const PetscInt *r,*ic,*ai = a->i,*aj = a->j,*xi;
543   PetscInt       prow,n = a->mbs,*ainew,*ajnew,jmax,*fill,nz,*im,*ajfill,*flev,*xitmp;
544   PetscInt       *dloc,idx,row,m,fm,nzf,nzi,reallocate = 0,dcount = 0;
545   PetscInt       incrlev,nnz,i,bs = A->rmap->bs,bs2 = a->bs2,levels,diagonal_fill,dd;
546   PetscBool      col_identity,row_identity,both_identity,flg;
547   PetscReal      f;
548 
549   PetscFunctionBegin;
550   ierr = MatMissingDiagonal_SeqBAIJ(A,&flg,&dd);CHKERRQ(ierr);
551   if (flg) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix A is missing diagonal entry in row %D",dd);
552 
553   f             = info->fill;
554   levels        = (PetscInt)info->levels;
555   diagonal_fill = (PetscInt)info->diagonal_fill;
556 
557   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
558 
559   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
560   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
561   both_identity = (PetscBool) (row_identity && col_identity);
562 
563   if (!levels && both_identity) {  /* special case copy the nonzero structure */
564     ierr = MatDuplicateNoCreate_SeqBAIJ(fact,A,MAT_DO_NOT_COPY_VALUES,PETSC_TRUE);CHKERRQ(ierr);
565     ierr = MatSeqBAIJSetNumericFactorization_inplace(fact,both_identity);CHKERRQ(ierr);
566 
567     fact->factortype = MAT_FACTOR_ILU;
568     b                = (Mat_SeqBAIJ*)fact->data;
569     b->row           = isrow;
570     b->col           = iscol;
571     ierr             = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
572     ierr             = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
573     b->icol          = isicol;
574     b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE;
575 
576     ierr         = PetscMalloc1((n+1)*bs,&b->solve_work);CHKERRQ(ierr);
577     PetscFunctionReturn(0);
578   }
579 
580   /* general case perform the symbolic factorization */
581   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
582   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
583 
584   /* get new row pointers */
585   ierr     = PetscMalloc1(n+1,&ainew);CHKERRQ(ierr);
586   ainew[0] = 0;
587   /* don't know how many column pointers are needed so estimate */
588   jmax = (PetscInt)(f*ai[n] + 1);
589   ierr = PetscMalloc1(jmax,&ajnew);CHKERRQ(ierr);
590   /* ajfill is level of fill for each fill entry */
591   ierr = PetscMalloc1(jmax,&ajfill);CHKERRQ(ierr);
592   /* fill is a linked list of nonzeros in active row */
593   ierr = PetscMalloc1(n+1,&fill);CHKERRQ(ierr);
594   /* im is level for each filled value */
595   ierr = PetscMalloc1(n+1,&im);CHKERRQ(ierr);
596   /* dloc is location of diagonal in factor */
597   ierr    = PetscMalloc1(n+1,&dloc);CHKERRQ(ierr);
598   dloc[0] = 0;
599   for (prow=0; prow<n; prow++) {
600 
601     /* copy prow into linked list */
602     nzf = nz = ai[r[prow]+1] - ai[r[prow]];
603     if (!nz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[prow],prow);
604     xi         = aj + ai[r[prow]];
605     fill[n]    = n;
606     fill[prow] = -1;   /* marker for diagonal entry */
607     while (nz--) {
608       fm  = n;
609       idx = ic[*xi++];
610       do {
611         m  = fm;
612         fm = fill[m];
613       } while (fm < idx);
614       fill[m]   = idx;
615       fill[idx] = fm;
616       im[idx]   = 0;
617     }
618 
619     /* make sure diagonal entry is included */
620     if (diagonal_fill && fill[prow] == -1) {
621       fm = n;
622       while (fill[fm] < prow) fm = fill[fm];
623       fill[prow] = fill[fm];    /* insert diagonal into linked list */
624       fill[fm]   = prow;
625       im[prow]   = 0;
626       nzf++;
627       dcount++;
628     }
629 
630     nzi = 0;
631     row = fill[n];
632     while (row < prow) {
633       incrlev = im[row] + 1;
634       nz      = dloc[row];
635       xi      = ajnew  + ainew[row] + nz + 1;
636       flev    = ajfill + ainew[row] + nz + 1;
637       nnz     = ainew[row+1] - ainew[row] - nz - 1;
638       fm      = row;
639       while (nnz-- > 0) {
640         idx = *xi++;
641         if (*flev + incrlev > levels) {
642           flev++;
643           continue;
644         }
645         do {
646           m  = fm;
647           fm = fill[m];
648         } while (fm < idx);
649         if (fm != idx) {
650           im[idx]   = *flev + incrlev;
651           fill[m]   = idx;
652           fill[idx] = fm;
653           fm        = idx;
654           nzf++;
655         } else if (im[idx] > *flev + incrlev) im[idx] = *flev+incrlev;
656         flev++;
657       }
658       row = fill[row];
659       nzi++;
660     }
661     /* copy new filled row into permanent storage */
662     ainew[prow+1] = ainew[prow] + nzf;
663     if (ainew[prow+1] > jmax) {
664 
665       /* estimate how much additional space we will need */
666       /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */
667       /* just double the memory each time */
668       PetscInt maxadd = jmax;
669       /* maxadd = (int)(((f*ai[n]+1)*(n-prow+5))/n); */
670       if (maxadd < nzf) maxadd = (n-prow)*(nzf+1);
671       jmax += maxadd;
672 
673       /* allocate a longer ajnew and ajfill */
674       ierr   = PetscMalloc1(jmax,&xitmp);CHKERRQ(ierr);
675       ierr   = PetscMemcpy(xitmp,ajnew,ainew[prow]*sizeof(PetscInt));CHKERRQ(ierr);
676       ierr   = PetscFree(ajnew);CHKERRQ(ierr);
677       ajnew  = xitmp;
678       ierr   = PetscMalloc1(jmax,&xitmp);CHKERRQ(ierr);
679       ierr   = PetscMemcpy(xitmp,ajfill,ainew[prow]*sizeof(PetscInt));CHKERRQ(ierr);
680       ierr   = PetscFree(ajfill);CHKERRQ(ierr);
681       ajfill = xitmp;
682       reallocate++;   /* count how many reallocations are needed */
683     }
684     xitmp      = ajnew + ainew[prow];
685     flev       = ajfill + ainew[prow];
686     dloc[prow] = nzi;
687     fm         = fill[n];
688     while (nzf--) {
689       *xitmp++ = fm;
690       *flev++  = im[fm];
691       fm       = fill[fm];
692     }
693     /* make sure row has diagonal entry */
694     if (ajnew[ainew[prow]+dloc[prow]] != prow) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Row %D has missing diagonal in factored matrix\n\
695                                                         try running with -pc_factor_nonzeros_along_diagonal or -pc_factor_diagonal_fill",prow);
696   }
697   ierr = PetscFree(ajfill);CHKERRQ(ierr);
698   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
699   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
700   ierr = PetscFree(fill);CHKERRQ(ierr);
701   ierr = PetscFree(im);CHKERRQ(ierr);
702 
703 #if defined(PETSC_USE_INFO)
704   {
705     PetscReal af = ((PetscReal)ainew[n])/((PetscReal)ai[n]);
706     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocate,(double)f,(double)af);CHKERRQ(ierr);
707     ierr = PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr);
708     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%g);\n",(double)af);CHKERRQ(ierr);
709     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
710     if (diagonal_fill) {
711       ierr = PetscInfo1(A,"Detected and replaced %D missing diagonals\n",dcount);CHKERRQ(ierr);
712     }
713   }
714 #endif
715 
716   /* put together the new matrix */
717   ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(fact,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr);
718   ierr = PetscLogObjectParent((PetscObject)fact,(PetscObject)isicol);CHKERRQ(ierr);
719   b    = (Mat_SeqBAIJ*)fact->data;
720 
721   b->free_a       = PETSC_TRUE;
722   b->free_ij      = PETSC_TRUE;
723   b->singlemalloc = PETSC_FALSE;
724 
725   ierr = PetscMalloc1(bs2*ainew[n],&b->a);CHKERRQ(ierr);
726 
727   b->j          = ajnew;
728   b->i          = ainew;
729   for (i=0; i<n; i++) dloc[i] += ainew[i];
730   b->diag          = dloc;
731   b->free_diag     = PETSC_TRUE;
732   b->ilen          = 0;
733   b->imax          = 0;
734   b->row           = isrow;
735   b->col           = iscol;
736   b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE;
737 
738   ierr    = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
739   ierr    = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
740   b->icol = isicol;
741   ierr    = PetscMalloc1(bs*n+bs,&b->solve_work);CHKERRQ(ierr);
742   /* In b structure:  Free imax, ilen, old a, old j.
743      Allocate dloc, solve_work, new a, new j */
744   ierr     = PetscLogObjectMemory((PetscObject)fact,(ainew[n]-n)*(sizeof(PetscInt))+bs2*ainew[n]*sizeof(PetscScalar));CHKERRQ(ierr);
745   b->maxnz = b->nz = ainew[n];
746 
747   fact->info.factor_mallocs    = reallocate;
748   fact->info.fill_ratio_given  = f;
749   fact->info.fill_ratio_needed = ((PetscReal)ainew[n])/((PetscReal)ai[prow]);
750 
751   ierr = MatSeqBAIJSetNumericFactorization_inplace(fact,both_identity);CHKERRQ(ierr);
752   PetscFunctionReturn(0);
753 }
754 
755 #undef __FUNCT__
756 #define __FUNCT__ "MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE"
757 PetscErrorCode MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE(Mat A)
758 {
759   /* Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; */
760   /* int i,*AJ=a->j,nz=a->nz; */
761 
762   PetscFunctionBegin;
763   /* Undo Column scaling */
764   /*    while (nz--) { */
765   /*      AJ[i] = AJ[i]/4; */
766   /*    } */
767   /* This should really invoke a push/pop logic, but we don't have that yet. */
768   A->ops->setunfactored = NULL;
769   PetscFunctionReturn(0);
770 }
771 
772 #undef __FUNCT__
773 #define __FUNCT__ "MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE_usj"
774 PetscErrorCode MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE_usj(Mat A)
775 {
776   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
777   PetscInt       *AJ=a->j,nz=a->nz;
778   unsigned short *aj=(unsigned short*)AJ;
779 
780   PetscFunctionBegin;
781   /* Is this really necessary? */
782   while (nz--) {
783     AJ[nz] = (int)((unsigned int)aj[nz]); /* First extend, then convert to signed. */
784   }
785   A->ops->setunfactored = NULL;
786   PetscFunctionReturn(0);
787 }
788 
789 
790