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