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