Actual source code: baijfact9.c
1: /*$Id: baijfact9.c,v 1.4 2001/03/23 23:22:07 balay Exp $*/
2: /*
3: Factorization code for BAIJ format.
4: */
5: #include src/mat/impls/baij/seq/baij.h
6: #include src/vec/vecimpl.h
7: #include src/inline/ilu.h
9: /* ------------------------------------------------------------*/
10: /*
11: Version for when blocks are 5 by 5
12: */
13: int MatLUFactorNumeric_SeqBAIJ_5(Mat A,Mat *B)
14: {
15: Mat C = *B;
16: Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data;
17: IS isrow = b->row,isicol = b->icol;
18: int *r,*ic,ierr,i,j,n = a->mbs,*bi = b->i,*bj = b->j;
19: int *ajtmpold,*ajtmp,nz,row;
20: int *diag_offset = b->diag,idx,*ai=a->i,*aj=a->j,*pj;
21: MatScalar *pv,*v,*rtmp,*pc,*w,*x;
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,*aa = a->a;
30: ISGetIndices(isrow,&r);
31: ISGetIndices(isicol,&ic);
32: PetscMalloc(25*(n+1)*sizeof(MatScalar),&rtmp);
34: for (i=0; i<n; i++) {
35: nz = bi[i+1] - bi[i];
36: ajtmp = bj + bi[i];
37: for (j=0; j<nz; j++) {
38: x = rtmp+25*ajtmp[j];
39: x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = x[9] = 0.0;
40: x[10] = x[11] = x[12] = x[13] = x[14] = x[15] = x[16] = x[17] = 0.0;
41: x[18] = x[19] = x[20] = x[21] = x[22] = x[23] = x[24] = 0.0;
42: }
43: /* load in initial (unfactored row) */
44: idx = r[i];
45: nz = ai[idx+1] - ai[idx];
46: ajtmpold = aj + ai[idx];
47: v = aa + 25*ai[idx];
48: for (j=0; j<nz; j++) {
49: x = rtmp+25*ic[ajtmpold[j]];
50: x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3];
51: x[4] = v[4]; x[5] = v[5]; x[6] = v[6]; x[7] = v[7]; x[8] = v[8];
52: x[9] = v[9]; x[10] = v[10]; x[11] = v[11]; x[12] = v[12]; x[13] = v[13];
53: x[14] = v[14]; x[15] = v[15]; x[16] = v[16]; x[17] = v[17];
54: x[18] = v[18]; x[19] = v[19]; x[20] = v[20]; x[21] = v[21];
55: x[22] = v[22]; x[23] = v[23]; x[24] = v[24];
56: v += 25;
57: }
58: row = *ajtmp++;
59: while (row < i) {
60: pc = rtmp + 25*row;
61: p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3];
62: p5 = pc[4]; p6 = pc[5]; p7 = pc[6]; p8 = pc[7]; p9 = pc[8];
63: p10 = pc[9]; p11 = pc[10]; p12 = pc[11]; p13 = pc[12]; p14 = pc[13];
64: p15 = pc[14]; p16 = pc[15]; p17 = pc[16]; p18 = pc[17]; p19 = pc[18];
65: p20 = pc[19]; p21 = pc[20]; p22 = pc[21]; p23 = pc[22]; p24 = pc[23];
66: p25 = pc[24];
67: if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 ||
68: p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0 || p10 != 0.0 ||
69: p11 != 0.0 || p12 != 0.0 || p13 != 0.0 || p14 != 0.0 || p15 != 0.0
70: || p16 != 0.0 || p17 != 0.0 || p18 != 0.0 || p19 != 0.0 ||
71: p20 != 0.0 || p21 != 0.0 || p22 != 0.0 || p23 != 0.0 ||
72: p24 != 0.0 || p25 != 0.0) {
73: pv = ba + 25*diag_offset[row];
74: pj = bj + diag_offset[row] + 1;
75: x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
76: x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8];
77: x10 = pv[9]; x11 = pv[10]; x12 = pv[11]; x13 = pv[12]; x14 = pv[13];
78: x15 = pv[14]; x16 = pv[15]; x17 = pv[16]; x18 = pv[17];
79: x19 = pv[18]; x20 = pv[19]; x21 = pv[20]; x22 = pv[21];
80: x23 = pv[22]; x24 = pv[23]; x25 = pv[24];
81: pc[0] = m1 = p1*x1 + p6*x2 + p11*x3 + p16*x4 + p21*x5;
82: pc[1] = m2 = p2*x1 + p7*x2 + p12*x3 + p17*x4 + p22*x5;
83: pc[2] = m3 = p3*x1 + p8*x2 + p13*x3 + p18*x4 + p23*x5;
84: pc[3] = m4 = p4*x1 + p9*x2 + p14*x3 + p19*x4 + p24*x5;
85: pc[4] = m5 = p5*x1 + p10*x2 + p15*x3 + p20*x4 + p25*x5;
87: pc[5] = m6 = p1*x6 + p6*x7 + p11*x8 + p16*x9 + p21*x10;
88: pc[6] = m7 = p2*x6 + p7*x7 + p12*x8 + p17*x9 + p22*x10;
89: pc[7] = m8 = p3*x6 + p8*x7 + p13*x8 + p18*x9 + p23*x10;
90: pc[8] = m9 = p4*x6 + p9*x7 + p14*x8 + p19*x9 + p24*x10;
91: pc[9] = m10 = p5*x6 + p10*x7 + p15*x8 + p20*x9 + p25*x10;
93: pc[10] = m11 = p1*x11 + p6*x12 + p11*x13 + p16*x14 + p21*x15;
94: pc[11] = m12 = p2*x11 + p7*x12 + p12*x13 + p17*x14 + p22*x15;
95: pc[12] = m13 = p3*x11 + p8*x12 + p13*x13 + p18*x14 + p23*x15;
96: pc[13] = m14 = p4*x11 + p9*x12 + p14*x13 + p19*x14 + p24*x15;
97: pc[14] = m15 = p5*x11 + p10*x12 + p15*x13 + p20*x14 + p25*x15;
99: pc[15] = m16 = p1*x16 + p6*x17 + p11*x18 + p16*x19 + p21*x20;
100: pc[16] = m17 = p2*x16 + p7*x17 + p12*x18 + p17*x19 + p22*x20;
101: pc[17] = m18 = p3*x16 + p8*x17 + p13*x18 + p18*x19 + p23*x20;
102: pc[18] = m19 = p4*x16 + p9*x17 + p14*x18 + p19*x19 + p24*x20;
103: pc[19] = m20 = p5*x16 + p10*x17 + p15*x18 + p20*x19 + p25*x20;
105: pc[20] = m21 = p1*x21 + p6*x22 + p11*x23 + p16*x24 + p21*x25;
106: pc[21] = m22 = p2*x21 + p7*x22 + p12*x23 + p17*x24 + p22*x25;
107: pc[22] = m23 = p3*x21 + p8*x22 + p13*x23 + p18*x24 + p23*x25;
108: pc[23] = m24 = p4*x21 + p9*x22 + p14*x23 + p19*x24 + p24*x25;
109: pc[24] = m25 = p5*x21 + p10*x22 + p15*x23 + p20*x24 + p25*x25;
111: nz = bi[row+1] - diag_offset[row] - 1;
112: pv += 25;
113: for (j=0; j<nz; j++) {
114: x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
115: x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8];
116: x10 = pv[9]; x11 = pv[10]; x12 = pv[11]; x13 = pv[12];
117: x14 = pv[13]; x15 = pv[14]; x16 = pv[15]; x17 = pv[16];
118: x18 = pv[17]; x19 = pv[18]; x20 = pv[19]; x21 = pv[20];
119: x22 = pv[21]; x23 = pv[22]; x24 = pv[23]; x25 = pv[24];
120: x = rtmp + 25*pj[j];
121: x[0] -= m1*x1 + m6*x2 + m11*x3 + m16*x4 + m21*x5;
122: x[1] -= m2*x1 + m7*x2 + m12*x3 + m17*x4 + m22*x5;
123: x[2] -= m3*x1 + m8*x2 + m13*x3 + m18*x4 + m23*x5;
124: x[3] -= m4*x1 + m9*x2 + m14*x3 + m19*x4 + m24*x5;
125: x[4] -= m5*x1 + m10*x2 + m15*x3 + m20*x4 + m25*x5;
127: x[5] -= m1*x6 + m6*x7 + m11*x8 + m16*x9 + m21*x10;
128: x[6] -= m2*x6 + m7*x7 + m12*x8 + m17*x9 + m22*x10;
129: x[7] -= m3*x6 + m8*x7 + m13*x8 + m18*x9 + m23*x10;
130: x[8] -= m4*x6 + m9*x7 + m14*x8 + m19*x9 + m24*x10;
131: x[9] -= m5*x6 + m10*x7 + m15*x8 + m20*x9 + m25*x10;
133: x[10] -= m1*x11 + m6*x12 + m11*x13 + m16*x14 + m21*x15;
134: x[11] -= m2*x11 + m7*x12 + m12*x13 + m17*x14 + m22*x15;
135: x[12] -= m3*x11 + m8*x12 + m13*x13 + m18*x14 + m23*x15;
136: x[13] -= m4*x11 + m9*x12 + m14*x13 + m19*x14 + m24*x15;
137: x[14] -= m5*x11 + m10*x12 + m15*x13 + m20*x14 + m25*x15;
139: x[15] -= m1*x16 + m6*x17 + m11*x18 + m16*x19 + m21*x20;
140: x[16] -= m2*x16 + m7*x17 + m12*x18 + m17*x19 + m22*x20;
141: x[17] -= m3*x16 + m8*x17 + m13*x18 + m18*x19 + m23*x20;
142: x[18] -= m4*x16 + m9*x17 + m14*x18 + m19*x19 + m24*x20;
143: x[19] -= m5*x16 + m10*x17 + m15*x18 + m20*x19 + m25*x20;
145: x[20] -= m1*x21 + m6*x22 + m11*x23 + m16*x24 + m21*x25;
146: x[21] -= m2*x21 + m7*x22 + m12*x23 + m17*x24 + m22*x25;
147: x[22] -= m3*x21 + m8*x22 + m13*x23 + m18*x24 + m23*x25;
148: x[23] -= m4*x21 + m9*x22 + m14*x23 + m19*x24 + m24*x25;
149: x[24] -= m5*x21 + m10*x22 + m15*x23 + m20*x24 + m25*x25;
151: pv += 25;
152: }
153: PetscLogFlops(250*nz+225);
154: }
155: row = *ajtmp++;
156: }
157: /* finished row so stick it into b->a */
158: pv = ba + 25*bi[i];
159: pj = bj + bi[i];
160: nz = bi[i+1] - bi[i];
161: for (j=0; j<nz; j++) {
162: x = rtmp+25*pj[j];
163: pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3];
164: pv[4] = x[4]; pv[5] = x[5]; pv[6] = x[6]; pv[7] = x[7]; pv[8] = x[8];
165: pv[9] = x[9]; pv[10] = x[10]; pv[11] = x[11]; pv[12] = x[12];
166: pv[13] = x[13]; pv[14] = x[14]; pv[15] = x[15]; pv[16] = x[16];
167: pv[17] = x[17]; pv[18] = x[18]; pv[19] = x[19]; pv[20] = x[20];
168: pv[21] = x[21]; pv[22] = x[22]; pv[23] = x[23]; pv[24] = x[24];
169: pv += 25;
170: }
171: /* invert diagonal block */
172: w = ba + 25*diag_offset[i];
173: Kernel_A_gets_inverse_A_5(w);
174: }
176: PetscFree(rtmp);
177: ISRestoreIndices(isicol,&ic);
178: ISRestoreIndices(isrow,&r);
179: C->factor = FACTOR_LU;
180: C->assembled = PETSC_TRUE;
181: PetscLogFlops(1.3333*125*b->mbs); /* from inverting diagonal blocks */
182: return(0);
183: }