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