Actual source code: baijfact7.c
petsc-3.13.6 2020-09-29
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 6 by 6
11: */
12: PetscErrorCode MatLUFactorNumeric_SeqBAIJ_6_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 *ajtmpold,*ajtmp,*diag_offset = b->diag,*r,*ic,*bi = b->i,*bj = b->j,*ai=a->i,*aj=a->j,*pj;
18: PetscInt nz,row,i,j,n = a->mbs,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 x26,x27,x28,x29,x30,x31,x32,x33,x34,x35,x36;
27: MatScalar m26,m27,m28,m29,m30,m31,m32,m33,m34,m35,m36;
28: MatScalar *ba = b->a,*aa = a->a;
29: PetscReal shift = info->shiftamount;
30: PetscBool allowzeropivot,zeropivotdetected;
33: allowzeropivot = PetscNot(A->erroriffailure);
34: ISGetIndices(isrow,&r);
35: ISGetIndices(isicol,&ic);
36: PetscMalloc1(36*(n+1),&rtmp);
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: x = rtmp+36*ajtmp[j];
43: x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = x[9] = 0.0;
44: x[10] = x[11] = x[12] = x[13] = x[14] = x[15] = x[16] = x[17] = 0.0;
45: x[18] = x[19] = x[20] = x[21] = x[22] = x[23] = x[24] = x[25] = 0.0;
46: x[26] = x[27] = x[28] = x[29] = x[30] = x[31] = x[32] = x[33] = 0.0;
47: x[34] = x[35] = 0.0;
48: }
49: /* load in initial (unfactored row) */
50: idx = r[i];
51: nz = ai[idx+1] - ai[idx];
52: ajtmpold = aj + ai[idx];
53: v = aa + 36*ai[idx];
54: for (j=0; j<nz; j++) {
55: x = rtmp+36*ic[ajtmpold[j]];
56: x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3];
57: x[4] = v[4]; x[5] = v[5]; x[6] = v[6]; x[7] = v[7];
58: x[8] = v[8]; x[9] = v[9]; x[10] = v[10]; x[11] = v[11];
59: x[12] = v[12]; x[13] = v[13]; x[14] = v[14]; x[15] = v[15];
60: x[16] = v[16]; x[17] = v[17]; x[18] = v[18]; x[19] = v[19];
61: x[20] = v[20]; x[21] = v[21]; x[22] = v[22]; x[23] = v[23];
62: x[24] = v[24]; x[25] = v[25]; x[26] = v[26]; x[27] = v[27];
63: x[28] = v[28]; x[29] = v[29]; x[30] = v[30]; x[31] = v[31];
64: x[32] = v[32]; x[33] = v[33]; x[34] = v[34]; x[35] = v[35];
65: v += 36;
66: }
67: row = *ajtmp++;
68: while (row < i) {
69: pc = rtmp + 36*row;
70: p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3];
71: p5 = pc[4]; p6 = pc[5]; p7 = pc[6]; p8 = pc[7];
72: p9 = pc[8]; p10 = pc[9]; p11 = pc[10]; p12 = pc[11];
73: p13 = pc[12]; p14 = pc[13]; p15 = pc[14]; p16 = pc[15];
74: p17 = pc[16]; p18 = pc[17]; p19 = pc[18]; p20 = pc[19];
75: p21 = pc[20]; p22 = pc[21]; p23 = pc[22]; p24 = pc[23];
76: p25 = pc[24]; p26 = pc[25]; p27 = pc[26]; p28 = pc[27];
77: p29 = pc[28]; p30 = pc[29]; p31 = pc[30]; p32 = pc[31];
78: p33 = pc[32]; p34 = pc[33]; p35 = pc[34]; p36 = pc[35];
79: if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 ||
80: p5 != 0.0 || p6 != 0.0 || p7 != 0.0 || p8 != 0.0 ||
81: p9 != 0.0 || p10 != 0.0 || p11 != 0.0 || p12 != 0.0 ||
82: p13 != 0.0 || p14 != 0.0 || p15 != 0.0 || p16 != 0.0 ||
83: p17 != 0.0 || p18 != 0.0 || p19 != 0.0 || p20 != 0.0 ||
84: p21 != 0.0 || p22 != 0.0 || p23 != 0.0 || p24 != 0.0 ||
85: p25 != 0.0 || p26 != 0.0 || p27 != 0.0 || p28 != 0.0 ||
86: p29 != 0.0 || p30 != 0.0 || p31 != 0.0 || p32 != 0.0 ||
87: p33 != 0.0 || p34 != 0.0 || p35 != 0.0 || p36 != 0.0) {
88: pv = ba + 36*diag_offset[row];
89: pj = bj + diag_offset[row] + 1;
90: x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
91: x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7];
92: x9 = pv[8]; x10 = pv[9]; x11 = pv[10]; x12 = pv[11];
93: x13 = pv[12]; x14 = pv[13]; x15 = pv[14]; x16 = pv[15];
94: x17 = pv[16]; x18 = pv[17]; x19 = pv[18]; x20 = pv[19];
95: x21 = pv[20]; x22 = pv[21]; x23 = pv[22]; x24 = pv[23];
96: x25 = pv[24]; x26 = pv[25]; x27 = pv[26]; x28 = pv[27];
97: x29 = pv[28]; x30 = pv[29]; x31 = pv[30]; x32 = pv[31];
98: x33 = pv[32]; x34 = pv[33]; x35 = pv[34]; x36 = pv[35];
99: pc[0] = m1 = p1*x1 + p7*x2 + p13*x3 + p19*x4 + p25*x5 + p31*x6;
100: pc[1] = m2 = p2*x1 + p8*x2 + p14*x3 + p20*x4 + p26*x5 + p32*x6;
101: pc[2] = m3 = p3*x1 + p9*x2 + p15*x3 + p21*x4 + p27*x5 + p33*x6;
102: pc[3] = m4 = p4*x1 + p10*x2 + p16*x3 + p22*x4 + p28*x5 + p34*x6;
103: pc[4] = m5 = p5*x1 + p11*x2 + p17*x3 + p23*x4 + p29*x5 + p35*x6;
104: pc[5] = m6 = p6*x1 + p12*x2 + p18*x3 + p24*x4 + p30*x5 + p36*x6;
106: pc[6] = m7 = p1*x7 + p7*x8 + p13*x9 + p19*x10 + p25*x11 + p31*x12;
107: pc[7] = m8 = p2*x7 + p8*x8 + p14*x9 + p20*x10 + p26*x11 + p32*x12;
108: pc[8] = m9 = p3*x7 + p9*x8 + p15*x9 + p21*x10 + p27*x11 + p33*x12;
109: pc[9] = m10 = p4*x7 + p10*x8 + p16*x9 + p22*x10 + p28*x11 + p34*x12;
110: pc[10] = m11 = p5*x7 + p11*x8 + p17*x9 + p23*x10 + p29*x11 + p35*x12;
111: pc[11] = m12 = p6*x7 + p12*x8 + p18*x9 + p24*x10 + p30*x11 + p36*x12;
113: pc[12] = m13 = p1*x13 + p7*x14 + p13*x15 + p19*x16 + p25*x17 + p31*x18;
114: pc[13] = m14 = p2*x13 + p8*x14 + p14*x15 + p20*x16 + p26*x17 + p32*x18;
115: pc[14] = m15 = p3*x13 + p9*x14 + p15*x15 + p21*x16 + p27*x17 + p33*x18;
116: pc[15] = m16 = p4*x13 + p10*x14 + p16*x15 + p22*x16 + p28*x17 + p34*x18;
117: pc[16] = m17 = p5*x13 + p11*x14 + p17*x15 + p23*x16 + p29*x17 + p35*x18;
118: pc[17] = m18 = p6*x13 + p12*x14 + p18*x15 + p24*x16 + p30*x17 + p36*x18;
120: pc[18] = m19 = p1*x19 + p7*x20 + p13*x21 + p19*x22 + p25*x23 + p31*x24;
121: pc[19] = m20 = p2*x19 + p8*x20 + p14*x21 + p20*x22 + p26*x23 + p32*x24;
122: pc[20] = m21 = p3*x19 + p9*x20 + p15*x21 + p21*x22 + p27*x23 + p33*x24;
123: pc[21] = m22 = p4*x19 + p10*x20 + p16*x21 + p22*x22 + p28*x23 + p34*x24;
124: pc[22] = m23 = p5*x19 + p11*x20 + p17*x21 + p23*x22 + p29*x23 + p35*x24;
125: pc[23] = m24 = p6*x19 + p12*x20 + p18*x21 + p24*x22 + p30*x23 + p36*x24;
127: pc[24] = m25 = p1*x25 + p7*x26 + p13*x27 + p19*x28 + p25*x29 + p31*x30;
128: pc[25] = m26 = p2*x25 + p8*x26 + p14*x27 + p20*x28 + p26*x29 + p32*x30;
129: pc[26] = m27 = p3*x25 + p9*x26 + p15*x27 + p21*x28 + p27*x29 + p33*x30;
130: pc[27] = m28 = p4*x25 + p10*x26 + p16*x27 + p22*x28 + p28*x29 + p34*x30;
131: pc[28] = m29 = p5*x25 + p11*x26 + p17*x27 + p23*x28 + p29*x29 + p35*x30;
132: pc[29] = m30 = p6*x25 + p12*x26 + p18*x27 + p24*x28 + p30*x29 + p36*x30;
134: pc[30] = m31 = p1*x31 + p7*x32 + p13*x33 + p19*x34 + p25*x35 + p31*x36;
135: pc[31] = m32 = p2*x31 + p8*x32 + p14*x33 + p20*x34 + p26*x35 + p32*x36;
136: pc[32] = m33 = p3*x31 + p9*x32 + p15*x33 + p21*x34 + p27*x35 + p33*x36;
137: pc[33] = m34 = p4*x31 + p10*x32 + p16*x33 + p22*x34 + p28*x35 + p34*x36;
138: pc[34] = m35 = p5*x31 + p11*x32 + p17*x33 + p23*x34 + p29*x35 + p35*x36;
139: pc[35] = m36 = p6*x31 + p12*x32 + p18*x33 + p24*x34 + p30*x35 + p36*x36;
141: nz = bi[row+1] - diag_offset[row] - 1;
142: pv += 36;
143: for (j=0; j<nz; j++) {
144: x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
145: x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7];
146: x9 = pv[8]; x10 = pv[9]; x11 = pv[10]; x12 = pv[11];
147: x13 = pv[12]; x14 = pv[13]; x15 = pv[14]; x16 = pv[15];
148: x17 = pv[16]; x18 = pv[17]; x19 = pv[18]; x20 = pv[19];
149: x21 = pv[20]; x22 = pv[21]; x23 = pv[22]; x24 = pv[23];
150: x25 = pv[24]; x26 = pv[25]; x27 = pv[26]; x28 = pv[27];
151: x29 = pv[28]; x30 = pv[29]; x31 = pv[30]; x32 = pv[31];
152: x33 = pv[32]; x34 = pv[33]; x35 = pv[34]; x36 = pv[35];
153: x = rtmp + 36*pj[j];
154: x[0] -= m1*x1 + m7*x2 + m13*x3 + m19*x4 + m25*x5 + m31*x6;
155: x[1] -= m2*x1 + m8*x2 + m14*x3 + m20*x4 + m26*x5 + m32*x6;
156: x[2] -= m3*x1 + m9*x2 + m15*x3 + m21*x4 + m27*x5 + m33*x6;
157: x[3] -= m4*x1 + m10*x2 + m16*x3 + m22*x4 + m28*x5 + m34*x6;
158: x[4] -= m5*x1 + m11*x2 + m17*x3 + m23*x4 + m29*x5 + m35*x6;
159: x[5] -= m6*x1 + m12*x2 + m18*x3 + m24*x4 + m30*x5 + m36*x6;
161: x[6] -= m1*x7 + m7*x8 + m13*x9 + m19*x10 + m25*x11 + m31*x12;
162: x[7] -= m2*x7 + m8*x8 + m14*x9 + m20*x10 + m26*x11 + m32*x12;
163: x[8] -= m3*x7 + m9*x8 + m15*x9 + m21*x10 + m27*x11 + m33*x12;
164: x[9] -= m4*x7 + m10*x8 + m16*x9 + m22*x10 + m28*x11 + m34*x12;
165: x[10] -= m5*x7 + m11*x8 + m17*x9 + m23*x10 + m29*x11 + m35*x12;
166: x[11] -= m6*x7 + m12*x8 + m18*x9 + m24*x10 + m30*x11 + m36*x12;
168: x[12] -= m1*x13 + m7*x14 + m13*x15 + m19*x16 + m25*x17 + m31*x18;
169: x[13] -= m2*x13 + m8*x14 + m14*x15 + m20*x16 + m26*x17 + m32*x18;
170: x[14] -= m3*x13 + m9*x14 + m15*x15 + m21*x16 + m27*x17 + m33*x18;
171: x[15] -= m4*x13 + m10*x14 + m16*x15 + m22*x16 + m28*x17 + m34*x18;
172: x[16] -= m5*x13 + m11*x14 + m17*x15 + m23*x16 + m29*x17 + m35*x18;
173: x[17] -= m6*x13 + m12*x14 + m18*x15 + m24*x16 + m30*x17 + m36*x18;
175: x[18] -= m1*x19 + m7*x20 + m13*x21 + m19*x22 + m25*x23 + m31*x24;
176: x[19] -= m2*x19 + m8*x20 + m14*x21 + m20*x22 + m26*x23 + m32*x24;
177: x[20] -= m3*x19 + m9*x20 + m15*x21 + m21*x22 + m27*x23 + m33*x24;
178: x[21] -= m4*x19 + m10*x20 + m16*x21 + m22*x22 + m28*x23 + m34*x24;
179: x[22] -= m5*x19 + m11*x20 + m17*x21 + m23*x22 + m29*x23 + m35*x24;
180: x[23] -= m6*x19 + m12*x20 + m18*x21 + m24*x22 + m30*x23 + m36*x24;
182: x[24] -= m1*x25 + m7*x26 + m13*x27 + m19*x28 + m25*x29 + m31*x30;
183: x[25] -= m2*x25 + m8*x26 + m14*x27 + m20*x28 + m26*x29 + m32*x30;
184: x[26] -= m3*x25 + m9*x26 + m15*x27 + m21*x28 + m27*x29 + m33*x30;
185: x[27] -= m4*x25 + m10*x26 + m16*x27 + m22*x28 + m28*x29 + m34*x30;
186: x[28] -= m5*x25 + m11*x26 + m17*x27 + m23*x28 + m29*x29 + m35*x30;
187: x[29] -= m6*x25 + m12*x26 + m18*x27 + m24*x28 + m30*x29 + m36*x30;
189: x[30] -= m1*x31 + m7*x32 + m13*x33 + m19*x34 + m25*x35 + m31*x36;
190: x[31] -= m2*x31 + m8*x32 + m14*x33 + m20*x34 + m26*x35 + m32*x36;
191: x[32] -= m3*x31 + m9*x32 + m15*x33 + m21*x34 + m27*x35 + m33*x36;
192: x[33] -= m4*x31 + m10*x32 + m16*x33 + m22*x34 + m28*x35 + m34*x36;
193: x[34] -= m5*x31 + m11*x32 + m17*x33 + m23*x34 + m29*x35 + m35*x36;
194: x[35] -= m6*x31 + m12*x32 + m18*x33 + m24*x34 + m30*x35 + m36*x36;
196: pv += 36;
197: }
198: PetscLogFlops(432.0*nz+396.0);
199: }
200: row = *ajtmp++;
201: }
202: /* finished row so stick it into b->a */
203: pv = ba + 36*bi[i];
204: pj = bj + bi[i];
205: nz = bi[i+1] - bi[i];
206: for (j=0; j<nz; j++) {
207: x = rtmp+36*pj[j];
208: pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3];
209: pv[4] = x[4]; pv[5] = x[5]; pv[6] = x[6]; pv[7] = x[7];
210: pv[8] = x[8]; pv[9] = x[9]; pv[10] = x[10]; pv[11] = x[11];
211: pv[12] = x[12]; pv[13] = x[13]; pv[14] = x[14]; pv[15] = x[15];
212: pv[16] = x[16]; pv[17] = x[17]; pv[18] = x[18]; pv[19] = x[19];
213: pv[20] = x[20]; pv[21] = x[21]; pv[22] = x[22]; pv[23] = x[23];
214: pv[24] = x[24]; pv[25] = x[25]; pv[26] = x[26]; pv[27] = x[27];
215: pv[28] = x[28]; pv[29] = x[29]; pv[30] = x[30]; pv[31] = x[31];
216: pv[32] = x[32]; pv[33] = x[33]; pv[34] = x[34]; pv[35] = x[35];
217: pv += 36;
218: }
219: /* invert diagonal block */
220: w = ba + 36*diag_offset[i];
221: PetscKernel_A_gets_inverse_A_6(w,shift,allowzeropivot,&zeropivotdetected);
222: if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
223: }
225: PetscFree(rtmp);
226: ISRestoreIndices(isicol,&ic);
227: ISRestoreIndices(isrow,&r);
229: C->ops->solve = MatSolve_SeqBAIJ_6_inplace;
230: C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_6_inplace;
231: C->assembled = PETSC_TRUE;
233: PetscLogFlops(1.333333333333*6*6*6*b->mbs); /* from inverting diagonal blocks */
234: return(0);
235: }
237: PetscErrorCode MatLUFactorNumeric_SeqBAIJ_6(Mat B,Mat A,const MatFactorInfo *info)
238: {
239: Mat C = B;
240: Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data;
241: IS isrow = b->row,isicol = b->icol;
243: const PetscInt *r,*ic;
244: PetscInt i,j,k,nz,nzL,row;
245: const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
246: const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag,*pj,bs2=a->bs2;
247: MatScalar *rtmp,*pc,*mwork,*v,*pv,*aa=a->a;
248: PetscInt flg;
249: PetscReal shift = info->shiftamount;
250: PetscBool allowzeropivot,zeropivotdetected;
253: allowzeropivot = PetscNot(A->erroriffailure);
254: ISGetIndices(isrow,&r);
255: ISGetIndices(isicol,&ic);
257: /* generate work space needed by the factorization */
258: PetscMalloc2(bs2*n,&rtmp,bs2,&mwork);
259: PetscArrayzero(rtmp,bs2*n);
261: for (i=0; i<n; i++) {
262: /* zero rtmp */
263: /* L part */
264: nz = bi[i+1] - bi[i];
265: bjtmp = bj + bi[i];
266: for (j=0; j<nz; j++) {
267: PetscArrayzero(rtmp+bs2*bjtmp[j],bs2);
268: }
270: /* U part */
271: nz = bdiag[i] - bdiag[i+1];
272: bjtmp = bj + bdiag[i+1]+1;
273: for (j=0; j<nz; j++) {
274: PetscArrayzero(rtmp+bs2*bjtmp[j],bs2);
275: }
277: /* load in initial (unfactored row) */
278: nz = ai[r[i]+1] - ai[r[i]];
279: ajtmp = aj + ai[r[i]];
280: v = aa + bs2*ai[r[i]];
281: for (j=0; j<nz; j++) {
282: PetscArraycpy(rtmp+bs2*ic[ajtmp[j]],v+bs2*j,bs2);
283: }
285: /* elimination */
286: bjtmp = bj + bi[i];
287: nzL = bi[i+1] - bi[i];
288: for (k=0; k < nzL; k++) {
289: row = bjtmp[k];
290: pc = rtmp + bs2*row;
291: for (flg=0,j=0; j<bs2; j++) {
292: if (pc[j]!=0.0) {
293: flg = 1;
294: break;
295: }
296: }
297: if (flg) {
298: pv = b->a + bs2*bdiag[row];
299: /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
300: PetscKernel_A_gets_A_times_B_6(pc,pv,mwork);
302: pj = b->j + bdiag[row+1]+1; /* begining of U(row,:) */
303: pv = b->a + bs2*(bdiag[row+1]+1);
304: nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries inU(row,:), excluding diag */
305: for (j=0; j<nz; j++) {
306: /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
307: /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
308: v = rtmp + bs2*pj[j];
309: PetscKernel_A_gets_A_minus_B_times_C_6(v,pc,pv);
310: pv += bs2;
311: }
312: PetscLogFlops(432.0*nz+396); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
313: }
314: }
316: /* finished row so stick it into b->a */
317: /* L part */
318: pv = b->a + bs2*bi[i];
319: pj = b->j + bi[i];
320: nz = bi[i+1] - bi[i];
321: for (j=0; j<nz; j++) {
322: PetscArraycpy(pv+bs2*j,rtmp+bs2*pj[j],bs2);
323: }
325: /* Mark diagonal and invert diagonal for simplier triangular solves */
326: pv = b->a + bs2*bdiag[i];
327: pj = b->j + bdiag[i];
328: PetscArraycpy(pv,rtmp+bs2*pj[0],bs2);
329: PetscKernel_A_gets_inverse_A_6(pv,shift,allowzeropivot,&zeropivotdetected);
330: if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
332: /* U part */
333: pv = b->a + bs2*(bdiag[i+1]+1);
334: pj = b->j + bdiag[i+1]+1;
335: nz = bdiag[i] - bdiag[i+1] - 1;
336: for (j=0; j<nz; j++) {
337: PetscArraycpy(pv+bs2*j,rtmp+bs2*pj[j],bs2);
338: }
339: }
341: PetscFree2(rtmp,mwork);
342: ISRestoreIndices(isicol,&ic);
343: ISRestoreIndices(isrow,&r);
345: C->ops->solve = MatSolve_SeqBAIJ_6;
346: C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_6;
347: C->assembled = PETSC_TRUE;
349: PetscLogFlops(1.333333333333*6*6*6*n); /* from inverting diagonal blocks */
350: return(0);
351: }
353: PetscErrorCode MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering_inplace(Mat C,Mat A,const MatFactorInfo *info)
354: {
355: Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)C->data;
357: PetscInt i,j,n = a->mbs,*bi = b->i,*bj = b->j;
358: PetscInt *ajtmpold,*ajtmp,nz,row;
359: PetscInt *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj;
360: MatScalar *pv,*v,*rtmp,*pc,*w,*x;
361: MatScalar x1,x2,x3,x4,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15;
362: MatScalar x16,x17,x18,x19,x20,x21,x22,x23,x24,x25;
363: MatScalar p1,p2,p3,p4,p5,p6,p7,p8,p9,p10,p11,p12,p13,p14,p15;
364: MatScalar p16,p17,p18,p19,p20,p21,p22,p23,p24,p25;
365: MatScalar m1,m2,m3,m4,m5,m6,m7,m8,m9,m10,m11,m12,m13,m14,m15;
366: MatScalar m16,m17,m18,m19,m20,m21,m22,m23,m24,m25;
367: MatScalar p26,p27,p28,p29,p30,p31,p32,p33,p34,p35,p36;
368: MatScalar x26,x27,x28,x29,x30,x31,x32,x33,x34,x35,x36;
369: MatScalar m26,m27,m28,m29,m30,m31,m32,m33,m34,m35,m36;
370: MatScalar *ba = b->a,*aa = a->a;
371: PetscReal shift = info->shiftamount;
372: PetscBool allowzeropivot,zeropivotdetected;
375: allowzeropivot = PetscNot(A->erroriffailure);
376: PetscMalloc1(36*(n+1),&rtmp);
377: for (i=0; i<n; i++) {
378: nz = bi[i+1] - bi[i];
379: ajtmp = bj + bi[i];
380: for (j=0; j<nz; j++) {
381: x = rtmp+36*ajtmp[j];
382: x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = x[9] = 0.0;
383: x[10] = x[11] = x[12] = x[13] = x[14] = x[15] = x[16] = x[17] = 0.0;
384: x[18] = x[19] = x[20] = x[21] = x[22] = x[23] = x[24] = x[25] = 0.0;
385: x[26] = x[27] = x[28] = x[29] = x[30] = x[31] = x[32] = x[33] = 0.0;
386: x[34] = x[35] = 0.0;
387: }
388: /* load in initial (unfactored row) */
389: nz = ai[i+1] - ai[i];
390: ajtmpold = aj + ai[i];
391: v = aa + 36*ai[i];
392: for (j=0; j<nz; j++) {
393: x = rtmp+36*ajtmpold[j];
394: x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3];
395: x[4] = v[4]; x[5] = v[5]; x[6] = v[6]; x[7] = v[7];
396: x[8] = v[8]; x[9] = v[9]; x[10] = v[10]; x[11] = v[11];
397: x[12] = v[12]; x[13] = v[13]; x[14] = v[14]; x[15] = v[15];
398: x[16] = v[16]; x[17] = v[17]; x[18] = v[18]; x[19] = v[19];
399: x[20] = v[20]; x[21] = v[21]; x[22] = v[22]; x[23] = v[23];
400: x[24] = v[24]; x[25] = v[25]; x[26] = v[26]; x[27] = v[27];
401: x[28] = v[28]; x[29] = v[29]; x[30] = v[30]; x[31] = v[31];
402: x[32] = v[32]; x[33] = v[33]; x[34] = v[34]; x[35] = v[35];
403: v += 36;
404: }
405: row = *ajtmp++;
406: while (row < i) {
407: pc = rtmp + 36*row;
408: p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3];
409: p5 = pc[4]; p6 = pc[5]; p7 = pc[6]; p8 = pc[7];
410: p9 = pc[8]; p10 = pc[9]; p11 = pc[10]; p12 = pc[11];
411: p13 = pc[12]; p14 = pc[13]; p15 = pc[14]; p16 = pc[15];
412: p17 = pc[16]; p18 = pc[17]; p19 = pc[18]; p20 = pc[19];
413: p21 = pc[20]; p22 = pc[21]; p23 = pc[22]; p24 = pc[23];
414: p25 = pc[24]; p26 = pc[25]; p27 = pc[26]; p28 = pc[27];
415: p29 = pc[28]; p30 = pc[29]; p31 = pc[30]; p32 = pc[31];
416: p33 = pc[32]; p34 = pc[33]; p35 = pc[34]; p36 = pc[35];
417: if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 ||
418: p5 != 0.0 || p6 != 0.0 || p7 != 0.0 || p8 != 0.0 ||
419: p9 != 0.0 || p10 != 0.0 || p11 != 0.0 || p12 != 0.0 ||
420: p13 != 0.0 || p14 != 0.0 || p15 != 0.0 || p16 != 0.0 ||
421: p17 != 0.0 || p18 != 0.0 || p19 != 0.0 || p20 != 0.0 ||
422: p21 != 0.0 || p22 != 0.0 || p23 != 0.0 || p24 != 0.0 ||
423: p25 != 0.0 || p26 != 0.0 || p27 != 0.0 || p28 != 0.0 ||
424: p29 != 0.0 || p30 != 0.0 || p31 != 0.0 || p32 != 0.0 ||
425: p33 != 0.0 || p34 != 0.0 || p35 != 0.0 || p36 != 0.0) {
426: pv = ba + 36*diag_offset[row];
427: pj = bj + diag_offset[row] + 1;
428: x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
429: x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7];
430: x9 = pv[8]; x10 = pv[9]; x11 = pv[10]; x12 = pv[11];
431: x13 = pv[12]; x14 = pv[13]; x15 = pv[14]; x16 = pv[15];
432: x17 = pv[16]; x18 = pv[17]; x19 = pv[18]; x20 = pv[19];
433: x21 = pv[20]; x22 = pv[21]; x23 = pv[22]; x24 = pv[23];
434: x25 = pv[24]; x26 = pv[25]; x27 = pv[26]; x28 = pv[27];
435: x29 = pv[28]; x30 = pv[29]; x31 = pv[30]; x32 = pv[31];
436: x33 = pv[32]; x34 = pv[33]; x35 = pv[34]; x36 = pv[35];
437: pc[0] = m1 = p1*x1 + p7*x2 + p13*x3 + p19*x4 + p25*x5 + p31*x6;
438: pc[1] = m2 = p2*x1 + p8*x2 + p14*x3 + p20*x4 + p26*x5 + p32*x6;
439: pc[2] = m3 = p3*x1 + p9*x2 + p15*x3 + p21*x4 + p27*x5 + p33*x6;
440: pc[3] = m4 = p4*x1 + p10*x2 + p16*x3 + p22*x4 + p28*x5 + p34*x6;
441: pc[4] = m5 = p5*x1 + p11*x2 + p17*x3 + p23*x4 + p29*x5 + p35*x6;
442: pc[5] = m6 = p6*x1 + p12*x2 + p18*x3 + p24*x4 + p30*x5 + p36*x6;
444: pc[6] = m7 = p1*x7 + p7*x8 + p13*x9 + p19*x10 + p25*x11 + p31*x12;
445: pc[7] = m8 = p2*x7 + p8*x8 + p14*x9 + p20*x10 + p26*x11 + p32*x12;
446: pc[8] = m9 = p3*x7 + p9*x8 + p15*x9 + p21*x10 + p27*x11 + p33*x12;
447: pc[9] = m10 = p4*x7 + p10*x8 + p16*x9 + p22*x10 + p28*x11 + p34*x12;
448: pc[10] = m11 = p5*x7 + p11*x8 + p17*x9 + p23*x10 + p29*x11 + p35*x12;
449: pc[11] = m12 = p6*x7 + p12*x8 + p18*x9 + p24*x10 + p30*x11 + p36*x12;
451: pc[12] = m13 = p1*x13 + p7*x14 + p13*x15 + p19*x16 + p25*x17 + p31*x18;
452: pc[13] = m14 = p2*x13 + p8*x14 + p14*x15 + p20*x16 + p26*x17 + p32*x18;
453: pc[14] = m15 = p3*x13 + p9*x14 + p15*x15 + p21*x16 + p27*x17 + p33*x18;
454: pc[15] = m16 = p4*x13 + p10*x14 + p16*x15 + p22*x16 + p28*x17 + p34*x18;
455: pc[16] = m17 = p5*x13 + p11*x14 + p17*x15 + p23*x16 + p29*x17 + p35*x18;
456: pc[17] = m18 = p6*x13 + p12*x14 + p18*x15 + p24*x16 + p30*x17 + p36*x18;
458: pc[18] = m19 = p1*x19 + p7*x20 + p13*x21 + p19*x22 + p25*x23 + p31*x24;
459: pc[19] = m20 = p2*x19 + p8*x20 + p14*x21 + p20*x22 + p26*x23 + p32*x24;
460: pc[20] = m21 = p3*x19 + p9*x20 + p15*x21 + p21*x22 + p27*x23 + p33*x24;
461: pc[21] = m22 = p4*x19 + p10*x20 + p16*x21 + p22*x22 + p28*x23 + p34*x24;
462: pc[22] = m23 = p5*x19 + p11*x20 + p17*x21 + p23*x22 + p29*x23 + p35*x24;
463: pc[23] = m24 = p6*x19 + p12*x20 + p18*x21 + p24*x22 + p30*x23 + p36*x24;
465: pc[24] = m25 = p1*x25 + p7*x26 + p13*x27 + p19*x28 + p25*x29 + p31*x30;
466: pc[25] = m26 = p2*x25 + p8*x26 + p14*x27 + p20*x28 + p26*x29 + p32*x30;
467: pc[26] = m27 = p3*x25 + p9*x26 + p15*x27 + p21*x28 + p27*x29 + p33*x30;
468: pc[27] = m28 = p4*x25 + p10*x26 + p16*x27 + p22*x28 + p28*x29 + p34*x30;
469: pc[28] = m29 = p5*x25 + p11*x26 + p17*x27 + p23*x28 + p29*x29 + p35*x30;
470: pc[29] = m30 = p6*x25 + p12*x26 + p18*x27 + p24*x28 + p30*x29 + p36*x30;
472: pc[30] = m31 = p1*x31 + p7*x32 + p13*x33 + p19*x34 + p25*x35 + p31*x36;
473: pc[31] = m32 = p2*x31 + p8*x32 + p14*x33 + p20*x34 + p26*x35 + p32*x36;
474: pc[32] = m33 = p3*x31 + p9*x32 + p15*x33 + p21*x34 + p27*x35 + p33*x36;
475: pc[33] = m34 = p4*x31 + p10*x32 + p16*x33 + p22*x34 + p28*x35 + p34*x36;
476: pc[34] = m35 = p5*x31 + p11*x32 + p17*x33 + p23*x34 + p29*x35 + p35*x36;
477: pc[35] = m36 = p6*x31 + p12*x32 + p18*x33 + p24*x34 + p30*x35 + p36*x36;
479: nz = bi[row+1] - diag_offset[row] - 1;
480: pv += 36;
481: for (j=0; j<nz; j++) {
482: x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
483: x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7];
484: x9 = pv[8]; x10 = pv[9]; x11 = pv[10]; x12 = pv[11];
485: x13 = pv[12]; x14 = pv[13]; x15 = pv[14]; x16 = pv[15];
486: x17 = pv[16]; x18 = pv[17]; x19 = pv[18]; x20 = pv[19];
487: x21 = pv[20]; x22 = pv[21]; x23 = pv[22]; x24 = pv[23];
488: x25 = pv[24]; x26 = pv[25]; x27 = pv[26]; x28 = pv[27];
489: x29 = pv[28]; x30 = pv[29]; x31 = pv[30]; x32 = pv[31];
490: x33 = pv[32]; x34 = pv[33]; x35 = pv[34]; x36 = pv[35];
491: x = rtmp + 36*pj[j];
492: x[0] -= m1*x1 + m7*x2 + m13*x3 + m19*x4 + m25*x5 + m31*x6;
493: x[1] -= m2*x1 + m8*x2 + m14*x3 + m20*x4 + m26*x5 + m32*x6;
494: x[2] -= m3*x1 + m9*x2 + m15*x3 + m21*x4 + m27*x5 + m33*x6;
495: x[3] -= m4*x1 + m10*x2 + m16*x3 + m22*x4 + m28*x5 + m34*x6;
496: x[4] -= m5*x1 + m11*x2 + m17*x3 + m23*x4 + m29*x5 + m35*x6;
497: x[5] -= m6*x1 + m12*x2 + m18*x3 + m24*x4 + m30*x5 + m36*x6;
499: x[6] -= m1*x7 + m7*x8 + m13*x9 + m19*x10 + m25*x11 + m31*x12;
500: x[7] -= m2*x7 + m8*x8 + m14*x9 + m20*x10 + m26*x11 + m32*x12;
501: x[8] -= m3*x7 + m9*x8 + m15*x9 + m21*x10 + m27*x11 + m33*x12;
502: x[9] -= m4*x7 + m10*x8 + m16*x9 + m22*x10 + m28*x11 + m34*x12;
503: x[10] -= m5*x7 + m11*x8 + m17*x9 + m23*x10 + m29*x11 + m35*x12;
504: x[11] -= m6*x7 + m12*x8 + m18*x9 + m24*x10 + m30*x11 + m36*x12;
506: x[12] -= m1*x13 + m7*x14 + m13*x15 + m19*x16 + m25*x17 + m31*x18;
507: x[13] -= m2*x13 + m8*x14 + m14*x15 + m20*x16 + m26*x17 + m32*x18;
508: x[14] -= m3*x13 + m9*x14 + m15*x15 + m21*x16 + m27*x17 + m33*x18;
509: x[15] -= m4*x13 + m10*x14 + m16*x15 + m22*x16 + m28*x17 + m34*x18;
510: x[16] -= m5*x13 + m11*x14 + m17*x15 + m23*x16 + m29*x17 + m35*x18;
511: x[17] -= m6*x13 + m12*x14 + m18*x15 + m24*x16 + m30*x17 + m36*x18;
513: x[18] -= m1*x19 + m7*x20 + m13*x21 + m19*x22 + m25*x23 + m31*x24;
514: x[19] -= m2*x19 + m8*x20 + m14*x21 + m20*x22 + m26*x23 + m32*x24;
515: x[20] -= m3*x19 + m9*x20 + m15*x21 + m21*x22 + m27*x23 + m33*x24;
516: x[21] -= m4*x19 + m10*x20 + m16*x21 + m22*x22 + m28*x23 + m34*x24;
517: x[22] -= m5*x19 + m11*x20 + m17*x21 + m23*x22 + m29*x23 + m35*x24;
518: x[23] -= m6*x19 + m12*x20 + m18*x21 + m24*x22 + m30*x23 + m36*x24;
520: x[24] -= m1*x25 + m7*x26 + m13*x27 + m19*x28 + m25*x29 + m31*x30;
521: x[25] -= m2*x25 + m8*x26 + m14*x27 + m20*x28 + m26*x29 + m32*x30;
522: x[26] -= m3*x25 + m9*x26 + m15*x27 + m21*x28 + m27*x29 + m33*x30;
523: x[27] -= m4*x25 + m10*x26 + m16*x27 + m22*x28 + m28*x29 + m34*x30;
524: x[28] -= m5*x25 + m11*x26 + m17*x27 + m23*x28 + m29*x29 + m35*x30;
525: x[29] -= m6*x25 + m12*x26 + m18*x27 + m24*x28 + m30*x29 + m36*x30;
527: x[30] -= m1*x31 + m7*x32 + m13*x33 + m19*x34 + m25*x35 + m31*x36;
528: x[31] -= m2*x31 + m8*x32 + m14*x33 + m20*x34 + m26*x35 + m32*x36;
529: x[32] -= m3*x31 + m9*x32 + m15*x33 + m21*x34 + m27*x35 + m33*x36;
530: x[33] -= m4*x31 + m10*x32 + m16*x33 + m22*x34 + m28*x35 + m34*x36;
531: x[34] -= m5*x31 + m11*x32 + m17*x33 + m23*x34 + m29*x35 + m35*x36;
532: x[35] -= m6*x31 + m12*x32 + m18*x33 + m24*x34 + m30*x35 + m36*x36;
534: pv += 36;
535: }
536: PetscLogFlops(432.0*nz+396.0);
537: }
538: row = *ajtmp++;
539: }
540: /* finished row so stick it into b->a */
541: pv = ba + 36*bi[i];
542: pj = bj + bi[i];
543: nz = bi[i+1] - bi[i];
544: for (j=0; j<nz; j++) {
545: x = rtmp+36*pj[j];
546: pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3];
547: pv[4] = x[4]; pv[5] = x[5]; pv[6] = x[6]; pv[7] = x[7];
548: pv[8] = x[8]; pv[9] = x[9]; pv[10] = x[10]; pv[11] = x[11];
549: pv[12] = x[12]; pv[13] = x[13]; pv[14] = x[14]; pv[15] = x[15];
550: pv[16] = x[16]; pv[17] = x[17]; pv[18] = x[18]; pv[19] = x[19];
551: pv[20] = x[20]; pv[21] = x[21]; pv[22] = x[22]; pv[23] = x[23];
552: pv[24] = x[24]; pv[25] = x[25]; pv[26] = x[26]; pv[27] = x[27];
553: pv[28] = x[28]; pv[29] = x[29]; pv[30] = x[30]; pv[31] = x[31];
554: pv[32] = x[32]; pv[33] = x[33]; pv[34] = x[34]; pv[35] = x[35];
555: pv += 36;
556: }
557: /* invert diagonal block */
558: w = ba + 36*diag_offset[i];
559: PetscKernel_A_gets_inverse_A_6(w,shift,allowzeropivot,&zeropivotdetected);
560: if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
561: }
563: PetscFree(rtmp);
565: C->ops->solve = MatSolve_SeqBAIJ_6_NaturalOrdering_inplace;
566: C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_6_NaturalOrdering_inplace;
567: C->assembled = PETSC_TRUE;
569: PetscLogFlops(1.333333333333*6*6*6*b->mbs); /* from inverting diagonal blocks */
570: return(0);
571: }
573: PetscErrorCode MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info)
574: {
575: Mat C =B;
576: Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data;
578: PetscInt i,j,k,nz,nzL,row;
579: const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
580: const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag,*pj,bs2=a->bs2;
581: MatScalar *rtmp,*pc,*mwork,*v,*pv,*aa=a->a;
582: PetscInt flg;
583: PetscReal shift = info->shiftamount;
584: PetscBool allowzeropivot,zeropivotdetected;
587: allowzeropivot = PetscNot(A->erroriffailure);
589: /* generate work space needed by the factorization */
590: PetscMalloc2(bs2*n,&rtmp,bs2,&mwork);
591: PetscArrayzero(rtmp,bs2*n);
593: for (i=0; i<n; i++) {
594: /* zero rtmp */
595: /* L part */
596: nz = bi[i+1] - bi[i];
597: bjtmp = bj + bi[i];
598: for (j=0; j<nz; j++) {
599: PetscArrayzero(rtmp+bs2*bjtmp[j],bs2);
600: }
602: /* U part */
603: nz = bdiag[i] - bdiag[i+1];
604: bjtmp = bj + bdiag[i+1]+1;
605: for (j=0; j<nz; j++) {
606: PetscArrayzero(rtmp+bs2*bjtmp[j],bs2);
607: }
609: /* load in initial (unfactored row) */
610: nz = ai[i+1] - ai[i];
611: ajtmp = aj + ai[i];
612: v = aa + bs2*ai[i];
613: for (j=0; j<nz; j++) {
614: PetscArraycpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2);
615: }
617: /* elimination */
618: bjtmp = bj + bi[i];
619: nzL = bi[i+1] - bi[i];
620: for (k=0; k < nzL; k++) {
621: row = bjtmp[k];
622: pc = rtmp + bs2*row;
623: for (flg=0,j=0; j<bs2; j++) {
624: if (pc[j]!=0.0) {
625: flg = 1;
626: break;
627: }
628: }
629: if (flg) {
630: pv = b->a + bs2*bdiag[row];
631: /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
632: PetscKernel_A_gets_A_times_B_6(pc,pv,mwork);
634: pj = b->j + bdiag[row+1]+1; /* begining of U(row,:) */
635: pv = b->a + bs2*(bdiag[row+1]+1);
636: nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries inU(row,:), excluding diag */
637: for (j=0; j<nz; j++) {
638: /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
639: /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
640: v = rtmp + bs2*pj[j];
641: PetscKernel_A_gets_A_minus_B_times_C_6(v,pc,pv);
642: pv += bs2;
643: }
644: PetscLogFlops(432.0*nz+396); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
645: }
646: }
648: /* finished row so stick it into b->a */
649: /* L part */
650: pv = b->a + bs2*bi[i];
651: pj = b->j + bi[i];
652: nz = bi[i+1] - bi[i];
653: for (j=0; j<nz; j++) {
654: PetscArraycpy(pv+bs2*j,rtmp+bs2*pj[j],bs2);
655: }
657: /* Mark diagonal and invert diagonal for simplier triangular solves */
658: pv = b->a + bs2*bdiag[i];
659: pj = b->j + bdiag[i];
660: PetscArraycpy(pv,rtmp+bs2*pj[0],bs2);
661: PetscKernel_A_gets_inverse_A_6(pv,shift,allowzeropivot,&zeropivotdetected);
662: if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
664: /* U part */
665: pv = b->a + bs2*(bdiag[i+1]+1);
666: pj = b->j + bdiag[i+1]+1;
667: nz = bdiag[i] - bdiag[i+1] - 1;
668: for (j=0; j<nz; j++) {
669: PetscArraycpy(pv+bs2*j,rtmp+bs2*pj[j],bs2);
670: }
671: }
672: PetscFree2(rtmp,mwork);
674: C->ops->solve = MatSolve_SeqBAIJ_6_NaturalOrdering;
675: C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_6_NaturalOrdering;
676: C->assembled = PETSC_TRUE;
678: PetscLogFlops(1.333333333333*6*6*6*n); /* from inverting diagonal blocks */
679: return(0);
680: }