Actual source code: baijfact.c
2: /*
3: Factorization code for BAIJ format.
4: */
5: #include <../src/mat/impls/baij/seq/baij.h>
6: #include <../src/mat/blockinvert.h>
10: PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2(Mat B,Mat A,const MatFactorInfo *info)
11: {
12: Mat C=B;
13: Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ *)C->data;
14: IS isrow = b->row,isicol = b->icol;
15: PetscErrorCode ierr;
16: const PetscInt *r,*ic;
17: PetscInt i,j,k,nz,nzL,row,*pj;
18: const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,bs2=a->bs2;
19: const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag;
20: MatScalar *rtmp,*pc,*mwork,*pv;
21: MatScalar *aa=a->a,*v;
22: PetscInt flg;
23: PetscReal shift = info->shiftamount;
26: ISGetIndices(isrow,&r);
27: ISGetIndices(isicol,&ic);
29: /* generate work space needed by the factorization */
30: PetscMalloc2(bs2*n,MatScalar,&rtmp,bs2,MatScalar,&mwork);
31: PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));
33: for (i=0; i<n; i++){
34: /* zero rtmp */
35: /* L part */
36: nz = bi[i+1] - bi[i];
37: bjtmp = bj + bi[i];
38: for (j=0; j<nz; j++){
39: PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));
40: }
42: /* U part */
43: nz = bdiag[i] - bdiag[i+1];
44: bjtmp = bj + bdiag[i+1]+1;
45: for (j=0; j<nz; j++){
46: PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));
47: }
48:
49: /* load in initial (unfactored row) */
50: nz = ai[r[i]+1] - ai[r[i]];
51: ajtmp = aj + ai[r[i]];
52: v = aa + bs2*ai[r[i]];
53: for (j=0; j<nz; j++) {
54: PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],v+bs2*j,bs2*sizeof(MatScalar));
55: }
57: /* elimination */
58: bjtmp = bj + bi[i];
59: nzL = bi[i+1] - bi[i];
60: for(k=0;k < nzL;k++) {
61: row = bjtmp[k];
62: pc = rtmp + bs2*row;
63: for (flg=0,j=0; j<bs2; j++) { if (pc[j] != (PetscScalar)0.0) { flg = 1; break; }}
64: if (flg) {
65: pv = b->a + bs2*bdiag[row];
66: /* Kernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
67: Kernel_A_gets_A_times_B_2(pc,pv,mwork);
68:
69: pj = b->j + bdiag[row+1]+1; /* begining of U(row,:) */
70: pv = b->a + bs2*(bdiag[row+1]+1);
71: nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries inU(row,:), excluding diag */
72: for (j=0; j<nz; j++) {
73: /* Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
74: /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
75: v = rtmp + 4*pj[j];
76: Kernel_A_gets_A_minus_B_times_C_2(v,pc,pv);
77: pv += 4;
78: }
79: PetscLogFlops(16*nz+12); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
80: }
81: }
83: /* finished row so stick it into b->a */
84: /* L part */
85: pv = b->a + bs2*bi[i] ;
86: pj = b->j + bi[i] ;
87: nz = bi[i+1] - bi[i];
88: for (j=0; j<nz; j++) {
89: PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));
90: }
92: /* Mark diagonal and invert diagonal for simplier triangular solves */
93: pv = b->a + bs2*bdiag[i];
94: pj = b->j + bdiag[i];
95: PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));
96: /* Kernel_A_gets_inverse_A(bs,pv,v_pivots,v_work); */
97: Kernel_A_gets_inverse_A_2(pv,shift);
98:
99: /* U part */
100: pv = b->a + bs2*(bdiag[i+1]+1);
101: pj = b->j + bdiag[i+1]+1;
102: nz = bdiag[i] - bdiag[i+1] - 1;
103: for (j=0; j<nz; j++){
104: PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));
105: }
106: }
108: PetscFree2(rtmp,mwork);
109: ISRestoreIndices(isicol,&ic);
110: ISRestoreIndices(isrow,&r);
111: C->ops->solve = MatSolve_SeqBAIJ_2;
112: C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2;
113:
114: C->assembled = PETSC_TRUE;
115: PetscLogFlops(1.333333333333*2*2*2*n); /* from inverting diagonal blocks */
116: return(0);
117: }
121: PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info)
122: {
123: Mat C=B;
124: Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ *)C->data;
125: PetscErrorCode ierr;
126: PetscInt i,j,k,nz,nzL,row,*pj;
127: const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,bs2=a->bs2;
128: const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag;
129: MatScalar *rtmp,*pc,*mwork,*pv;
130: MatScalar *aa=a->a,*v;
131: PetscInt flg;
132: PetscReal shift = info->shiftamount;
135: /* generate work space needed by the factorization */
136: PetscMalloc2(bs2*n,MatScalar,&rtmp,bs2,MatScalar,&mwork);
137: PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));
139: for (i=0; i<n; i++){
140: /* zero rtmp */
141: /* L part */
142: nz = bi[i+1] - bi[i];
143: bjtmp = bj + bi[i];
144: for (j=0; j<nz; j++){
145: PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));
146: }
148: /* U part */
149: nz = bdiag[i] - bdiag[i+1];
150: bjtmp = bj + bdiag[i+1]+1;
151: for (j=0; j<nz; j++){
152: PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));
153: }
154:
155: /* load in initial (unfactored row) */
156: nz = ai[i+1] - ai[i];
157: ajtmp = aj + ai[i];
158: v = aa + bs2*ai[i];
159: for (j=0; j<nz; j++) {
160: PetscMemcpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2*sizeof(MatScalar));
161: }
163: /* elimination */
164: bjtmp = bj + bi[i];
165: nzL = bi[i+1] - bi[i];
166: for(k=0;k < nzL;k++) {
167: row = bjtmp[k];
168: pc = rtmp + bs2*row;
169: for (flg=0,j=0; j<bs2; j++) { if (pc[j]!=(PetscScalar)0.0) { flg = 1; break; }}
170: if (flg) {
171: pv = b->a + bs2*bdiag[row];
172: /* Kernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
173: Kernel_A_gets_A_times_B_2(pc,pv,mwork);
174:
175: pj = b->j + bdiag[row+1]+1; /* beginning of U(row,:) */
176: pv = b->a + bs2*(bdiag[row+1]+1);
177: nz = bdiag[row]-bdiag[row+1] - 1; /* num of entries in U(row,:) excluding diag */
178: for (j=0; j<nz; j++) {
179: /* Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
180: /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
181: v = rtmp + 4*pj[j];
182: Kernel_A_gets_A_minus_B_times_C_2(v,pc,pv);
183: pv += 4;
184: }
185: PetscLogFlops(16*nz+12); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
186: }
187: }
189: /* finished row so stick it into b->a */
190: /* L part */
191: pv = b->a + bs2*bi[i] ;
192: pj = b->j + bi[i] ;
193: nz = bi[i+1] - bi[i];
194: for (j=0; j<nz; j++) {
195: PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));
196: }
198: /* Mark diagonal and invert diagonal for simplier triangular solves */
199: pv = b->a + bs2*bdiag[i];
200: pj = b->j + bdiag[i];
201: PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));
202: /* Kernel_A_gets_inverse_A(bs,pv,v_pivots,v_work); */
203: Kernel_A_gets_inverse_A_2(pv,shift);
204:
205: /* U part */
206: /*
207: pv = b->a + bs2*bi[2*n-i];
208: pj = b->j + bi[2*n-i];
209: nz = bi[2*n-i+1] - bi[2*n-i] - 1;
210: */
211: pv = b->a + bs2*(bdiag[i+1]+1);
212: pj = b->j + bdiag[i+1]+1;
213: nz = bdiag[i] - bdiag[i+1] - 1;
214: for (j=0; j<nz; j++){
215: PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));
216: }
217: }
218: PetscFree2(rtmp,mwork);
220: C->ops->solve = MatSolve_SeqBAIJ_2_NaturalOrdering;
221: C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_NaturalOrdering;
222: C->assembled = PETSC_TRUE;
223: PetscLogFlops(1.333333333333*2*2*2*n); /* from inverting diagonal blocks */
224: return(0);
225: }
229: PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_inplace(Mat B,Mat A,const MatFactorInfo *info)
230: {
231: Mat C = B;
232: Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data;
233: IS isrow = b->row,isicol = b->icol;
235: const PetscInt *r,*ic;
236: PetscInt i,j,n = a->mbs,*bi = b->i,*bj = b->j;
237: PetscInt *ajtmpold,*ajtmp,nz,row;
238: PetscInt *diag_offset=b->diag,idx,*ai=a->i,*aj=a->j,*pj;
239: MatScalar *pv,*v,*rtmp,m1,m2,m3,m4,*pc,*w,*x,x1,x2,x3,x4;
240: MatScalar p1,p2,p3,p4;
241: MatScalar *ba = b->a,*aa = a->a;
242: PetscReal shift = info->shiftamount;
245: ISGetIndices(isrow,&r);
246: ISGetIndices(isicol,&ic);
247: PetscMalloc(4*(n+1)*sizeof(MatScalar),&rtmp);
249: for (i=0; i<n; i++) {
250: nz = bi[i+1] - bi[i];
251: ajtmp = bj + bi[i];
252: for (j=0; j<nz; j++) {
253: x = rtmp+4*ajtmp[j]; x[0] = x[1] = x[2] = x[3] = 0.0;
254: }
255: /* load in initial (unfactored row) */
256: idx = r[i];
257: nz = ai[idx+1] - ai[idx];
258: ajtmpold = aj + ai[idx];
259: v = aa + 4*ai[idx];
260: for (j=0; j<nz; j++) {
261: x = rtmp+4*ic[ajtmpold[j]];
262: x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3];
263: v += 4;
264: }
265: row = *ajtmp++;
266: while (row < i) {
267: pc = rtmp + 4*row;
268: p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3];
269: if (p1 != (PetscScalar)0.0 || p2 != (PetscScalar)0.0 || p3 != (PetscScalar)0.0 || p4 != (PetscScalar)0.0) {
270: pv = ba + 4*diag_offset[row];
271: pj = bj + diag_offset[row] + 1;
272: x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
273: pc[0] = m1 = p1*x1 + p3*x2;
274: pc[1] = m2 = p2*x1 + p4*x2;
275: pc[2] = m3 = p1*x3 + p3*x4;
276: pc[3] = m4 = p2*x3 + p4*x4;
277: nz = bi[row+1] - diag_offset[row] - 1;
278: pv += 4;
279: for (j=0; j<nz; j++) {
280: x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
281: x = rtmp + 4*pj[j];
282: x[0] -= m1*x1 + m3*x2;
283: x[1] -= m2*x1 + m4*x2;
284: x[2] -= m1*x3 + m3*x4;
285: x[3] -= m2*x3 + m4*x4;
286: pv += 4;
287: }
288: PetscLogFlops(16.0*nz+12.0);
289: }
290: row = *ajtmp++;
291: }
292: /* finished row so stick it into b->a */
293: pv = ba + 4*bi[i];
294: pj = bj + bi[i];
295: nz = bi[i+1] - bi[i];
296: for (j=0; j<nz; j++) {
297: x = rtmp+4*pj[j];
298: pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3];
299: pv += 4;
300: }
301: /* invert diagonal block */
302: w = ba + 4*diag_offset[i];
303: Kernel_A_gets_inverse_A_2(w,shift);
304: }
306: PetscFree(rtmp);
307: ISRestoreIndices(isicol,&ic);
308: ISRestoreIndices(isrow,&r);
309: C->ops->solve = MatSolve_SeqBAIJ_2_inplace;
310: C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_inplace;
311: C->assembled = PETSC_TRUE;
312: PetscLogFlops(1.333333333333*8*b->mbs); /* from inverting diagonal blocks */
313: return(0);
314: }
315: /*
316: Version for when blocks are 2 by 2 Using natural ordering
317: */
320: PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_inplace(Mat C,Mat A,const MatFactorInfo *info)
321: {
322: Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data;
324: PetscInt i,j,n = a->mbs,*bi = b->i,*bj = b->j;
325: PetscInt *ajtmpold,*ajtmp,nz,row;
326: PetscInt *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj;
327: MatScalar *pv,*v,*rtmp,*pc,*w,*x;
328: MatScalar p1,p2,p3,p4,m1,m2,m3,m4,x1,x2,x3,x4;
329: MatScalar *ba = b->a,*aa = a->a;
330: PetscReal shift = info->shiftamount;
333: PetscMalloc(4*(n+1)*sizeof(MatScalar),&rtmp);
334: for (i=0; i<n; i++) {
335: nz = bi[i+1] - bi[i];
336: ajtmp = bj + bi[i];
337: for (j=0; j<nz; j++) {
338: x = rtmp+4*ajtmp[j];
339: x[0] = x[1] = x[2] = x[3] = 0.0;
340: }
341: /* load in initial (unfactored row) */
342: nz = ai[i+1] - ai[i];
343: ajtmpold = aj + ai[i];
344: v = aa + 4*ai[i];
345: for (j=0; j<nz; j++) {
346: x = rtmp+4*ajtmpold[j];
347: x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3];
348: v += 4;
349: }
350: row = *ajtmp++;
351: while (row < i) {
352: pc = rtmp + 4*row;
353: p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3];
354: if (p1 != (PetscScalar)0.0 || p2 != (PetscScalar)0.0 || p3 != (PetscScalar)0.0 || p4 != (PetscScalar)0.0) {
355: pv = ba + 4*diag_offset[row];
356: pj = bj + diag_offset[row] + 1;
357: x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
358: pc[0] = m1 = p1*x1 + p3*x2;
359: pc[1] = m2 = p2*x1 + p4*x2;
360: pc[2] = m3 = p1*x3 + p3*x4;
361: pc[3] = m4 = p2*x3 + p4*x4;
362: nz = bi[row+1] - diag_offset[row] - 1;
363: pv += 4;
364: for (j=0; j<nz; j++) {
365: x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
366: x = rtmp + 4*pj[j];
367: x[0] -= m1*x1 + m3*x2;
368: x[1] -= m2*x1 + m4*x2;
369: x[2] -= m1*x3 + m3*x4;
370: x[3] -= m2*x3 + m4*x4;
371: pv += 4;
372: }
373: PetscLogFlops(16.0*nz+12.0);
374: }
375: row = *ajtmp++;
376: }
377: /* finished row so stick it into b->a */
378: pv = ba + 4*bi[i];
379: pj = bj + bi[i];
380: nz = bi[i+1] - bi[i];
381: for (j=0; j<nz; j++) {
382: x = rtmp+4*pj[j];
383: pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3];
384: /*
385: printf(" col %d:",pj[j]);
386: PetscInt j1;
387: for (j1=0; j1<4; j1++) printf(" %g,",*(pv+j1));
388: printf("\n");
389: */
390: pv += 4;
391: }
392: /* invert diagonal block */
393: w = ba + 4*diag_offset[i];
394: /*
395: printf(" \n%d -th: diag: ",i);
396: for (j=0; j<4; j++){
397: printf(" %g,",w[j]);
398: }
399: printf("\n----------------------------\n");
400: */
401: Kernel_A_gets_inverse_A_2(w,shift);
402: }
404: PetscFree(rtmp);
405: C->ops->solve = MatSolve_SeqBAIJ_2_NaturalOrdering_inplace;
406: C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_NaturalOrdering_inplace;
407: C->assembled = PETSC_TRUE;
408: PetscLogFlops(1.333333333333*8*b->mbs); /* from inverting diagonal blocks */
409: return(0);
410: }
412: /* ----------------------------------------------------------- */
413: /*
414: Version for when blocks are 1 by 1.
415: */
418: PetscErrorCode MatLUFactorNumeric_SeqBAIJ_1(Mat B,Mat A,const MatFactorInfo *info)
419: {
420: Mat C=B;
421: Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ *)C->data;
422: IS isrow = b->row,isicol = b->icol;
423: PetscErrorCode ierr;
424: const PetscInt *r,*ic,*ics;
425: const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bdiag=b->diag;
426: PetscInt i,j,k,nz,nzL,row,*pj;
427: const PetscInt *ajtmp,*bjtmp;
428: MatScalar *rtmp,*pc,multiplier,*pv;
429: const MatScalar *aa=a->a,*v;
430: PetscBool row_identity,col_identity;
431: FactorShiftCtx sctx;
432: const PetscInt *ddiag;
433: PetscReal rs;
434: MatScalar d;
437: /* MatPivotSetUp(): initialize shift context sctx */
438: PetscMemzero(&sctx,sizeof(FactorShiftCtx));
440: if (info->shifttype == (PetscReal) MAT_SHIFT_POSITIVE_DEFINITE) { /* set sctx.shift_top=max{rs} */
441: ddiag = a->diag;
442: sctx.shift_top = info->zeropivot;
443: for (i=0; i<n; i++) {
444: /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
445: d = (aa)[ddiag[i]];
446: rs = -PetscAbsScalar(d) - PetscRealPart(d);
447: v = aa+ai[i];
448: nz = ai[i+1] - ai[i];
449: for (j=0; j<nz; j++)
450: rs += PetscAbsScalar(v[j]);
451: if (rs>sctx.shift_top) sctx.shift_top = rs;
452: }
453: sctx.shift_top *= 1.1;
454: sctx.nshift_max = 5;
455: sctx.shift_lo = 0.;
456: sctx.shift_hi = 1.;
457: }
459: ISGetIndices(isrow,&r);
460: ISGetIndices(isicol,&ic);
461: PetscMalloc((n+1)*sizeof(MatScalar),&rtmp);
462: ics = ic;
464: do {
465: sctx.newshift = PETSC_FALSE;
466: for (i=0; i<n; i++){
467: /* zero rtmp */
468: /* L part */
469: nz = bi[i+1] - bi[i];
470: bjtmp = bj + bi[i];
471: for (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0;
473: /* U part */
474: nz = bdiag[i]-bdiag[i+1];
475: bjtmp = bj + bdiag[i+1]+1;
476: for (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0;
477:
478: /* load in initial (unfactored row) */
479: nz = ai[r[i]+1] - ai[r[i]];
480: ajtmp = aj + ai[r[i]];
481: v = aa + ai[r[i]];
482: for (j=0; j<nz; j++) {
483: rtmp[ics[ajtmp[j]]] = v[j];
484: }
485: /* ZeropivotApply() */
486: rtmp[i] += sctx.shift_amount; /* shift the diagonal of the matrix */
487:
488: /* elimination */
489: bjtmp = bj + bi[i];
490: row = *bjtmp++;
491: nzL = bi[i+1] - bi[i];
492: for(k=0; k < nzL;k++) {
493: pc = rtmp + row;
494: if (*pc != (PetscScalar)0.0) {
495: pv = b->a + bdiag[row];
496: multiplier = *pc * (*pv);
497: *pc = multiplier;
498: pj = b->j + bdiag[row+1]+1; /* beginning of U(row,:) */
499: pv = b->a + bdiag[row+1]+1;
500: nz = bdiag[row]-bdiag[row+1]-1; /* num of entries in U(row,:) excluding diag */
501: for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j];
502: PetscLogFlops(2.0*nz);
503: }
504: row = *bjtmp++;
505: }
507: /* finished row so stick it into b->a */
508: rs = 0.0;
509: /* L part */
510: pv = b->a + bi[i] ;
511: pj = b->j + bi[i] ;
512: nz = bi[i+1] - bi[i];
513: for (j=0; j<nz; j++) {
514: pv[j] = rtmp[pj[j]]; rs += PetscAbsScalar(pv[j]);
515: }
517: /* U part */
518: pv = b->a + bdiag[i+1]+1;
519: pj = b->j + bdiag[i+1]+1;
520: nz = bdiag[i] - bdiag[i+1]-1;
521: for (j=0; j<nz; j++) {
522: pv[j] = rtmp[pj[j]]; rs += PetscAbsScalar(pv[j]);
523: }
525: sctx.rs = rs;
526: sctx.pv = rtmp[i];
527: MatPivotCheck(A,info,&sctx,i);
528: if(sctx.newshift) break; /* break for-loop */
529: rtmp[i] = sctx.pv; /* sctx.pv might be updated in the case of MAT_SHIFT_INBLOCKS */
531: /* Mark diagonal and invert diagonal for simplier triangular solves */
532: pv = b->a + bdiag[i];
533: *pv = (PetscScalar)1.0/rtmp[i];
535: } /* endof for (i=0; i<n; i++){ */
537: /* MatPivotRefine() */
538: if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE && !sctx.newshift && sctx.shift_fraction>0 && sctx.nshift<sctx.nshift_max){
539: /*
540: * if no shift in this attempt & shifting & started shifting & can refine,
541: * then try lower shift
542: */
543: sctx.shift_hi = sctx.shift_fraction;
544: sctx.shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.;
545: sctx.shift_amount = sctx.shift_fraction * sctx.shift_top;
546: sctx.newshift = PETSC_TRUE;
547: sctx.nshift++;
548: }
549: } while (sctx.newshift);
551: PetscFree(rtmp);
552: ISRestoreIndices(isicol,&ic);
553: ISRestoreIndices(isrow,&r);
554:
555: ISIdentity(isrow,&row_identity);
556: ISIdentity(isicol,&col_identity);
557: if (row_identity && col_identity) {
558: C->ops->solve = MatSolve_SeqBAIJ_1_NaturalOrdering;
559: C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1_NaturalOrdering;
560: } else {
561: C->ops->solve = MatSolve_SeqBAIJ_1;
562: C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1;
563: }
564: C->assembled = PETSC_TRUE;
565: PetscLogFlops(C->cmap->n);
567: /* MatShiftView(A,info,&sctx) */
568: if (sctx.nshift){
569: if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
570: PetscInfo4(A,"number of shift_pd tries %D, shift_amount %G, diagonal shifted up by %e fraction top_value %e\n",sctx.nshift,sctx.shift_amount,sctx.shift_fraction,sctx.shift_top);
571: } else if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
572: PetscInfo2(A,"number of shift_nz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);
573: } else if (info->shifttype == (PetscReal)MAT_SHIFT_INBLOCKS){
574: PetscInfo2(A,"number of shift_inblocks applied %D, each shift_amount %G\n",sctx.nshift,info->shiftamount);
575: }
576: }
577: return(0);
578: }
582: PetscErrorCode MatLUFactorNumeric_SeqBAIJ_1_inplace(Mat C,Mat A,const MatFactorInfo *info)
583: {
584: Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data;
585: IS isrow = b->row,isicol = b->icol;
587: const PetscInt *r,*ic;
588: PetscInt i,j,n = a->mbs,*bi = b->i,*bj = b->j;
589: PetscInt *ajtmpold,*ajtmp,nz,row,*ai = a->i,*aj = a->j;
590: PetscInt *diag_offset = b->diag,diag,*pj;
591: MatScalar *pv,*v,*rtmp,multiplier,*pc;
592: MatScalar *ba = b->a,*aa = a->a;
593: PetscBool row_identity, col_identity;
596: ISGetIndices(isrow,&r);
597: ISGetIndices(isicol,&ic);
598: PetscMalloc((n+1)*sizeof(MatScalar),&rtmp);
600: for (i=0; i<n; i++) {
601: nz = bi[i+1] - bi[i];
602: ajtmp = bj + bi[i];
603: for (j=0; j<nz; j++) rtmp[ajtmp[j]] = 0.0;
605: /* load in initial (unfactored row) */
606: nz = ai[r[i]+1] - ai[r[i]];
607: ajtmpold = aj + ai[r[i]];
608: v = aa + ai[r[i]];
609: for (j=0; j<nz; j++) rtmp[ic[ajtmpold[j]]] = v[j];
611: row = *ajtmp++;
612: while (row < i) {
613: pc = rtmp + row;
614: if (*pc != 0.0) {
615: pv = ba + diag_offset[row];
616: pj = bj + diag_offset[row] + 1;
617: multiplier = *pc * *pv++;
618: *pc = multiplier;
619: nz = bi[row+1] - diag_offset[row] - 1;
620: for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j];
621: PetscLogFlops(1.0+2.0*nz);
622: }
623: row = *ajtmp++;
624: }
625: /* finished row so stick it into b->a */
626: pv = ba + bi[i];
627: pj = bj + bi[i];
628: nz = bi[i+1] - bi[i];
629: for (j=0; j<nz; j++) {pv[j] = rtmp[pj[j]];}
630: diag = diag_offset[i] - bi[i];
631: /* check pivot entry for current row */
632: if (pv[diag] == 0.0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot: row in original ordering %D in permuted ordering %D",r[i],i);
633: pv[diag] = 1.0/pv[diag];
634: }
636: PetscFree(rtmp);
637: ISRestoreIndices(isicol,&ic);
638: ISRestoreIndices(isrow,&r);
639: ISIdentity(isrow,&row_identity);
640: ISIdentity(isicol,&col_identity);
641: if (row_identity && col_identity) {
642: C->ops->solve = MatSolve_SeqBAIJ_1_NaturalOrdering_inplace;
643: C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1_NaturalOrdering_inplace;
644: } else {
645: C->ops->solve = MatSolve_SeqBAIJ_1_inplace;
646: C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1_inplace;
647: }
648: C->assembled = PETSC_TRUE;
649: PetscLogFlops(C->cmap->n);
650: return(0);
651: }
656: PetscErrorCode MatGetFactor_seqbaij_petsc(Mat A,MatFactorType ftype,Mat *B)
657: {
658: PetscInt n = A->rmap->n;
659: PetscErrorCode ierr;
662: MatCreate(((PetscObject)A)->comm,B);
663: MatSetSizes(*B,n,n,n,n);
664: if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU || ftype == MAT_FACTOR_ILUDT) {
665: MatSetType(*B,MATSEQBAIJ);
666: (*B)->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqBAIJ;
667: (*B)->ops->ilufactorsymbolic = MatILUFactorSymbolic_SeqBAIJ;
668: } else if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) {
669: MatSetType(*B,MATSEQSBAIJ);
670: MatSeqSBAIJSetPreallocation(*B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);
671: (*B)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqBAIJ;
672: (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqBAIJ;
673: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
674: (*B)->factortype = ftype;
675: return(0);
676: }
682: PetscErrorCode MatGetFactorAvailable_seqbaij_petsc(Mat A,MatFactorType ftype,PetscBool *flg)
683: {
685: *flg = PETSC_TRUE;
686: return(0);
687: }
690: /* ----------------------------------------------------------- */
693: PetscErrorCode MatLUFactor_SeqBAIJ(Mat A,IS row,IS col,const MatFactorInfo *info)
694: {
696: Mat C;
699: MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_LU,&C);
700: MatLUFactorSymbolic(C,A,row,col,info);
701: MatLUFactorNumeric(C,A,info);
702: A->ops->solve = C->ops->solve;
703: A->ops->solvetranspose = C->ops->solvetranspose;
704: MatHeaderMerge(A,C);
705: PetscLogObjectParent(A,((Mat_SeqBAIJ*)(A->data))->icol);
706: return(0);
707: }
709: #include <../src/mat/impls/sbaij/seq/sbaij.h>
712: PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N(Mat C,Mat A,const MatFactorInfo *info)
713: {
715: Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data;
716: Mat_SeqSBAIJ *b=(Mat_SeqSBAIJ*)C->data;
717: IS ip=b->row;
718: const PetscInt *rip;
719: PetscInt i,j,mbs=a->mbs,bs=A->rmap->bs,*bi=b->i,*bj=b->j,*bcol;
720: PetscInt *ai=a->i,*aj=a->j;
721: PetscInt k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
722: MatScalar *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
723: PetscReal rs;
724: FactorShiftCtx sctx;
727: if (bs > 1) {/* convert A to a SBAIJ matrix and apply Cholesky factorization from it */
728: if (!a->sbaijMat){
729: MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);
730: }
731: (a->sbaijMat)->ops->choleskyfactornumeric(C,a->sbaijMat,info);
732: MatDestroy(&a->sbaijMat);
733: return(0);
734: }
735:
736: /* MatPivotSetUp(): initialize shift context sctx */
737: PetscMemzero(&sctx,sizeof(FactorShiftCtx));
739: ISGetIndices(ip,&rip);
740: PetscMalloc3(mbs,MatScalar,&rtmp,mbs,PetscInt,&il,mbs,PetscInt,&jl);
742: sctx.shift_amount = 0.;
743: sctx.nshift = 0;
744: do {
745: sctx.newshift = PETSC_FALSE;
746: for (i=0; i<mbs; i++) {
747: rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
748: }
749:
750: for (k = 0; k<mbs; k++){
751: bval = ba + bi[k];
752: /* initialize k-th row by the perm[k]-th row of A */
753: jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
754: for (j = jmin; j < jmax; j++){
755: col = rip[aj[j]];
756: if (col >= k){ /* only take upper triangular entry */
757: rtmp[col] = aa[j];
758: *bval++ = 0.0; /* for in-place factorization */
759: }
760: }
761:
762: /* shift the diagonal of the matrix */
763: if (sctx.nshift) rtmp[k] += sctx.shift_amount;
765: /* modify k-th row by adding in those rows i with U(i,k)!=0 */
766: dk = rtmp[k];
767: i = jl[k]; /* first row to be added to k_th row */
769: while (i < k){
770: nexti = jl[i]; /* next row to be added to k_th row */
772: /* compute multiplier, update diag(k) and U(i,k) */
773: ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
774: uikdi = - ba[ili]*ba[bi[i]]; /* diagonal(k) */
775: dk += uikdi*ba[ili];
776: ba[ili] = uikdi; /* -U(i,k) */
778: /* add multiple of row i to k-th row */
779: jmin = ili + 1; jmax = bi[i+1];
780: if (jmin < jmax){
781: for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
782: /* update il and jl for row i */
783: il[i] = jmin;
784: j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
785: }
786: i = nexti;
787: }
789: /* shift the diagonals when zero pivot is detected */
790: /* compute rs=sum of abs(off-diagonal) */
791: rs = 0.0;
792: jmin = bi[k]+1;
793: nz = bi[k+1] - jmin;
794: if (nz){
795: bcol = bj + jmin;
796: while (nz--){
797: rs += PetscAbsScalar(rtmp[*bcol]);
798: bcol++;
799: }
800: }
802: sctx.rs = rs;
803: sctx.pv = dk;
804: MatPivotCheck(A,info,&sctx,k);
805: if (sctx.newshift) break;
806: dk = sctx.pv;
808: /* copy data into U(k,:) */
809: ba[bi[k]] = 1.0/dk; /* U(k,k) */
810: jmin = bi[k]+1; jmax = bi[k+1];
811: if (jmin < jmax) {
812: for (j=jmin; j<jmax; j++){
813: col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
814: }
815: /* add the k-th row into il and jl */
816: il[k] = jmin;
817: i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
818: }
819: }
820: } while (sctx.newshift);
821: PetscFree3(rtmp,il,jl);
823: ISRestoreIndices(ip,&rip);
824: C->assembled = PETSC_TRUE;
825: C->preallocated = PETSC_TRUE;
826: PetscLogFlops(C->rmap->N);
827: if (sctx.nshift){
828: if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
829: PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);
830: } else if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
831: PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);
832: }
833: }
834: return(0);
835: }
839: PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info)
840: {
841: Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data;
842: Mat_SeqSBAIJ *b=(Mat_SeqSBAIJ*)C->data;
844: PetscInt i,j,am=a->mbs;
845: PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
846: PetscInt k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz;
847: MatScalar *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval;
848: PetscReal rs;
849: FactorShiftCtx sctx;
852: /* MatPivotSetUp(): initialize shift context sctx */
853: PetscMemzero(&sctx,sizeof(FactorShiftCtx));
854:
855: PetscMalloc3(am,MatScalar,&rtmp,am,PetscInt,&il,am,PetscInt,&jl);
857: do {
858: sctx.newshift = PETSC_FALSE;
859: for (i=0; i<am; i++) {
860: rtmp[i] = 0.0; jl[i] = am; il[0] = 0;
861: }
863: for (k = 0; k<am; k++){
864: /* initialize k-th row with elements nonzero in row perm(k) of A */
865: nz = ai[k+1] - ai[k];
866: acol = aj + ai[k];
867: aval = aa + ai[k];
868: bval = ba + bi[k];
869: while (nz -- ){
870: if (*acol < k) { /* skip lower triangular entries */
871: acol++; aval++;
872: } else {
873: rtmp[*acol++] = *aval++;
874: *bval++ = 0.0; /* for in-place factorization */
875: }
876: }
877:
878: /* shift the diagonal of the matrix */
879: if (sctx.nshift) rtmp[k] += sctx.shift_amount;
880:
881: /* modify k-th row by adding in those rows i with U(i,k)!=0 */
882: dk = rtmp[k];
883: i = jl[k]; /* first row to be added to k_th row */
885: while (i < k){
886: nexti = jl[i]; /* next row to be added to k_th row */
887: /* compute multiplier, update D(k) and U(i,k) */
888: ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
889: uikdi = - ba[ili]*ba[bi[i]];
890: dk += uikdi*ba[ili];
891: ba[ili] = uikdi; /* -U(i,k) */
893: /* add multiple of row i to k-th row ... */
894: jmin = ili + 1;
895: nz = bi[i+1] - jmin;
896: if (nz > 0){
897: bcol = bj + jmin;
898: bval = ba + jmin;
899: while (nz --) rtmp[*bcol++] += uikdi*(*bval++);
900: /* update il and jl for i-th row */
901: il[i] = jmin;
902: j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
903: }
904: i = nexti;
905: }
907: /* shift the diagonals when zero pivot is detected */
908: /* compute rs=sum of abs(off-diagonal) */
909: rs = 0.0;
910: jmin = bi[k]+1;
911: nz = bi[k+1] - jmin;
912: if (nz){
913: bcol = bj + jmin;
914: while (nz--){
915: rs += PetscAbsScalar(rtmp[*bcol]);
916: bcol++;
917: }
918: }
920: sctx.rs = rs;
921: sctx.pv = dk;
922: MatPivotCheck(A,info,&sctx,k);
923: if (sctx.newshift) break; /* sctx.shift_amount is updated */
924: dk = sctx.pv;
926: /* copy data into U(k,:) */
927: ba[bi[k]] = 1.0/dk;
928: jmin = bi[k]+1;
929: nz = bi[k+1] - jmin;
930: if (nz){
931: bcol = bj + jmin;
932: bval = ba + jmin;
933: while (nz--){
934: *bval++ = rtmp[*bcol];
935: rtmp[*bcol++] = 0.0;
936: }
937: /* add k-th row into il and jl */
938: il[k] = jmin;
939: i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
940: }
941: }
942: } while (sctx.newshift);
943: PetscFree3(rtmp,il,jl);
944:
945: C->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
946: C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
947: C->assembled = PETSC_TRUE;
948: C->preallocated = PETSC_TRUE;
949: PetscLogFlops(C->rmap->N);
950: if (sctx.nshift){
951: if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
952: PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);
953: } else if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
954: PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);
955: }
956: }
957: return(0);
958: }
960: #include <petscbt.h>
961: #include <../src/mat/utils/freespace.h>
964: PetscErrorCode MatICCFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
965: {
966: Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
967: Mat_SeqSBAIJ *b;
968: Mat B;
969: PetscErrorCode ierr;
970: PetscBool perm_identity;
971: PetscInt reallocs=0,i,*ai=a->i,*aj=a->j,am=a->mbs,bs=A->rmap->bs,*ui;
972: const PetscInt *rip;
973: PetscInt jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
974: PetscInt nlnk,*lnk,*lnk_lvl=PETSC_NULL,ncols,ncols_upper,*cols,*cols_lvl,*uj,**uj_ptr,**uj_lvl_ptr;
975: PetscReal fill=info->fill,levels=info->levels;
976: PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
977: PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
978: PetscBT lnkbt;
981: if (bs > 1){
982: if (!a->sbaijMat){
983: MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);
984: }
985: (fact)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqSBAIJ; /* undue the change made in MatGetFactor_seqbaij_petsc */
986: MatICCFactorSymbolic(fact,a->sbaijMat,perm,info);
987: return(0);
988: }
990: ISIdentity(perm,&perm_identity);
991: ISGetIndices(perm,&rip);
993: /* special case that simply copies fill pattern */
994: if (!levels && perm_identity) {
995: MatMarkDiagonal_SeqBAIJ(A);
996: PetscMalloc((am+1)*sizeof(PetscInt),&ui);
997: for (i=0; i<am; i++) {
998: ui[i] = ai[i+1] - a->diag[i]; /* ui: rowlengths - changes when !perm_identity */
999: }
1000: B = fact;
1001: MatSeqSBAIJSetPreallocation(B,1,0,ui);
1004: b = (Mat_SeqSBAIJ*)B->data;
1005: uj = b->j;
1006: for (i=0; i<am; i++) {
1007: aj = a->j + a->diag[i];
1008: for (j=0; j<ui[i]; j++){
1009: *uj++ = *aj++;
1010: }
1011: b->ilen[i] = ui[i];
1012: }
1013: PetscFree(ui);
1014: B->factortype = MAT_FACTOR_NONE;
1015: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
1016: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
1017: B->factortype = MAT_FACTOR_ICC;
1019: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
1020: return(0);
1021: }
1023: /* initialization */
1024: PetscMalloc((am+1)*sizeof(PetscInt),&ui);
1025: ui[0] = 0;
1026: PetscMalloc((2*am+1)*sizeof(PetscInt),&cols_lvl);
1028: /* jl: linked list for storing indices of the pivot rows
1029: il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
1030: PetscMalloc4(am,PetscInt*,&uj_ptr,am,PetscInt*,&uj_lvl_ptr,am,PetscInt,&il,am,PetscInt,&jl);
1031: for (i=0; i<am; i++){
1032: jl[i] = am; il[i] = 0;
1033: }
1035: /* create and initialize a linked list for storing column indices of the active row k */
1036: nlnk = am + 1;
1037: PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);
1039: /* initial FreeSpace size is fill*(ai[am]+am)/2 */
1040: PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+am)/2),&free_space);
1041: current_space = free_space;
1042: PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+am)/2),&free_space_lvl);
1043: current_space_lvl = free_space_lvl;
1045: for (k=0; k<am; k++){ /* for each active row k */
1046: /* initialize lnk by the column indices of row rip[k] of A */
1047: nzk = 0;
1048: ncols = ai[rip[k]+1] - ai[rip[k]];
1049: ncols_upper = 0;
1050: cols = cols_lvl + am;
1051: for (j=0; j<ncols; j++){
1052: i = rip[*(aj + ai[rip[k]] + j)];
1053: if (i >= k){ /* only take upper triangular entry */
1054: cols[ncols_upper] = i;
1055: cols_lvl[ncols_upper] = -1; /* initialize level for nonzero entries */
1056: ncols_upper++;
1057: }
1058: }
1059: PetscIncompleteLLAdd(ncols_upper,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);
1060: nzk += nlnk;
1062: /* update lnk by computing fill-in for each pivot row to be merged in */
1063: prow = jl[k]; /* 1st pivot row */
1064:
1065: while (prow < k){
1066: nextprow = jl[prow];
1067:
1068: /* merge prow into k-th row */
1069: jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */
1070: jmax = ui[prow+1];
1071: ncols = jmax-jmin;
1072: i = jmin - ui[prow];
1073: cols = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
1074: for (j=0; j<ncols; j++) cols_lvl[j] = *(uj_lvl_ptr[prow] + i + j);
1075: PetscIncompleteLLAddSorted(ncols,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);
1076: nzk += nlnk;
1078: /* update il and jl for prow */
1079: if (jmin < jmax){
1080: il[prow] = jmin;
1081: j = *cols; jl[prow] = jl[j]; jl[j] = prow;
1082: }
1083: prow = nextprow;
1084: }
1086: /* if free space is not available, make more free space */
1087: if (current_space->local_remaining<nzk) {
1088: i = am - k + 1; /* num of unfactored rows */
1089: i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
1090: PetscFreeSpaceGet(i,¤t_space);
1091: PetscFreeSpaceGet(i,¤t_space_lvl);
1092: reallocs++;
1093: }
1095: /* copy data into free_space and free_space_lvl, then initialize lnk */
1096: PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);
1098: /* add the k-th row into il and jl */
1099: if (nzk-1 > 0){
1100: i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
1101: jl[k] = jl[i]; jl[i] = k;
1102: il[k] = ui[k] + 1;
1103: }
1104: uj_ptr[k] = current_space->array;
1105: uj_lvl_ptr[k] = current_space_lvl->array;
1107: current_space->array += nzk;
1108: current_space->local_used += nzk;
1109: current_space->local_remaining -= nzk;
1111: current_space_lvl->array += nzk;
1112: current_space_lvl->local_used += nzk;
1113: current_space_lvl->local_remaining -= nzk;
1115: ui[k+1] = ui[k] + nzk;
1116: }
1118: ISRestoreIndices(perm,&rip);
1119: PetscFree4(uj_ptr,uj_lvl_ptr,il,jl);
1120: PetscFree(cols_lvl);
1122: /* copy free_space into uj and free free_space; set uj in new datastructure; */
1123: PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);
1124: PetscFreeSpaceContiguous(&free_space,uj);
1125: PetscIncompleteLLDestroy(lnk,lnkbt);
1126: PetscFreeSpaceDestroy(free_space_lvl);
1128: /* put together the new matrix in MATSEQSBAIJ format */
1129: B = fact;
1130: MatSeqSBAIJSetPreallocation(B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);
1132: b = (Mat_SeqSBAIJ*)B->data;
1133: b->singlemalloc = PETSC_FALSE;
1134: b->free_a = PETSC_TRUE;
1135: b->free_ij = PETSC_TRUE;
1136: PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);
1137: b->j = uj;
1138: b->i = ui;
1139: b->diag = 0;
1140: b->ilen = 0;
1141: b->imax = 0;
1142: b->row = perm;
1143: b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
1144: PetscObjectReference((PetscObject)perm);
1145: b->icol = perm;
1146: PetscObjectReference((PetscObject)perm);
1147: PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);
1148: PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));
1149: b->maxnz = b->nz = ui[am];
1150:
1151: B->info.factor_mallocs = reallocs;
1152: B->info.fill_ratio_given = fill;
1153: if (ai[am] != 0.) {
1154: /* nonzeros in lower triangular part of A (includign diagonals)= (ai[am]+am)/2 */
1155: B->info.fill_ratio_needed = ((PetscReal)2*ui[am])/(ai[am]+am);
1156: } else {
1157: B->info.fill_ratio_needed = 0.0;
1158: }
1159: #if defined(PETSC_USE_INFO)
1160: if (ai[am] != 0) {
1161: PetscReal af = B->info.fill_ratio_needed;
1162: PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);
1163: PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);
1164: PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);
1165: } else {
1166: PetscInfo(A,"Empty matrix.\n");
1167: }
1168: #endif
1169: if (perm_identity){
1170: B->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
1171: B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
1172: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
1173: } else {
1174: (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
1175: }
1176: return(0);
1177: }
1181: PetscErrorCode MatCholeskyFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
1182: {
1183: Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
1184: Mat_SeqSBAIJ *b;
1185: Mat B;
1186: PetscErrorCode ierr;
1187: PetscBool perm_identity;
1188: PetscReal fill = info->fill;
1189: const PetscInt *rip;
1190: PetscInt i,mbs=a->mbs,bs=A->rmap->bs,*ai=a->i,*aj=a->j,reallocs=0,prow;
1191: PetscInt *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
1192: PetscInt nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr;
1193: PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
1194: PetscBT lnkbt;
1197: if (bs > 1) { /* convert to seqsbaij */
1198: if (!a->sbaijMat){
1199: MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);
1200: }
1201: (fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqSBAIJ; /* undue the change made in MatGetFactor_seqbaij_petsc */
1202: MatCholeskyFactorSymbolic(fact,a->sbaijMat,perm,info);
1203: return(0);
1204: }
1206: /* check whether perm is the identity mapping */
1207: ISIdentity(perm,&perm_identity);
1208: if (!perm_identity) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix reordering is not supported");
1209: ISGetIndices(perm,&rip);
1211: /* initialization */
1212: PetscMalloc((mbs+1)*sizeof(PetscInt),&ui);
1213: ui[0] = 0;
1215: /* jl: linked list for storing indices of the pivot rows
1216: il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
1217: PetscMalloc4(mbs,PetscInt*,&ui_ptr,mbs,PetscInt,&il,mbs,PetscInt,&jl,mbs,PetscInt,&cols);
1218: for (i=0; i<mbs; i++){
1219: jl[i] = mbs; il[i] = 0;
1220: }
1222: /* create and initialize a linked list for storing column indices of the active row k */
1223: nlnk = mbs + 1;
1224: PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);
1226: /* initial FreeSpace size is fill* (ai[mbs]+mbs)/2 */
1227: PetscFreeSpaceGet((PetscInt)(fill*(ai[mbs]+mbs)/2),&free_space);
1228: current_space = free_space;
1230: for (k=0; k<mbs; k++){ /* for each active row k */
1231: /* initialize lnk by the column indices of row rip[k] of A */
1232: nzk = 0;
1233: ncols = ai[rip[k]+1] - ai[rip[k]];
1234: ncols_upper = 0;
1235: for (j=0; j<ncols; j++){
1236: i = rip[*(aj + ai[rip[k]] + j)];
1237: if (i >= k){ /* only take upper triangular entry */
1238: cols[ncols_upper] = i;
1239: ncols_upper++;
1240: }
1241: }
1242: PetscLLAdd(ncols_upper,cols,mbs,nlnk,lnk,lnkbt);
1243: nzk += nlnk;
1245: /* update lnk by computing fill-in for each pivot row to be merged in */
1246: prow = jl[k]; /* 1st pivot row */
1247:
1248: while (prow < k){
1249: nextprow = jl[prow];
1250: /* merge prow into k-th row */
1251: jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
1252: jmax = ui[prow+1];
1253: ncols = jmax-jmin;
1254: uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
1255: PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);
1256: nzk += nlnk;
1258: /* update il and jl for prow */
1259: if (jmin < jmax){
1260: il[prow] = jmin;
1261: j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
1262: }
1263: prow = nextprow;
1264: }
1266: /* if free space is not available, make more free space */
1267: if (current_space->local_remaining<nzk) {
1268: i = mbs - k + 1; /* num of unfactored rows */
1269: i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
1270: PetscFreeSpaceGet(i,¤t_space);
1271: reallocs++;
1272: }
1274: /* copy data into free space, then initialize lnk */
1275: PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);
1277: /* add the k-th row into il and jl */
1278: if (nzk-1 > 0){
1279: i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
1280: jl[k] = jl[i]; jl[i] = k;
1281: il[k] = ui[k] + 1;
1282: }
1283: ui_ptr[k] = current_space->array;
1284: current_space->array += nzk;
1285: current_space->local_used += nzk;
1286: current_space->local_remaining -= nzk;
1288: ui[k+1] = ui[k] + nzk;
1289: }
1291: ISRestoreIndices(perm,&rip);
1292: PetscFree4(ui_ptr,il,jl,cols);
1294: /* copy free_space into uj and free free_space; set uj in new datastructure; */
1295: PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&uj);
1296: PetscFreeSpaceContiguous(&free_space,uj);
1297: PetscLLDestroy(lnk,lnkbt);
1299: /* put together the new matrix in MATSEQSBAIJ format */
1300: B = fact;
1301: MatSeqSBAIJSetPreallocation(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);
1303: b = (Mat_SeqSBAIJ*)B->data;
1304: b->singlemalloc = PETSC_FALSE;
1305: b->free_a = PETSC_TRUE;
1306: b->free_ij = PETSC_TRUE;
1307: PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&b->a);
1308: b->j = uj;
1309: b->i = ui;
1310: b->diag = 0;
1311: b->ilen = 0;
1312: b->imax = 0;
1313: b->row = perm;
1314: b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
1315: PetscObjectReference((PetscObject)perm);
1316: b->icol = perm;
1317: PetscObjectReference((PetscObject)perm);
1318: PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);
1319: PetscLogObjectMemory(B,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));
1320: b->maxnz = b->nz = ui[mbs];
1321:
1322: B->info.factor_mallocs = reallocs;
1323: B->info.fill_ratio_given = fill;
1324: if (ai[mbs] != 0.) {
1325: /* nonzeros in lower triangular part of A = (ai[mbs]+mbs)/2 */
1326: B->info.fill_ratio_needed = ((PetscReal)2*ui[mbs])/(ai[mbs]+mbs);
1327: } else {
1328: B->info.fill_ratio_needed = 0.0;
1329: }
1330: #if defined(PETSC_USE_INFO)
1331: if (ai[mbs] != 0.) {
1332: PetscReal af = B->info.fill_ratio_needed;
1333: PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);
1334: PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);
1335: PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);
1336: } else {
1337: PetscInfo(A,"Empty matrix.\n");
1338: }
1339: #endif
1340: if (perm_identity){
1341: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
1342: } else {
1343: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
1344: }
1345: return(0);
1346: }
1350: PetscErrorCode MatSolve_SeqBAIJ_N_NaturalOrdering(Mat A,Vec bb,Vec xx)
1351: {
1352: Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data;
1354: const PetscInt *ai=a->i,*aj=a->j,*adiag=a->diag,*vi;
1355: PetscInt i,k,n=a->mbs;
1356: PetscInt nz,bs=A->rmap->bs,bs2=a->bs2;
1357: MatScalar *aa=a->a,*v;
1358: PetscScalar *x,*b,*s,*t,*ls;
1361: VecGetArray(bb,&b);
1362: VecGetArray(xx,&x);
1363: t = a->solve_work;
1365: /* forward solve the lower triangular */
1366: PetscMemcpy(t,b,bs*sizeof(PetscScalar)); /* copy 1st block of b to t */
1368: for (i=1; i<n; i++) {
1369: v = aa + bs2*ai[i];
1370: vi = aj + ai[i];
1371: nz = ai[i+1] - ai[i];
1372: s = t + bs*i;
1373: PetscMemcpy(s,b+bs*i,bs*sizeof(PetscScalar)); /* copy i_th block of b to t */
1374: for(k=0;k<nz;k++){
1375: Kernel_v_gets_v_minus_A_times_w(bs,s,v,t+bs*vi[k]);
1376: v += bs2;
1377: }
1378: }
1379:
1380: /* backward solve the upper triangular */
1381: ls = a->solve_work + A->cmap->n;
1382: for (i=n-1; i>=0; i--){
1383: v = aa + bs2*(adiag[i+1]+1);
1384: vi = aj + adiag[i+1]+1;
1385: nz = adiag[i] - adiag[i+1]-1;
1386: PetscMemcpy(ls,t+i*bs,bs*sizeof(PetscScalar));
1387: for(k=0;k<nz;k++){
1388: Kernel_v_gets_v_minus_A_times_w(bs,ls,v,t+bs*vi[k]);
1389: v += bs2;
1390: }
1391: Kernel_w_gets_A_times_v(bs,ls,aa+bs2*adiag[i],t+i*bs); /* *inv(diagonal[i]) */
1392: PetscMemcpy(x+i*bs,t+i*bs,bs*sizeof(PetscScalar));
1393: }
1394:
1395: VecRestoreArray(bb,&b);
1396: VecRestoreArray(xx,&x);
1397: PetscLogFlops(2.0*(a->bs2)*(a->nz) - A->rmap->bs*A->cmap->n);
1398: return(0);
1399: }
1403: PetscErrorCode MatSolve_SeqBAIJ_N(Mat A,Vec bb,Vec xx)
1404: {
1405: Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data;
1406: IS iscol=a->col,isrow=a->row;
1408: const PetscInt *r,*c,*rout,*cout,*ai=a->i,*aj=a->j,*adiag=a->diag,*vi;
1409: PetscInt i,m,n=a->mbs;
1410: PetscInt nz,bs=A->rmap->bs,bs2=a->bs2;
1411: MatScalar *aa=a->a,*v;
1412: PetscScalar *x,*b,*s,*t,*ls;
1415: VecGetArray(bb,&b);
1416: VecGetArray(xx,&x);
1417: t = a->solve_work;
1419: ISGetIndices(isrow,&rout); r = rout;
1420: ISGetIndices(iscol,&cout); c = cout;
1422: /* forward solve the lower triangular */
1423: PetscMemcpy(t,b+bs*r[0],bs*sizeof(PetscScalar));
1424: for (i=1; i<n; i++) {
1425: v = aa + bs2*ai[i];
1426: vi = aj + ai[i];
1427: nz = ai[i+1] - ai[i];
1428: s = t + bs*i;
1429: PetscMemcpy(s,b+bs*r[i],bs*sizeof(PetscScalar));
1430: for(m=0;m<nz;m++){
1431: Kernel_v_gets_v_minus_A_times_w(bs,s,v,t+bs*vi[m]);
1432: v += bs2;
1433: }
1434: }
1436: /* backward solve the upper triangular */
1437: ls = a->solve_work + A->cmap->n;
1438: for (i=n-1; i>=0; i--){
1439: v = aa + bs2*(adiag[i+1]+1);
1440: vi = aj + adiag[i+1]+1;
1441: nz = adiag[i] - adiag[i+1] - 1;
1442: PetscMemcpy(ls,t+i*bs,bs*sizeof(PetscScalar));
1443: for(m=0;m<nz;m++){
1444: Kernel_v_gets_v_minus_A_times_w(bs,ls,v,t+bs*vi[m]);
1445: v += bs2;
1446: }
1447: Kernel_w_gets_A_times_v(bs,ls,v,t+i*bs); /* *inv(diagonal[i]) */
1448: PetscMemcpy(x + bs*c[i],t+i*bs,bs*sizeof(PetscScalar));
1449: }
1450: ISRestoreIndices(isrow,&rout);
1451: ISRestoreIndices(iscol,&cout);
1452: VecRestoreArray(bb,&b);
1453: VecRestoreArray(xx,&x);
1454: PetscLogFlops(2.0*(a->bs2)*(a->nz) - A->rmap->bs*A->cmap->n);
1455: return(0);
1456: }
1460: PetscErrorCode BlockAbs_private(PetscInt nbs,PetscInt bs2,PetscScalar *blockarray,PetscReal *absarray)
1461: {
1462: PetscErrorCode ierr;
1463: PetscInt i,j;
1465: PetscMemzero(absarray,(nbs+1)*sizeof(PetscReal));
1466: for (i=0; i<nbs; i++){
1467: for (j=0; j<bs2; j++){
1468: if (absarray[i] < PetscAbsScalar(blockarray[i*nbs+j])) absarray[i] = PetscAbsScalar(blockarray[i*nbs+j]);
1469: }
1470: }
1471: return(0);
1472: }
1476: /*
1477: This needs to be renamed and called by the regular MatILUFactor_SeqBAIJ when drop tolerance is used
1478: */
1479: PetscErrorCode MatILUDTFactor_SeqBAIJ(Mat A,IS isrow,IS iscol,const MatFactorInfo *info,Mat *fact)
1480: {
1481: Mat B = *fact;
1482: Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data,*b;
1483: IS isicol;
1484: PetscErrorCode ierr;
1485: const PetscInt *r,*ic;
1486: PetscInt i,mbs=a->mbs,bs=A->rmap->bs,bs2=a->bs2,*ai=a->i,*aj=a->j,*ajtmp,*adiag;
1487: PetscInt *bi,*bj,*bdiag;
1488:
1489: PetscInt row,nzi,nzi_bl,nzi_bu,*im,dtcount,nzi_al,nzi_au;
1490: PetscInt nlnk,*lnk;
1491: PetscBT lnkbt;
1492: PetscBool row_identity,icol_identity;
1493: MatScalar *aatmp,*pv,*batmp,*ba,*rtmp,*pc,*multiplier,*vtmp;
1494: PetscInt j,nz,*pj,*bjtmp,k,ncut,*jtmp;
1495:
1496: PetscReal dt=info->dt; /* shift=info->shiftamount; */
1497: PetscInt nnz_max;
1498: PetscBool missing;
1499: PetscReal *vtmp_abs;
1500: MatScalar *v_work;
1501: PetscInt *v_pivots;
1504: /* ------- symbolic factorization, can be reused ---------*/
1505: MatMissingDiagonal(A,&missing,&i);
1506: if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i);
1507: adiag=a->diag;
1509: ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);
1511: /* bdiag is location of diagonal in factor */
1512: PetscMalloc((mbs+1)*sizeof(PetscInt),&bdiag);
1514: /* allocate row pointers bi */
1515: PetscMalloc((2*mbs+2)*sizeof(PetscInt),&bi);
1517: /* allocate bj and ba; max num of nonzero entries is (ai[n]+2*n*dtcount+2) */
1518: dtcount = (PetscInt)info->dtcount;
1519: if (dtcount > mbs-1) dtcount = mbs-1;
1520: nnz_max = ai[mbs]+2*mbs*dtcount +2;
1521: /* printf("MatILUDTFactor_SeqBAIJ, bs %d, ai[mbs] %d, nnz_max %d, dtcount %d\n",bs,ai[mbs],nnz_max,dtcount); */
1522: PetscMalloc(nnz_max*sizeof(PetscInt),&bj);
1523: nnz_max = nnz_max*bs2;
1524: PetscMalloc(nnz_max*sizeof(MatScalar),&ba);
1526: /* put together the new matrix */
1527: MatSeqBAIJSetPreallocation_SeqBAIJ(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);
1528: PetscLogObjectParent(B,isicol);
1529: b = (Mat_SeqBAIJ*)(B)->data;
1530: b->free_a = PETSC_TRUE;
1531: b->free_ij = PETSC_TRUE;
1532: b->singlemalloc = PETSC_FALSE;
1533: b->a = ba;
1534: b->j = bj;
1535: b->i = bi;
1536: b->diag = bdiag;
1537: b->ilen = 0;
1538: b->imax = 0;
1539: b->row = isrow;
1540: b->col = iscol;
1541: PetscObjectReference((PetscObject)isrow);
1542: PetscObjectReference((PetscObject)iscol);
1543: b->icol = isicol;
1544: PetscMalloc((bs*(mbs+1))*sizeof(PetscScalar),&b->solve_work);
1546: PetscLogObjectMemory(B,nnz_max*(sizeof(PetscInt)+sizeof(MatScalar)));
1547: b->maxnz = nnz_max/bs2;
1549: (B)->factortype = MAT_FACTOR_ILUDT;
1550: (B)->info.factor_mallocs = 0;
1551: (B)->info.fill_ratio_given = ((PetscReal)nnz_max)/((PetscReal)(ai[mbs]*bs2));
1552: CHKMEMQ;
1553: /* ------- end of symbolic factorization ---------*/
1554: ISGetIndices(isrow,&r);
1555: ISGetIndices(isicol,&ic);
1557: /* linked list for storing column indices of the active row */
1558: nlnk = mbs + 1;
1559: PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);
1561: /* im: used by PetscLLAddSortedLU(); jtmp: working array for column indices of active row */
1562: PetscMalloc2(mbs,PetscInt,&im,mbs,PetscInt,&jtmp);
1563: /* rtmp, vtmp: working arrays for sparse and contiguous row entries of active row */
1564: PetscMalloc2(mbs*bs2,MatScalar,&rtmp,mbs*bs2,MatScalar,&vtmp);
1565: PetscMalloc((mbs+1)*sizeof(PetscReal),&vtmp_abs);
1566: PetscMalloc3(bs,MatScalar,&v_work,bs2,MatScalar,&multiplier,bs,PetscInt,&v_pivots);
1568: bi[0] = 0;
1569: bdiag[0] = (nnz_max/bs2)-1; /* location of diagonal in factor B */
1570: bi[2*mbs+1] = bdiag[0]+1; /* endof bj and ba array */
1571: for (i=0; i<mbs; i++) {
1572: /* copy initial fill into linked list */
1573: nzi = 0; /* nonzeros for active row i */
1574: nzi = ai[r[i]+1] - ai[r[i]];
1575: if (!nzi) 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);
1576: nzi_al = adiag[r[i]] - ai[r[i]];
1577: nzi_au = ai[r[i]+1] - adiag[r[i]] -1;
1578: /* printf("row %d, nzi_al/au %d %d\n",i,nzi_al,nzi_au); */
1579:
1580: /* load in initial unfactored row */
1581: ajtmp = aj + ai[r[i]];
1582: PetscLLAddPerm(nzi,ajtmp,ic,mbs,nlnk,lnk,lnkbt);
1583: PetscMemzero(rtmp,mbs*bs2*sizeof(PetscScalar));
1584: aatmp = a->a + bs2*ai[r[i]];
1585: for (j=0; j<nzi; j++) {
1586: PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],aatmp+bs2*j,bs2*sizeof(MatScalar));
1587: }
1588:
1589: /* add pivot rows into linked list */
1590: row = lnk[mbs];
1591: while (row < i) {
1592: nzi_bl = bi[row+1] - bi[row] + 1;
1593: bjtmp = bj + bdiag[row+1]+1; /* points to 1st column next to the diagonal in U */
1594: PetscLLAddSortedLU(bjtmp,row,nlnk,lnk,lnkbt,i,nzi_bl,im);
1595: nzi += nlnk;
1596: row = lnk[row];
1597: }
1598:
1599: /* copy data from lnk into jtmp, then initialize lnk */
1600: PetscLLClean(mbs,mbs,nzi,lnk,jtmp,lnkbt);
1602: /* numerical factorization */
1603: bjtmp = jtmp;
1604: row = *bjtmp++; /* 1st pivot row */
1606: while (row < i) {
1607: pc = rtmp + bs2*row;
1608: pv = ba + bs2*bdiag[row]; /* inv(diag) of the pivot row */
1609: Kernel_A_gets_A_times_B(bs,pc,pv,multiplier); /* pc= multiplier = pc*inv(diag[row]) */
1610: BlockAbs_private(1,bs2,pc,vtmp_abs);
1611: if (vtmp_abs[0] > dt){ /* apply tolerance dropping rule */
1612: pj = bj + bdiag[row+1] + 1; /* point to 1st entry of U(row,:) */
1613: pv = ba + bs2*(bdiag[row+1] + 1);
1614: nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:), excluding diagonal */
1615: for (j=0; j<nz; j++){
1616: Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j);
1617: }
1618: /* PetscLogFlops(bslog*(nz+1.0)-bs); */
1619: }
1620: row = *bjtmp++;
1621: }
1623: /* copy sparse rtmp into contiguous vtmp; separate L and U part */
1624: nzi_bl = 0; j = 0;
1625: while (jtmp[j] < i){ /* L-part. Note: jtmp is sorted */
1626: PetscMemcpy(vtmp+bs2*j,rtmp+bs2*jtmp[j],bs2*sizeof(MatScalar));
1627: nzi_bl++; j++;
1628: }
1629: nzi_bu = nzi - nzi_bl -1;
1630: /* printf("nzi %d, nzi_bl %d, nzi_bu %d\n",nzi,nzi_bl,nzi_bu); */
1632: while (j < nzi){ /* U-part */
1633: PetscMemcpy(vtmp+bs2*j,rtmp+bs2*jtmp[j],bs2*sizeof(MatScalar));
1634: /*
1635: printf(" col %d: ",jtmp[j]);
1636: for (j1=0; j1<bs2; j1++) printf(" %g",*(vtmp+bs2*j+j1));
1637: printf(" \n");
1638: */
1639: j++;
1640: }
1642: BlockAbs_private(nzi,bs2,vtmp,vtmp_abs);
1643: /*
1644: printf(" row %d, nzi %d, vtmp_abs\n",i,nzi);
1645: for (j1=0; j1<nzi; j1++) printf(" (%d %g),",jtmp[j1],vtmp_abs[j1]);
1646: printf(" \n");
1647: */
1648: bjtmp = bj + bi[i];
1649: batmp = ba + bs2*bi[i];
1650: /* apply level dropping rule to L part */
1651: ncut = nzi_al + dtcount;
1652: if (ncut < nzi_bl){
1653: PetscSortSplitReal(ncut,nzi_bl,vtmp_abs,jtmp);
1654: PetscSortIntWithScalarArray(ncut,jtmp,vtmp);
1655: } else {
1656: ncut = nzi_bl;
1657: }
1658: for (j=0; j<ncut; j++){
1659: bjtmp[j] = jtmp[j];
1660: PetscMemcpy(batmp+bs2*j,rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));
1661: /*
1662: printf(" col %d: ",bjtmp[j]);
1663: for (j1=0; j1<bs2; j1++) printf(" %g,",*(batmp+bs2*j+j1));
1664: printf("\n");
1665: */
1666: }
1667: bi[i+1] = bi[i] + ncut;
1668: nzi = ncut + 1;
1669:
1670: /* apply level dropping rule to U part */
1671: ncut = nzi_au + dtcount;
1672: if (ncut < nzi_bu){
1673: PetscSortSplitReal(ncut,nzi_bu,vtmp_abs+nzi_bl+1,jtmp+nzi_bl+1);
1674: PetscSortIntWithScalarArray(ncut,jtmp+nzi_bl+1,vtmp+nzi_bl+1);
1675: } else {
1676: ncut = nzi_bu;
1677: }
1678: nzi += ncut;
1679:
1680: /* mark bdiagonal */
1681: bdiag[i+1] = bdiag[i] - (ncut + 1);
1682: bi[2*mbs - i] = bi[2*mbs - i +1] - (ncut + 1);
1683:
1684: bjtmp = bj + bdiag[i];
1685: batmp = ba + bs2*bdiag[i];
1686: PetscMemcpy(batmp,rtmp+bs2*i,bs2*sizeof(MatScalar));
1687: *bjtmp = i;
1688: /*
1689: printf(" diag %d: ",*bjtmp);
1690: for (j=0; j<bs2; j++){
1691: printf(" %g,",batmp[j]);
1692: }
1693: printf("\n");
1694: */
1695: bjtmp = bj + bdiag[i+1]+1;
1696: batmp = ba + (bdiag[i+1]+1)*bs2;
1697:
1698: for (k=0; k<ncut; k++){
1699: bjtmp[k] = jtmp[nzi_bl+1+k];
1700: PetscMemcpy(batmp+bs2*k,rtmp+bs2*bjtmp[k],bs2*sizeof(MatScalar));
1701: /*
1702: printf(" col %d:",bjtmp[k]);
1703: for (j1=0; j1<bs2; j1++) printf(" %g,",*(batmp+bs2*k+j1));
1704: printf("\n");
1705: */
1706: }
1707:
1708: im[i] = nzi; /* used by PetscLLAddSortedLU() */
1709:
1710: /* invert diagonal block for simplier triangular solves - add shift??? */
1711: batmp = ba + bs2*bdiag[i];
1712: Kernel_A_gets_inverse_A(bs,batmp,v_pivots,v_work);
1713: } /* for (i=0; i<mbs; i++) */
1714: PetscFree3(v_work,multiplier,v_pivots);
1716: /* printf("end of L %d, beginning of U %d\n",bi[mbs],bdiag[mbs]); */
1717: if (bi[mbs] >= bdiag[mbs]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"end of L array %d cannot >= the beginning of U array %d",bi[mbs],bdiag[mbs]);
1719: ISRestoreIndices(isrow,&r);
1720: ISRestoreIndices(isicol,&ic);
1722: PetscLLDestroy(lnk,lnkbt);
1724: PetscFree2(im,jtmp);
1725: PetscFree2(rtmp,vtmp);
1726:
1727: PetscLogFlops(bs2*B->cmap->n);
1728: b->maxnz = b->nz = bi[mbs] + bdiag[0] - bdiag[mbs];
1730: ISIdentity(isrow,&row_identity);
1731: ISIdentity(isicol,&icol_identity);
1732: if (row_identity && icol_identity) {
1733: B->ops->solve = MatSolve_SeqBAIJ_N_NaturalOrdering;
1734: } else {
1735: B->ops->solve = MatSolve_SeqBAIJ_N;
1736: }
1737:
1738: B->ops->solveadd = 0;
1739: B->ops->solvetranspose = 0;
1740: B->ops->solvetransposeadd = 0;
1741: B->ops->matsolve = 0;
1742: B->assembled = PETSC_TRUE;
1743: B->preallocated = PETSC_TRUE;
1744: return(0);
1745: }