Actual source code: baijfact13.c
2: /*
3: Factorization code for BAIJ format.
4: */
5: #include <../src/mat/impls/baij/seq/baij.h>
6: #include <petsc/private/kernels/blockinvert.h>
8: /*
9: Version for when blocks are 3 by 3
10: */
11: PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3_inplace(Mat C,Mat A,const MatFactorInfo *info)
12: {
13: Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)C->data;
14: IS isrow = b->row,isicol = b->icol;
15: const PetscInt *r,*ic;
16: PetscInt i,j,n = a->mbs,*bi = b->i,*bj = b->j;
17: PetscInt *ajtmpold,*ajtmp,nz,row,*ai=a->i,*aj=a->j;
18: PetscInt *diag_offset = b->diag,idx,*pj;
19: MatScalar *pv,*v,*rtmp,*pc,*w,*x;
20: MatScalar p1,p2,p3,p4,m1,m2,m3,m4,m5,m6,m7,m8,m9,x1,x2,x3,x4;
21: MatScalar p5,p6,p7,p8,p9,x5,x6,x7,x8,x9;
22: MatScalar *ba = b->a,*aa = a->a;
23: PetscReal shift = info->shiftamount;
24: PetscBool allowzeropivot,zeropivotdetected;
26: ISGetIndices(isrow,&r);
27: ISGetIndices(isicol,&ic);
28: PetscMalloc1(9*(n+1),&rtmp);
29: allowzeropivot = PetscNot(A->erroriffailure);
31: for (i=0; i<n; i++) {
32: nz = bi[i+1] - bi[i];
33: ajtmp = bj + bi[i];
34: for (j=0; j<nz; j++) {
35: x = rtmp + 9*ajtmp[j];
36: x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = 0.0;
37: }
38: /* load in initial (unfactored row) */
39: idx = r[i];
40: nz = ai[idx+1] - ai[idx];
41: ajtmpold = aj + ai[idx];
42: v = aa + 9*ai[idx];
43: for (j=0; j<nz; j++) {
44: x = rtmp + 9*ic[ajtmpold[j]];
45: x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3];
46: x[4] = v[4]; x[5] = v[5]; x[6] = v[6]; x[7] = v[7]; x[8] = v[8];
47: v += 9;
48: }
49: row = *ajtmp++;
50: while (row < i) {
51: pc = rtmp + 9*row;
52: p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3];
53: p5 = pc[4]; p6 = pc[5]; p7 = pc[6]; p8 = pc[7]; p9 = pc[8];
54: if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 ||
55: p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0) {
56: pv = ba + 9*diag_offset[row];
57: pj = bj + diag_offset[row] + 1;
58: x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
59: x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8];
60: pc[0] = m1 = p1*x1 + p4*x2 + p7*x3;
61: pc[1] = m2 = p2*x1 + p5*x2 + p8*x3;
62: pc[2] = m3 = p3*x1 + p6*x2 + p9*x3;
64: pc[3] = m4 = p1*x4 + p4*x5 + p7*x6;
65: pc[4] = m5 = p2*x4 + p5*x5 + p8*x6;
66: pc[5] = m6 = p3*x4 + p6*x5 + p9*x6;
68: pc[6] = m7 = p1*x7 + p4*x8 + p7*x9;
69: pc[7] = m8 = p2*x7 + p5*x8 + p8*x9;
70: pc[8] = m9 = p3*x7 + p6*x8 + p9*x9;
71: nz = bi[row+1] - diag_offset[row] - 1;
72: pv += 9;
73: for (j=0; j<nz; j++) {
74: x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
75: x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8];
76: x = rtmp + 9*pj[j];
77: x[0] -= m1*x1 + m4*x2 + m7*x3;
78: x[1] -= m2*x1 + m5*x2 + m8*x3;
79: x[2] -= m3*x1 + m6*x2 + m9*x3;
81: x[3] -= m1*x4 + m4*x5 + m7*x6;
82: x[4] -= m2*x4 + m5*x5 + m8*x6;
83: x[5] -= m3*x4 + m6*x5 + m9*x6;
85: x[6] -= m1*x7 + m4*x8 + m7*x9;
86: x[7] -= m2*x7 + m5*x8 + m8*x9;
87: x[8] -= m3*x7 + m6*x8 + m9*x9;
88: pv += 9;
89: }
90: PetscLogFlops(54.0*nz+36.0);
91: }
92: row = *ajtmp++;
93: }
94: /* finished row so stick it into b->a */
95: pv = ba + 9*bi[i];
96: pj = bj + bi[i];
97: nz = bi[i+1] - bi[i];
98: for (j=0; j<nz; j++) {
99: x = rtmp + 9*pj[j];
100: pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3];
101: pv[4] = x[4]; pv[5] = x[5]; pv[6] = x[6]; pv[7] = x[7]; pv[8] = x[8];
102: pv += 9;
103: }
104: /* invert diagonal block */
105: w = ba + 9*diag_offset[i];
106: PetscKernel_A_gets_inverse_A_3(w,shift,allowzeropivot,&zeropivotdetected);
107: if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
108: }
110: PetscFree(rtmp);
111: ISRestoreIndices(isicol,&ic);
112: ISRestoreIndices(isrow,&r);
114: C->ops->solve = MatSolve_SeqBAIJ_3_inplace;
115: C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3_inplace;
116: C->assembled = PETSC_TRUE;
118: PetscLogFlops(1.333333333333*3*3*3*b->mbs); /* from inverting diagonal blocks */
119: return 0;
120: }
122: /* MatLUFactorNumeric_SeqBAIJ_3 -
123: copied from MatLUFactorNumeric_SeqBAIJ_N_inplace() and manually re-implemented
124: PetscKernel_A_gets_A_times_B()
125: PetscKernel_A_gets_A_minus_B_times_C()
126: PetscKernel_A_gets_inverse_A()
127: */
128: PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3(Mat B,Mat A,const MatFactorInfo *info)
129: {
130: Mat C =B;
131: Mat_SeqBAIJ *a =(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data;
132: IS isrow = b->row,isicol = b->icol;
133: const PetscInt *r,*ic;
134: PetscInt i,j,k,nz,nzL,row;
135: const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
136: const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag,*pj,bs2=a->bs2;
137: MatScalar *rtmp,*pc,*mwork,*v,*pv,*aa=a->a;
138: PetscInt flg;
139: PetscReal shift = info->shiftamount;
140: PetscBool allowzeropivot,zeropivotdetected;
142: ISGetIndices(isrow,&r);
143: ISGetIndices(isicol,&ic);
144: allowzeropivot = PetscNot(A->erroriffailure);
146: /* generate work space needed by the factorization */
147: PetscMalloc2(bs2*n,&rtmp,bs2,&mwork);
148: PetscArrayzero(rtmp,bs2*n);
150: for (i=0; i<n; i++) {
151: /* zero rtmp */
152: /* L part */
153: nz = bi[i+1] - bi[i];
154: bjtmp = bj + bi[i];
155: for (j=0; j<nz; j++) {
156: PetscArrayzero(rtmp+bs2*bjtmp[j],bs2);
157: }
159: /* U part */
160: nz = bdiag[i] - bdiag[i+1];
161: bjtmp = bj + bdiag[i+1]+1;
162: for (j=0; j<nz; j++) {
163: PetscArrayzero(rtmp+bs2*bjtmp[j],bs2);
164: }
166: /* load in initial (unfactored row) */
167: nz = ai[r[i]+1] - ai[r[i]];
168: ajtmp = aj + ai[r[i]];
169: v = aa + bs2*ai[r[i]];
170: for (j=0; j<nz; j++) {
171: PetscArraycpy(rtmp+bs2*ic[ajtmp[j]],v+bs2*j,bs2);
172: }
174: /* elimination */
175: bjtmp = bj + bi[i];
176: nzL = bi[i+1] - bi[i];
177: for (k = 0; k < nzL; k++) {
178: row = bjtmp[k];
179: pc = rtmp + bs2*row;
180: for (flg=0,j=0; j<bs2; j++) {
181: if (pc[j]!=0.0) {
182: flg = 1;
183: break;
184: }
185: }
186: if (flg) {
187: pv = b->a + bs2*bdiag[row];
188: /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
189: PetscKernel_A_gets_A_times_B_3(pc,pv,mwork);
191: pj = b->j + bdiag[row+1] + 1; /* beginning of U(row,:) */
192: pv = b->a + bs2*(bdiag[row+1]+1);
193: nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:) excluding diag */
194: for (j=0; j<nz; j++) {
195: /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
196: /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
197: v = rtmp + bs2*pj[j];
198: PetscKernel_A_gets_A_minus_B_times_C_3(v,pc,pv);
199: pv += bs2;
200: }
201: PetscLogFlops(54.0*nz+45); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
202: }
203: }
205: /* finished row so stick it into b->a */
206: /* L part */
207: pv = b->a + bs2*bi[i];
208: pj = b->j + bi[i];
209: nz = bi[i+1] - bi[i];
210: for (j=0; j<nz; j++) {
211: PetscArraycpy(pv+bs2*j,rtmp+bs2*pj[j],bs2);
212: }
214: /* Mark diagonal and invert diagonal for simpler triangular solves */
215: pv = b->a + bs2*bdiag[i];
216: pj = b->j + bdiag[i];
217: PetscArraycpy(pv,rtmp+bs2*pj[0],bs2);
218: PetscKernel_A_gets_inverse_A_3(pv,shift,allowzeropivot,&zeropivotdetected);
219: if (zeropivotdetected) B->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
221: /* U part */
222: pj = b->j + bdiag[i+1] + 1;
223: pv = b->a + bs2*(bdiag[i+1]+1);
224: nz = bdiag[i] - bdiag[i+1] - 1;
225: for (j=0; j<nz; j++) {
226: PetscArraycpy(pv+bs2*j,rtmp+bs2*pj[j],bs2);
227: }
228: }
230: PetscFree2(rtmp,mwork);
231: ISRestoreIndices(isicol,&ic);
232: ISRestoreIndices(isrow,&r);
234: C->ops->solve = MatSolve_SeqBAIJ_3;
235: C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3;
236: C->assembled = PETSC_TRUE;
238: PetscLogFlops(1.333333333333*3*3*3*n); /* from inverting diagonal blocks */
239: return 0;
240: }
242: PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering_inplace(Mat C,Mat A,const MatFactorInfo *info)
243: {
244: Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)C->data;
245: PetscInt i,j,n = a->mbs,*bi = b->i,*bj = b->j;
246: PetscInt *ajtmpold,*ajtmp,nz,row;
247: PetscInt *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj;
248: MatScalar *pv,*v,*rtmp,*pc,*w,*x;
249: MatScalar p1,p2,p3,p4,m1,m2,m3,m4,m5,m6,m7,m8,m9,x1,x2,x3,x4;
250: MatScalar p5,p6,p7,p8,p9,x5,x6,x7,x8,x9;
251: MatScalar *ba = b->a,*aa = a->a;
252: PetscReal shift = info->shiftamount;
253: PetscBool allowzeropivot,zeropivotdetected;
255: PetscMalloc1(9*(n+1),&rtmp);
256: allowzeropivot = PetscNot(A->erroriffailure);
258: for (i=0; i<n; i++) {
259: nz = bi[i+1] - bi[i];
260: ajtmp = bj + bi[i];
261: for (j=0; j<nz; j++) {
262: x = rtmp+9*ajtmp[j];
263: x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = 0.0;
264: }
265: /* load in initial (unfactored row) */
266: nz = ai[i+1] - ai[i];
267: ajtmpold = aj + ai[i];
268: v = aa + 9*ai[i];
269: for (j=0; j<nz; j++) {
270: x = rtmp+9*ajtmpold[j];
271: x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3];
272: x[4] = v[4]; x[5] = v[5]; x[6] = v[6]; x[7] = v[7]; x[8] = v[8];
273: v += 9;
274: }
275: row = *ajtmp++;
276: while (row < i) {
277: pc = rtmp + 9*row;
278: p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3];
279: p5 = pc[4]; p6 = pc[5]; p7 = pc[6]; p8 = pc[7]; p9 = pc[8];
280: if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 ||
281: p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0) {
282: pv = ba + 9*diag_offset[row];
283: pj = bj + diag_offset[row] + 1;
284: x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
285: x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8];
286: pc[0] = m1 = p1*x1 + p4*x2 + p7*x3;
287: pc[1] = m2 = p2*x1 + p5*x2 + p8*x3;
288: pc[2] = m3 = p3*x1 + p6*x2 + p9*x3;
290: pc[3] = m4 = p1*x4 + p4*x5 + p7*x6;
291: pc[4] = m5 = p2*x4 + p5*x5 + p8*x6;
292: pc[5] = m6 = p3*x4 + p6*x5 + p9*x6;
294: pc[6] = m7 = p1*x7 + p4*x8 + p7*x9;
295: pc[7] = m8 = p2*x7 + p5*x8 + p8*x9;
296: pc[8] = m9 = p3*x7 + p6*x8 + p9*x9;
298: nz = bi[row+1] - diag_offset[row] - 1;
299: pv += 9;
300: for (j=0; j<nz; j++) {
301: x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
302: x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8];
303: x = rtmp + 9*pj[j];
304: x[0] -= m1*x1 + m4*x2 + m7*x3;
305: x[1] -= m2*x1 + m5*x2 + m8*x3;
306: x[2] -= m3*x1 + m6*x2 + m9*x3;
308: x[3] -= m1*x4 + m4*x5 + m7*x6;
309: x[4] -= m2*x4 + m5*x5 + m8*x6;
310: x[5] -= m3*x4 + m6*x5 + m9*x6;
312: x[6] -= m1*x7 + m4*x8 + m7*x9;
313: x[7] -= m2*x7 + m5*x8 + m8*x9;
314: x[8] -= m3*x7 + m6*x8 + m9*x9;
315: pv += 9;
316: }
317: PetscLogFlops(54.0*nz+36.0);
318: }
319: row = *ajtmp++;
320: }
321: /* finished row so stick it into b->a */
322: pv = ba + 9*bi[i];
323: pj = bj + bi[i];
324: nz = bi[i+1] - bi[i];
325: for (j=0; j<nz; j++) {
326: x = rtmp+9*pj[j];
327: pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3];
328: pv[4] = x[4]; pv[5] = x[5]; pv[6] = x[6]; pv[7] = x[7]; pv[8] = x[8];
329: pv += 9;
330: }
331: /* invert diagonal block */
332: w = ba + 9*diag_offset[i];
333: PetscKernel_A_gets_inverse_A_3(w,shift,allowzeropivot,&zeropivotdetected);
334: if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
335: }
337: PetscFree(rtmp);
339: C->ops->solve = MatSolve_SeqBAIJ_3_NaturalOrdering_inplace;
340: C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3_NaturalOrdering_inplace;
341: C->assembled = PETSC_TRUE;
343: PetscLogFlops(1.333333333333*3*3*3*b->mbs); /* from inverting diagonal blocks */
344: return 0;
345: }
347: /*
348: MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering -
349: copied from MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_inplace()
350: */
351: PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info)
352: {
353: Mat C =B;
354: Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data;
355: PetscInt i,j,k,nz,nzL,row;
356: const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
357: const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag,*pj,bs2=a->bs2;
358: MatScalar *rtmp,*pc,*mwork,*v,*pv,*aa=a->a;
359: PetscInt flg;
360: PetscReal shift = info->shiftamount;
361: PetscBool allowzeropivot,zeropivotdetected;
363: allowzeropivot = PetscNot(A->erroriffailure);
365: /* generate work space needed by the factorization */
366: PetscMalloc2(bs2*n,&rtmp,bs2,&mwork);
367: PetscArrayzero(rtmp,bs2*n);
369: for (i=0; i<n; i++) {
370: /* zero rtmp */
371: /* L part */
372: nz = bi[i+1] - bi[i];
373: bjtmp = bj + bi[i];
374: for (j=0; j<nz; j++) {
375: PetscArrayzero(rtmp+bs2*bjtmp[j],bs2);
376: }
378: /* U part */
379: nz = bdiag[i] - bdiag[i+1];
380: bjtmp = bj + bdiag[i+1] + 1;
381: for (j=0; j<nz; j++) {
382: PetscArrayzero(rtmp+bs2*bjtmp[j],bs2);
383: }
385: /* load in initial (unfactored row) */
386: nz = ai[i+1] - ai[i];
387: ajtmp = aj + ai[i];
388: v = aa + bs2*ai[i];
389: for (j=0; j<nz; j++) {
390: PetscArraycpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2);
391: }
393: /* elimination */
394: bjtmp = bj + bi[i];
395: nzL = bi[i+1] - bi[i];
396: for (k=0; k<nzL; k++) {
397: row = bjtmp[k];
398: pc = rtmp + bs2*row;
399: for (flg=0,j=0; j<bs2; j++) {
400: if (pc[j]!=0.0) {
401: flg = 1;
402: break;
403: }
404: }
405: if (flg) {
406: pv = b->a + bs2*bdiag[row];
407: /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
408: PetscKernel_A_gets_A_times_B_3(pc,pv,mwork);
410: pj = b->j + bdiag[row+1]+1; /* beginning of U(row,:) */
411: pv = b->a + bs2*(bdiag[row+1]+1);
412: nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:) excluding diag */
413: for (j=0; j<nz; j++) {
414: /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
415: /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
416: v = rtmp + bs2*pj[j];
417: PetscKernel_A_gets_A_minus_B_times_C_3(v,pc,pv);
418: pv += bs2;
419: }
420: PetscLogFlops(54.0*nz+45); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
421: }
422: }
424: /* finished row so stick it into b->a */
425: /* L part */
426: pv = b->a + bs2*bi[i];
427: pj = b->j + bi[i];
428: nz = bi[i+1] - bi[i];
429: for (j=0; j<nz; j++) {
430: PetscArraycpy(pv+bs2*j,rtmp+bs2*pj[j],bs2);
431: }
433: /* Mark diagonal and invert diagonal for simpler triangular solves */
434: pv = b->a + bs2*bdiag[i];
435: pj = b->j + bdiag[i];
436: PetscArraycpy(pv,rtmp+bs2*pj[0],bs2);
437: PetscKernel_A_gets_inverse_A_3(pv,shift,allowzeropivot,&zeropivotdetected);
438: if (zeropivotdetected) B->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
440: /* U part */
441: pv = b->a + bs2*(bdiag[i+1]+1);
442: pj = b->j + bdiag[i+1]+1;
443: nz = bdiag[i] - bdiag[i+1] - 1;
444: for (j=0; j<nz; j++) {
445: PetscArraycpy(pv+bs2*j,rtmp+bs2*pj[j],bs2);
446: }
447: }
448: PetscFree2(rtmp,mwork);
450: C->ops->solve = MatSolve_SeqBAIJ_3_NaturalOrdering;
451: C->ops->forwardsolve = MatForwardSolve_SeqBAIJ_3_NaturalOrdering;
452: C->ops->backwardsolve = MatBackwardSolve_SeqBAIJ_3_NaturalOrdering;
453: C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3_NaturalOrdering;
454: C->assembled = PETSC_TRUE;
456: PetscLogFlops(1.333333333333*3*3*3*n); /* from inverting diagonal blocks */
457: return 0;
458: }