xref: /petsc/src/mat/impls/baij/seq/baijfact2.c (revision 96e086a24b9cdc04c368e2851f46791dcaa837d7)
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 /*
15    This is not much faster than MatLUFactorNumeric_SeqBAIJ_N() but the solve is faster at least sometimes
16 */
17 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_15_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info)
18 {
19   Mat             C =B;
20   Mat_SeqBAIJ     *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data;
21   PetscErrorCode  ierr;
22   PetscInt        i,j,k,ipvt[15];
23   const PetscInt  n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ajtmp,*bjtmp,*bdiag=b->diag,*pj;
24   PetscInt        nz,nzL,row;
25   MatScalar       *rtmp,*pc,*mwork,*pv,*vv,work[225];
26   const MatScalar *v,*aa=a->a;
27   PetscInt        bs2 = a->bs2,bs=A->rmap->bs,flg;
28   PetscInt        sol_ver;
29   PetscBool       allowzeropivot,zeropivotdetected;
30 
31   PetscFunctionBegin;
32   allowzeropivot = PetscNot(A->erroriffailure);
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     ierr = PetscKernel_A_gets_inverse_A_15(pv,ipvt,work,info->shiftamount,allowzeropivot,&zeropivotdetected);CHKERRQ(ierr);
106     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
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 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_N(Mat B,Mat A,const MatFactorInfo *info)
128 {
129   Mat            C     =B;
130   Mat_SeqBAIJ    *a    =(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data;
131   IS             isrow = b->row,isicol = b->icol;
132   PetscErrorCode ierr;
133   const PetscInt *r,*ic;
134   PetscInt       i,j,k,n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
135   PetscInt       *ajtmp,*bjtmp,nz,nzL,row,*bdiag=b->diag,*pj;
136   MatScalar      *rtmp,*pc,*mwork,*v,*pv,*aa=a->a;
137   PetscInt       bs=A->rmap->bs,bs2 = a->bs2,*v_pivots,flg;
138   MatScalar      *v_work;
139   PetscBool      col_identity,row_identity,both_identity;
140   PetscBool      allowzeropivot,zeropivotdetected;
141 
142   PetscFunctionBegin;
143   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
144   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
145   allowzeropivot = PetscNot(A->erroriffailure);
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     ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr);
215 
216     ierr = PetscKernel_A_gets_inverse_A(bs,pv,v_pivots,v_work,allowzeropivot,&zeropivotdetected);CHKERRQ(ierr);
217     if (zeropivotdetected) B->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
218 
219     /* U part */
220     pv = b->a + bs2*(bdiag[i+1]+1);
221     pj = b->j + bdiag[i+1]+1;
222     nz = bdiag[i] - bdiag[i+1] - 1;
223     for (j=0; j<nz; j++) {
224       ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr);
225     }
226   }
227 
228   ierr = PetscFree(rtmp);CHKERRQ(ierr);
229   ierr = PetscFree3(v_work,mwork,v_pivots);CHKERRQ(ierr);
230   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
231   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
232 
233   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
234   ierr = ISIdentity(isicol,&col_identity);CHKERRQ(ierr);
235 
236   both_identity = (PetscBool) (row_identity && col_identity);
237   if (both_identity) {
238     switch (bs) {
239     case  9:
240       C->ops->solve = MatSolve_SeqBAIJ_9_NaturalOrdering;
241       break;
242     case 11:
243       C->ops->solve = MatSolve_SeqBAIJ_11_NaturalOrdering;
244       break;
245     case 12:
246       C->ops->solve = MatSolve_SeqBAIJ_12_NaturalOrdering;
247       break;
248     case 13:
249       C->ops->solve = MatSolve_SeqBAIJ_13_NaturalOrdering;
250       break;
251     case 14:
252       C->ops->solve = MatSolve_SeqBAIJ_14_NaturalOrdering;
253       break;
254     default:
255       C->ops->solve = MatSolve_SeqBAIJ_N_NaturalOrdering;
256       break;
257     }
258   } else {
259     C->ops->solve = MatSolve_SeqBAIJ_N;
260   }
261   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_N;
262 
263   C->assembled = PETSC_TRUE;
264 
265   ierr = PetscLogFlops(1.333333333333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
266   PetscFunctionReturn(0);
267 }
268 
269 /*
270    ilu(0) with natural ordering under new data structure.
271    See MatILUFactorSymbolic_SeqAIJ_ilu0() for detailed description
272    because this code is almost identical to MatILUFactorSymbolic_SeqAIJ_ilu0_inplace().
273 */
274 
275 PetscErrorCode MatILUFactorSymbolic_SeqBAIJ_ilu0(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
276 {
277 
278   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b;
279   PetscErrorCode ierr;
280   PetscInt       n=a->mbs,*ai=a->i,*aj,*adiag=a->diag,bs2 = a->bs2;
281   PetscInt       i,j,nz,*bi,*bj,*bdiag,bi_temp;
282 
283   PetscFunctionBegin;
284   ierr = MatDuplicateNoCreate_SeqBAIJ(fact,A,MAT_DO_NOT_COPY_VALUES,PETSC_FALSE);CHKERRQ(ierr);
285   b    = (Mat_SeqBAIJ*)(fact)->data;
286 
287   /* allocate matrix arrays for new data structure */
288   ierr = PetscMalloc3(bs2*ai[n]+1,&b->a,ai[n]+1,&b->j,n+1,&b->i);CHKERRQ(ierr);
289   ierr = PetscLogObjectMemory((PetscObject)fact,ai[n]*(bs2*sizeof(PetscScalar)+sizeof(PetscInt))+(n+1)*sizeof(PetscInt));CHKERRQ(ierr);
290 
291   b->singlemalloc    = PETSC_TRUE;
292   b->free_a          = PETSC_TRUE;
293   b->free_ij         = PETSC_TRUE;
294   fact->preallocated = PETSC_TRUE;
295   fact->assembled    = PETSC_TRUE;
296   if (!b->diag) {
297     ierr = PetscMalloc1(n+1,&b->diag);CHKERRQ(ierr);
298     ierr = PetscLogObjectMemory((PetscObject)fact,(n+1)*sizeof(PetscInt));CHKERRQ(ierr);
299   }
300   bdiag = b->diag;
301 
302   if (n > 0) {
303     ierr = PetscMemzero(b->a,bs2*ai[n]*sizeof(MatScalar));CHKERRQ(ierr);
304   }
305 
306   /* set bi and bj with new data structure */
307   bi = b->i;
308   bj = b->j;
309 
310   /* L part */
311   bi[0] = 0;
312   for (i=0; i<n; i++) {
313     nz      = adiag[i] - ai[i];
314     bi[i+1] = bi[i] + nz;
315     aj      = a->j + ai[i];
316     for (j=0; j<nz; j++) {
317       *bj = aj[j]; bj++;
318     }
319   }
320 
321   /* U part */
322   bi_temp  = bi[n];
323   bdiag[n] = bi[n]-1;
324   for (i=n-1; i>=0; i--) {
325     nz      = ai[i+1] - adiag[i] - 1;
326     bi_temp = bi_temp + nz + 1;
327     aj      = a->j + adiag[i] + 1;
328     for (j=0; j<nz; j++) {
329       *bj = aj[j]; bj++;
330     }
331     /* diag[i] */
332     *bj      = i; bj++;
333     bdiag[i] = bi_temp - 1;
334   }
335   PetscFunctionReturn(0);
336 }
337 
338 PetscErrorCode MatILUFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
339 {
340   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data,*b;
341   IS                 isicol;
342   PetscErrorCode     ierr;
343   const PetscInt     *r,*ic;
344   PetscInt           n=a->mbs,*ai=a->i,*aj=a->j,d;
345   PetscInt           *bi,*cols,nnz,*cols_lvl;
346   PetscInt           *bdiag,prow,fm,nzbd,reallocs=0,dcount=0;
347   PetscInt           i,levels,diagonal_fill;
348   PetscBool          col_identity,row_identity,both_identity;
349   PetscReal          f;
350   PetscInt           nlnk,*lnk,*lnk_lvl=NULL;
351   PetscBT            lnkbt;
352   PetscInt           nzi,*bj,**bj_ptr,**bjlvl_ptr;
353   PetscFreeSpaceList free_space    =NULL,current_space=NULL;
354   PetscFreeSpaceList free_space_lvl=NULL,current_space_lvl=NULL;
355   PetscBool          missing;
356   PetscInt           bs=A->rmap->bs,bs2=a->bs2;
357 
358   PetscFunctionBegin;
359   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);
360   if (bs>1) {  /* check shifttype */
361     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");
362   }
363 
364   ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr);
365   if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
366 
367   f             = info->fill;
368   levels        = (PetscInt)info->levels;
369   diagonal_fill = (PetscInt)info->diagonal_fill;
370 
371   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
372 
373   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
374   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
375 
376   both_identity = (PetscBool) (row_identity && col_identity);
377 
378   if (!levels && both_identity) {
379     /* special case: ilu(0) with natural ordering */
380     ierr = MatILUFactorSymbolic_SeqBAIJ_ilu0(fact,A,isrow,iscol,info);CHKERRQ(ierr);
381     ierr = MatSeqBAIJSetNumericFactorization(fact,both_identity);CHKERRQ(ierr);
382 
383     fact->factortype               = MAT_FACTOR_ILU;
384     (fact)->info.factor_mallocs    = 0;
385     (fact)->info.fill_ratio_given  = info->fill;
386     (fact)->info.fill_ratio_needed = 1.0;
387 
388     b                = (Mat_SeqBAIJ*)(fact)->data;
389     b->row           = isrow;
390     b->col           = iscol;
391     b->icol          = isicol;
392     ierr             = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
393     ierr             = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
394     b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE;
395 
396     ierr = PetscMalloc1((n+1)*bs,&b->solve_work);CHKERRQ(ierr);
397     PetscFunctionReturn(0);
398   }
399 
400   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
401   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
402 
403   /* get new row pointers */
404   ierr  = PetscMalloc1(n+1,&bi);CHKERRQ(ierr);
405   bi[0] = 0;
406   /* bdiag is location of diagonal in factor */
407   ierr     = PetscMalloc1(n+1,&bdiag);CHKERRQ(ierr);
408   bdiag[0] = 0;
409 
410   ierr = PetscMalloc2(n,&bj_ptr,n,&bjlvl_ptr);CHKERRQ(ierr);
411 
412   /* create a linked list for storing column indices of the active row */
413   nlnk = n + 1;
414   ierr = PetscIncompleteLLCreate(n,n,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
415 
416   /* initial FreeSpace size is f*(ai[n]+1) */
417   ierr              = PetscFreeSpaceGet(PetscRealIntMultTruncate(f,ai[n]+1),&free_space);CHKERRQ(ierr);
418   current_space     = free_space;
419   ierr              = PetscFreeSpaceGet(PetscRealIntMultTruncate(f,ai[n]+1),&free_space_lvl);CHKERRQ(ierr);
420   current_space_lvl = free_space_lvl;
421 
422   for (i=0; i<n; i++) {
423     nzi = 0;
424     /* copy current row into linked list */
425     nnz = ai[r[i]+1] - ai[r[i]];
426     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);
427     cols   = aj + ai[r[i]];
428     lnk[i] = -1; /* marker to indicate if diagonal exists */
429     ierr   = PetscIncompleteLLInit(nnz,cols,n,ic,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
430     nzi   += nlnk;
431 
432     /* make sure diagonal entry is included */
433     if (diagonal_fill && lnk[i] == -1) {
434       fm = n;
435       while (lnk[fm] < i) fm = lnk[fm];
436       lnk[i]     = lnk[fm]; /* insert diagonal into linked list */
437       lnk[fm]    = i;
438       lnk_lvl[i] = 0;
439       nzi++; dcount++;
440     }
441 
442     /* add pivot rows into the active row */
443     nzbd = 0;
444     prow = lnk[n];
445     while (prow < i) {
446       nnz      = bdiag[prow];
447       cols     = bj_ptr[prow] + nnz + 1;
448       cols_lvl = bjlvl_ptr[prow] + nnz + 1;
449       nnz      = bi[prow+1] - bi[prow] - nnz - 1;
450 
451       ierr = PetscILULLAddSorted(nnz,cols,levels,cols_lvl,prow,nlnk,lnk,lnk_lvl,lnkbt,prow);CHKERRQ(ierr);
452       nzi += nlnk;
453       prow = lnk[prow];
454       nzbd++;
455     }
456     bdiag[i] = nzbd;
457     bi[i+1]  = bi[i] + nzi;
458 
459     /* if free space is not available, make more free space */
460     if (current_space->local_remaining<nzi) {
461       nnz  = PetscIntMultTruncate(2,PetscIntMultTruncate(nzi,(n - i))); /* estimated and max additional space needed */
462       ierr = PetscFreeSpaceGet(nnz,&current_space);CHKERRQ(ierr);
463       ierr = PetscFreeSpaceGet(nnz,&current_space_lvl);CHKERRQ(ierr);
464       reallocs++;
465     }
466 
467     /* copy data into free_space and free_space_lvl, then initialize lnk */
468     ierr = PetscIncompleteLLClean(n,n,nzi,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
469 
470     bj_ptr[i]    = current_space->array;
471     bjlvl_ptr[i] = current_space_lvl->array;
472 
473     /* make sure the active row i has diagonal entry */
474     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);
475 
476     current_space->array           += nzi;
477     current_space->local_used      += nzi;
478     current_space->local_remaining -= nzi;
479 
480     current_space_lvl->array           += nzi;
481     current_space_lvl->local_used      += nzi;
482     current_space_lvl->local_remaining -= nzi;
483   }
484 
485   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
486   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
487 
488   /* copy free_space into bj and free free_space; set bi, bj, bdiag in new datastructure; */
489   ierr = PetscMalloc1(bi[n]+1,&bj);CHKERRQ(ierr);
490   ierr = PetscFreeSpaceContiguous_LU(&free_space,bj,n,bi,bdiag);CHKERRQ(ierr);
491 
492   ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
493   ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
494   ierr = PetscFree2(bj_ptr,bjlvl_ptr);CHKERRQ(ierr);
495 
496 #if defined(PETSC_USE_INFO)
497   {
498     PetscReal af = ((PetscReal)(bdiag[0]+1))/((PetscReal)ai[n]);
499     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)f,(double)af);CHKERRQ(ierr);
500     ierr = PetscInfo1(A,"Run with -[sub_]pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr);
501     ierr = PetscInfo1(A,"PCFactorSetFill([sub]pc,%g);\n",(double)af);CHKERRQ(ierr);
502     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
503     if (diagonal_fill) {
504       ierr = PetscInfo1(A,"Detected and replaced %D missing diagonals\n",dcount);CHKERRQ(ierr);
505     }
506   }
507 #endif
508 
509   /* put together the new matrix */
510   ierr = MatSeqBAIJSetPreallocation(fact,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr);
511   ierr = PetscLogObjectParent((PetscObject)fact,(PetscObject)isicol);CHKERRQ(ierr);
512 
513   b               = (Mat_SeqBAIJ*)(fact)->data;
514   b->free_a       = PETSC_TRUE;
515   b->free_ij      = PETSC_TRUE;
516   b->singlemalloc = PETSC_FALSE;
517 
518   ierr = PetscMalloc1(bs2*(bdiag[0]+1),&b->a);CHKERRQ(ierr);
519 
520   b->j          = bj;
521   b->i          = bi;
522   b->diag       = bdiag;
523   b->free_diag  = PETSC_TRUE;
524   b->ilen       = 0;
525   b->imax       = 0;
526   b->row        = isrow;
527   b->col        = iscol;
528   ierr          = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
529   ierr          = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
530   b->icol       = isicol;
531 
532   ierr = PetscMalloc1(bs*n+bs,&b->solve_work);CHKERRQ(ierr);
533   /* In b structure:  Free imax, ilen, old a, old j.
534      Allocate bdiag, solve_work, new a, new j */
535   ierr     = PetscLogObjectMemory((PetscObject)fact,(bdiag[0]+1) * (sizeof(PetscInt)+bs2*sizeof(PetscScalar)));CHKERRQ(ierr);
536   b->maxnz = b->nz = bdiag[0]+1;
537 
538   fact->info.factor_mallocs    = reallocs;
539   fact->info.fill_ratio_given  = f;
540   fact->info.fill_ratio_needed = ((PetscReal)(bdiag[0]+1))/((PetscReal)ai[n]);
541 
542   ierr = MatSeqBAIJSetNumericFactorization(fact,both_identity);CHKERRQ(ierr);
543   PetscFunctionReturn(0);
544 }
545 
546 /*
547      This code is virtually identical to MatILUFactorSymbolic_SeqAIJ
548    except that the data structure of Mat_SeqAIJ is slightly different.
549    Not a good example of code reuse.
550 */
551 PetscErrorCode MatILUFactorSymbolic_SeqBAIJ_inplace(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
552 {
553   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b;
554   IS             isicol;
555   PetscErrorCode ierr;
556   const PetscInt *r,*ic,*ai = a->i,*aj = a->j,*xi;
557   PetscInt       prow,n = a->mbs,*ainew,*ajnew,jmax,*fill,nz,*im,*ajfill,*flev,*xitmp;
558   PetscInt       *dloc,idx,row,m,fm,nzf,nzi,reallocate = 0,dcount = 0;
559   PetscInt       incrlev,nnz,i,bs = A->rmap->bs,bs2 = a->bs2,levels,diagonal_fill,dd;
560   PetscBool      col_identity,row_identity,both_identity,flg;
561   PetscReal      f;
562 
563   PetscFunctionBegin;
564   ierr = MatMissingDiagonal_SeqBAIJ(A,&flg,&dd);CHKERRQ(ierr);
565   if (flg) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix A is missing diagonal entry in row %D",dd);
566 
567   f             = info->fill;
568   levels        = (PetscInt)info->levels;
569   diagonal_fill = (PetscInt)info->diagonal_fill;
570 
571   ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr);
572 
573   ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr);
574   ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr);
575   both_identity = (PetscBool) (row_identity && col_identity);
576 
577   if (!levels && both_identity) {  /* special case copy the nonzero structure */
578     ierr = MatDuplicateNoCreate_SeqBAIJ(fact,A,MAT_DO_NOT_COPY_VALUES,PETSC_TRUE);CHKERRQ(ierr);
579     ierr = MatSeqBAIJSetNumericFactorization_inplace(fact,both_identity);CHKERRQ(ierr);
580 
581     fact->factortype = MAT_FACTOR_ILU;
582     b                = (Mat_SeqBAIJ*)fact->data;
583     b->row           = isrow;
584     b->col           = iscol;
585     ierr             = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
586     ierr             = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
587     b->icol          = isicol;
588     b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE;
589 
590     ierr         = PetscMalloc1((n+1)*bs,&b->solve_work);CHKERRQ(ierr);
591     PetscFunctionReturn(0);
592   }
593 
594   /* general case perform the symbolic factorization */
595   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
596   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
597 
598   /* get new row pointers */
599   ierr     = PetscMalloc1(n+1,&ainew);CHKERRQ(ierr);
600   ainew[0] = 0;
601   /* don't know how many column pointers are needed so estimate */
602   jmax = (PetscInt)(f*ai[n] + 1);
603   ierr = PetscMalloc1(jmax,&ajnew);CHKERRQ(ierr);
604   /* ajfill is level of fill for each fill entry */
605   ierr = PetscMalloc1(jmax,&ajfill);CHKERRQ(ierr);
606   /* fill is a linked list of nonzeros in active row */
607   ierr = PetscMalloc1(n+1,&fill);CHKERRQ(ierr);
608   /* im is level for each filled value */
609   ierr = PetscMalloc1(n+1,&im);CHKERRQ(ierr);
610   /* dloc is location of diagonal in factor */
611   ierr    = PetscMalloc1(n+1,&dloc);CHKERRQ(ierr);
612   dloc[0] = 0;
613   for (prow=0; prow<n; prow++) {
614 
615     /* copy prow into linked list */
616     nzf = nz = ai[r[prow]+1] - ai[r[prow]];
617     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);
618     xi         = aj + ai[r[prow]];
619     fill[n]    = n;
620     fill[prow] = -1;   /* marker for diagonal entry */
621     while (nz--) {
622       fm  = n;
623       idx = ic[*xi++];
624       do {
625         m  = fm;
626         fm = fill[m];
627       } while (fm < idx);
628       fill[m]   = idx;
629       fill[idx] = fm;
630       im[idx]   = 0;
631     }
632 
633     /* make sure diagonal entry is included */
634     if (diagonal_fill && fill[prow] == -1) {
635       fm = n;
636       while (fill[fm] < prow) fm = fill[fm];
637       fill[prow] = fill[fm];    /* insert diagonal into linked list */
638       fill[fm]   = prow;
639       im[prow]   = 0;
640       nzf++;
641       dcount++;
642     }
643 
644     nzi = 0;
645     row = fill[n];
646     while (row < prow) {
647       incrlev = im[row] + 1;
648       nz      = dloc[row];
649       xi      = ajnew  + ainew[row] + nz + 1;
650       flev    = ajfill + ainew[row] + nz + 1;
651       nnz     = ainew[row+1] - ainew[row] - nz - 1;
652       fm      = row;
653       while (nnz-- > 0) {
654         idx = *xi++;
655         if (*flev + incrlev > levels) {
656           flev++;
657           continue;
658         }
659         do {
660           m  = fm;
661           fm = fill[m];
662         } while (fm < idx);
663         if (fm != idx) {
664           im[idx]   = *flev + incrlev;
665           fill[m]   = idx;
666           fill[idx] = fm;
667           fm        = idx;
668           nzf++;
669         } else if (im[idx] > *flev + incrlev) im[idx] = *flev+incrlev;
670         flev++;
671       }
672       row = fill[row];
673       nzi++;
674     }
675     /* copy new filled row into permanent storage */
676     ainew[prow+1] = ainew[prow] + nzf;
677     if (ainew[prow+1] > jmax) {
678 
679       /* estimate how much additional space we will need */
680       /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */
681       /* just double the memory each time */
682       PetscInt maxadd = jmax;
683       /* maxadd = (int)(((f*ai[n]+1)*(n-prow+5))/n); */
684       if (maxadd < nzf) maxadd = (n-prow)*(nzf+1);
685       jmax += maxadd;
686 
687       /* allocate a longer ajnew and ajfill */
688       ierr   = PetscMalloc1(jmax,&xitmp);CHKERRQ(ierr);
689       ierr   = PetscMemcpy(xitmp,ajnew,ainew[prow]*sizeof(PetscInt));CHKERRQ(ierr);
690       ierr   = PetscFree(ajnew);CHKERRQ(ierr);
691       ajnew  = xitmp;
692       ierr   = PetscMalloc1(jmax,&xitmp);CHKERRQ(ierr);
693       ierr   = PetscMemcpy(xitmp,ajfill,ainew[prow]*sizeof(PetscInt));CHKERRQ(ierr);
694       ierr   = PetscFree(ajfill);CHKERRQ(ierr);
695       ajfill = xitmp;
696       reallocate++;   /* count how many reallocations are needed */
697     }
698     xitmp      = ajnew + ainew[prow];
699     flev       = ajfill + ainew[prow];
700     dloc[prow] = nzi;
701     fm         = fill[n];
702     while (nzf--) {
703       *xitmp++ = fm;
704       *flev++  = im[fm];
705       fm       = fill[fm];
706     }
707     /* make sure row has diagonal entry */
708     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\
709                                                         try running with -pc_factor_nonzeros_along_diagonal or -pc_factor_diagonal_fill",prow);
710   }
711   ierr = PetscFree(ajfill);CHKERRQ(ierr);
712   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
713   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
714   ierr = PetscFree(fill);CHKERRQ(ierr);
715   ierr = PetscFree(im);CHKERRQ(ierr);
716 
717 #if defined(PETSC_USE_INFO)
718   {
719     PetscReal af = ((PetscReal)ainew[n])/((PetscReal)ai[n]);
720     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocate,(double)f,(double)af);CHKERRQ(ierr);
721     ierr = PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr);
722     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%g);\n",(double)af);CHKERRQ(ierr);
723     ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr);
724     if (diagonal_fill) {
725       ierr = PetscInfo1(A,"Detected and replaced %D missing diagonals\n",dcount);CHKERRQ(ierr);
726     }
727   }
728 #endif
729 
730   /* put together the new matrix */
731   ierr = MatSeqBAIJSetPreallocation(fact,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr);
732   ierr = PetscLogObjectParent((PetscObject)fact,(PetscObject)isicol);CHKERRQ(ierr);
733   b    = (Mat_SeqBAIJ*)fact->data;
734 
735   b->free_a       = PETSC_TRUE;
736   b->free_ij      = PETSC_TRUE;
737   b->singlemalloc = PETSC_FALSE;
738 
739   ierr = PetscMalloc1(bs2*ainew[n],&b->a);CHKERRQ(ierr);
740 
741   b->j          = ajnew;
742   b->i          = ainew;
743   for (i=0; i<n; i++) dloc[i] += ainew[i];
744   b->diag          = dloc;
745   b->free_diag     = PETSC_TRUE;
746   b->ilen          = 0;
747   b->imax          = 0;
748   b->row           = isrow;
749   b->col           = iscol;
750   b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE;
751 
752   ierr    = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr);
753   ierr    = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr);
754   b->icol = isicol;
755   ierr    = PetscMalloc1(bs*n+bs,&b->solve_work);CHKERRQ(ierr);
756   /* In b structure:  Free imax, ilen, old a, old j.
757      Allocate dloc, solve_work, new a, new j */
758   ierr     = PetscLogObjectMemory((PetscObject)fact,(ainew[n]-n)*(sizeof(PetscInt))+bs2*ainew[n]*sizeof(PetscScalar));CHKERRQ(ierr);
759   b->maxnz = b->nz = ainew[n];
760 
761   fact->info.factor_mallocs    = reallocate;
762   fact->info.fill_ratio_given  = f;
763   fact->info.fill_ratio_needed = ((PetscReal)ainew[n])/((PetscReal)ai[prow]);
764 
765   ierr = MatSeqBAIJSetNumericFactorization_inplace(fact,both_identity);CHKERRQ(ierr);
766   PetscFunctionReturn(0);
767 }
768 
769 PetscErrorCode MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE(Mat A)
770 {
771   /* Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; */
772   /* int i,*AJ=a->j,nz=a->nz; */
773 
774   PetscFunctionBegin;
775   /* Undo Column scaling */
776   /*    while (nz--) { */
777   /*      AJ[i] = AJ[i]/4; */
778   /*    } */
779   /* This should really invoke a push/pop logic, but we don't have that yet. */
780   A->ops->setunfactored = NULL;
781   PetscFunctionReturn(0);
782 }
783 
784 PetscErrorCode MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE_usj(Mat A)
785 {
786   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
787   PetscInt       *AJ=a->j,nz=a->nz;
788   unsigned short *aj=(unsigned short*)AJ;
789 
790   PetscFunctionBegin;
791   /* Is this really necessary? */
792   while (nz--) {
793     AJ[nz] = (int)((unsigned int)aj[nz]); /* First extend, then convert to signed. */
794   }
795   A->ops->setunfactored = NULL;
796   PetscFunctionReturn(0);
797 }
798 
799 
800