Actual source code: baijfact5.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>
7: /*
8: Version for when blocks are 7 by 7
9: */
12: PetscErrorCode MatLUFactorNumeric_SeqBAIJ_7_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;
17: const PetscInt *r,*ic,*bi = b->i,*bj = b->j,*ajtmp,*diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj,*ajtmpold;
18: PetscInt i,j,n = a->mbs,nz,row,idx;
19: MatScalar *pv,*v,*rtmp,*pc,*w,*x;
20: MatScalar p1,p2,p3,p4,m1,m2,m3,m4,m5,m6,m7,m8,m9,x1,x2,x3,x4;
21: MatScalar p5,p6,p7,p8,p9,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15,x16;
22: MatScalar x17,x18,x19,x20,x21,x22,x23,x24,x25,p10,p11,p12,p13,p14;
23: MatScalar p15,p16,p17,p18,p19,p20,p21,p22,p23,p24,p25,m10,m11,m12;
24: MatScalar m13,m14,m15,m16,m17,m18,m19,m20,m21,m22,m23,m24,m25;
25: MatScalar p26,p27,p28,p29,p30,p31,p32,p33,p34,p35,p36;
26: MatScalar p37,p38,p39,p40,p41,p42,p43,p44,p45,p46,p47,p48,p49;
27: MatScalar x26,x27,x28,x29,x30,x31,x32,x33,x34,x35,x36;
28: MatScalar x37,x38,x39,x40,x41,x42,x43,x44,x45,x46,x47,x48,x49;
29: MatScalar m26,m27,m28,m29,m30,m31,m32,m33,m34,m35,m36;
30: MatScalar m37,m38,m39,m40,m41,m42,m43,m44,m45,m46,m47,m48,m49;
31: MatScalar *ba = b->a,*aa = a->a;
32: PetscReal shift = info->shiftamount;
33: PetscBool allowzeropivot,zeropivotdetected;
36: allowzeropivot = PetscNot(A->erroriffailure);
37: ISGetIndices(isrow,&r);
38: ISGetIndices(isicol,&ic);
39: PetscMalloc1(49*(n+1),&rtmp);
41: for (i=0; i<n; i++) {
42: nz = bi[i+1] - bi[i];
43: ajtmp = bj + bi[i];
44: for (j=0; j<nz; j++) {
45: x = rtmp+49*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] = x[25] = 0.0;
49: x[26] = x[27] = x[28] = x[29] = x[30] = x[31] = x[32] = x[33] = 0.0;
50: x[34] = x[35] = x[36] = x[37] = x[38] = x[39] = x[40] = x[41] = 0.0;
51: x[42] = x[43] = x[44] = x[45] = x[46] = x[47] = x[48] = 0.0;
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 + 49*ai[idx];
58: for (j=0; j<nz; j++) {
59: x = rtmp+49*ic[ajtmpold[j]];
60: x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3];
61: x[4] = v[4]; x[5] = v[5]; x[6] = v[6]; x[7] = v[7];
62: x[8] = v[8]; x[9] = v[9]; x[10] = v[10]; x[11] = v[11];
63: x[12] = v[12]; x[13] = v[13]; x[14] = v[14]; x[15] = v[15];
64: x[16] = v[16]; x[17] = v[17]; x[18] = v[18]; x[19] = v[19];
65: x[20] = v[20]; x[21] = v[21]; x[22] = v[22]; x[23] = v[23];
66: x[24] = v[24]; x[25] = v[25]; x[26] = v[26]; x[27] = v[27];
67: x[28] = v[28]; x[29] = v[29]; x[30] = v[30]; x[31] = v[31];
68: x[32] = v[32]; x[33] = v[33]; x[34] = v[34]; x[35] = v[35];
69: x[36] = v[36]; x[37] = v[37]; x[38] = v[38]; x[39] = v[39];
70: x[40] = v[40]; x[41] = v[41]; x[42] = v[42]; x[43] = v[43];
71: x[44] = v[44]; x[45] = v[45]; x[46] = v[46]; x[47] = v[47];
72: x[48] = v[48];
73: v += 49;
74: }
75: row = *ajtmp++;
76: while (row < i) {
77: pc = rtmp + 49*row;
78: p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3];
79: p5 = pc[4]; p6 = pc[5]; p7 = pc[6]; p8 = pc[7];
80: p9 = pc[8]; p10 = pc[9]; p11 = pc[10]; p12 = pc[11];
81: p13 = pc[12]; p14 = pc[13]; p15 = pc[14]; p16 = pc[15];
82: p17 = pc[16]; p18 = pc[17]; p19 = pc[18]; p20 = pc[19];
83: p21 = pc[20]; p22 = pc[21]; p23 = pc[22]; p24 = pc[23];
84: p25 = pc[24]; p26 = pc[25]; p27 = pc[26]; p28 = pc[27];
85: p29 = pc[28]; p30 = pc[29]; p31 = pc[30]; p32 = pc[31];
86: p33 = pc[32]; p34 = pc[33]; p35 = pc[34]; p36 = pc[35];
87: p37 = pc[36]; p38 = pc[37]; p39 = pc[38]; p40 = pc[39];
88: p41 = pc[40]; p42 = pc[41]; p43 = pc[42]; p44 = pc[43];
89: p45 = pc[44]; p46 = pc[45]; p47 = pc[46]; p48 = pc[47];
90: p49 = pc[48];
91: if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 ||
92: p5 != 0.0 || p6 != 0.0 || p7 != 0.0 || p8 != 0.0 ||
93: p9 != 0.0 || p10 != 0.0 || p11 != 0.0 || p12 != 0.0 ||
94: p13 != 0.0 || p14 != 0.0 || p15 != 0.0 || p16 != 0.0 ||
95: p17 != 0.0 || p18 != 0.0 || p19 != 0.0 || p20 != 0.0 ||
96: p21 != 0.0 || p22 != 0.0 || p23 != 0.0 || p24 != 0.0 ||
97: p25 != 0.0 || p26 != 0.0 || p27 != 0.0 || p28 != 0.0 ||
98: p29 != 0.0 || p30 != 0.0 || p31 != 0.0 || p32 != 0.0 ||
99: p33 != 0.0 || p34 != 0.0 || p35 != 0.0 || p36 != 0.0 ||
100: p37 != 0.0 || p38 != 0.0 || p39 != 0.0 || p40 != 0.0 ||
101: p41 != 0.0 || p42 != 0.0 || p43 != 0.0 || p44 != 0.0 ||
102: p45 != 0.0 || p46 != 0.0 || p47 != 0.0 || p48 != 0.0 ||
103: p49 != 0.0) {
104: pv = ba + 49*diag_offset[row];
105: pj = bj + diag_offset[row] + 1;
106: x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
107: x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7];
108: x9 = pv[8]; x10 = pv[9]; x11 = pv[10]; x12 = pv[11];
109: x13 = pv[12]; x14 = pv[13]; x15 = pv[14]; x16 = pv[15];
110: x17 = pv[16]; x18 = pv[17]; x19 = pv[18]; x20 = pv[19];
111: x21 = pv[20]; x22 = pv[21]; x23 = pv[22]; x24 = pv[23];
112: x25 = pv[24]; x26 = pv[25]; x27 = pv[26]; x28 = pv[27];
113: x29 = pv[28]; x30 = pv[29]; x31 = pv[30]; x32 = pv[31];
114: x33 = pv[32]; x34 = pv[33]; x35 = pv[34]; x36 = pv[35];
115: x37 = pv[36]; x38 = pv[37]; x39 = pv[38]; x40 = pv[39];
116: x41 = pv[40]; x42 = pv[41]; x43 = pv[42]; x44 = pv[43];
117: x45 = pv[44]; x46 = pv[45]; x47 = pv[46]; x48 = pv[47];
118: x49 = pv[48];
119: pc[0] = m1 = p1*x1 + p8*x2 + p15*x3 + p22*x4 + p29*x5 + p36*x6 + p43*x7;
120: pc[1] = m2 = p2*x1 + p9*x2 + p16*x3 + p23*x4 + p30*x5 + p37*x6 + p44*x7;
121: pc[2] = m3 = p3*x1 + p10*x2 + p17*x3 + p24*x4 + p31*x5 + p38*x6 + p45*x7;
122: pc[3] = m4 = p4*x1 + p11*x2 + p18*x3 + p25*x4 + p32*x5 + p39*x6 + p46*x7;
123: pc[4] = m5 = p5*x1 + p12*x2 + p19*x3 + p26*x4 + p33*x5 + p40*x6 + p47*x7;
124: pc[5] = m6 = p6*x1 + p13*x2 + p20*x3 + p27*x4 + p34*x5 + p41*x6 + p48*x7;
125: pc[6] = m7 = p7*x1 + p14*x2 + p21*x3 + p28*x4 + p35*x5 + p42*x6 + p49*x7;
127: pc[7] = m8 = p1*x8 + p8*x9 + p15*x10 + p22*x11 + p29*x12 + p36*x13 + p43*x14;
128: pc[8] = m9 = p2*x8 + p9*x9 + p16*x10 + p23*x11 + p30*x12 + p37*x13 + p44*x14;
129: pc[9] = m10 = p3*x8 + p10*x9 + p17*x10 + p24*x11 + p31*x12 + p38*x13 + p45*x14;
130: pc[10] = m11 = p4*x8 + p11*x9 + p18*x10 + p25*x11 + p32*x12 + p39*x13 + p46*x14;
131: pc[11] = m12 = p5*x8 + p12*x9 + p19*x10 + p26*x11 + p33*x12 + p40*x13 + p47*x14;
132: pc[12] = m13 = p6*x8 + p13*x9 + p20*x10 + p27*x11 + p34*x12 + p41*x13 + p48*x14;
133: pc[13] = m14 = p7*x8 + p14*x9 + p21*x10 + p28*x11 + p35*x12 + p42*x13 + p49*x14;
135: pc[14] = m15 = p1*x15 + p8*x16 + p15*x17 + p22*x18 + p29*x19 + p36*x20 + p43*x21;
136: pc[15] = m16 = p2*x15 + p9*x16 + p16*x17 + p23*x18 + p30*x19 + p37*x20 + p44*x21;
137: pc[16] = m17 = p3*x15 + p10*x16 + p17*x17 + p24*x18 + p31*x19 + p38*x20 + p45*x21;
138: pc[17] = m18 = p4*x15 + p11*x16 + p18*x17 + p25*x18 + p32*x19 + p39*x20 + p46*x21;
139: pc[18] = m19 = p5*x15 + p12*x16 + p19*x17 + p26*x18 + p33*x19 + p40*x20 + p47*x21;
140: pc[19] = m20 = p6*x15 + p13*x16 + p20*x17 + p27*x18 + p34*x19 + p41*x20 + p48*x21;
141: pc[20] = m21 = p7*x15 + p14*x16 + p21*x17 + p28*x18 + p35*x19 + p42*x20 + p49*x21;
143: pc[21] = m22 = p1*x22 + p8*x23 + p15*x24 + p22*x25 + p29*x26 + p36*x27 + p43*x28;
144: pc[22] = m23 = p2*x22 + p9*x23 + p16*x24 + p23*x25 + p30*x26 + p37*x27 + p44*x28;
145: pc[23] = m24 = p3*x22 + p10*x23 + p17*x24 + p24*x25 + p31*x26 + p38*x27 + p45*x28;
146: pc[24] = m25 = p4*x22 + p11*x23 + p18*x24 + p25*x25 + p32*x26 + p39*x27 + p46*x28;
147: pc[25] = m26 = p5*x22 + p12*x23 + p19*x24 + p26*x25 + p33*x26 + p40*x27 + p47*x28;
148: pc[26] = m27 = p6*x22 + p13*x23 + p20*x24 + p27*x25 + p34*x26 + p41*x27 + p48*x28;
149: pc[27] = m28 = p7*x22 + p14*x23 + p21*x24 + p28*x25 + p35*x26 + p42*x27 + p49*x28;
151: pc[28] = m29 = p1*x29 + p8*x30 + p15*x31 + p22*x32 + p29*x33 + p36*x34 + p43*x35;
152: pc[29] = m30 = p2*x29 + p9*x30 + p16*x31 + p23*x32 + p30*x33 + p37*x34 + p44*x35;
153: pc[30] = m31 = p3*x29 + p10*x30 + p17*x31 + p24*x32 + p31*x33 + p38*x34 + p45*x35;
154: pc[31] = m32 = p4*x29 + p11*x30 + p18*x31 + p25*x32 + p32*x33 + p39*x34 + p46*x35;
155: pc[32] = m33 = p5*x29 + p12*x30 + p19*x31 + p26*x32 + p33*x33 + p40*x34 + p47*x35;
156: pc[33] = m34 = p6*x29 + p13*x30 + p20*x31 + p27*x32 + p34*x33 + p41*x34 + p48*x35;
157: pc[34] = m35 = p7*x29 + p14*x30 + p21*x31 + p28*x32 + p35*x33 + p42*x34 + p49*x35;
159: pc[35] = m36 = p1*x36 + p8*x37 + p15*x38 + p22*x39 + p29*x40 + p36*x41 + p43*x42;
160: pc[36] = m37 = p2*x36 + p9*x37 + p16*x38 + p23*x39 + p30*x40 + p37*x41 + p44*x42;
161: pc[37] = m38 = p3*x36 + p10*x37 + p17*x38 + p24*x39 + p31*x40 + p38*x41 + p45*x42;
162: pc[38] = m39 = p4*x36 + p11*x37 + p18*x38 + p25*x39 + p32*x40 + p39*x41 + p46*x42;
163: pc[39] = m40 = p5*x36 + p12*x37 + p19*x38 + p26*x39 + p33*x40 + p40*x41 + p47*x42;
164: pc[40] = m41 = p6*x36 + p13*x37 + p20*x38 + p27*x39 + p34*x40 + p41*x41 + p48*x42;
165: pc[41] = m42 = p7*x36 + p14*x37 + p21*x38 + p28*x39 + p35*x40 + p42*x41 + p49*x42;
167: pc[42] = m43 = p1*x43 + p8*x44 + p15*x45 + p22*x46 + p29*x47 + p36*x48 + p43*x49;
168: pc[43] = m44 = p2*x43 + p9*x44 + p16*x45 + p23*x46 + p30*x47 + p37*x48 + p44*x49;
169: pc[44] = m45 = p3*x43 + p10*x44 + p17*x45 + p24*x46 + p31*x47 + p38*x48 + p45*x49;
170: pc[45] = m46 = p4*x43 + p11*x44 + p18*x45 + p25*x46 + p32*x47 + p39*x48 + p46*x49;
171: pc[46] = m47 = p5*x43 + p12*x44 + p19*x45 + p26*x46 + p33*x47 + p40*x48 + p47*x49;
172: pc[47] = m48 = p6*x43 + p13*x44 + p20*x45 + p27*x46 + p34*x47 + p41*x48 + p48*x49;
173: pc[48] = m49 = p7*x43 + p14*x44 + p21*x45 + p28*x46 + p35*x47 + p42*x48 + p49*x49;
175: nz = bi[row+1] - diag_offset[row] - 1;
176: pv += 49;
177: for (j=0; j<nz; j++) {
178: x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
179: x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7];
180: x9 = pv[8]; x10 = pv[9]; x11 = pv[10]; x12 = pv[11];
181: x13 = pv[12]; x14 = pv[13]; x15 = pv[14]; x16 = pv[15];
182: x17 = pv[16]; x18 = pv[17]; x19 = pv[18]; x20 = pv[19];
183: x21 = pv[20]; x22 = pv[21]; x23 = pv[22]; x24 = pv[23];
184: x25 = pv[24]; x26 = pv[25]; x27 = pv[26]; x28 = pv[27];
185: x29 = pv[28]; x30 = pv[29]; x31 = pv[30]; x32 = pv[31];
186: x33 = pv[32]; x34 = pv[33]; x35 = pv[34]; x36 = pv[35];
187: x37 = pv[36]; x38 = pv[37]; x39 = pv[38]; x40 = pv[39];
188: x41 = pv[40]; x42 = pv[41]; x43 = pv[42]; x44 = pv[43];
189: x45 = pv[44]; x46 = pv[45]; x47 = pv[46]; x48 = pv[47];
190: x49 = pv[48];
191: x = rtmp + 49*pj[j];
192: x[0] -= m1*x1 + m8*x2 + m15*x3 + m22*x4 + m29*x5 + m36*x6 + m43*x7;
193: x[1] -= m2*x1 + m9*x2 + m16*x3 + m23*x4 + m30*x5 + m37*x6 + m44*x7;
194: x[2] -= m3*x1 + m10*x2 + m17*x3 + m24*x4 + m31*x5 + m38*x6 + m45*x7;
195: x[3] -= m4*x1 + m11*x2 + m18*x3 + m25*x4 + m32*x5 + m39*x6 + m46*x7;
196: x[4] -= m5*x1 + m12*x2 + m19*x3 + m26*x4 + m33*x5 + m40*x6 + m47*x7;
197: x[5] -= m6*x1 + m13*x2 + m20*x3 + m27*x4 + m34*x5 + m41*x6 + m48*x7;
198: x[6] -= m7*x1 + m14*x2 + m21*x3 + m28*x4 + m35*x5 + m42*x6 + m49*x7;
200: x[7] -= m1*x8 + m8*x9 + m15*x10 + m22*x11 + m29*x12 + m36*x13 + m43*x14;
201: x[8] -= m2*x8 + m9*x9 + m16*x10 + m23*x11 + m30*x12 + m37*x13 + m44*x14;
202: x[9] -= m3*x8 + m10*x9 + m17*x10 + m24*x11 + m31*x12 + m38*x13 + m45*x14;
203: x[10] -= m4*x8 + m11*x9 + m18*x10 + m25*x11 + m32*x12 + m39*x13 + m46*x14;
204: x[11] -= m5*x8 + m12*x9 + m19*x10 + m26*x11 + m33*x12 + m40*x13 + m47*x14;
205: x[12] -= m6*x8 + m13*x9 + m20*x10 + m27*x11 + m34*x12 + m41*x13 + m48*x14;
206: x[13] -= m7*x8 + m14*x9 + m21*x10 + m28*x11 + m35*x12 + m42*x13 + m49*x14;
208: x[14] -= m1*x15 + m8*x16 + m15*x17 + m22*x18 + m29*x19 + m36*x20 + m43*x21;
209: x[15] -= m2*x15 + m9*x16 + m16*x17 + m23*x18 + m30*x19 + m37*x20 + m44*x21;
210: x[16] -= m3*x15 + m10*x16 + m17*x17 + m24*x18 + m31*x19 + m38*x20 + m45*x21;
211: x[17] -= m4*x15 + m11*x16 + m18*x17 + m25*x18 + m32*x19 + m39*x20 + m46*x21;
212: x[18] -= m5*x15 + m12*x16 + m19*x17 + m26*x18 + m33*x19 + m40*x20 + m47*x21;
213: x[19] -= m6*x15 + m13*x16 + m20*x17 + m27*x18 + m34*x19 + m41*x20 + m48*x21;
214: x[20] -= m7*x15 + m14*x16 + m21*x17 + m28*x18 + m35*x19 + m42*x20 + m49*x21;
216: x[21] -= m1*x22 + m8*x23 + m15*x24 + m22*x25 + m29*x26 + m36*x27 + m43*x28;
217: x[22] -= m2*x22 + m9*x23 + m16*x24 + m23*x25 + m30*x26 + m37*x27 + m44*x28;
218: x[23] -= m3*x22 + m10*x23 + m17*x24 + m24*x25 + m31*x26 + m38*x27 + m45*x28;
219: x[24] -= m4*x22 + m11*x23 + m18*x24 + m25*x25 + m32*x26 + m39*x27 + m46*x28;
220: x[25] -= m5*x22 + m12*x23 + m19*x24 + m26*x25 + m33*x26 + m40*x27 + m47*x28;
221: x[26] -= m6*x22 + m13*x23 + m20*x24 + m27*x25 + m34*x26 + m41*x27 + m48*x28;
222: x[27] -= m7*x22 + m14*x23 + m21*x24 + m28*x25 + m35*x26 + m42*x27 + m49*x28;
224: x[28] -= m1*x29 + m8*x30 + m15*x31 + m22*x32 + m29*x33 + m36*x34 + m43*x35;
225: x[29] -= m2*x29 + m9*x30 + m16*x31 + m23*x32 + m30*x33 + m37*x34 + m44*x35;
226: x[30] -= m3*x29 + m10*x30 + m17*x31 + m24*x32 + m31*x33 + m38*x34 + m45*x35;
227: x[31] -= m4*x29 + m11*x30 + m18*x31 + m25*x32 + m32*x33 + m39*x34 + m46*x35;
228: x[32] -= m5*x29 + m12*x30 + m19*x31 + m26*x32 + m33*x33 + m40*x34 + m47*x35;
229: x[33] -= m6*x29 + m13*x30 + m20*x31 + m27*x32 + m34*x33 + m41*x34 + m48*x35;
230: x[34] -= m7*x29 + m14*x30 + m21*x31 + m28*x32 + m35*x33 + m42*x34 + m49*x35;
232: x[35] -= m1*x36 + m8*x37 + m15*x38 + m22*x39 + m29*x40 + m36*x41 + m43*x42;
233: x[36] -= m2*x36 + m9*x37 + m16*x38 + m23*x39 + m30*x40 + m37*x41 + m44*x42;
234: x[37] -= m3*x36 + m10*x37 + m17*x38 + m24*x39 + m31*x40 + m38*x41 + m45*x42;
235: x[38] -= m4*x36 + m11*x37 + m18*x38 + m25*x39 + m32*x40 + m39*x41 + m46*x42;
236: x[39] -= m5*x36 + m12*x37 + m19*x38 + m26*x39 + m33*x40 + m40*x41 + m47*x42;
237: x[40] -= m6*x36 + m13*x37 + m20*x38 + m27*x39 + m34*x40 + m41*x41 + m48*x42;
238: x[41] -= m7*x36 + m14*x37 + m21*x38 + m28*x39 + m35*x40 + m42*x41 + m49*x42;
240: x[42] -= m1*x43 + m8*x44 + m15*x45 + m22*x46 + m29*x47 + m36*x48 + m43*x49;
241: x[43] -= m2*x43 + m9*x44 + m16*x45 + m23*x46 + m30*x47 + m37*x48 + m44*x49;
242: x[44] -= m3*x43 + m10*x44 + m17*x45 + m24*x46 + m31*x47 + m38*x48 + m45*x49;
243: x[45] -= m4*x43 + m11*x44 + m18*x45 + m25*x46 + m32*x47 + m39*x48 + m46*x49;
244: x[46] -= m5*x43 + m12*x44 + m19*x45 + m26*x46 + m33*x47 + m40*x48 + m47*x49;
245: x[47] -= m6*x43 + m13*x44 + m20*x45 + m27*x46 + m34*x47 + m41*x48 + m48*x49;
246: x[48] -= m7*x43 + m14*x44 + m21*x45 + m28*x46 + m35*x47 + m42*x48 + m49*x49;
247: pv += 49;
248: }
249: PetscLogFlops(686.0*nz+637.0);
250: }
251: row = *ajtmp++;
252: }
253: /* finished row so stick it into b->a */
254: pv = ba + 49*bi[i];
255: pj = bj + bi[i];
256: nz = bi[i+1] - bi[i];
257: for (j=0; j<nz; j++) {
258: x = rtmp+49*pj[j];
259: pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3];
260: pv[4] = x[4]; pv[5] = x[5]; pv[6] = x[6]; pv[7] = x[7];
261: pv[8] = x[8]; pv[9] = x[9]; pv[10] = x[10]; pv[11] = x[11];
262: pv[12] = x[12]; pv[13] = x[13]; pv[14] = x[14]; pv[15] = x[15];
263: pv[16] = x[16]; pv[17] = x[17]; pv[18] = x[18]; pv[19] = x[19];
264: pv[20] = x[20]; pv[21] = x[21]; pv[22] = x[22]; pv[23] = x[23];
265: pv[24] = x[24]; pv[25] = x[25]; pv[26] = x[26]; pv[27] = x[27];
266: pv[28] = x[28]; pv[29] = x[29]; pv[30] = x[30]; pv[31] = x[31];
267: pv[32] = x[32]; pv[33] = x[33]; pv[34] = x[34]; pv[35] = x[35];
268: pv[36] = x[36]; pv[37] = x[37]; pv[38] = x[38]; pv[39] = x[39];
269: pv[40] = x[40]; pv[41] = x[41]; pv[42] = x[42]; pv[43] = x[43];
270: pv[44] = x[44]; pv[45] = x[45]; pv[46] = x[46]; pv[47] = x[47];
271: pv[48] = x[48];
272: pv += 49;
273: }
274: /* invert diagonal block */
275: w = ba + 49*diag_offset[i];
276: PetscKernel_A_gets_inverse_A_7(w,shift,allowzeropivot,&zeropivotdetected);
277: if (zeropivotdetected) C->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
278: }
280: PetscFree(rtmp);
281: ISRestoreIndices(isicol,&ic);
282: ISRestoreIndices(isrow,&r);
284: C->ops->solve = MatSolve_SeqBAIJ_7_inplace;
285: C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_7_inplace;
286: C->assembled = PETSC_TRUE;
288: PetscLogFlops(1.333333333333*7*7*7*b->mbs); /* from inverting diagonal blocks */
289: return(0);
290: }
295: PetscErrorCode MatLUFactorNumeric_SeqBAIJ_7(Mat B,Mat A,const MatFactorInfo *info)
296: {
297: Mat C =B;
298: Mat_SeqBAIJ *a =(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data;
299: IS isrow = b->row,isicol = b->icol;
301: const PetscInt *r,*ic;
302: PetscInt i,j,k,nz,nzL,row;
303: const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
304: const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag,*pj,bs2=a->bs2;
305: MatScalar *rtmp,*pc,*mwork,*v,*pv,*aa=a->a;
306: PetscInt flg;
307: PetscReal shift = info->shiftamount;
308: PetscBool allowzeropivot,zeropivotdetected;
311: allowzeropivot = PetscNot(A->erroriffailure);
312: ISGetIndices(isrow,&r);
313: ISGetIndices(isicol,&ic);
315: /* generate work space needed by the factorization */
316: PetscMalloc2(bs2*n,&rtmp,bs2,&mwork);
317: PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));
319: for (i=0; i<n; i++) {
320: /* zero rtmp */
321: /* L part */
322: nz = bi[i+1] - bi[i];
323: bjtmp = bj + bi[i];
324: for (j=0; j<nz; j++) {
325: PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));
326: }
328: /* U part */
329: nz = bdiag[i] - bdiag[i+1];
330: bjtmp = bj + bdiag[i+1]+1;
331: for (j=0; j<nz; j++) {
332: PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));
333: }
335: /* load in initial (unfactored row) */
336: nz = ai[r[i]+1] - ai[r[i]];
337: ajtmp = aj + ai[r[i]];
338: v = aa + bs2*ai[r[i]];
339: for (j=0; j<nz; j++) {
340: PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],v+bs2*j,bs2*sizeof(MatScalar));
341: }
343: /* elimination */
344: bjtmp = bj + bi[i];
345: nzL = bi[i+1] - bi[i];
346: for (k=0; k < nzL; k++) {
347: row = bjtmp[k];
348: pc = rtmp + bs2*row;
349: for (flg=0,j=0; j<bs2; j++) {
350: if (pc[j]!=0.0) {
351: flg = 1;
352: break;
353: }
354: }
355: if (flg) {
356: pv = b->a + bs2*bdiag[row];
357: /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
358: PetscKernel_A_gets_A_times_B_7(pc,pv,mwork);
360: pj = b->j + bdiag[row+1]+1; /* begining of U(row,:) */
361: pv = b->a + bs2*(bdiag[row+1]+1);
362: nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries inU(row,:), excluding diag */
363: for (j=0; j<nz; j++) {
364: /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
365: /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
366: v = rtmp + bs2*pj[j];
367: PetscKernel_A_gets_A_minus_B_times_C_7(v,pc,pv);
368: pv += bs2;
369: }
370: PetscLogFlops(686*nz+637); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
371: }
372: }
374: /* finished row so stick it into b->a */
375: /* L part */
376: pv = b->a + bs2*bi[i];
377: pj = b->j + bi[i];
378: nz = bi[i+1] - bi[i];
379: for (j=0; j<nz; j++) {
380: PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));
381: }
383: /* Mark diagonal and invert diagonal for simplier triangular solves */
384: pv = b->a + bs2*bdiag[i];
385: pj = b->j + bdiag[i];
386: PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));
387: PetscKernel_A_gets_inverse_A_7(pv,shift,allowzeropivot,&zeropivotdetected);
388: if (zeropivotdetected) C->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
390: /* U part */
391: pv = b->a + bs2*(bdiag[i+1]+1);
392: pj = b->j + bdiag[i+1]+1;
393: nz = bdiag[i] - bdiag[i+1] - 1;
394: for (j=0; j<nz; j++) {
395: PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));
396: }
397: }
399: PetscFree2(rtmp,mwork);
400: ISRestoreIndices(isicol,&ic);
401: ISRestoreIndices(isrow,&r);
403: C->ops->solve = MatSolve_SeqBAIJ_7;
404: C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_7;
405: C->assembled = PETSC_TRUE;
407: PetscLogFlops(1.333333333333*7*7*7*n); /* from inverting diagonal blocks */
408: return(0);
409: }
413: PetscErrorCode MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering_inplace(Mat C,Mat A,const MatFactorInfo *info)
414: {
415: Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)C->data;
417: PetscInt i,j,n = a->mbs,*bi = b->i,*bj = b->j;
418: PetscInt *ajtmpold,*ajtmp,nz,row;
419: PetscInt *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj;
420: MatScalar *pv,*v,*rtmp,*pc,*w,*x;
421: MatScalar x1,x2,x3,x4,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15;
422: MatScalar x16,x17,x18,x19,x20,x21,x22,x23,x24,x25;
423: MatScalar p1,p2,p3,p4,p5,p6,p7,p8,p9,p10,p11,p12,p13,p14,p15;
424: MatScalar p16,p17,p18,p19,p20,p21,p22,p23,p24,p25;
425: MatScalar m1,m2,m3,m4,m5,m6,m7,m8,m9,m10,m11,m12,m13,m14,m15;
426: MatScalar m16,m17,m18,m19,m20,m21,m22,m23,m24,m25;
427: MatScalar p26,p27,p28,p29,p30,p31,p32,p33,p34,p35,p36;
428: MatScalar p37,p38,p39,p40,p41,p42,p43,p44,p45,p46,p47,p48,p49;
429: MatScalar x26,x27,x28,x29,x30,x31,x32,x33,x34,x35,x36;
430: MatScalar x37,x38,x39,x40,x41,x42,x43,x44,x45,x46,x47,x48,x49;
431: MatScalar m26,m27,m28,m29,m30,m31,m32,m33,m34,m35,m36;
432: MatScalar m37,m38,m39,m40,m41,m42,m43,m44,m45,m46,m47,m48,m49;
433: MatScalar *ba = b->a,*aa = a->a;
434: PetscReal shift = info->shiftamount;
435: PetscBool allowzeropivot,zeropivotdetected;
438: allowzeropivot = PetscNot(A->erroriffailure);
439: PetscMalloc1(49*(n+1),&rtmp);
440: for (i=0; i<n; i++) {
441: nz = bi[i+1] - bi[i];
442: ajtmp = bj + bi[i];
443: for (j=0; j<nz; j++) {
444: x = rtmp+49*ajtmp[j];
445: x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = x[9] = 0.0;
446: x[10] = x[11] = x[12] = x[13] = x[14] = x[15] = x[16] = x[17] = 0.0;
447: x[18] = x[19] = x[20] = x[21] = x[22] = x[23] = x[24] = x[25] = 0.0;
448: x[26] = x[27] = x[28] = x[29] = x[30] = x[31] = x[32] = x[33] = 0.0;
449: x[34] = x[35] = x[36] = x[37] = x[38] = x[39] = x[40] = x[41] = 0.0;
450: x[42] = x[43] = x[44] = x[45] = x[46] = x[47] = x[48] = 0.0;
451: }
452: /* load in initial (unfactored row) */
453: nz = ai[i+1] - ai[i];
454: ajtmpold = aj + ai[i];
455: v = aa + 49*ai[i];
456: for (j=0; j<nz; j++) {
457: x = rtmp+49*ajtmpold[j];
458: x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3];
459: x[4] = v[4]; x[5] = v[5]; x[6] = v[6]; x[7] = v[7];
460: x[8] = v[8]; x[9] = v[9]; x[10] = v[10]; x[11] = v[11];
461: x[12] = v[12]; x[13] = v[13]; x[14] = v[14]; x[15] = v[15];
462: x[16] = v[16]; x[17] = v[17]; x[18] = v[18]; x[19] = v[19];
463: x[20] = v[20]; x[21] = v[21]; x[22] = v[22]; x[23] = v[23];
464: x[24] = v[24]; x[25] = v[25]; x[26] = v[26]; x[27] = v[27];
465: x[28] = v[28]; x[29] = v[29]; x[30] = v[30]; x[31] = v[31];
466: x[32] = v[32]; x[33] = v[33]; x[34] = v[34]; x[35] = v[35];
467: x[36] = v[36]; x[37] = v[37]; x[38] = v[38]; x[39] = v[39];
468: x[40] = v[40]; x[41] = v[41]; x[42] = v[42]; x[43] = v[43];
469: x[44] = v[44]; x[45] = v[45]; x[46] = v[46]; x[47] = v[47];
470: x[48] = v[48];
471: v += 49;
472: }
473: row = *ajtmp++;
474: while (row < i) {
475: pc = rtmp + 49*row;
476: p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3];
477: p5 = pc[4]; p6 = pc[5]; p7 = pc[6]; p8 = pc[7];
478: p9 = pc[8]; p10 = pc[9]; p11 = pc[10]; p12 = pc[11];
479: p13 = pc[12]; p14 = pc[13]; p15 = pc[14]; p16 = pc[15];
480: p17 = pc[16]; p18 = pc[17]; p19 = pc[18]; p20 = pc[19];
481: p21 = pc[20]; p22 = pc[21]; p23 = pc[22]; p24 = pc[23];
482: p25 = pc[24]; p26 = pc[25]; p27 = pc[26]; p28 = pc[27];
483: p29 = pc[28]; p30 = pc[29]; p31 = pc[30]; p32 = pc[31];
484: p33 = pc[32]; p34 = pc[33]; p35 = pc[34]; p36 = pc[35];
485: p37 = pc[36]; p38 = pc[37]; p39 = pc[38]; p40 = pc[39];
486: p41 = pc[40]; p42 = pc[41]; p43 = pc[42]; p44 = pc[43];
487: p45 = pc[44]; p46 = pc[45]; p47 = pc[46]; p48 = pc[47];
488: p49 = pc[48];
489: if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 ||
490: p5 != 0.0 || p6 != 0.0 || p7 != 0.0 || p8 != 0.0 ||
491: p9 != 0.0 || p10 != 0.0 || p11 != 0.0 || p12 != 0.0 ||
492: p13 != 0.0 || p14 != 0.0 || p15 != 0.0 || p16 != 0.0 ||
493: p17 != 0.0 || p18 != 0.0 || p19 != 0.0 || p20 != 0.0 ||
494: p21 != 0.0 || p22 != 0.0 || p23 != 0.0 || p24 != 0.0 ||
495: p25 != 0.0 || p26 != 0.0 || p27 != 0.0 || p28 != 0.0 ||
496: p29 != 0.0 || p30 != 0.0 || p31 != 0.0 || p32 != 0.0 ||
497: p33 != 0.0 || p34 != 0.0 || p35 != 0.0 || p36 != 0.0 ||
498: p37 != 0.0 || p38 != 0.0 || p39 != 0.0 || p40 != 0.0 ||
499: p41 != 0.0 || p42 != 0.0 || p43 != 0.0 || p44 != 0.0 ||
500: p45 != 0.0 || p46 != 0.0 || p47 != 0.0 || p48 != 0.0 ||
501: p49 != 0.0) {
502: pv = ba + 49*diag_offset[row];
503: pj = bj + diag_offset[row] + 1;
504: x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
505: x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7];
506: x9 = pv[8]; x10 = pv[9]; x11 = pv[10]; x12 = pv[11];
507: x13 = pv[12]; x14 = pv[13]; x15 = pv[14]; x16 = pv[15];
508: x17 = pv[16]; x18 = pv[17]; x19 = pv[18]; x20 = pv[19];
509: x21 = pv[20]; x22 = pv[21]; x23 = pv[22]; x24 = pv[23];
510: x25 = pv[24]; x26 = pv[25]; x27 = pv[26]; x28 = pv[27];
511: x29 = pv[28]; x30 = pv[29]; x31 = pv[30]; x32 = pv[31];
512: x33 = pv[32]; x34 = pv[33]; x35 = pv[34]; x36 = pv[35];
513: x37 = pv[36]; x38 = pv[37]; x39 = pv[38]; x40 = pv[39];
514: x41 = pv[40]; x42 = pv[41]; x43 = pv[42]; x44 = pv[43];
515: x45 = pv[44]; x46 = pv[45]; x47 = pv[46]; x48 = pv[47];
516: x49 = pv[48];
517: pc[0] = m1 = p1*x1 + p8*x2 + p15*x3 + p22*x4 + p29*x5 + p36*x6 + p43*x7;
518: pc[1] = m2 = p2*x1 + p9*x2 + p16*x3 + p23*x4 + p30*x5 + p37*x6 + p44*x7;
519: pc[2] = m3 = p3*x1 + p10*x2 + p17*x3 + p24*x4 + p31*x5 + p38*x6 + p45*x7;
520: pc[3] = m4 = p4*x1 + p11*x2 + p18*x3 + p25*x4 + p32*x5 + p39*x6 + p46*x7;
521: pc[4] = m5 = p5*x1 + p12*x2 + p19*x3 + p26*x4 + p33*x5 + p40*x6 + p47*x7;
522: pc[5] = m6 = p6*x1 + p13*x2 + p20*x3 + p27*x4 + p34*x5 + p41*x6 + p48*x7;
523: pc[6] = m7 = p7*x1 + p14*x2 + p21*x3 + p28*x4 + p35*x5 + p42*x6 + p49*x7;
525: pc[7] = m8 = p1*x8 + p8*x9 + p15*x10 + p22*x11 + p29*x12 + p36*x13 + p43*x14;
526: pc[8] = m9 = p2*x8 + p9*x9 + p16*x10 + p23*x11 + p30*x12 + p37*x13 + p44*x14;
527: pc[9] = m10 = p3*x8 + p10*x9 + p17*x10 + p24*x11 + p31*x12 + p38*x13 + p45*x14;
528: pc[10] = m11 = p4*x8 + p11*x9 + p18*x10 + p25*x11 + p32*x12 + p39*x13 + p46*x14;
529: pc[11] = m12 = p5*x8 + p12*x9 + p19*x10 + p26*x11 + p33*x12 + p40*x13 + p47*x14;
530: pc[12] = m13 = p6*x8 + p13*x9 + p20*x10 + p27*x11 + p34*x12 + p41*x13 + p48*x14;
531: pc[13] = m14 = p7*x8 + p14*x9 + p21*x10 + p28*x11 + p35*x12 + p42*x13 + p49*x14;
533: pc[14] = m15 = p1*x15 + p8*x16 + p15*x17 + p22*x18 + p29*x19 + p36*x20 + p43*x21;
534: pc[15] = m16 = p2*x15 + p9*x16 + p16*x17 + p23*x18 + p30*x19 + p37*x20 + p44*x21;
535: pc[16] = m17 = p3*x15 + p10*x16 + p17*x17 + p24*x18 + p31*x19 + p38*x20 + p45*x21;
536: pc[17] = m18 = p4*x15 + p11*x16 + p18*x17 + p25*x18 + p32*x19 + p39*x20 + p46*x21;
537: pc[18] = m19 = p5*x15 + p12*x16 + p19*x17 + p26*x18 + p33*x19 + p40*x20 + p47*x21;
538: pc[19] = m20 = p6*x15 + p13*x16 + p20*x17 + p27*x18 + p34*x19 + p41*x20 + p48*x21;
539: pc[20] = m21 = p7*x15 + p14*x16 + p21*x17 + p28*x18 + p35*x19 + p42*x20 + p49*x21;
541: pc[21] = m22 = p1*x22 + p8*x23 + p15*x24 + p22*x25 + p29*x26 + p36*x27 + p43*x28;
542: pc[22] = m23 = p2*x22 + p9*x23 + p16*x24 + p23*x25 + p30*x26 + p37*x27 + p44*x28;
543: pc[23] = m24 = p3*x22 + p10*x23 + p17*x24 + p24*x25 + p31*x26 + p38*x27 + p45*x28;
544: pc[24] = m25 = p4*x22 + p11*x23 + p18*x24 + p25*x25 + p32*x26 + p39*x27 + p46*x28;
545: pc[25] = m26 = p5*x22 + p12*x23 + p19*x24 + p26*x25 + p33*x26 + p40*x27 + p47*x28;
546: pc[26] = m27 = p6*x22 + p13*x23 + p20*x24 + p27*x25 + p34*x26 + p41*x27 + p48*x28;
547: pc[27] = m28 = p7*x22 + p14*x23 + p21*x24 + p28*x25 + p35*x26 + p42*x27 + p49*x28;
549: pc[28] = m29 = p1*x29 + p8*x30 + p15*x31 + p22*x32 + p29*x33 + p36*x34 + p43*x35;
550: pc[29] = m30 = p2*x29 + p9*x30 + p16*x31 + p23*x32 + p30*x33 + p37*x34 + p44*x35;
551: pc[30] = m31 = p3*x29 + p10*x30 + p17*x31 + p24*x32 + p31*x33 + p38*x34 + p45*x35;
552: pc[31] = m32 = p4*x29 + p11*x30 + p18*x31 + p25*x32 + p32*x33 + p39*x34 + p46*x35;
553: pc[32] = m33 = p5*x29 + p12*x30 + p19*x31 + p26*x32 + p33*x33 + p40*x34 + p47*x35;
554: pc[33] = m34 = p6*x29 + p13*x30 + p20*x31 + p27*x32 + p34*x33 + p41*x34 + p48*x35;
555: pc[34] = m35 = p7*x29 + p14*x30 + p21*x31 + p28*x32 + p35*x33 + p42*x34 + p49*x35;
557: pc[35] = m36 = p1*x36 + p8*x37 + p15*x38 + p22*x39 + p29*x40 + p36*x41 + p43*x42;
558: pc[36] = m37 = p2*x36 + p9*x37 + p16*x38 + p23*x39 + p30*x40 + p37*x41 + p44*x42;
559: pc[37] = m38 = p3*x36 + p10*x37 + p17*x38 + p24*x39 + p31*x40 + p38*x41 + p45*x42;
560: pc[38] = m39 = p4*x36 + p11*x37 + p18*x38 + p25*x39 + p32*x40 + p39*x41 + p46*x42;
561: pc[39] = m40 = p5*x36 + p12*x37 + p19*x38 + p26*x39 + p33*x40 + p40*x41 + p47*x42;
562: pc[40] = m41 = p6*x36 + p13*x37 + p20*x38 + p27*x39 + p34*x40 + p41*x41 + p48*x42;
563: pc[41] = m42 = p7*x36 + p14*x37 + p21*x38 + p28*x39 + p35*x40 + p42*x41 + p49*x42;
565: pc[42] = m43 = p1*x43 + p8*x44 + p15*x45 + p22*x46 + p29*x47 + p36*x48 + p43*x49;
566: pc[43] = m44 = p2*x43 + p9*x44 + p16*x45 + p23*x46 + p30*x47 + p37*x48 + p44*x49;
567: pc[44] = m45 = p3*x43 + p10*x44 + p17*x45 + p24*x46 + p31*x47 + p38*x48 + p45*x49;
568: pc[45] = m46 = p4*x43 + p11*x44 + p18*x45 + p25*x46 + p32*x47 + p39*x48 + p46*x49;
569: pc[46] = m47 = p5*x43 + p12*x44 + p19*x45 + p26*x46 + p33*x47 + p40*x48 + p47*x49;
570: pc[47] = m48 = p6*x43 + p13*x44 + p20*x45 + p27*x46 + p34*x47 + p41*x48 + p48*x49;
571: pc[48] = m49 = p7*x43 + p14*x44 + p21*x45 + p28*x46 + p35*x47 + p42*x48 + p49*x49;
573: nz = bi[row+1] - diag_offset[row] - 1;
574: pv += 49;
575: for (j=0; j<nz; j++) {
576: x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
577: x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7];
578: x9 = pv[8]; x10 = pv[9]; x11 = pv[10]; x12 = pv[11];
579: x13 = pv[12]; x14 = pv[13]; x15 = pv[14]; x16 = pv[15];
580: x17 = pv[16]; x18 = pv[17]; x19 = pv[18]; x20 = pv[19];
581: x21 = pv[20]; x22 = pv[21]; x23 = pv[22]; x24 = pv[23];
582: x25 = pv[24]; x26 = pv[25]; x27 = pv[26]; x28 = pv[27];
583: x29 = pv[28]; x30 = pv[29]; x31 = pv[30]; x32 = pv[31];
584: x33 = pv[32]; x34 = pv[33]; x35 = pv[34]; x36 = pv[35];
585: x37 = pv[36]; x38 = pv[37]; x39 = pv[38]; x40 = pv[39];
586: x41 = pv[40]; x42 = pv[41]; x43 = pv[42]; x44 = pv[43];
587: x45 = pv[44]; x46 = pv[45]; x47 = pv[46]; x48 = pv[47];
588: x49 = pv[48];
589: x = rtmp + 49*pj[j];
590: x[0] -= m1*x1 + m8*x2 + m15*x3 + m22*x4 + m29*x5 + m36*x6 + m43*x7;
591: x[1] -= m2*x1 + m9*x2 + m16*x3 + m23*x4 + m30*x5 + m37*x6 + m44*x7;
592: x[2] -= m3*x1 + m10*x2 + m17*x3 + m24*x4 + m31*x5 + m38*x6 + m45*x7;
593: x[3] -= m4*x1 + m11*x2 + m18*x3 + m25*x4 + m32*x5 + m39*x6 + m46*x7;
594: x[4] -= m5*x1 + m12*x2 + m19*x3 + m26*x4 + m33*x5 + m40*x6 + m47*x7;
595: x[5] -= m6*x1 + m13*x2 + m20*x3 + m27*x4 + m34*x5 + m41*x6 + m48*x7;
596: x[6] -= m7*x1 + m14*x2 + m21*x3 + m28*x4 + m35*x5 + m42*x6 + m49*x7;
598: x[7] -= m1*x8 + m8*x9 + m15*x10 + m22*x11 + m29*x12 + m36*x13 + m43*x14;
599: x[8] -= m2*x8 + m9*x9 + m16*x10 + m23*x11 + m30*x12 + m37*x13 + m44*x14;
600: x[9] -= m3*x8 + m10*x9 + m17*x10 + m24*x11 + m31*x12 + m38*x13 + m45*x14;
601: x[10] -= m4*x8 + m11*x9 + m18*x10 + m25*x11 + m32*x12 + m39*x13 + m46*x14;
602: x[11] -= m5*x8 + m12*x9 + m19*x10 + m26*x11 + m33*x12 + m40*x13 + m47*x14;
603: x[12] -= m6*x8 + m13*x9 + m20*x10 + m27*x11 + m34*x12 + m41*x13 + m48*x14;
604: x[13] -= m7*x8 + m14*x9 + m21*x10 + m28*x11 + m35*x12 + m42*x13 + m49*x14;
606: x[14] -= m1*x15 + m8*x16 + m15*x17 + m22*x18 + m29*x19 + m36*x20 + m43*x21;
607: x[15] -= m2*x15 + m9*x16 + m16*x17 + m23*x18 + m30*x19 + m37*x20 + m44*x21;
608: x[16] -= m3*x15 + m10*x16 + m17*x17 + m24*x18 + m31*x19 + m38*x20 + m45*x21;
609: x[17] -= m4*x15 + m11*x16 + m18*x17 + m25*x18 + m32*x19 + m39*x20 + m46*x21;
610: x[18] -= m5*x15 + m12*x16 + m19*x17 + m26*x18 + m33*x19 + m40*x20 + m47*x21;
611: x[19] -= m6*x15 + m13*x16 + m20*x17 + m27*x18 + m34*x19 + m41*x20 + m48*x21;
612: x[20] -= m7*x15 + m14*x16 + m21*x17 + m28*x18 + m35*x19 + m42*x20 + m49*x21;
614: x[21] -= m1*x22 + m8*x23 + m15*x24 + m22*x25 + m29*x26 + m36*x27 + m43*x28;
615: x[22] -= m2*x22 + m9*x23 + m16*x24 + m23*x25 + m30*x26 + m37*x27 + m44*x28;
616: x[23] -= m3*x22 + m10*x23 + m17*x24 + m24*x25 + m31*x26 + m38*x27 + m45*x28;
617: x[24] -= m4*x22 + m11*x23 + m18*x24 + m25*x25 + m32*x26 + m39*x27 + m46*x28;
618: x[25] -= m5*x22 + m12*x23 + m19*x24 + m26*x25 + m33*x26 + m40*x27 + m47*x28;
619: x[26] -= m6*x22 + m13*x23 + m20*x24 + m27*x25 + m34*x26 + m41*x27 + m48*x28;
620: x[27] -= m7*x22 + m14*x23 + m21*x24 + m28*x25 + m35*x26 + m42*x27 + m49*x28;
622: x[28] -= m1*x29 + m8*x30 + m15*x31 + m22*x32 + m29*x33 + m36*x34 + m43*x35;
623: x[29] -= m2*x29 + m9*x30 + m16*x31 + m23*x32 + m30*x33 + m37*x34 + m44*x35;
624: x[30] -= m3*x29 + m10*x30 + m17*x31 + m24*x32 + m31*x33 + m38*x34 + m45*x35;
625: x[31] -= m4*x29 + m11*x30 + m18*x31 + m25*x32 + m32*x33 + m39*x34 + m46*x35;
626: x[32] -= m5*x29 + m12*x30 + m19*x31 + m26*x32 + m33*x33 + m40*x34 + m47*x35;
627: x[33] -= m6*x29 + m13*x30 + m20*x31 + m27*x32 + m34*x33 + m41*x34 + m48*x35;
628: x[34] -= m7*x29 + m14*x30 + m21*x31 + m28*x32 + m35*x33 + m42*x34 + m49*x35;
630: x[35] -= m1*x36 + m8*x37 + m15*x38 + m22*x39 + m29*x40 + m36*x41 + m43*x42;
631: x[36] -= m2*x36 + m9*x37 + m16*x38 + m23*x39 + m30*x40 + m37*x41 + m44*x42;
632: x[37] -= m3*x36 + m10*x37 + m17*x38 + m24*x39 + m31*x40 + m38*x41 + m45*x42;
633: x[38] -= m4*x36 + m11*x37 + m18*x38 + m25*x39 + m32*x40 + m39*x41 + m46*x42;
634: x[39] -= m5*x36 + m12*x37 + m19*x38 + m26*x39 + m33*x40 + m40*x41 + m47*x42;
635: x[40] -= m6*x36 + m13*x37 + m20*x38 + m27*x39 + m34*x40 + m41*x41 + m48*x42;
636: x[41] -= m7*x36 + m14*x37 + m21*x38 + m28*x39 + m35*x40 + m42*x41 + m49*x42;
638: x[42] -= m1*x43 + m8*x44 + m15*x45 + m22*x46 + m29*x47 + m36*x48 + m43*x49;
639: x[43] -= m2*x43 + m9*x44 + m16*x45 + m23*x46 + m30*x47 + m37*x48 + m44*x49;
640: x[44] -= m3*x43 + m10*x44 + m17*x45 + m24*x46 + m31*x47 + m38*x48 + m45*x49;
641: x[45] -= m4*x43 + m11*x44 + m18*x45 + m25*x46 + m32*x47 + m39*x48 + m46*x49;
642: x[46] -= m5*x43 + m12*x44 + m19*x45 + m26*x46 + m33*x47 + m40*x48 + m47*x49;
643: x[47] -= m6*x43 + m13*x44 + m20*x45 + m27*x46 + m34*x47 + m41*x48 + m48*x49;
644: x[48] -= m7*x43 + m14*x44 + m21*x45 + m28*x46 + m35*x47 + m42*x48 + m49*x49;
645: pv += 49;
646: }
647: PetscLogFlops(686.0*nz+637.0);
648: }
649: row = *ajtmp++;
650: }
651: /* finished row so stick it into b->a */
652: pv = ba + 49*bi[i];
653: pj = bj + bi[i];
654: nz = bi[i+1] - bi[i];
655: for (j=0; j<nz; j++) {
656: x = rtmp+49*pj[j];
657: pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3];
658: pv[4] = x[4]; pv[5] = x[5]; pv[6] = x[6]; pv[7] = x[7];
659: pv[8] = x[8]; pv[9] = x[9]; pv[10] = x[10]; pv[11] = x[11];
660: pv[12] = x[12]; pv[13] = x[13]; pv[14] = x[14]; pv[15] = x[15];
661: pv[16] = x[16]; pv[17] = x[17]; pv[18] = x[18]; pv[19] = x[19];
662: pv[20] = x[20]; pv[21] = x[21]; pv[22] = x[22]; pv[23] = x[23];
663: pv[24] = x[24]; pv[25] = x[25]; pv[26] = x[26]; pv[27] = x[27];
664: pv[28] = x[28]; pv[29] = x[29]; pv[30] = x[30]; pv[31] = x[31];
665: pv[32] = x[32]; pv[33] = x[33]; pv[34] = x[34]; pv[35] = x[35];
666: pv[36] = x[36]; pv[37] = x[37]; pv[38] = x[38]; pv[39] = x[39];
667: pv[40] = x[40]; pv[41] = x[41]; pv[42] = x[42]; pv[43] = x[43];
668: pv[44] = x[44]; pv[45] = x[45]; pv[46] = x[46]; pv[47] = x[47];
669: pv[48] = x[48];
670: pv += 49;
671: }
672: /* invert diagonal block */
673: w = ba + 49*diag_offset[i];
674: PetscKernel_A_gets_inverse_A_7(w,shift,allowzeropivot,&zeropivotdetected);
675: if (zeropivotdetected) C->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
676: }
678: PetscFree(rtmp);
680: C->ops->solve = MatSolve_SeqBAIJ_7_NaturalOrdering_inplace;
681: C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_7_NaturalOrdering_inplace;
682: C->assembled = PETSC_TRUE;
684: PetscLogFlops(1.333333333333*7*7*7*b->mbs); /* from inverting diagonal blocks */
685: return(0);
686: }
690: PetscErrorCode MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info)
691: {
692: Mat C =B;
693: Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data;
695: PetscInt i,j,k,nz,nzL,row;
696: const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
697: const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag,*pj,bs2=a->bs2;
698: MatScalar *rtmp,*pc,*mwork,*v,*pv,*aa=a->a;
699: PetscInt flg;
700: PetscReal shift = info->shiftamount;
701: PetscBool allowzeropivot,zeropivotdetected;
704: allowzeropivot = PetscNot(A->erroriffailure);
706: /* generate work space needed by the factorization */
707: PetscMalloc2(bs2*n,&rtmp,bs2,&mwork);
708: PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));
710: for (i=0; i<n; i++) {
711: /* zero rtmp */
712: /* L part */
713: nz = bi[i+1] - bi[i];
714: bjtmp = bj + bi[i];
715: for (j=0; j<nz; j++) {
716: PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));
717: }
719: /* U part */
720: nz = bdiag[i] - bdiag[i+1];
721: bjtmp = bj + bdiag[i+1]+1;
722: for (j=0; j<nz; j++) {
723: PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));
724: }
726: /* load in initial (unfactored row) */
727: nz = ai[i+1] - ai[i];
728: ajtmp = aj + ai[i];
729: v = aa + bs2*ai[i];
730: for (j=0; j<nz; j++) {
731: PetscMemcpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2*sizeof(MatScalar));
732: }
734: /* elimination */
735: bjtmp = bj + bi[i];
736: nzL = bi[i+1] - bi[i];
737: for (k=0; k < nzL; k++) {
738: row = bjtmp[k];
739: pc = rtmp + bs2*row;
740: for (flg=0,j=0; j<bs2; j++) {
741: if (pc[j]!=0.0) {
742: flg = 1;
743: break;
744: }
745: }
746: if (flg) {
747: pv = b->a + bs2*bdiag[row];
748: /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
749: PetscKernel_A_gets_A_times_B_7(pc,pv,mwork);
751: pj = b->j + bdiag[row+1]+1; /* begining of U(row,:) */
752: pv = b->a + bs2*(bdiag[row+1]+1);
753: nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries inU(row,:), excluding diag */
754: for (j=0; j<nz; j++) {
755: /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
756: /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
757: v = rtmp + bs2*pj[j];
758: PetscKernel_A_gets_A_minus_B_times_C_7(v,pc,pv);
759: pv += bs2;
760: }
761: PetscLogFlops(686*nz+637); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
762: }
763: }
765: /* finished row so stick it into b->a */
766: /* L part */
767: pv = b->a + bs2*bi[i];
768: pj = b->j + bi[i];
769: nz = bi[i+1] - bi[i];
770: for (j=0; j<nz; j++) {
771: PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));
772: }
774: /* Mark diagonal and invert diagonal for simplier triangular solves */
775: pv = b->a + bs2*bdiag[i];
776: pj = b->j + bdiag[i];
777: PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));
778: PetscKernel_A_gets_inverse_A_7(pv,shift,allowzeropivot,&zeropivotdetected);
779: if (zeropivotdetected) C->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
781: /* U part */
782: pv = b->a + bs2*(bdiag[i+1]+1);
783: pj = b->j + bdiag[i+1]+1;
784: nz = bdiag[i] - bdiag[i+1] - 1;
785: for (j=0; j<nz; j++) {
786: PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));
787: }
788: }
789: PetscFree2(rtmp,mwork);
791: C->ops->solve = MatSolve_SeqBAIJ_7_NaturalOrdering;
792: C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_7_NaturalOrdering;
793: C->assembled = PETSC_TRUE;
795: PetscLogFlops(1.333333333333*7*7*7*n); /* from inverting diagonal blocks */
796: return(0);
797: }