Actual source code: baijfact9.c
petsc-3.9.4 2018-09-11
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: /*
10: Version for when blocks are 5 by 5
11: */
12: PetscErrorCode MatLUFactorNumeric_SeqBAIJ_5_inplace(Mat C,Mat A,const MatFactorInfo *info)
13: {
14: Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)C->data;
15: IS isrow = b->row,isicol = b->icol;
16: PetscErrorCode ierr;
17: const PetscInt *r,*ic,*bi = b->i,*bj = b->j,*ajtmpold,*ajtmp;
18: PetscInt i,j,n = a->mbs,nz,row,idx,ipvt[5];
19: const PetscInt *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj;
20: MatScalar *w,*pv,*rtmp,*x,*pc;
21: const MatScalar *v,*aa = a->a;
22: MatScalar p1,p2,p3,p4,m1,m2,m3,m4,m5,m6,m7,m8,m9,x1,x2,x3,x4;
23: MatScalar p5,p6,p7,p8,p9,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15,x16;
24: MatScalar x17,x18,x19,x20,x21,x22,x23,x24,x25,p10,p11,p12,p13,p14;
25: MatScalar p15,p16,p17,p18,p19,p20,p21,p22,p23,p24,p25,m10,m11,m12;
26: MatScalar m13,m14,m15,m16,m17,m18,m19,m20,m21,m22,m23,m24,m25;
27: MatScalar *ba = b->a,work[25];
28: PetscReal shift = info->shiftamount;
29: PetscBool allowzeropivot,zeropivotdetected;
32: allowzeropivot = PetscNot(A->erroriffailure);
33: ISGetIndices(isrow,&r);
34: ISGetIndices(isicol,&ic);
35: PetscMalloc1(25*(n+1),&rtmp);
37: #define PETSC_USE_MEMZERO 1
38: #define PETSC_USE_MEMCPY 1
40: for (i=0; i<n; i++) {
41: nz = bi[i+1] - bi[i];
42: ajtmp = bj + bi[i];
43: for (j=0; j<nz; j++) {
44: #if defined(PETSC_USE_MEMZERO)
45: PetscMemzero(rtmp+25*ajtmp[j],25*sizeof(PetscScalar));
46: #else
47: x = rtmp+25*ajtmp[j];
48: x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = x[9] = 0.0;
49: x[10] = x[11] = x[12] = x[13] = x[14] = x[15] = x[16] = x[17] = 0.0;
50: x[18] = x[19] = x[20] = x[21] = x[22] = x[23] = x[24] = 0.0;
51: #endif
52: }
53: /* load in initial (unfactored row) */
54: idx = r[i];
55: nz = ai[idx+1] - ai[idx];
56: ajtmpold = aj + ai[idx];
57: v = aa + 25*ai[idx];
58: for (j=0; j<nz; j++) {
59: #if defined(PETSC_USE_MEMCPY)
60: PetscMemcpy(rtmp+25*ic[ajtmpold[j]],v,25*sizeof(PetscScalar));
61: #else
62: x = rtmp+25*ic[ajtmpold[j]];
63: x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3];
64: x[4] = v[4]; x[5] = v[5]; x[6] = v[6]; x[7] = v[7]; x[8] = v[8];
65: x[9] = v[9]; x[10] = v[10]; x[11] = v[11]; x[12] = v[12]; x[13] = v[13];
66: x[14] = v[14]; x[15] = v[15]; x[16] = v[16]; x[17] = v[17];
67: x[18] = v[18]; x[19] = v[19]; x[20] = v[20]; x[21] = v[21];
68: x[22] = v[22]; x[23] = v[23]; x[24] = v[24];
69: #endif
70: v += 25;
71: }
72: row = *ajtmp++;
73: while (row < i) {
74: pc = rtmp + 25*row;
75: p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3];
76: p5 = pc[4]; p6 = pc[5]; p7 = pc[6]; p8 = pc[7]; p9 = pc[8];
77: p10 = pc[9]; p11 = pc[10]; p12 = pc[11]; p13 = pc[12]; p14 = pc[13];
78: p15 = pc[14]; p16 = pc[15]; p17 = pc[16]; p18 = pc[17]; p19 = pc[18];
79: p20 = pc[19]; p21 = pc[20]; p22 = pc[21]; p23 = pc[22]; p24 = pc[23];
80: p25 = pc[24];
81: if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 ||
82: p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0 || p10 != 0.0 ||
83: p11 != 0.0 || p12 != 0.0 || p13 != 0.0 || p14 != 0.0 || p15 != 0.0
84: || p16 != 0.0 || p17 != 0.0 || p18 != 0.0 || p19 != 0.0 ||
85: p20 != 0.0 || p21 != 0.0 || p22 != 0.0 || p23 != 0.0 ||
86: p24 != 0.0 || p25 != 0.0) {
87: pv = ba + 25*diag_offset[row];
88: pj = bj + diag_offset[row] + 1;
89: x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
90: x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8];
91: x10 = pv[9]; x11 = pv[10]; x12 = pv[11]; x13 = pv[12]; x14 = pv[13];
92: x15 = pv[14]; x16 = pv[15]; x17 = pv[16]; x18 = pv[17];
93: x19 = pv[18]; x20 = pv[19]; x21 = pv[20]; x22 = pv[21];
94: x23 = pv[22]; x24 = pv[23]; x25 = pv[24];
95: pc[0] = m1 = p1*x1 + p6*x2 + p11*x3 + p16*x4 + p21*x5;
96: pc[1] = m2 = p2*x1 + p7*x2 + p12*x3 + p17*x4 + p22*x5;
97: pc[2] = m3 = p3*x1 + p8*x2 + p13*x3 + p18*x4 + p23*x5;
98: pc[3] = m4 = p4*x1 + p9*x2 + p14*x3 + p19*x4 + p24*x5;
99: pc[4] = m5 = p5*x1 + p10*x2 + p15*x3 + p20*x4 + p25*x5;
101: pc[5] = m6 = p1*x6 + p6*x7 + p11*x8 + p16*x9 + p21*x10;
102: pc[6] = m7 = p2*x6 + p7*x7 + p12*x8 + p17*x9 + p22*x10;
103: pc[7] = m8 = p3*x6 + p8*x7 + p13*x8 + p18*x9 + p23*x10;
104: pc[8] = m9 = p4*x6 + p9*x7 + p14*x8 + p19*x9 + p24*x10;
105: pc[9] = m10 = p5*x6 + p10*x7 + p15*x8 + p20*x9 + p25*x10;
107: pc[10] = m11 = p1*x11 + p6*x12 + p11*x13 + p16*x14 + p21*x15;
108: pc[11] = m12 = p2*x11 + p7*x12 + p12*x13 + p17*x14 + p22*x15;
109: pc[12] = m13 = p3*x11 + p8*x12 + p13*x13 + p18*x14 + p23*x15;
110: pc[13] = m14 = p4*x11 + p9*x12 + p14*x13 + p19*x14 + p24*x15;
111: pc[14] = m15 = p5*x11 + p10*x12 + p15*x13 + p20*x14 + p25*x15;
113: pc[15] = m16 = p1*x16 + p6*x17 + p11*x18 + p16*x19 + p21*x20;
114: pc[16] = m17 = p2*x16 + p7*x17 + p12*x18 + p17*x19 + p22*x20;
115: pc[17] = m18 = p3*x16 + p8*x17 + p13*x18 + p18*x19 + p23*x20;
116: pc[18] = m19 = p4*x16 + p9*x17 + p14*x18 + p19*x19 + p24*x20;
117: pc[19] = m20 = p5*x16 + p10*x17 + p15*x18 + p20*x19 + p25*x20;
119: pc[20] = m21 = p1*x21 + p6*x22 + p11*x23 + p16*x24 + p21*x25;
120: pc[21] = m22 = p2*x21 + p7*x22 + p12*x23 + p17*x24 + p22*x25;
121: pc[22] = m23 = p3*x21 + p8*x22 + p13*x23 + p18*x24 + p23*x25;
122: pc[23] = m24 = p4*x21 + p9*x22 + p14*x23 + p19*x24 + p24*x25;
123: pc[24] = m25 = p5*x21 + p10*x22 + p15*x23 + p20*x24 + p25*x25;
125: nz = bi[row+1] - diag_offset[row] - 1;
126: pv += 25;
127: for (j=0; j<nz; j++) {
128: x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
129: x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8];
130: x10 = pv[9]; x11 = pv[10]; x12 = pv[11]; x13 = pv[12];
131: x14 = pv[13]; x15 = pv[14]; x16 = pv[15]; x17 = pv[16];
132: x18 = pv[17]; x19 = pv[18]; x20 = pv[19]; x21 = pv[20];
133: x22 = pv[21]; x23 = pv[22]; x24 = pv[23]; x25 = pv[24];
134: x = rtmp + 25*pj[j];
135: x[0] -= m1*x1 + m6*x2 + m11*x3 + m16*x4 + m21*x5;
136: x[1] -= m2*x1 + m7*x2 + m12*x3 + m17*x4 + m22*x5;
137: x[2] -= m3*x1 + m8*x2 + m13*x3 + m18*x4 + m23*x5;
138: x[3] -= m4*x1 + m9*x2 + m14*x3 + m19*x4 + m24*x5;
139: x[4] -= m5*x1 + m10*x2 + m15*x3 + m20*x4 + m25*x5;
141: x[5] -= m1*x6 + m6*x7 + m11*x8 + m16*x9 + m21*x10;
142: x[6] -= m2*x6 + m7*x7 + m12*x8 + m17*x9 + m22*x10;
143: x[7] -= m3*x6 + m8*x7 + m13*x8 + m18*x9 + m23*x10;
144: x[8] -= m4*x6 + m9*x7 + m14*x8 + m19*x9 + m24*x10;
145: x[9] -= m5*x6 + m10*x7 + m15*x8 + m20*x9 + m25*x10;
147: x[10] -= m1*x11 + m6*x12 + m11*x13 + m16*x14 + m21*x15;
148: x[11] -= m2*x11 + m7*x12 + m12*x13 + m17*x14 + m22*x15;
149: x[12] -= m3*x11 + m8*x12 + m13*x13 + m18*x14 + m23*x15;
150: x[13] -= m4*x11 + m9*x12 + m14*x13 + m19*x14 + m24*x15;
151: x[14] -= m5*x11 + m10*x12 + m15*x13 + m20*x14 + m25*x15;
153: x[15] -= m1*x16 + m6*x17 + m11*x18 + m16*x19 + m21*x20;
154: x[16] -= m2*x16 + m7*x17 + m12*x18 + m17*x19 + m22*x20;
155: x[17] -= m3*x16 + m8*x17 + m13*x18 + m18*x19 + m23*x20;
156: x[18] -= m4*x16 + m9*x17 + m14*x18 + m19*x19 + m24*x20;
157: x[19] -= m5*x16 + m10*x17 + m15*x18 + m20*x19 + m25*x20;
159: x[20] -= m1*x21 + m6*x22 + m11*x23 + m16*x24 + m21*x25;
160: x[21] -= m2*x21 + m7*x22 + m12*x23 + m17*x24 + m22*x25;
161: x[22] -= m3*x21 + m8*x22 + m13*x23 + m18*x24 + m23*x25;
162: x[23] -= m4*x21 + m9*x22 + m14*x23 + m19*x24 + m24*x25;
163: x[24] -= m5*x21 + m10*x22 + m15*x23 + m20*x24 + m25*x25;
165: pv += 25;
166: }
167: PetscLogFlops(250.0*nz+225.0);
168: }
169: row = *ajtmp++;
170: }
171: /* finished row so stick it into b->a */
172: pv = ba + 25*bi[i];
173: pj = bj + bi[i];
174: nz = bi[i+1] - bi[i];
175: for (j=0; j<nz; j++) {
176: #if defined(PETSC_USE_MEMCPY)
177: PetscMemcpy(pv,rtmp+25*pj[j],25*sizeof(PetscScalar));
178: #else
179: x = rtmp+25*pj[j];
180: pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3];
181: pv[4] = x[4]; pv[5] = x[5]; pv[6] = x[6]; pv[7] = x[7]; pv[8] = x[8];
182: pv[9] = x[9]; pv[10] = x[10]; pv[11] = x[11]; pv[12] = x[12];
183: pv[13] = x[13]; pv[14] = x[14]; pv[15] = x[15]; pv[16] = x[16];
184: pv[17] = x[17]; pv[18] = x[18]; pv[19] = x[19]; pv[20] = x[20];
185: pv[21] = x[21]; pv[22] = x[22]; pv[23] = x[23]; pv[24] = x[24];
186: #endif
187: pv += 25;
188: }
189: /* invert diagonal block */
190: w = ba + 25*diag_offset[i];
191: PetscKernel_A_gets_inverse_A_5(w,ipvt,work,shift,allowzeropivot,&zeropivotdetected);
192: if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
193: }
195: PetscFree(rtmp);
196: ISRestoreIndices(isicol,&ic);
197: ISRestoreIndices(isrow,&r);
199: C->ops->solve = MatSolve_SeqBAIJ_5_inplace;
200: C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_5_inplace;
201: C->assembled = PETSC_TRUE;
203: PetscLogFlops(1.333333333333*5*5*5*b->mbs); /* from inverting diagonal blocks */
204: return(0);
205: }
207: /* MatLUFactorNumeric_SeqBAIJ_5 -
208: copied from MatLUFactorNumeric_SeqBAIJ_N_inplace() and manually re-implemented
209: PetscKernel_A_gets_A_times_B()
210: PetscKernel_A_gets_A_minus_B_times_C()
211: PetscKernel_A_gets_inverse_A()
212: */
214: PetscErrorCode MatLUFactorNumeric_SeqBAIJ_5(Mat B,Mat A,const MatFactorInfo *info)
215: {
216: Mat C =B;
217: Mat_SeqBAIJ *a =(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data;
218: IS isrow = b->row,isicol = b->icol;
220: const PetscInt *r,*ic;
221: PetscInt i,j,k,nz,nzL,row;
222: const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
223: const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag,*pj,bs2=a->bs2;
224: MatScalar *rtmp,*pc,*mwork,*v,*pv,*aa=a->a,work[25];
225: PetscInt flg,ipvt[5];
226: PetscReal shift = info->shiftamount;
227: PetscBool allowzeropivot,zeropivotdetected;
230: allowzeropivot = PetscNot(A->erroriffailure);
231: ISGetIndices(isrow,&r);
232: ISGetIndices(isicol,&ic);
234: /* generate work space needed by the factorization */
235: PetscMalloc2(bs2*n,&rtmp,bs2,&mwork);
236: PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));
238: for (i=0; i<n; i++) {
239: /* zero rtmp */
240: /* L part */
241: nz = bi[i+1] - bi[i];
242: bjtmp = bj + bi[i];
243: for (j=0; j<nz; j++) {
244: PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));
245: }
247: /* U part */
248: nz = bdiag[i] - bdiag[i+1];
249: bjtmp = bj + bdiag[i+1]+1;
250: for (j=0; j<nz; j++) {
251: PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));
252: }
254: /* load in initial (unfactored row) */
255: nz = ai[r[i]+1] - ai[r[i]];
256: ajtmp = aj + ai[r[i]];
257: v = aa + bs2*ai[r[i]];
258: for (j=0; j<nz; j++) {
259: PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],v+bs2*j,bs2*sizeof(MatScalar));
260: }
262: /* elimination */
263: bjtmp = bj + bi[i];
264: nzL = bi[i+1] - bi[i];
265: for (k=0; k < nzL; k++) {
266: row = bjtmp[k];
267: pc = rtmp + bs2*row;
268: for (flg=0,j=0; j<bs2; j++) {
269: if (pc[j]!=0.0) {
270: flg = 1;
271: break;
272: }
273: }
274: if (flg) {
275: pv = b->a + bs2*bdiag[row];
276: /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
277: PetscKernel_A_gets_A_times_B_5(pc,pv,mwork);
279: pj = b->j + bdiag[row+1]+1; /* begining of U(row,:) */
280: pv = b->a + bs2*(bdiag[row+1]+1);
281: nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries inU(row,:), excluding diag */
282: for (j=0; j<nz; j++) {
283: /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
284: /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
285: v = rtmp + bs2*pj[j];
286: PetscKernel_A_gets_A_minus_B_times_C_5(v,pc,pv);
287: pv += bs2;
288: }
289: PetscLogFlops(250*nz+225); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
290: }
291: }
293: /* finished row so stick it into b->a */
294: /* L part */
295: pv = b->a + bs2*bi[i];
296: pj = b->j + bi[i];
297: nz = bi[i+1] - bi[i];
298: for (j=0; j<nz; j++) {
299: PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));
300: }
302: /* Mark diagonal and invert diagonal for simplier triangular solves */
303: pv = b->a + bs2*bdiag[i];
304: pj = b->j + bdiag[i];
305: PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));
306: PetscKernel_A_gets_inverse_A_5(pv,ipvt,work,shift,allowzeropivot,&zeropivotdetected);
307: if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
309: /* U part */
310: pv = b->a + bs2*(bdiag[i+1]+1);
311: pj = b->j + bdiag[i+1]+1;
312: nz = bdiag[i] - bdiag[i+1] - 1;
313: for (j=0; j<nz; j++) {
314: PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));
315: }
316: }
318: PetscFree2(rtmp,mwork);
319: ISRestoreIndices(isicol,&ic);
320: ISRestoreIndices(isrow,&r);
322: C->ops->solve = MatSolve_SeqBAIJ_5;
323: C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_5;
324: C->assembled = PETSC_TRUE;
326: PetscLogFlops(1.333333333333*5*5*5*n); /* from inverting diagonal blocks */
327: return(0);
328: }
330: /*
331: Version for when blocks are 5 by 5 Using natural ordering
332: */
333: PetscErrorCode MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering_inplace(Mat C,Mat A,const MatFactorInfo *info)
334: {
335: Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)C->data;
337: PetscInt i,j,n = a->mbs,*bi = b->i,*bj = b->j,ipvt[5];
338: PetscInt *ajtmpold,*ajtmp,nz,row;
339: PetscInt *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj;
340: MatScalar *pv,*v,*rtmp,*pc,*w,*x;
341: MatScalar x1,x2,x3,x4,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15;
342: MatScalar x16,x17,x18,x19,x20,x21,x22,x23,x24,x25;
343: MatScalar p1,p2,p3,p4,p5,p6,p7,p8,p9,p10,p11,p12,p13,p14,p15;
344: MatScalar p16,p17,p18,p19,p20,p21,p22,p23,p24,p25;
345: MatScalar m1,m2,m3,m4,m5,m6,m7,m8,m9,m10,m11,m12,m13,m14,m15;
346: MatScalar m16,m17,m18,m19,m20,m21,m22,m23,m24,m25;
347: MatScalar *ba = b->a,*aa = a->a,work[25];
348: PetscReal shift = info->shiftamount;
349: PetscBool allowzeropivot,zeropivotdetected;
352: allowzeropivot = PetscNot(A->erroriffailure);
353: PetscMalloc1(25*(n+1),&rtmp);
354: for (i=0; i<n; i++) {
355: nz = bi[i+1] - bi[i];
356: ajtmp = bj + bi[i];
357: for (j=0; j<nz; j++) {
358: x = rtmp+25*ajtmp[j];
359: x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = x[9] = 0.0;
360: x[10] = x[11] = x[12] = x[13] = x[14] = x[15] = 0.0;
361: x[16] = x[17] = x[18] = x[19] = x[20] = x[21] = x[22] = x[23] = x[24] = 0.0;
362: }
363: /* load in initial (unfactored row) */
364: nz = ai[i+1] - ai[i];
365: ajtmpold = aj + ai[i];
366: v = aa + 25*ai[i];
367: for (j=0; j<nz; j++) {
368: x = rtmp+25*ajtmpold[j];
369: x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3];
370: x[4] = v[4]; x[5] = v[5]; x[6] = v[6]; x[7] = v[7]; x[8] = v[8];
371: x[9] = v[9]; x[10] = v[10]; x[11] = v[11]; x[12] = v[12]; x[13] = v[13];
372: x[14] = v[14]; x[15] = v[15]; x[16] = v[16]; x[17] = v[17]; x[18] = v[18];
373: x[19] = v[19]; x[20] = v[20]; x[21] = v[21]; x[22] = v[22]; x[23] = v[23];
374: x[24] = v[24];
375: v += 25;
376: }
377: row = *ajtmp++;
378: while (row < i) {
379: pc = rtmp + 25*row;
380: p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3];
381: p5 = pc[4]; p6 = pc[5]; p7 = pc[6]; p8 = pc[7]; p9 = pc[8];
382: p10 = pc[9]; p11 = pc[10]; p12 = pc[11]; p13 = pc[12]; p14 = pc[13];
383: p15 = pc[14]; p16 = pc[15]; p17 = pc[16]; p18 = pc[17];
384: p19 = pc[18]; p20 = pc[19]; p21 = pc[20]; p22 = pc[21]; p23 = pc[22];
385: p24 = pc[23]; p25 = pc[24];
386: if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 ||
387: p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0 || p10 != 0.0 ||
388: p11 != 0.0 || p12 != 0.0 || p13 != 0.0 || p14 != 0.0 || p15 != 0.0
389: || p16 != 0.0 || p17 != 0.0 || p18 != 0.0 || p19 != 0.0 || p20 != 0.0
390: || p21 != 0.0 || p22 != 0.0 || p23 != 0.0 || p24 != 0.0 || p25 != 0.0) {
391: pv = ba + 25*diag_offset[row];
392: pj = bj + diag_offset[row] + 1;
393: x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
394: x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8];
395: x10 = pv[9]; x11 = pv[10]; x12 = pv[11]; x13 = pv[12]; x14 = pv[13];
396: x15 = pv[14]; x16 = pv[15]; x17 = pv[16]; x18 = pv[17]; x19 = pv[18];
397: x20 = pv[19]; x21 = pv[20]; x22 = pv[21]; x23 = pv[22]; x24 = pv[23];
398: x25 = pv[24];
399: pc[0] = m1 = p1*x1 + p6*x2 + p11*x3 + p16*x4 + p21*x5;
400: pc[1] = m2 = p2*x1 + p7*x2 + p12*x3 + p17*x4 + p22*x5;
401: pc[2] = m3 = p3*x1 + p8*x2 + p13*x3 + p18*x4 + p23*x5;
402: pc[3] = m4 = p4*x1 + p9*x2 + p14*x3 + p19*x4 + p24*x5;
403: pc[4] = m5 = p5*x1 + p10*x2 + p15*x3 + p20*x4 + p25*x5;
405: pc[5] = m6 = p1*x6 + p6*x7 + p11*x8 + p16*x9 + p21*x10;
406: pc[6] = m7 = p2*x6 + p7*x7 + p12*x8 + p17*x9 + p22*x10;
407: pc[7] = m8 = p3*x6 + p8*x7 + p13*x8 + p18*x9 + p23*x10;
408: pc[8] = m9 = p4*x6 + p9*x7 + p14*x8 + p19*x9 + p24*x10;
409: pc[9] = m10 = p5*x6 + p10*x7 + p15*x8 + p20*x9 + p25*x10;
411: pc[10] = m11 = p1*x11 + p6*x12 + p11*x13 + p16*x14 + p21*x15;
412: pc[11] = m12 = p2*x11 + p7*x12 + p12*x13 + p17*x14 + p22*x15;
413: pc[12] = m13 = p3*x11 + p8*x12 + p13*x13 + p18*x14 + p23*x15;
414: pc[13] = m14 = p4*x11 + p9*x12 + p14*x13 + p19*x14 + p24*x15;
415: pc[14] = m15 = p5*x11 + p10*x12 + p15*x13 + p20*x14 + p25*x15;
417: pc[15] = m16 = p1*x16 + p6*x17 + p11*x18 + p16*x19 + p21*x20;
418: pc[16] = m17 = p2*x16 + p7*x17 + p12*x18 + p17*x19 + p22*x20;
419: pc[17] = m18 = p3*x16 + p8*x17 + p13*x18 + p18*x19 + p23*x20;
420: pc[18] = m19 = p4*x16 + p9*x17 + p14*x18 + p19*x19 + p24*x20;
421: pc[19] = m20 = p5*x16 + p10*x17 + p15*x18 + p20*x19 + p25*x20;
423: pc[20] = m21 = p1*x21 + p6*x22 + p11*x23 + p16*x24 + p21*x25;
424: pc[21] = m22 = p2*x21 + p7*x22 + p12*x23 + p17*x24 + p22*x25;
425: pc[22] = m23 = p3*x21 + p8*x22 + p13*x23 + p18*x24 + p23*x25;
426: pc[23] = m24 = p4*x21 + p9*x22 + p14*x23 + p19*x24 + p24*x25;
427: pc[24] = m25 = p5*x21 + p10*x22 + p15*x23 + p20*x24 + p25*x25;
429: nz = bi[row+1] - diag_offset[row] - 1;
430: pv += 25;
431: for (j=0; j<nz; j++) {
432: x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
433: x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8];
434: x10 = pv[9]; x11 = pv[10]; x12 = pv[11]; x13 = pv[12];
435: x14 = pv[13]; x15 = pv[14]; x16 = pv[15]; x17 = pv[16]; x18 = pv[17];
436: x19 = pv[18]; x20 = pv[19]; x21 = pv[20]; x22 = pv[21]; x23 = pv[22];
437: x24 = pv[23]; x25 = pv[24];
438: x = rtmp + 25*pj[j];
439: x[0] -= m1*x1 + m6*x2 + m11*x3 + m16*x4 + m21*x5;
440: x[1] -= m2*x1 + m7*x2 + m12*x3 + m17*x4 + m22*x5;
441: x[2] -= m3*x1 + m8*x2 + m13*x3 + m18*x4 + m23*x5;
442: x[3] -= m4*x1 + m9*x2 + m14*x3 + m19*x4 + m24*x5;
443: x[4] -= m5*x1 + m10*x2 + m15*x3 + m20*x4 + m25*x5;
445: x[5] -= m1*x6 + m6*x7 + m11*x8 + m16*x9 + m21*x10;
446: x[6] -= m2*x6 + m7*x7 + m12*x8 + m17*x9 + m22*x10;
447: x[7] -= m3*x6 + m8*x7 + m13*x8 + m18*x9 + m23*x10;
448: x[8] -= m4*x6 + m9*x7 + m14*x8 + m19*x9 + m24*x10;
449: x[9] -= m5*x6 + m10*x7 + m15*x8 + m20*x9 + m25*x10;
451: x[10] -= m1*x11 + m6*x12 + m11*x13 + m16*x14 + m21*x15;
452: x[11] -= m2*x11 + m7*x12 + m12*x13 + m17*x14 + m22*x15;
453: x[12] -= m3*x11 + m8*x12 + m13*x13 + m18*x14 + m23*x15;
454: x[13] -= m4*x11 + m9*x12 + m14*x13 + m19*x14 + m24*x15;
455: x[14] -= m5*x11 + m10*x12 + m15*x13 + m20*x14 + m25*x15;
457: x[15] -= m1*x16 + m6*x17 + m11*x18 + m16*x19 + m21*x20;
458: x[16] -= m2*x16 + m7*x17 + m12*x18 + m17*x19 + m22*x20;
459: x[17] -= m3*x16 + m8*x17 + m13*x18 + m18*x19 + m23*x20;
460: x[18] -= m4*x16 + m9*x17 + m14*x18 + m19*x19 + m24*x20;
461: x[19] -= m5*x16 + m10*x17 + m15*x18 + m20*x19 + m25*x20;
463: x[20] -= m1*x21 + m6*x22 + m11*x23 + m16*x24 + m21*x25;
464: x[21] -= m2*x21 + m7*x22 + m12*x23 + m17*x24 + m22*x25;
465: x[22] -= m3*x21 + m8*x22 + m13*x23 + m18*x24 + m23*x25;
466: x[23] -= m4*x21 + m9*x22 + m14*x23 + m19*x24 + m24*x25;
467: x[24] -= m5*x21 + m10*x22 + m15*x23 + m20*x24 + m25*x25;
468: pv += 25;
469: }
470: PetscLogFlops(250.0*nz+225.0);
471: }
472: row = *ajtmp++;
473: }
474: /* finished row so stick it into b->a */
475: pv = ba + 25*bi[i];
476: pj = bj + bi[i];
477: nz = bi[i+1] - bi[i];
478: for (j=0; j<nz; j++) {
479: x = rtmp+25*pj[j];
480: pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3];
481: pv[4] = x[4]; pv[5] = x[5]; pv[6] = x[6]; pv[7] = x[7]; pv[8] = x[8];
482: pv[9] = x[9]; pv[10] = x[10]; pv[11] = x[11]; pv[12] = x[12];
483: pv[13] = x[13]; pv[14] = x[14]; pv[15] = x[15]; pv[16] = x[16]; pv[17] = x[17];
484: pv[18] = x[18]; pv[19] = x[19]; pv[20] = x[20]; pv[21] = x[21]; pv[22] = x[22];
485: pv[23] = x[23]; pv[24] = x[24];
486: pv += 25;
487: }
488: /* invert diagonal block */
489: w = ba + 25*diag_offset[i];
490: PetscKernel_A_gets_inverse_A_5(w,ipvt,work,shift,allowzeropivot,&zeropivotdetected);
491: if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
492: }
494: PetscFree(rtmp);
496: C->ops->solve = MatSolve_SeqBAIJ_5_NaturalOrdering_inplace;
497: C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_5_NaturalOrdering_inplace;
498: C->assembled = PETSC_TRUE;
500: PetscLogFlops(1.333333333333*5*5*5*b->mbs); /* from inverting diagonal blocks */
501: return(0);
502: }
504: PetscErrorCode MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info)
505: {
506: Mat C =B;
507: Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data;
509: PetscInt i,j,k,nz,nzL,row;
510: const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
511: const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag,*pj,bs2=a->bs2;
512: MatScalar *rtmp,*pc,*mwork,*v,*vv,*pv,*aa=a->a,work[25];
513: PetscInt flg,ipvt[5];
514: PetscReal shift = info->shiftamount;
515: PetscBool allowzeropivot,zeropivotdetected;
518: allowzeropivot = PetscNot(A->erroriffailure);
520: /* generate work space needed by the factorization */
521: PetscMalloc2(bs2*n,&rtmp,bs2,&mwork);
522: PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));
524: for (i=0; i<n; i++) {
525: /* zero rtmp */
526: /* L part */
527: nz = bi[i+1] - bi[i];
528: bjtmp = bj + bi[i];
529: for (j=0; j<nz; j++) {
530: PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));
531: }
533: /* U part */
534: nz = bdiag[i] - bdiag[i+1];
535: bjtmp = bj + bdiag[i+1]+1;
536: for (j=0; j<nz; j++) {
537: PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));
538: }
540: /* load in initial (unfactored row) */
541: nz = ai[i+1] - ai[i];
542: ajtmp = aj + ai[i];
543: v = aa + bs2*ai[i];
544: for (j=0; j<nz; j++) {
545: PetscMemcpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2*sizeof(MatScalar));
546: }
548: /* elimination */
549: bjtmp = bj + bi[i];
550: nzL = bi[i+1] - bi[i];
551: for (k=0; k < nzL; k++) {
552: row = bjtmp[k];
553: pc = rtmp + bs2*row;
554: for (flg=0,j=0; j<bs2; j++) {
555: if (pc[j]!=0.0) {
556: flg = 1;
557: break;
558: }
559: }
560: if (flg) {
561: pv = b->a + bs2*bdiag[row];
562: /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
563: PetscKernel_A_gets_A_times_B_5(pc,pv,mwork);
565: pj = b->j + bdiag[row+1]+1; /* begining of U(row,:) */
566: pv = b->a + bs2*(bdiag[row+1]+1);
567: nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries inU(row,:), excluding diag */
568: for (j=0; j<nz; j++) {
569: /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
570: /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
571: vv = rtmp + bs2*pj[j];
572: PetscKernel_A_gets_A_minus_B_times_C_5(vv,pc,pv);
573: pv += bs2;
574: }
575: PetscLogFlops(250*nz+225); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
576: }
577: }
579: /* finished row so stick it into b->a */
580: /* L part */
581: pv = b->a + bs2*bi[i];
582: pj = b->j + bi[i];
583: nz = bi[i+1] - bi[i];
584: for (j=0; j<nz; j++) {
585: PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));
586: }
588: /* Mark diagonal and invert diagonal for simplier triangular solves */
589: pv = b->a + bs2*bdiag[i];
590: pj = b->j + bdiag[i];
591: PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));
592: PetscKernel_A_gets_inverse_A_5(pv,ipvt,work,shift,allowzeropivot,&zeropivotdetected);
593: if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
595: /* U part */
596: pv = b->a + bs2*(bdiag[i+1]+1);
597: pj = b->j + bdiag[i+1]+1;
598: nz = bdiag[i] - bdiag[i+1] - 1;
599: for (j=0; j<nz; j++) {
600: PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));
601: }
602: }
603: PetscFree2(rtmp,mwork);
605: C->ops->solve = MatSolve_SeqBAIJ_5_NaturalOrdering;
606: C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_5_NaturalOrdering;
607: C->assembled = PETSC_TRUE;
609: PetscLogFlops(1.333333333333*5*5*5*n); /* from inverting diagonal blocks */
610: return(0);
611: }