Actual source code: baijfact7.c
1: #define PETSCMAT_DLL
3: /*
4: Factorization code for BAIJ format.
5: */
6: #include ../src/mat/impls/baij/seq/baij.h
7: #include ../src/inline/ilu.h
9: /* ------------------------------------------------------------*/
10: /*
11: Version for when blocks are 6 by 6
12: */
15: PetscErrorCode MatLUFactorNumeric_SeqBAIJ_6(Mat C,Mat A,const MatFactorInfo *info)
16: {
17: Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data;
18: IS isrow = b->row,isicol = b->icol;
20: const PetscInt *ajtmpold,*ajtmp,*diag_offset = b->diag,*r,*ic,*bi = b->i,*bj = b->j,*ai=a->i,*aj=a->j,*pj;
21: PetscInt nz,row,i,j,n = a->mbs,idx;
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,x10,x11,x12,x13,x14,x15,x16;
25: MatScalar x17,x18,x19,x20,x21,x22,x23,x24,x25,p10,p11,p12,p13,p14;
26: MatScalar p15,p16,p17,p18,p19,p20,p21,p22,p23,p24,p25,m10,m11,m12;
27: MatScalar m13,m14,m15,m16,m17,m18,m19,m20,m21,m22,m23,m24,m25;
28: MatScalar p26,p27,p28,p29,p30,p31,p32,p33,p34,p35,p36;
29: MatScalar x26,x27,x28,x29,x30,x31,x32,x33,x34,x35,x36;
30: MatScalar m26,m27,m28,m29,m30,m31,m32,m33,m34,m35,m36;
31: MatScalar *ba = b->a,*aa = a->a;
32: PetscReal shift = info->shiftinblocks;
35: ISGetIndices(isrow,&r);
36: ISGetIndices(isicol,&ic);
37: PetscMalloc(36*(n+1)*sizeof(MatScalar),&rtmp);
39: for (i=0; i<n; i++) {
40: nz = bi[i+1] - bi[i];
41: ajtmp = bj + bi[i];
42: for (j=0; j<nz; j++) {
43: x = rtmp+36*ajtmp[j];
44: x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = x[9] = 0.0;
45: x[10] = x[11] = x[12] = x[13] = x[14] = x[15] = x[16] = x[17] = 0.0;
46: x[18] = x[19] = x[20] = x[21] = x[22] = x[23] = x[24] = x[25] = 0.0 ;
47: x[26] = x[27] = x[28] = x[29] = x[30] = x[31] = x[32] = x[33] = 0.0 ;
48: x[34] = x[35] = 0.0 ;
49: }
50: /* load in initial (unfactored row) */
51: idx = r[i];
52: nz = ai[idx+1] - ai[idx];
53: ajtmpold = aj + ai[idx];
54: v = aa + 36*ai[idx];
55: for (j=0; j<nz; j++) {
56: x = rtmp+36*ic[ajtmpold[j]];
57: x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3];
58: x[4] = v[4]; x[5] = v[5]; x[6] = v[6]; x[7] = v[7];
59: x[8] = v[8]; x[9] = v[9]; x[10] = v[10]; x[11] = v[11];
60: x[12] = v[12]; x[13] = v[13]; x[14] = v[14]; x[15] = v[15];
61: x[16] = v[16]; x[17] = v[17]; x[18] = v[18]; x[19] = v[19];
62: x[20] = v[20]; x[21] = v[21]; x[22] = v[22]; x[23] = v[23];
63: x[24] = v[24]; x[25] = v[25]; x[26] = v[26]; x[27] = v[27];
64: x[28] = v[28]; x[29] = v[29]; x[30] = v[30]; x[31] = v[31];
65: x[32] = v[32]; x[33] = v[33]; x[34] = v[34]; x[35] = v[35];
66: v += 36;
67: }
68: row = *ajtmp++;
69: while (row < i) {
70: pc = rtmp + 36*row;
71: p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3];
72: p5 = pc[4]; p6 = pc[5]; p7 = pc[6]; p8 = pc[7];
73: p9 = pc[8]; p10 = pc[9]; p11 = pc[10]; p12 = pc[11];
74: p13 = pc[12]; p14 = pc[13]; p15 = pc[14]; p16 = pc[15];
75: p17 = pc[16]; p18 = pc[17]; p19 = pc[18]; p20 = pc[19];
76: p21 = pc[20]; p22 = pc[21]; p23 = pc[22]; p24 = pc[23];
77: p25 = pc[24]; p26 = pc[25]; p27 = pc[26]; p28 = pc[27];
78: p29 = pc[28]; p30 = pc[29]; p31 = pc[30]; p32 = pc[31];
79: p33 = pc[32]; p34 = pc[33]; p35 = pc[34]; p36 = pc[35];
80: if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 ||
81: p5 != 0.0 || p6 != 0.0 || p7 != 0.0 || p8 != 0.0 ||
82: p9 != 0.0 || p10 != 0.0 || p11 != 0.0 || p12 != 0.0 ||
83: p13 != 0.0 || p14 != 0.0 || p15 != 0.0 || p16 != 0.0 ||
84: p17 != 0.0 || p18 != 0.0 || p19 != 0.0 || p20 != 0.0 ||
85: p21 != 0.0 || p22 != 0.0 || p23 != 0.0 || p24 != 0.0 ||
86: p25 != 0.0 || p26 != 0.0 || p27 != 0.0 || p28 != 0.0 ||
87: p29 != 0.0 || p30 != 0.0 || p31 != 0.0 || p32 != 0.0 ||
88: p33 != 0.0 || p34 != 0.0 || p35 != 0.0 || p36 != 0.0) {
89: pv = ba + 36*diag_offset[row];
90: pj = bj + diag_offset[row] + 1;
91: x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
92: x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7];
93: x9 = pv[8]; x10 = pv[9]; x11 = pv[10]; x12 = pv[11];
94: x13 = pv[12]; x14 = pv[13]; x15 = pv[14]; x16 = pv[15];
95: x17 = pv[16]; x18 = pv[17]; x19 = pv[18]; x20 = pv[19];
96: x21 = pv[20]; x22 = pv[21]; x23 = pv[22]; x24 = pv[23];
97: x25 = pv[24]; x26 = pv[25]; x27 = pv[26]; x28 = pv[27];
98: x29 = pv[28]; x30 = pv[29]; x31 = pv[30]; x32 = pv[31];
99: x33 = pv[32]; x34 = pv[33]; x35 = pv[34]; x36 = pv[35];
100: pc[0] = m1 = p1*x1 + p7*x2 + p13*x3 + p19*x4 + p25*x5 + p31*x6;
101: pc[1] = m2 = p2*x1 + p8*x2 + p14*x3 + p20*x4 + p26*x5 + p32*x6;
102: pc[2] = m3 = p3*x1 + p9*x2 + p15*x3 + p21*x4 + p27*x5 + p33*x6;
103: pc[3] = m4 = p4*x1 + p10*x2 + p16*x3 + p22*x4 + p28*x5 + p34*x6;
104: pc[4] = m5 = p5*x1 + p11*x2 + p17*x3 + p23*x4 + p29*x5 + p35*x6;
105: pc[5] = m6 = p6*x1 + p12*x2 + p18*x3 + p24*x4 + p30*x5 + p36*x6;
107: pc[6] = m7 = p1*x7 + p7*x8 + p13*x9 + p19*x10 + p25*x11 + p31*x12;
108: pc[7] = m8 = p2*x7 + p8*x8 + p14*x9 + p20*x10 + p26*x11 + p32*x12;
109: pc[8] = m9 = p3*x7 + p9*x8 + p15*x9 + p21*x10 + p27*x11 + p33*x12;
110: pc[9] = m10 = p4*x7 + p10*x8 + p16*x9 + p22*x10 + p28*x11 + p34*x12;
111: pc[10] = m11 = p5*x7 + p11*x8 + p17*x9 + p23*x10 + p29*x11 + p35*x12;
112: pc[11] = m12 = p6*x7 + p12*x8 + p18*x9 + p24*x10 + p30*x11 + p36*x12;
114: pc[12] = m13 = p1*x13 + p7*x14 + p13*x15 + p19*x16 + p25*x17 + p31*x18;
115: pc[13] = m14 = p2*x13 + p8*x14 + p14*x15 + p20*x16 + p26*x17 + p32*x18;
116: pc[14] = m15 = p3*x13 + p9*x14 + p15*x15 + p21*x16 + p27*x17 + p33*x18;
117: pc[15] = m16 = p4*x13 + p10*x14 + p16*x15 + p22*x16 + p28*x17 + p34*x18;
118: pc[16] = m17 = p5*x13 + p11*x14 + p17*x15 + p23*x16 + p29*x17 + p35*x18;
119: pc[17] = m18 = p6*x13 + p12*x14 + p18*x15 + p24*x16 + p30*x17 + p36*x18;
121: pc[18] = m19 = p1*x19 + p7*x20 + p13*x21 + p19*x22 + p25*x23 + p31*x24;
122: pc[19] = m20 = p2*x19 + p8*x20 + p14*x21 + p20*x22 + p26*x23 + p32*x24;
123: pc[20] = m21 = p3*x19 + p9*x20 + p15*x21 + p21*x22 + p27*x23 + p33*x24;
124: pc[21] = m22 = p4*x19 + p10*x20 + p16*x21 + p22*x22 + p28*x23 + p34*x24;
125: pc[22] = m23 = p5*x19 + p11*x20 + p17*x21 + p23*x22 + p29*x23 + p35*x24;
126: pc[23] = m24 = p6*x19 + p12*x20 + p18*x21 + p24*x22 + p30*x23 + p36*x24;
128: pc[24] = m25 = p1*x25 + p7*x26 + p13*x27 + p19*x28 + p25*x29 + p31*x30;
129: pc[25] = m26 = p2*x25 + p8*x26 + p14*x27 + p20*x28 + p26*x29 + p32*x30;
130: pc[26] = m27 = p3*x25 + p9*x26 + p15*x27 + p21*x28 + p27*x29 + p33*x30;
131: pc[27] = m28 = p4*x25 + p10*x26 + p16*x27 + p22*x28 + p28*x29 + p34*x30;
132: pc[28] = m29 = p5*x25 + p11*x26 + p17*x27 + p23*x28 + p29*x29 + p35*x30;
133: pc[29] = m30 = p6*x25 + p12*x26 + p18*x27 + p24*x28 + p30*x29 + p36*x30;
135: pc[30] = m31 = p1*x31 + p7*x32 + p13*x33 + p19*x34 + p25*x35 + p31*x36;
136: pc[31] = m32 = p2*x31 + p8*x32 + p14*x33 + p20*x34 + p26*x35 + p32*x36;
137: pc[32] = m33 = p3*x31 + p9*x32 + p15*x33 + p21*x34 + p27*x35 + p33*x36;
138: pc[33] = m34 = p4*x31 + p10*x32 + p16*x33 + p22*x34 + p28*x35 + p34*x36;
139: pc[34] = m35 = p5*x31 + p11*x32 + p17*x33 + p23*x34 + p29*x35 + p35*x36;
140: pc[35] = m36 = p6*x31 + p12*x32 + p18*x33 + p24*x34 + p30*x35 + p36*x36;
142: nz = bi[row+1] - diag_offset[row] - 1;
143: pv += 36;
144: for (j=0; j<nz; j++) {
145: x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
146: x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7];
147: x9 = pv[8]; x10 = pv[9]; x11 = pv[10]; x12 = pv[11];
148: x13 = pv[12]; x14 = pv[13]; x15 = pv[14]; x16 = pv[15];
149: x17 = pv[16]; x18 = pv[17]; x19 = pv[18]; x20 = pv[19];
150: x21 = pv[20]; x22 = pv[21]; x23 = pv[22]; x24 = pv[23];
151: x25 = pv[24]; x26 = pv[25]; x27 = pv[26]; x28 = pv[27];
152: x29 = pv[28]; x30 = pv[29]; x31 = pv[30]; x32 = pv[31];
153: x33 = pv[32]; x34 = pv[33]; x35 = pv[34]; x36 = pv[35];
154: x = rtmp + 36*pj[j];
155: x[0] -= m1*x1 + m7*x2 + m13*x3 + m19*x4 + m25*x5 + m31*x6;
156: x[1] -= m2*x1 + m8*x2 + m14*x3 + m20*x4 + m26*x5 + m32*x6;
157: x[2] -= m3*x1 + m9*x2 + m15*x3 + m21*x4 + m27*x5 + m33*x6;
158: x[3] -= m4*x1 + m10*x2 + m16*x3 + m22*x4 + m28*x5 + m34*x6;
159: x[4] -= m5*x1 + m11*x2 + m17*x3 + m23*x4 + m29*x5 + m35*x6;
160: x[5] -= m6*x1 + m12*x2 + m18*x3 + m24*x4 + m30*x5 + m36*x6;
162: x[6] -= m1*x7 + m7*x8 + m13*x9 + m19*x10 + m25*x11 + m31*x12;
163: x[7] -= m2*x7 + m8*x8 + m14*x9 + m20*x10 + m26*x11 + m32*x12;
164: x[8] -= m3*x7 + m9*x8 + m15*x9 + m21*x10 + m27*x11 + m33*x12;
165: x[9] -= m4*x7 + m10*x8 + m16*x9 + m22*x10 + m28*x11 + m34*x12;
166: x[10] -= m5*x7 + m11*x8 + m17*x9 + m23*x10 + m29*x11 + m35*x12;
167: x[11] -= m6*x7 + m12*x8 + m18*x9 + m24*x10 + m30*x11 + m36*x12;
169: x[12] -= m1*x13 + m7*x14 + m13*x15 + m19*x16 + m25*x17 + m31*x18;
170: x[13] -= m2*x13 + m8*x14 + m14*x15 + m20*x16 + m26*x17 + m32*x18;
171: x[14] -= m3*x13 + m9*x14 + m15*x15 + m21*x16 + m27*x17 + m33*x18;
172: x[15] -= m4*x13 + m10*x14 + m16*x15 + m22*x16 + m28*x17 + m34*x18;
173: x[16] -= m5*x13 + m11*x14 + m17*x15 + m23*x16 + m29*x17 + m35*x18;
174: x[17] -= m6*x13 + m12*x14 + m18*x15 + m24*x16 + m30*x17 + m36*x18;
176: x[18] -= m1*x19 + m7*x20 + m13*x21 + m19*x22 + m25*x23 + m31*x24;
177: x[19] -= m2*x19 + m8*x20 + m14*x21 + m20*x22 + m26*x23 + m32*x24;
178: x[20] -= m3*x19 + m9*x20 + m15*x21 + m21*x22 + m27*x23 + m33*x24;
179: x[21] -= m4*x19 + m10*x20 + m16*x21 + m22*x22 + m28*x23 + m34*x24;
180: x[22] -= m5*x19 + m11*x20 + m17*x21 + m23*x22 + m29*x23 + m35*x24;
181: x[23] -= m6*x19 + m12*x20 + m18*x21 + m24*x22 + m30*x23 + m36*x24;
183: x[24] -= m1*x25 + m7*x26 + m13*x27 + m19*x28 + m25*x29 + m31*x30;
184: x[25] -= m2*x25 + m8*x26 + m14*x27 + m20*x28 + m26*x29 + m32*x30;
185: x[26] -= m3*x25 + m9*x26 + m15*x27 + m21*x28 + m27*x29 + m33*x30;
186: x[27] -= m4*x25 + m10*x26 + m16*x27 + m22*x28 + m28*x29 + m34*x30;
187: x[28] -= m5*x25 + m11*x26 + m17*x27 + m23*x28 + m29*x29 + m35*x30;
188: x[29] -= m6*x25 + m12*x26 + m18*x27 + m24*x28 + m30*x29 + m36*x30;
190: x[30] -= m1*x31 + m7*x32 + m13*x33 + m19*x34 + m25*x35 + m31*x36;
191: x[31] -= m2*x31 + m8*x32 + m14*x33 + m20*x34 + m26*x35 + m32*x36;
192: x[32] -= m3*x31 + m9*x32 + m15*x33 + m21*x34 + m27*x35 + m33*x36;
193: x[33] -= m4*x31 + m10*x32 + m16*x33 + m22*x34 + m28*x35 + m34*x36;
194: x[34] -= m5*x31 + m11*x32 + m17*x33 + m23*x34 + m29*x35 + m35*x36;
195: x[35] -= m6*x31 + m12*x32 + m18*x33 + m24*x34 + m30*x35 + m36*x36;
197: pv += 36;
198: }
199: PetscLogFlops(432*nz+396);
200: }
201: row = *ajtmp++;
202: }
203: /* finished row so stick it into b->a */
204: pv = ba + 36*bi[i];
205: pj = bj + bi[i];
206: nz = bi[i+1] - bi[i];
207: for (j=0; j<nz; j++) {
208: x = rtmp+36*pj[j];
209: pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3];
210: pv[4] = x[4]; pv[5] = x[5]; pv[6] = x[6]; pv[7] = x[7];
211: pv[8] = x[8]; pv[9] = x[9]; pv[10] = x[10]; pv[11] = x[11];
212: pv[12] = x[12]; pv[13] = x[13]; pv[14] = x[14]; pv[15] = x[15];
213: pv[16] = x[16]; pv[17] = x[17]; pv[18] = x[18]; pv[19] = x[19];
214: pv[20] = x[20]; pv[21] = x[21]; pv[22] = x[22]; pv[23] = x[23];
215: pv[24] = x[24]; pv[25] = x[25]; pv[26] = x[26]; pv[27] = x[27];
216: pv[28] = x[28]; pv[29] = x[29]; pv[30] = x[30]; pv[31] = x[31];
217: pv[32] = x[32]; pv[33] = x[33]; pv[34] = x[34]; pv[35] = x[35];
218: pv += 36;
219: }
220: /* invert diagonal block */
221: w = ba + 36*diag_offset[i];
222: Kernel_A_gets_inverse_A_6(w,shift);
223: }
225: PetscFree(rtmp);
226: ISRestoreIndices(isicol,&ic);
227: ISRestoreIndices(isrow,&r);
228: C->ops->solve = MatSolve_SeqBAIJ_6;
229: C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_6;
230: C->assembled = PETSC_TRUE;
231: PetscLogFlops(1.3333*216*b->mbs); /* from inverting diagonal blocks */
232: return(0);
233: }