xref: /petsc/src/mat/impls/baij/seq/baijfact.c (revision ce0a2cd1da0658c2b28aad1be2e2c8e41567bece)
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/inline/ilu.h"
8 
9 /* ------------------------------------------------------------*/
10 /*
11       Version for when blocks are 2 by 2
12 */
13 #undef __FUNCT__
14 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2"
15 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2(Mat A,MatFactorInfo *info,Mat *B)
16 {
17   Mat            C = *B;
18   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data;
19   IS             isrow = b->row,isicol = b->icol;
20   PetscErrorCode ierr;
21   PetscInt       *r,*ic,i,j,n = a->mbs,*bi = b->i,*bj = b->j;
22   PetscInt       *ajtmpold,*ajtmp,nz,row;
23   PetscInt       *diag_offset=b->diag,idx,*ai=a->i,*aj=a->j,*pj;
24   MatScalar      *pv,*v,*rtmp,m1,m2,m3,m4,*pc,*w,*x,x1,x2,x3,x4;
25   MatScalar      p1,p2,p3,p4;
26   MatScalar      *ba = b->a,*aa = a->a;
27   PetscReal      shift = info->shiftinblocks;
28 
29   PetscFunctionBegin;
30   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
31   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
32   ierr  = PetscMalloc(4*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
33 
34   for (i=0; i<n; i++) {
35     nz    = bi[i+1] - bi[i];
36     ajtmp = bj + bi[i];
37     for  (j=0; j<nz; j++) {
38       x = rtmp+4*ajtmp[j]; x[0] = x[1] = x[2] = x[3] = 0.0;
39     }
40     /* load in initial (unfactored row) */
41     idx      = r[i];
42     nz       = ai[idx+1] - ai[idx];
43     ajtmpold = aj + ai[idx];
44     v        = aa + 4*ai[idx];
45     for (j=0; j<nz; j++) {
46       x    = rtmp+4*ic[ajtmpold[j]];
47       x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3];
48       v    += 4;
49     }
50     row = *ajtmp++;
51     while (row < i) {
52       pc = rtmp + 4*row;
53       p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3];
54       if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0) {
55         pv = ba + 4*diag_offset[row];
56         pj = bj + diag_offset[row] + 1;
57         x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
58         pc[0] = m1 = p1*x1 + p3*x2;
59         pc[1] = m2 = p2*x1 + p4*x2;
60         pc[2] = m3 = p1*x3 + p3*x4;
61         pc[3] = m4 = p2*x3 + p4*x4;
62         nz = bi[row+1] - diag_offset[row] - 1;
63         pv += 4;
64         for (j=0; j<nz; j++) {
65           x1   = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
66           x    = rtmp + 4*pj[j];
67           x[0] -= m1*x1 + m3*x2;
68           x[1] -= m2*x1 + m4*x2;
69           x[2] -= m1*x3 + m3*x4;
70           x[3] -= m2*x3 + m4*x4;
71           pv   += 4;
72         }
73         ierr = PetscLogFlops(16*nz+12);CHKERRQ(ierr);
74       }
75       row = *ajtmp++;
76     }
77     /* finished row so stick it into b->a */
78     pv = ba + 4*bi[i];
79     pj = bj + bi[i];
80     nz = bi[i+1] - bi[i];
81     for (j=0; j<nz; j++) {
82       x     = rtmp+4*pj[j];
83       pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3];
84       pv   += 4;
85     }
86     /* invert diagonal block */
87     w = ba + 4*diag_offset[i];
88     ierr = Kernel_A_gets_inverse_A_2(w,shift);CHKERRQ(ierr);
89   }
90 
91   ierr = PetscFree(rtmp);CHKERRQ(ierr);
92   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
93   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
94   C->factor = FACTOR_LU;
95   C->assembled = PETSC_TRUE;
96   ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
97   PetscFunctionReturn(0);
98 }
99 /*
100       Version for when blocks are 2 by 2 Using natural ordering
101 */
102 #undef __FUNCT__
103 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering"
104 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering(Mat A,MatFactorInfo *info,Mat *B)
105 {
106   Mat            C = *B;
107   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data;
108   PetscErrorCode ierr;
109   PetscInt       i,j,n = a->mbs,*bi = b->i,*bj = b->j;
110   PetscInt       *ajtmpold,*ajtmp,nz,row;
111   PetscInt       *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj;
112   MatScalar      *pv,*v,*rtmp,*pc,*w,*x;
113   MatScalar      p1,p2,p3,p4,m1,m2,m3,m4,x1,x2,x3,x4;
114   MatScalar      *ba = b->a,*aa = a->a;
115   PetscReal      shift = info->shiftinblocks;
116 
117   PetscFunctionBegin;
118   ierr = PetscMalloc(4*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
119 
120   for (i=0; i<n; i++) {
121     nz    = bi[i+1] - bi[i];
122     ajtmp = bj + bi[i];
123     for  (j=0; j<nz; j++) {
124       x = rtmp+4*ajtmp[j];
125       x[0]  = x[1]  = x[2]  = x[3]  = 0.0;
126     }
127     /* load in initial (unfactored row) */
128     nz       = ai[i+1] - ai[i];
129     ajtmpold = aj + ai[i];
130     v        = aa + 4*ai[i];
131     for (j=0; j<nz; j++) {
132       x    = rtmp+4*ajtmpold[j];
133       x[0]  = v[0];  x[1]  = v[1];  x[2]  = v[2];  x[3]  = v[3];
134       v    += 4;
135     }
136     row = *ajtmp++;
137     while (row < i) {
138       pc  = rtmp + 4*row;
139       p1  = pc[0];  p2  = pc[1];  p3  = pc[2];  p4  = pc[3];
140       if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0) {
141         pv = ba + 4*diag_offset[row];
142         pj = bj + diag_offset[row] + 1;
143         x1  = pv[0];  x2  = pv[1];  x3  = pv[2];  x4  = pv[3];
144         pc[0] = m1 = p1*x1 + p3*x2;
145         pc[1] = m2 = p2*x1 + p4*x2;
146         pc[2] = m3 = p1*x3 + p3*x4;
147         pc[3] = m4 = p2*x3 + p4*x4;
148         nz = bi[row+1] - diag_offset[row] - 1;
149         pv += 4;
150         for (j=0; j<nz; j++) {
151           x1   = pv[0];  x2  = pv[1];   x3 = pv[2];  x4  = pv[3];
152           x    = rtmp + 4*pj[j];
153           x[0] -= m1*x1 + m3*x2;
154           x[1] -= m2*x1 + m4*x2;
155           x[2] -= m1*x3 + m3*x4;
156           x[3] -= m2*x3 + m4*x4;
157           pv   += 4;
158         }
159         ierr = PetscLogFlops(16*nz+12);CHKERRQ(ierr);
160       }
161       row = *ajtmp++;
162     }
163     /* finished row so stick it into b->a */
164     pv = ba + 4*bi[i];
165     pj = bj + bi[i];
166     nz = bi[i+1] - bi[i];
167     for (j=0; j<nz; j++) {
168       x      = rtmp+4*pj[j];
169       pv[0]  = x[0];  pv[1]  = x[1];  pv[2]  = x[2];  pv[3]  = x[3];
170       pv   += 4;
171     }
172     /* invert diagonal block */
173     w = ba + 4*diag_offset[i];
174     ierr = Kernel_A_gets_inverse_A_2(w,shift);CHKERRQ(ierr);
175   }
176 
177   ierr = PetscFree(rtmp);CHKERRQ(ierr);
178   C->factor    = FACTOR_LU;
179   C->assembled = PETSC_TRUE;
180   ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
181   PetscFunctionReturn(0);
182 }
183 
184 /* ----------------------------------------------------------- */
185 /*
186      Version for when blocks are 1 by 1.
187 */
188 #undef __FUNCT__
189 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_1"
190 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_1(Mat A,MatFactorInfo *info,Mat *B)
191 {
192   Mat            C = *B;
193   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data;
194   IS             isrow = b->row,isicol = b->icol;
195   PetscErrorCode ierr;
196   PetscInt       *r,*ic,i,j,n = a->mbs,*bi = b->i,*bj = b->j;
197   PetscInt       *ajtmpold,*ajtmp,nz,row,*ai = a->i,*aj = a->j;
198   PetscInt       *diag_offset = b->diag,diag,*pj;
199   MatScalar      *pv,*v,*rtmp,multiplier,*pc;
200   MatScalar      *ba = b->a,*aa = a->a;
201 
202   PetscFunctionBegin;
203   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
204   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
205   ierr  = PetscMalloc((n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
206 
207   for (i=0; i<n; i++) {
208     nz    = bi[i+1] - bi[i];
209     ajtmp = bj + bi[i];
210     for  (j=0; j<nz; j++) rtmp[ajtmp[j]] = 0.0;
211 
212     /* load in initial (unfactored row) */
213     nz       = ai[r[i]+1] - ai[r[i]];
214     ajtmpold = aj + ai[r[i]];
215     v        = aa + ai[r[i]];
216     for (j=0; j<nz; j++) rtmp[ic[ajtmpold[j]]] =  v[j];
217 
218     row = *ajtmp++;
219     while (row < i) {
220       pc = rtmp + row;
221       if (*pc != 0.0) {
222         pv         = ba + diag_offset[row];
223         pj         = bj + diag_offset[row] + 1;
224         multiplier = *pc * *pv++;
225         *pc        = multiplier;
226         nz         = bi[row+1] - diag_offset[row] - 1;
227         for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j];
228         ierr = PetscLogFlops(1+2*nz);CHKERRQ(ierr);
229       }
230       row = *ajtmp++;
231     }
232     /* finished row so stick it into b->a */
233     pv = ba + bi[i];
234     pj = bj + bi[i];
235     nz = bi[i+1] - bi[i];
236     for (j=0; j<nz; j++) {pv[j] = rtmp[pj[j]];}
237     diag = diag_offset[i] - bi[i];
238     /* check pivot entry for current row */
239     if (pv[diag] == 0.0) {
240       SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot");
241     }
242     pv[diag] = 1.0/pv[diag];
243   }
244 
245   ierr = PetscFree(rtmp);CHKERRQ(ierr);
246   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
247   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
248   C->factor    = FACTOR_LU;
249   C->assembled = PETSC_TRUE;
250   ierr = PetscLogFlops(C->cmap.n);CHKERRQ(ierr);
251   PetscFunctionReturn(0);
252 }
253 
254 
255 /* ----------------------------------------------------------- */
256 #undef __FUNCT__
257 #define __FUNCT__ "MatLUFactor_SeqBAIJ"
258 PetscErrorCode MatLUFactor_SeqBAIJ(Mat A,IS row,IS col,MatFactorInfo *info)
259 {
260   PetscErrorCode ierr;
261   Mat            C;
262 
263   PetscFunctionBegin;
264   ierr = MatLUFactorSymbolic(A,row,col,info,&C);CHKERRQ(ierr);
265   ierr = MatLUFactorNumeric(A,info,&C);CHKERRQ(ierr);
266   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
267   ierr = PetscLogObjectParent(A,((Mat_SeqBAIJ*)(A->data))->icol);CHKERRQ(ierr);
268   PetscFunctionReturn(0);
269 }
270 
271 #include "src/mat/impls/sbaij/seq/sbaij.h"
272 #undef __FUNCT__
273 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N"
274 PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N(Mat A,MatFactorInfo *info,Mat *B)
275 {
276   PetscErrorCode ierr;
277   Mat            C = *B;
278   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data;
279   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
280   IS             ip=b->row;
281   PetscInt       *rip,i,j,mbs=a->mbs,bs=A->rmap.bs,*bi=b->i,*bj=b->j,*bcol;
282   PetscInt       *ai=a->i,*aj=a->j;
283   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
284   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
285   PetscReal      zeropivot,rs,shiftnz;
286   PetscReal      shiftpd;
287   ChShift_Ctx    sctx;
288   PetscInt       newshift;
289 
290   PetscFunctionBegin;
291   if (bs > 1) {
292     if (!a->sbaijMat){
293       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
294     }
295     ierr = (a->sbaijMat)->ops->choleskyfactornumeric(a->sbaijMat,info,B);CHKERRQ(ierr);
296     ierr = MatDestroy(a->sbaijMat);CHKERRQ(ierr);
297     a->sbaijMat = PETSC_NULL;
298     PetscFunctionReturn(0);
299   }
300 
301   /* initialization */
302   shiftnz   = info->shiftnz;
303   shiftpd   = info->shiftpd;
304   zeropivot = info->zeropivot;
305 
306   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
307   nz   = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar);
308   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
309   jl   = il + mbs;
310   rtmp = (MatScalar*)(jl + mbs);
311 
312   sctx.shift_amount = 0;
313   sctx.nshift       = 0;
314   do {
315     sctx.chshift = PETSC_FALSE;
316     for (i=0; i<mbs; i++) {
317       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
318     }
319 
320     for (k = 0; k<mbs; k++){
321       bval = ba + bi[k];
322       /* initialize k-th row by the perm[k]-th row of A */
323       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
324       for (j = jmin; j < jmax; j++){
325         col = rip[aj[j]];
326         if (col >= k){ /* only take upper triangular entry */
327           rtmp[col] = aa[j];
328           *bval++  = 0.0; /* for in-place factorization */
329         }
330       }
331 
332       /* shift the diagonal of the matrix */
333       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
334 
335       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
336       dk = rtmp[k];
337       i = jl[k]; /* first row to be added to k_th row  */
338 
339       while (i < k){
340         nexti = jl[i]; /* next row to be added to k_th row */
341 
342         /* compute multiplier, update diag(k) and U(i,k) */
343         ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
344         uikdi = - ba[ili]*ba[bi[i]];  /* diagonal(k) */
345         dk += uikdi*ba[ili];
346         ba[ili] = uikdi; /* -U(i,k) */
347 
348         /* add multiple of row i to k-th row */
349         jmin = ili + 1; jmax = bi[i+1];
350         if (jmin < jmax){
351           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
352           /* update il and jl for row i */
353           il[i] = jmin;
354           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
355         }
356         i = nexti;
357       }
358 
359       /* shift the diagonals when zero pivot is detected */
360       /* compute rs=sum of abs(off-diagonal) */
361       rs   = 0.0;
362       jmin = bi[k]+1;
363       nz   = bi[k+1] - jmin;
364       if (nz){
365         bcol = bj + jmin;
366         while (nz--){
367           rs += PetscAbsScalar(rtmp[*bcol]);
368           bcol++;
369         }
370       }
371 
372       sctx.rs = rs;
373       sctx.pv = dk;
374       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
375       if (newshift == 1) break;
376 
377       /* copy data into U(k,:) */
378       ba[bi[k]] = 1.0/dk; /* U(k,k) */
379       jmin = bi[k]+1; jmax = bi[k+1];
380       if (jmin < jmax) {
381         for (j=jmin; j<jmax; j++){
382           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
383         }
384         /* add the k-th row into il and jl */
385         il[k] = jmin;
386         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
387       }
388     }
389   } while (sctx.chshift);
390   ierr = PetscFree(il);CHKERRQ(ierr);
391 
392   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
393   C->factor       = FACTOR_CHOLESKY;
394   C->assembled    = PETSC_TRUE;
395   C->preallocated = PETSC_TRUE;
396   ierr = PetscLogFlops(C->rmap.N);CHKERRQ(ierr);
397   if (sctx.nshift){
398     if (shiftnz) {
399       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
400     } else if (shiftpd) {
401       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
402     }
403   }
404   PetscFunctionReturn(0);
405 }
406 
407 #undef __FUNCT__
408 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering"
409 PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering(Mat A,MatFactorInfo *info,Mat *fact)
410 {
411   Mat            C = *fact;
412   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data;
413   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
414   PetscErrorCode ierr;
415   PetscInt       i,j,am=a->mbs;
416   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
417   PetscInt       k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz;
418   MatScalar      *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval;
419   PetscReal      zeropivot,rs,shiftnz;
420   PetscReal      shiftpd;
421   ChShift_Ctx    sctx;
422   PetscInt       newshift;
423 
424   PetscFunctionBegin;
425   /* initialization */
426   shiftnz   = info->shiftnz;
427   shiftpd   = info->shiftpd;
428   zeropivot = info->zeropivot;
429 
430   nz   = (2*am+1)*sizeof(PetscInt)+am*sizeof(MatScalar);
431   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
432   jl   = il + am;
433   rtmp = (MatScalar*)(jl + am);
434 
435   sctx.shift_amount = 0;
436   sctx.nshift       = 0;
437   do {
438     sctx.chshift = PETSC_FALSE;
439     for (i=0; i<am; i++) {
440       rtmp[i] = 0.0; jl[i] = am; il[0] = 0;
441     }
442 
443     for (k = 0; k<am; k++){
444     /* initialize k-th row with elements nonzero in row perm(k) of A */
445       nz   = ai[k+1] - ai[k];
446       acol = aj + ai[k];
447       aval = aa + ai[k];
448       bval = ba + bi[k];
449       while (nz -- ){
450         if (*acol < k) { /* skip lower triangular entries */
451           acol++; aval++;
452         } else {
453           rtmp[*acol++] = *aval++;
454           *bval++       = 0.0; /* for in-place factorization */
455         }
456       }
457 
458       /* shift the diagonal of the matrix */
459       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
460 
461       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
462       dk = rtmp[k];
463       i  = jl[k]; /* first row to be added to k_th row  */
464 
465       while (i < k){
466         nexti = jl[i]; /* next row to be added to k_th row */
467         /* compute multiplier, update D(k) and U(i,k) */
468         ili   = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
469         uikdi = - ba[ili]*ba[bi[i]];
470         dk   += uikdi*ba[ili];
471         ba[ili] = uikdi; /* -U(i,k) */
472 
473         /* add multiple of row i to k-th row ... */
474         jmin = ili + 1;
475         nz   = bi[i+1] - jmin;
476         if (nz > 0){
477           bcol = bj + jmin;
478           bval = ba + jmin;
479           while (nz --) rtmp[*bcol++] += uikdi*(*bval++);
480           /* update il and jl for i-th row */
481           il[i] = jmin;
482           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
483         }
484         i = nexti;
485       }
486 
487       /* shift the diagonals when zero pivot is detected */
488       /* compute rs=sum of abs(off-diagonal) */
489       rs   = 0.0;
490       jmin = bi[k]+1;
491       nz   = bi[k+1] - jmin;
492       if (nz){
493         bcol = bj + jmin;
494         while (nz--){
495           rs += PetscAbsScalar(rtmp[*bcol]);
496           bcol++;
497         }
498       }
499 
500       sctx.rs = rs;
501       sctx.pv = dk;
502       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
503       if (newshift == 1) break;    /* sctx.shift_amount is updated */
504 
505       /* copy data into U(k,:) */
506       ba[bi[k]] = 1.0/dk;
507       jmin      = bi[k]+1;
508       nz        = bi[k+1] - jmin;
509       if (nz){
510         bcol = bj + jmin;
511         bval = ba + jmin;
512         while (nz--){
513           *bval++       = rtmp[*bcol];
514           rtmp[*bcol++] = 0.0;
515         }
516         /* add k-th row into il and jl */
517         il[k] = jmin;
518         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
519       }
520     }
521   } while (sctx.chshift);
522   ierr = PetscFree(il);CHKERRQ(ierr);
523 
524   C->factor       = FACTOR_CHOLESKY;
525   C->assembled    = PETSC_TRUE;
526   C->preallocated = PETSC_TRUE;
527   ierr = PetscLogFlops(C->rmap.N);CHKERRQ(ierr);
528     if (sctx.nshift){
529     if (shiftnz) {
530       ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
531     } else if (shiftpd) {
532       ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
533     }
534   }
535   PetscFunctionReturn(0);
536 }
537 
538 #include "petscbt.h"
539 #include "src/mat/utils/freespace.h"
540 #undef __FUNCT__
541 #define __FUNCT__ "MatICCFactorSymbolic_SeqBAIJ"
542 PetscErrorCode MatICCFactorSymbolic_SeqBAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
543 {
544   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data;
545   Mat_SeqSBAIJ       *b;
546   Mat                B;
547   PetscErrorCode     ierr;
548   PetscTruth         perm_identity;
549   PetscInt           reallocs=0,*rip,i,*ai=a->i,*aj=a->j,am=a->mbs,bs=A->rmap.bs,*ui;
550   PetscInt           jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
551   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL,ncols,ncols_upper,*cols,*cols_lvl,*uj,**uj_ptr,**uj_lvl_ptr;
552   PetscReal          fill=info->fill,levels=info->levels;
553   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
554   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
555   PetscBT            lnkbt;
556 
557   PetscFunctionBegin;
558   if (bs > 1){
559     if (!a->sbaijMat){
560       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
561     }
562     ierr = MatICCFactorSymbolic(a->sbaijMat,perm,info,fact);CHKERRQ(ierr);
563     B = *fact;
564     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
565     PetscFunctionReturn(0);
566   }
567 
568   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
569   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
570 
571   /* special case that simply copies fill pattern */
572   if (!levels && perm_identity) {
573     ierr = MatMarkDiagonal_SeqBAIJ(A);CHKERRQ(ierr);
574     ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
575     for (i=0; i<am; i++) {
576       ui[i] = ai[i+1] - a->diag[i]; /* ui: rowlengths - changes when !perm_identity */
577     }
578     ierr = MatCreate(PETSC_COMM_SELF,fact);CHKERRQ(ierr);
579     ierr = MatSetSizes(*fact,am,am,am,am);CHKERRQ(ierr);
580     B = *fact;
581     ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
582     ierr = MatSeqSBAIJSetPreallocation(B,1,0,ui);CHKERRQ(ierr);
583 
584     b  = (Mat_SeqSBAIJ*)B->data;
585     uj = b->j;
586     for (i=0; i<am; i++) {
587       aj = a->j + a->diag[i];
588       for (j=0; j<ui[i]; j++){
589         *uj++ = *aj++;
590       }
591       b->ilen[i] = ui[i];
592     }
593     ierr = PetscFree(ui);CHKERRQ(ierr);
594     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
595     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
596 
597     B->ops->solve                 = MatSolve_SeqSBAIJ_1_NaturalOrdering;
598     B->ops->solvetranspose        = MatSolve_SeqSBAIJ_1_NaturalOrdering;
599     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
600     PetscFunctionReturn(0);
601   }
602 
603   /* initialization */
604   ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
605   ui[0] = 0;
606   ierr  = PetscMalloc((2*am+1)*sizeof(PetscInt),&cols_lvl);CHKERRQ(ierr);
607 
608   /* jl: linked list for storing indices of the pivot rows
609      il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
610   ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt*),&jl);CHKERRQ(ierr);
611   il         = jl + am;
612   uj_ptr     = (PetscInt**)(il + am);
613   uj_lvl_ptr = (PetscInt**)(uj_ptr + am);
614   for (i=0; i<am; i++){
615     jl[i] = am; il[i] = 0;
616   }
617 
618   /* create and initialize a linked list for storing column indices of the active row k */
619   nlnk = am + 1;
620   ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
621 
622   /* initial FreeSpace size is fill*(ai[am]+1) */
623   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
624   current_space = free_space;
625   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr);
626   current_space_lvl = free_space_lvl;
627 
628   for (k=0; k<am; k++){  /* for each active row k */
629     /* initialize lnk by the column indices of row rip[k] of A */
630     nzk   = 0;
631     ncols = ai[rip[k]+1] - ai[rip[k]];
632     ncols_upper = 0;
633     cols        = cols_lvl + am;
634     for (j=0; j<ncols; j++){
635       i = rip[*(aj + ai[rip[k]] + j)];
636       if (i >= k){ /* only take upper triangular entry */
637         cols[ncols_upper] = i;
638         cols_lvl[ncols_upper] = -1;  /* initialize level for nonzero entries */
639         ncols_upper++;
640       }
641     }
642     ierr = PetscIncompleteLLAdd(ncols_upper,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
643     nzk += nlnk;
644 
645     /* update lnk by computing fill-in for each pivot row to be merged in */
646     prow = jl[k]; /* 1st pivot row */
647 
648     while (prow < k){
649       nextprow = jl[prow];
650 
651       /* merge prow into k-th row */
652       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
653       jmax = ui[prow+1];
654       ncols = jmax-jmin;
655       i     = jmin - ui[prow];
656       cols = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
657       for (j=0; j<ncols; j++) cols_lvl[j] = *(uj_lvl_ptr[prow] + i + j);
658       ierr = PetscIncompleteLLAddSorted(ncols,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
659       nzk += nlnk;
660 
661       /* update il and jl for prow */
662       if (jmin < jmax){
663         il[prow] = jmin;
664         j = *cols; jl[prow] = jl[j]; jl[j] = prow;
665       }
666       prow = nextprow;
667     }
668 
669     /* if free space is not available, make more free space */
670     if (current_space->local_remaining<nzk) {
671       i = am - k + 1; /* num of unfactored rows */
672       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
673       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
674       ierr = PetscFreeSpaceGet(i,&current_space_lvl);CHKERRQ(ierr);
675       reallocs++;
676     }
677 
678     /* copy data into free_space and free_space_lvl, then initialize lnk */
679     ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
680 
681     /* add the k-th row into il and jl */
682     if (nzk-1 > 0){
683       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
684       jl[k] = jl[i]; jl[i] = k;
685       il[k] = ui[k] + 1;
686     }
687     uj_ptr[k]     = current_space->array;
688     uj_lvl_ptr[k] = current_space_lvl->array;
689 
690     current_space->array           += nzk;
691     current_space->local_used      += nzk;
692     current_space->local_remaining -= nzk;
693 
694     current_space_lvl->array           += nzk;
695     current_space_lvl->local_used      += nzk;
696     current_space_lvl->local_remaining -= nzk;
697 
698     ui[k+1] = ui[k] + nzk;
699   }
700 
701 #if defined(PETSC_USE_INFO)
702   if (ai[am] != 0) {
703     PetscReal af = ((PetscReal)(2*ui[am]-am))/((PetscReal)ai[am]);
704     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
705     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
706     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
707   } else {
708     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
709   }
710 #endif
711 
712   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
713   ierr = PetscFree(jl);CHKERRQ(ierr);
714   ierr = PetscFree(cols_lvl);CHKERRQ(ierr);
715 
716   /* destroy list of free space and other temporary array(s) */
717   ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
718   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
719   ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
720   ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
721 
722   /* put together the new matrix in MATSEQSBAIJ format */
723   ierr = MatCreate(PETSC_COMM_SELF,fact);CHKERRQ(ierr);
724   ierr = MatSetSizes(*fact,am,am,am,am);CHKERRQ(ierr);
725   B = *fact;
726   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
727   ierr = MatSeqSBAIJSetPreallocation(B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
728 
729   b = (Mat_SeqSBAIJ*)B->data;
730   b->singlemalloc = PETSC_FALSE;
731   b->free_a       = PETSC_TRUE;
732   b->free_ij       = PETSC_TRUE;
733   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
734   b->j    = uj;
735   b->i    = ui;
736   b->diag = 0;
737   b->ilen = 0;
738   b->imax = 0;
739   b->row  = perm;
740   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
741   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
742   b->icol = perm;
743   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
744   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
745   ierr    = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
746   b->maxnz = b->nz = ui[am];
747 
748   B->factor                 = FACTOR_CHOLESKY;
749   B->info.factor_mallocs    = reallocs;
750   B->info.fill_ratio_given  = fill;
751   if (ai[am] != 0) {
752     B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
753   } else {
754     B->info.fill_ratio_needed = 0.0;
755   }
756   if (perm_identity){
757     B->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
758     B->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
759     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
760   } else {
761     (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
762   }
763   PetscFunctionReturn(0);
764 }
765 
766 #undef __FUNCT__
767 #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqBAIJ"
768 PetscErrorCode MatCholeskyFactorSymbolic_SeqBAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
769 {
770   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data;
771   Mat_SeqSBAIJ       *b;
772   Mat                B;
773   PetscErrorCode     ierr;
774   PetscTruth         perm_identity;
775   PetscReal          fill = info->fill;
776   PetscInt           *rip,i,mbs=a->mbs,bs=A->rmap.bs,*ai=a->i,*aj=a->j,reallocs=0,prow;
777   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
778   PetscInt           nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr;
779   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
780   PetscBT            lnkbt;
781 
782   PetscFunctionBegin;
783   if (bs > 1) { /* convert to seqsbaij */
784     if (!a->sbaijMat){
785       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
786     }
787     ierr = MatCholeskyFactorSymbolic(a->sbaijMat,perm,info,fact);CHKERRQ(ierr);
788     B    = *fact;
789     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
790     PetscFunctionReturn(0);
791   }
792 
793   /* check whether perm is the identity mapping */
794   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
795   if (!perm_identity) SETERRQ(PETSC_ERR_SUP,"Matrix reordering is not supported");
796   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
797 
798   /* initialization */
799   ierr  = PetscMalloc((mbs+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
800   ui[0] = 0;
801 
802   /* jl: linked list for storing indices of the pivot rows
803      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
804   ierr = PetscMalloc((3*mbs+1)*sizeof(PetscInt)+mbs*sizeof(PetscInt*),&jl);CHKERRQ(ierr);
805   il     = jl + mbs;
806   cols   = il + mbs;
807   ui_ptr = (PetscInt**)(cols + mbs);
808   for (i=0; i<mbs; i++){
809     jl[i] = mbs; il[i] = 0;
810   }
811 
812   /* create and initialize a linked list for storing column indices of the active row k */
813   nlnk = mbs + 1;
814   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
815 
816   /* initial FreeSpace size is fill*(ai[mbs]+1) */
817   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[mbs]+1)),&free_space);CHKERRQ(ierr);
818   current_space = free_space;
819 
820   for (k=0; k<mbs; k++){  /* for each active row k */
821     /* initialize lnk by the column indices of row rip[k] of A */
822     nzk   = 0;
823     ncols = ai[rip[k]+1] - ai[rip[k]];
824     ncols_upper = 0;
825     for (j=0; j<ncols; j++){
826       i = rip[*(aj + ai[rip[k]] + j)];
827       if (i >= k){ /* only take upper triangular entry */
828         cols[ncols_upper] = i;
829         ncols_upper++;
830       }
831     }
832     ierr = PetscLLAdd(ncols_upper,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
833     nzk += nlnk;
834 
835     /* update lnk by computing fill-in for each pivot row to be merged in */
836     prow = jl[k]; /* 1st pivot row */
837 
838     while (prow < k){
839       nextprow = jl[prow];
840       /* merge prow into k-th row */
841       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
842       jmax = ui[prow+1];
843       ncols = jmax-jmin;
844       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
845       ierr = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
846       nzk += nlnk;
847 
848       /* update il and jl for prow */
849       if (jmin < jmax){
850         il[prow] = jmin;
851         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
852       }
853       prow = nextprow;
854     }
855 
856     /* if free space is not available, make more free space */
857     if (current_space->local_remaining<nzk) {
858       i = mbs - k + 1; /* num of unfactored rows */
859       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
860       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
861       reallocs++;
862     }
863 
864     /* copy data into free space, then initialize lnk */
865     ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
866 
867     /* add the k-th row into il and jl */
868     if (nzk-1 > 0){
869       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
870       jl[k] = jl[i]; jl[i] = k;
871       il[k] = ui[k] + 1;
872     }
873     ui_ptr[k] = current_space->array;
874     current_space->array           += nzk;
875     current_space->local_used      += nzk;
876     current_space->local_remaining -= nzk;
877 
878     ui[k+1] = ui[k] + nzk;
879   }
880 
881 #if defined(PETSC_USE_INFO)
882   if (ai[mbs] != 0) {
883     PetscReal af = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
884     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
885     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
886     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
887   } else {
888     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
889   }
890 #endif
891 
892   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
893   ierr = PetscFree(jl);CHKERRQ(ierr);
894 
895   /* destroy list of free space and other temporary array(s) */
896   ierr = PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
897   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
898   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
899 
900   /* put together the new matrix in MATSEQSBAIJ format */
901   ierr = MatCreate(PETSC_COMM_SELF,fact);CHKERRQ(ierr);
902   ierr = MatSetSizes(*fact,mbs,mbs,mbs,mbs);CHKERRQ(ierr);
903   B    = *fact;
904   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
905   ierr = MatSeqSBAIJSetPreallocation(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
906 
907   b = (Mat_SeqSBAIJ*)B->data;
908   b->singlemalloc = PETSC_FALSE;
909   b->free_a       = PETSC_TRUE;
910   b->free_ij      = PETSC_TRUE;
911   ierr = PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
912   b->j    = uj;
913   b->i    = ui;
914   b->diag = 0;
915   b->ilen = 0;
916   b->imax = 0;
917   b->row  = perm;
918   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
919   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
920   b->icol = perm;
921   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
922   ierr    = PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
923   ierr    = PetscLogObjectMemory(B,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
924   b->maxnz = b->nz = ui[mbs];
925 
926   B->factor                 = FACTOR_CHOLESKY;
927   B->info.factor_mallocs    = reallocs;
928   B->info.fill_ratio_given  = fill;
929   if (ai[mbs] != 0) {
930     B->info.fill_ratio_needed = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
931   } else {
932     B->info.fill_ratio_needed = 0.0;
933   }
934   if (perm_identity){
935     B->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
936     B->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
937     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
938   } else {
939     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
940   }
941   PetscFunctionReturn(0);
942 }
943