183287d42SBarry Smith /* 283287d42SBarry Smith Factorization code for BAIJ format. 383287d42SBarry Smith */ 4c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h> 5af0996ceSBarry Smith #include <petsc/private/kernels/blockinvert.h> 683287d42SBarry Smith 783287d42SBarry Smith /* 883287d42SBarry Smith Version for when blocks are 6 by 6 983287d42SBarry Smith */ 10*421480d9SBarry Smith PetscErrorCode MatILUFactorNumeric_SeqBAIJ_6_inplace(Mat C, Mat A, const MatFactorInfo *info) 11d71ae5a4SJacob Faibussowitsch { 1283287d42SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data; 1383287d42SBarry Smith IS isrow = b->row, isicol = b->icol; 145d0c19d7SBarry Smith const PetscInt *ajtmpold, *ajtmp, *diag_offset = b->diag, *r, *ic, *bi = b->i, *bj = b->j, *ai = a->i, *aj = a->j, *pj; 155d0c19d7SBarry Smith PetscInt nz, row, i, j, n = a->mbs, idx; 1683287d42SBarry Smith MatScalar *pv, *v, *rtmp, *pc, *w, *x; 1783287d42SBarry Smith MatScalar p1, p2, p3, p4, m1, m2, m3, m4, m5, m6, m7, m8, m9, x1, x2, x3, x4; 1883287d42SBarry Smith MatScalar p5, p6, p7, p8, p9, x5, x6, x7, x8, x9, x10, x11, x12, x13, x14, x15, x16; 1983287d42SBarry Smith MatScalar x17, x18, x19, x20, x21, x22, x23, x24, x25, p10, p11, p12, p13, p14; 2083287d42SBarry Smith MatScalar p15, p16, p17, p18, p19, p20, p21, p22, p23, p24, p25, m10, m11, m12; 2183287d42SBarry Smith MatScalar m13, m14, m15, m16, m17, m18, m19, m20, m21, m22, m23, m24, m25; 2283287d42SBarry Smith MatScalar p26, p27, p28, p29, p30, p31, p32, p33, p34, p35, p36; 2383287d42SBarry Smith MatScalar x26, x27, x28, x29, x30, x31, x32, x33, x34, x35, x36; 2483287d42SBarry Smith MatScalar m26, m27, m28, m29, m30, m31, m32, m33, m34, m35, m36; 2583287d42SBarry Smith MatScalar *ba = b->a, *aa = a->a; 26182b8fbaSHong Zhang PetscReal shift = info->shiftamount; 27a455e926SHong Zhang PetscBool allowzeropivot, zeropivotdetected; 2883287d42SBarry Smith 2983287d42SBarry Smith PetscFunctionBegin; 30*421480d9SBarry Smith /* Since A is C and C is labeled as a factored matrix we need to lie to MatGetDiagonalMarkers_SeqBAIJ() to get it to compute the diagonals */ 31*421480d9SBarry Smith A->factortype = MAT_FACTOR_NONE; 32*421480d9SBarry Smith PetscCall(MatGetDiagonalMarkers_SeqBAIJ(A, &diag_offset, NULL)); 33*421480d9SBarry Smith A->factortype = MAT_FACTOR_ILU; 340164db54SHong Zhang allowzeropivot = PetscNot(A->erroriffailure); 359566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isrow, &r)); 369566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isicol, &ic)); 379566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(36 * (n + 1), &rtmp)); 3883287d42SBarry Smith 3983287d42SBarry Smith for (i = 0; i < n; i++) { 4083287d42SBarry Smith nz = bi[i + 1] - bi[i]; 4183287d42SBarry Smith ajtmp = bj + bi[i]; 4283287d42SBarry Smith for (j = 0; j < nz; j++) { 4383287d42SBarry Smith x = rtmp + 36 * ajtmp[j]; 4483287d42SBarry Smith x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = x[9] = 0.0; 4583287d42SBarry Smith x[10] = x[11] = x[12] = x[13] = x[14] = x[15] = x[16] = x[17] = 0.0; 4683287d42SBarry Smith x[18] = x[19] = x[20] = x[21] = x[22] = x[23] = x[24] = x[25] = 0.0; 4783287d42SBarry Smith x[26] = x[27] = x[28] = x[29] = x[30] = x[31] = x[32] = x[33] = 0.0; 4883287d42SBarry Smith x[34] = x[35] = 0.0; 4983287d42SBarry Smith } 5083287d42SBarry Smith /* load in initial (unfactored row) */ 5183287d42SBarry Smith idx = r[i]; 5283287d42SBarry Smith nz = ai[idx + 1] - ai[idx]; 5383287d42SBarry Smith ajtmpold = aj + ai[idx]; 5483287d42SBarry Smith v = aa + 36 * ai[idx]; 5583287d42SBarry Smith for (j = 0; j < nz; j++) { 5683287d42SBarry Smith x = rtmp + 36 * ic[ajtmpold[j]]; 579371c9d4SSatish Balay x[0] = v[0]; 589371c9d4SSatish Balay x[1] = v[1]; 599371c9d4SSatish Balay x[2] = v[2]; 609371c9d4SSatish Balay x[3] = v[3]; 619371c9d4SSatish Balay x[4] = v[4]; 629371c9d4SSatish Balay x[5] = v[5]; 639371c9d4SSatish Balay x[6] = v[6]; 649371c9d4SSatish Balay x[7] = v[7]; 659371c9d4SSatish Balay x[8] = v[8]; 669371c9d4SSatish Balay x[9] = v[9]; 679371c9d4SSatish Balay x[10] = v[10]; 689371c9d4SSatish Balay x[11] = v[11]; 699371c9d4SSatish Balay x[12] = v[12]; 709371c9d4SSatish Balay x[13] = v[13]; 719371c9d4SSatish Balay x[14] = v[14]; 729371c9d4SSatish Balay x[15] = v[15]; 739371c9d4SSatish Balay x[16] = v[16]; 749371c9d4SSatish Balay x[17] = v[17]; 759371c9d4SSatish Balay x[18] = v[18]; 769371c9d4SSatish Balay x[19] = v[19]; 779371c9d4SSatish Balay x[20] = v[20]; 789371c9d4SSatish Balay x[21] = v[21]; 799371c9d4SSatish Balay x[22] = v[22]; 809371c9d4SSatish Balay x[23] = v[23]; 819371c9d4SSatish Balay x[24] = v[24]; 829371c9d4SSatish Balay x[25] = v[25]; 839371c9d4SSatish Balay x[26] = v[26]; 849371c9d4SSatish Balay x[27] = v[27]; 859371c9d4SSatish Balay x[28] = v[28]; 869371c9d4SSatish Balay x[29] = v[29]; 879371c9d4SSatish Balay x[30] = v[30]; 889371c9d4SSatish Balay x[31] = v[31]; 899371c9d4SSatish Balay x[32] = v[32]; 909371c9d4SSatish Balay x[33] = v[33]; 919371c9d4SSatish Balay x[34] = v[34]; 929371c9d4SSatish Balay x[35] = v[35]; 9383287d42SBarry Smith v += 36; 9483287d42SBarry Smith } 9583287d42SBarry Smith row = *ajtmp++; 9683287d42SBarry Smith while (row < i) { 9783287d42SBarry Smith pc = rtmp + 36 * row; 989371c9d4SSatish Balay p1 = pc[0]; 999371c9d4SSatish Balay p2 = pc[1]; 1009371c9d4SSatish Balay p3 = pc[2]; 1019371c9d4SSatish Balay p4 = pc[3]; 1029371c9d4SSatish Balay p5 = pc[4]; 1039371c9d4SSatish Balay p6 = pc[5]; 1049371c9d4SSatish Balay p7 = pc[6]; 1059371c9d4SSatish Balay p8 = pc[7]; 1069371c9d4SSatish Balay p9 = pc[8]; 1079371c9d4SSatish Balay p10 = pc[9]; 1089371c9d4SSatish Balay p11 = pc[10]; 1099371c9d4SSatish Balay p12 = pc[11]; 1109371c9d4SSatish Balay p13 = pc[12]; 1119371c9d4SSatish Balay p14 = pc[13]; 1129371c9d4SSatish Balay p15 = pc[14]; 1139371c9d4SSatish Balay p16 = pc[15]; 1149371c9d4SSatish Balay p17 = pc[16]; 1159371c9d4SSatish Balay p18 = pc[17]; 1169371c9d4SSatish Balay p19 = pc[18]; 1179371c9d4SSatish Balay p20 = pc[19]; 1189371c9d4SSatish Balay p21 = pc[20]; 1199371c9d4SSatish Balay p22 = pc[21]; 1209371c9d4SSatish Balay p23 = pc[22]; 1219371c9d4SSatish Balay p24 = pc[23]; 1229371c9d4SSatish Balay p25 = pc[24]; 1239371c9d4SSatish Balay p26 = pc[25]; 1249371c9d4SSatish Balay p27 = pc[26]; 1259371c9d4SSatish Balay p28 = pc[27]; 1269371c9d4SSatish Balay p29 = pc[28]; 1279371c9d4SSatish Balay p30 = pc[29]; 1289371c9d4SSatish Balay p31 = pc[30]; 1299371c9d4SSatish Balay p32 = pc[31]; 1309371c9d4SSatish Balay p33 = pc[32]; 1319371c9d4SSatish Balay p34 = pc[33]; 1329371c9d4SSatish Balay p35 = pc[34]; 1339371c9d4SSatish Balay p36 = pc[35]; 1349371c9d4SSatish Balay if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 || p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0 || p10 != 0.0 || p11 != 0.0 || p12 != 0.0 || p13 != 0.0 || p14 != 0.0 || p15 != 0.0 || p16 != 0.0 || p17 != 0.0 || p18 != 0.0 || p19 != 0.0 || p20 != 0.0 || p21 != 0.0 || p22 != 0.0 || p23 != 0.0 || p24 != 0.0 || p25 != 0.0 || p26 != 0.0 || p27 != 0.0 || p28 != 0.0 || p29 != 0.0 || p30 != 0.0 || p31 != 0.0 || p32 != 0.0 || p33 != 0.0 || p34 != 0.0 || p35 != 0.0 || p36 != 0.0) { 13583287d42SBarry Smith pv = ba + 36 * diag_offset[row]; 13683287d42SBarry Smith pj = bj + diag_offset[row] + 1; 1379371c9d4SSatish Balay x1 = pv[0]; 1389371c9d4SSatish Balay x2 = pv[1]; 1399371c9d4SSatish Balay x3 = pv[2]; 1409371c9d4SSatish Balay x4 = pv[3]; 1419371c9d4SSatish Balay x5 = pv[4]; 1429371c9d4SSatish Balay x6 = pv[5]; 1439371c9d4SSatish Balay x7 = pv[6]; 1449371c9d4SSatish Balay x8 = pv[7]; 1459371c9d4SSatish Balay x9 = pv[8]; 1469371c9d4SSatish Balay x10 = pv[9]; 1479371c9d4SSatish Balay x11 = pv[10]; 1489371c9d4SSatish Balay x12 = pv[11]; 1499371c9d4SSatish Balay x13 = pv[12]; 1509371c9d4SSatish Balay x14 = pv[13]; 1519371c9d4SSatish Balay x15 = pv[14]; 1529371c9d4SSatish Balay x16 = pv[15]; 1539371c9d4SSatish Balay x17 = pv[16]; 1549371c9d4SSatish Balay x18 = pv[17]; 1559371c9d4SSatish Balay x19 = pv[18]; 1569371c9d4SSatish Balay x20 = pv[19]; 1579371c9d4SSatish Balay x21 = pv[20]; 1589371c9d4SSatish Balay x22 = pv[21]; 1599371c9d4SSatish Balay x23 = pv[22]; 1609371c9d4SSatish Balay x24 = pv[23]; 1619371c9d4SSatish Balay x25 = pv[24]; 1629371c9d4SSatish Balay x26 = pv[25]; 1639371c9d4SSatish Balay x27 = pv[26]; 1649371c9d4SSatish Balay x28 = pv[27]; 1659371c9d4SSatish Balay x29 = pv[28]; 1669371c9d4SSatish Balay x30 = pv[29]; 1679371c9d4SSatish Balay x31 = pv[30]; 1689371c9d4SSatish Balay x32 = pv[31]; 1699371c9d4SSatish Balay x33 = pv[32]; 1709371c9d4SSatish Balay x34 = pv[33]; 1719371c9d4SSatish Balay x35 = pv[34]; 1729371c9d4SSatish Balay x36 = pv[35]; 17383287d42SBarry Smith pc[0] = m1 = p1 * x1 + p7 * x2 + p13 * x3 + p19 * x4 + p25 * x5 + p31 * x6; 17483287d42SBarry Smith pc[1] = m2 = p2 * x1 + p8 * x2 + p14 * x3 + p20 * x4 + p26 * x5 + p32 * x6; 17583287d42SBarry Smith pc[2] = m3 = p3 * x1 + p9 * x2 + p15 * x3 + p21 * x4 + p27 * x5 + p33 * x6; 17683287d42SBarry Smith pc[3] = m4 = p4 * x1 + p10 * x2 + p16 * x3 + p22 * x4 + p28 * x5 + p34 * x6; 17783287d42SBarry Smith pc[4] = m5 = p5 * x1 + p11 * x2 + p17 * x3 + p23 * x4 + p29 * x5 + p35 * x6; 17883287d42SBarry Smith pc[5] = m6 = p6 * x1 + p12 * x2 + p18 * x3 + p24 * x4 + p30 * x5 + p36 * x6; 17983287d42SBarry Smith 18083287d42SBarry Smith pc[6] = m7 = p1 * x7 + p7 * x8 + p13 * x9 + p19 * x10 + p25 * x11 + p31 * x12; 18183287d42SBarry Smith pc[7] = m8 = p2 * x7 + p8 * x8 + p14 * x9 + p20 * x10 + p26 * x11 + p32 * x12; 18283287d42SBarry Smith pc[8] = m9 = p3 * x7 + p9 * x8 + p15 * x9 + p21 * x10 + p27 * x11 + p33 * x12; 18383287d42SBarry Smith pc[9] = m10 = p4 * x7 + p10 * x8 + p16 * x9 + p22 * x10 + p28 * x11 + p34 * x12; 18483287d42SBarry Smith pc[10] = m11 = p5 * x7 + p11 * x8 + p17 * x9 + p23 * x10 + p29 * x11 + p35 * x12; 18583287d42SBarry Smith pc[11] = m12 = p6 * x7 + p12 * x8 + p18 * x9 + p24 * x10 + p30 * x11 + p36 * x12; 18683287d42SBarry Smith 18783287d42SBarry Smith pc[12] = m13 = p1 * x13 + p7 * x14 + p13 * x15 + p19 * x16 + p25 * x17 + p31 * x18; 18883287d42SBarry Smith pc[13] = m14 = p2 * x13 + p8 * x14 + p14 * x15 + p20 * x16 + p26 * x17 + p32 * x18; 18983287d42SBarry Smith pc[14] = m15 = p3 * x13 + p9 * x14 + p15 * x15 + p21 * x16 + p27 * x17 + p33 * x18; 19083287d42SBarry Smith pc[15] = m16 = p4 * x13 + p10 * x14 + p16 * x15 + p22 * x16 + p28 * x17 + p34 * x18; 19183287d42SBarry Smith pc[16] = m17 = p5 * x13 + p11 * x14 + p17 * x15 + p23 * x16 + p29 * x17 + p35 * x18; 19283287d42SBarry Smith pc[17] = m18 = p6 * x13 + p12 * x14 + p18 * x15 + p24 * x16 + p30 * x17 + p36 * x18; 19383287d42SBarry Smith 19483287d42SBarry Smith pc[18] = m19 = p1 * x19 + p7 * x20 + p13 * x21 + p19 * x22 + p25 * x23 + p31 * x24; 19583287d42SBarry Smith pc[19] = m20 = p2 * x19 + p8 * x20 + p14 * x21 + p20 * x22 + p26 * x23 + p32 * x24; 19683287d42SBarry Smith pc[20] = m21 = p3 * x19 + p9 * x20 + p15 * x21 + p21 * x22 + p27 * x23 + p33 * x24; 19783287d42SBarry Smith pc[21] = m22 = p4 * x19 + p10 * x20 + p16 * x21 + p22 * x22 + p28 * x23 + p34 * x24; 19883287d42SBarry Smith pc[22] = m23 = p5 * x19 + p11 * x20 + p17 * x21 + p23 * x22 + p29 * x23 + p35 * x24; 19983287d42SBarry Smith pc[23] = m24 = p6 * x19 + p12 * x20 + p18 * x21 + p24 * x22 + p30 * x23 + p36 * x24; 20083287d42SBarry Smith 20183287d42SBarry Smith pc[24] = m25 = p1 * x25 + p7 * x26 + p13 * x27 + p19 * x28 + p25 * x29 + p31 * x30; 20283287d42SBarry Smith pc[25] = m26 = p2 * x25 + p8 * x26 + p14 * x27 + p20 * x28 + p26 * x29 + p32 * x30; 20383287d42SBarry Smith pc[26] = m27 = p3 * x25 + p9 * x26 + p15 * x27 + p21 * x28 + p27 * x29 + p33 * x30; 20483287d42SBarry Smith pc[27] = m28 = p4 * x25 + p10 * x26 + p16 * x27 + p22 * x28 + p28 * x29 + p34 * x30; 20583287d42SBarry Smith pc[28] = m29 = p5 * x25 + p11 * x26 + p17 * x27 + p23 * x28 + p29 * x29 + p35 * x30; 20683287d42SBarry Smith pc[29] = m30 = p6 * x25 + p12 * x26 + p18 * x27 + p24 * x28 + p30 * x29 + p36 * x30; 20783287d42SBarry Smith 20883287d42SBarry Smith pc[30] = m31 = p1 * x31 + p7 * x32 + p13 * x33 + p19 * x34 + p25 * x35 + p31 * x36; 20983287d42SBarry Smith pc[31] = m32 = p2 * x31 + p8 * x32 + p14 * x33 + p20 * x34 + p26 * x35 + p32 * x36; 21083287d42SBarry Smith pc[32] = m33 = p3 * x31 + p9 * x32 + p15 * x33 + p21 * x34 + p27 * x35 + p33 * x36; 21183287d42SBarry Smith pc[33] = m34 = p4 * x31 + p10 * x32 + p16 * x33 + p22 * x34 + p28 * x35 + p34 * x36; 21283287d42SBarry Smith pc[34] = m35 = p5 * x31 + p11 * x32 + p17 * x33 + p23 * x34 + p29 * x35 + p35 * x36; 21383287d42SBarry Smith pc[35] = m36 = p6 * x31 + p12 * x32 + p18 * x33 + p24 * x34 + p30 * x35 + p36 * x36; 21483287d42SBarry Smith 21583287d42SBarry Smith nz = bi[row + 1] - diag_offset[row] - 1; 21683287d42SBarry Smith pv += 36; 21783287d42SBarry Smith for (j = 0; j < nz; j++) { 2189371c9d4SSatish Balay x1 = pv[0]; 2199371c9d4SSatish Balay x2 = pv[1]; 2209371c9d4SSatish Balay x3 = pv[2]; 2219371c9d4SSatish Balay x4 = pv[3]; 2229371c9d4SSatish Balay x5 = pv[4]; 2239371c9d4SSatish Balay x6 = pv[5]; 2249371c9d4SSatish Balay x7 = pv[6]; 2259371c9d4SSatish Balay x8 = pv[7]; 2269371c9d4SSatish Balay x9 = pv[8]; 2279371c9d4SSatish Balay x10 = pv[9]; 2289371c9d4SSatish Balay x11 = pv[10]; 2299371c9d4SSatish Balay x12 = pv[11]; 2309371c9d4SSatish Balay x13 = pv[12]; 2319371c9d4SSatish Balay x14 = pv[13]; 2329371c9d4SSatish Balay x15 = pv[14]; 2339371c9d4SSatish Balay x16 = pv[15]; 2349371c9d4SSatish Balay x17 = pv[16]; 2359371c9d4SSatish Balay x18 = pv[17]; 2369371c9d4SSatish Balay x19 = pv[18]; 2379371c9d4SSatish Balay x20 = pv[19]; 2389371c9d4SSatish Balay x21 = pv[20]; 2399371c9d4SSatish Balay x22 = pv[21]; 2409371c9d4SSatish Balay x23 = pv[22]; 2419371c9d4SSatish Balay x24 = pv[23]; 2429371c9d4SSatish Balay x25 = pv[24]; 2439371c9d4SSatish Balay x26 = pv[25]; 2449371c9d4SSatish Balay x27 = pv[26]; 2459371c9d4SSatish Balay x28 = pv[27]; 2469371c9d4SSatish Balay x29 = pv[28]; 2479371c9d4SSatish Balay x30 = pv[29]; 2489371c9d4SSatish Balay x31 = pv[30]; 2499371c9d4SSatish Balay x32 = pv[31]; 2509371c9d4SSatish Balay x33 = pv[32]; 2519371c9d4SSatish Balay x34 = pv[33]; 2529371c9d4SSatish Balay x35 = pv[34]; 2539371c9d4SSatish Balay x36 = pv[35]; 25483287d42SBarry Smith x = rtmp + 36 * pj[j]; 25583287d42SBarry Smith x[0] -= m1 * x1 + m7 * x2 + m13 * x3 + m19 * x4 + m25 * x5 + m31 * x6; 25683287d42SBarry Smith x[1] -= m2 * x1 + m8 * x2 + m14 * x3 + m20 * x4 + m26 * x5 + m32 * x6; 25783287d42SBarry Smith x[2] -= m3 * x1 + m9 * x2 + m15 * x3 + m21 * x4 + m27 * x5 + m33 * x6; 25883287d42SBarry Smith x[3] -= m4 * x1 + m10 * x2 + m16 * x3 + m22 * x4 + m28 * x5 + m34 * x6; 25983287d42SBarry Smith x[4] -= m5 * x1 + m11 * x2 + m17 * x3 + m23 * x4 + m29 * x5 + m35 * x6; 26083287d42SBarry Smith x[5] -= m6 * x1 + m12 * x2 + m18 * x3 + m24 * x4 + m30 * x5 + m36 * x6; 26183287d42SBarry Smith 26283287d42SBarry Smith x[6] -= m1 * x7 + m7 * x8 + m13 * x9 + m19 * x10 + m25 * x11 + m31 * x12; 26383287d42SBarry Smith x[7] -= m2 * x7 + m8 * x8 + m14 * x9 + m20 * x10 + m26 * x11 + m32 * x12; 26483287d42SBarry Smith x[8] -= m3 * x7 + m9 * x8 + m15 * x9 + m21 * x10 + m27 * x11 + m33 * x12; 26583287d42SBarry Smith x[9] -= m4 * x7 + m10 * x8 + m16 * x9 + m22 * x10 + m28 * x11 + m34 * x12; 26683287d42SBarry Smith x[10] -= m5 * x7 + m11 * x8 + m17 * x9 + m23 * x10 + m29 * x11 + m35 * x12; 26783287d42SBarry Smith x[11] -= m6 * x7 + m12 * x8 + m18 * x9 + m24 * x10 + m30 * x11 + m36 * x12; 26883287d42SBarry Smith 26983287d42SBarry Smith x[12] -= m1 * x13 + m7 * x14 + m13 * x15 + m19 * x16 + m25 * x17 + m31 * x18; 27083287d42SBarry Smith x[13] -= m2 * x13 + m8 * x14 + m14 * x15 + m20 * x16 + m26 * x17 + m32 * x18; 27183287d42SBarry Smith x[14] -= m3 * x13 + m9 * x14 + m15 * x15 + m21 * x16 + m27 * x17 + m33 * x18; 27283287d42SBarry Smith x[15] -= m4 * x13 + m10 * x14 + m16 * x15 + m22 * x16 + m28 * x17 + m34 * x18; 27383287d42SBarry Smith x[16] -= m5 * x13 + m11 * x14 + m17 * x15 + m23 * x16 + m29 * x17 + m35 * x18; 27483287d42SBarry Smith x[17] -= m6 * x13 + m12 * x14 + m18 * x15 + m24 * x16 + m30 * x17 + m36 * x18; 27583287d42SBarry Smith 27683287d42SBarry Smith x[18] -= m1 * x19 + m7 * x20 + m13 * x21 + m19 * x22 + m25 * x23 + m31 * x24; 27783287d42SBarry Smith x[19] -= m2 * x19 + m8 * x20 + m14 * x21 + m20 * x22 + m26 * x23 + m32 * x24; 27883287d42SBarry Smith x[20] -= m3 * x19 + m9 * x20 + m15 * x21 + m21 * x22 + m27 * x23 + m33 * x24; 27983287d42SBarry Smith x[21] -= m4 * x19 + m10 * x20 + m16 * x21 + m22 * x22 + m28 * x23 + m34 * x24; 28083287d42SBarry Smith x[22] -= m5 * x19 + m11 * x20 + m17 * x21 + m23 * x22 + m29 * x23 + m35 * x24; 28183287d42SBarry Smith x[23] -= m6 * x19 + m12 * x20 + m18 * x21 + m24 * x22 + m30 * x23 + m36 * x24; 28283287d42SBarry Smith 28383287d42SBarry Smith x[24] -= m1 * x25 + m7 * x26 + m13 * x27 + m19 * x28 + m25 * x29 + m31 * x30; 28483287d42SBarry Smith x[25] -= m2 * x25 + m8 * x26 + m14 * x27 + m20 * x28 + m26 * x29 + m32 * x30; 28583287d42SBarry Smith x[26] -= m3 * x25 + m9 * x26 + m15 * x27 + m21 * x28 + m27 * x29 + m33 * x30; 28683287d42SBarry Smith x[27] -= m4 * x25 + m10 * x26 + m16 * x27 + m22 * x28 + m28 * x29 + m34 * x30; 28783287d42SBarry Smith x[28] -= m5 * x25 + m11 * x26 + m17 * x27 + m23 * x28 + m29 * x29 + m35 * x30; 28883287d42SBarry Smith x[29] -= m6 * x25 + m12 * x26 + m18 * x27 + m24 * x28 + m30 * x29 + m36 * x30; 28983287d42SBarry Smith 29083287d42SBarry Smith x[30] -= m1 * x31 + m7 * x32 + m13 * x33 + m19 * x34 + m25 * x35 + m31 * x36; 29183287d42SBarry Smith x[31] -= m2 * x31 + m8 * x32 + m14 * x33 + m20 * x34 + m26 * x35 + m32 * x36; 29283287d42SBarry Smith x[32] -= m3 * x31 + m9 * x32 + m15 * x33 + m21 * x34 + m27 * x35 + m33 * x36; 29383287d42SBarry Smith x[33] -= m4 * x31 + m10 * x32 + m16 * x33 + m22 * x34 + m28 * x35 + m34 * x36; 29483287d42SBarry Smith x[34] -= m5 * x31 + m11 * x32 + m17 * x33 + m23 * x34 + m29 * x35 + m35 * x36; 29583287d42SBarry Smith x[35] -= m6 * x31 + m12 * x32 + m18 * x33 + m24 * x34 + m30 * x35 + m36 * x36; 29683287d42SBarry Smith 29783287d42SBarry Smith pv += 36; 29883287d42SBarry Smith } 2999566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(432.0 * nz + 396.0)); 30083287d42SBarry Smith } 30183287d42SBarry Smith row = *ajtmp++; 30283287d42SBarry Smith } 30383287d42SBarry Smith /* finished row so stick it into b->a */ 30483287d42SBarry Smith pv = ba + 36 * bi[i]; 30583287d42SBarry Smith pj = bj + bi[i]; 30683287d42SBarry Smith nz = bi[i + 1] - bi[i]; 30783287d42SBarry Smith for (j = 0; j < nz; j++) { 30883287d42SBarry Smith x = rtmp + 36 * pj[j]; 3099371c9d4SSatish Balay pv[0] = x[0]; 3109371c9d4SSatish Balay pv[1] = x[1]; 3119371c9d4SSatish Balay pv[2] = x[2]; 3129371c9d4SSatish Balay pv[3] = x[3]; 3139371c9d4SSatish Balay pv[4] = x[4]; 3149371c9d4SSatish Balay pv[5] = x[5]; 3159371c9d4SSatish Balay pv[6] = x[6]; 3169371c9d4SSatish Balay pv[7] = x[7]; 3179371c9d4SSatish Balay pv[8] = x[8]; 3189371c9d4SSatish Balay pv[9] = x[9]; 3199371c9d4SSatish Balay pv[10] = x[10]; 3209371c9d4SSatish Balay pv[11] = x[11]; 3219371c9d4SSatish Balay pv[12] = x[12]; 3229371c9d4SSatish Balay pv[13] = x[13]; 3239371c9d4SSatish Balay pv[14] = x[14]; 3249371c9d4SSatish Balay pv[15] = x[15]; 3259371c9d4SSatish Balay pv[16] = x[16]; 3269371c9d4SSatish Balay pv[17] = x[17]; 3279371c9d4SSatish Balay pv[18] = x[18]; 3289371c9d4SSatish Balay pv[19] = x[19]; 3299371c9d4SSatish Balay pv[20] = x[20]; 3309371c9d4SSatish Balay pv[21] = x[21]; 3319371c9d4SSatish Balay pv[22] = x[22]; 3329371c9d4SSatish Balay pv[23] = x[23]; 3339371c9d4SSatish Balay pv[24] = x[24]; 3349371c9d4SSatish Balay pv[25] = x[25]; 3359371c9d4SSatish Balay pv[26] = x[26]; 3369371c9d4SSatish Balay pv[27] = x[27]; 3379371c9d4SSatish Balay pv[28] = x[28]; 3389371c9d4SSatish Balay pv[29] = x[29]; 3399371c9d4SSatish Balay pv[30] = x[30]; 3409371c9d4SSatish Balay pv[31] = x[31]; 3419371c9d4SSatish Balay pv[32] = x[32]; 3429371c9d4SSatish Balay pv[33] = x[33]; 3439371c9d4SSatish Balay pv[34] = x[34]; 3449371c9d4SSatish Balay pv[35] = x[35]; 34583287d42SBarry Smith pv += 36; 34683287d42SBarry Smith } 34783287d42SBarry Smith /* invert diagonal block */ 34883287d42SBarry Smith w = ba + 36 * diag_offset[i]; 3499566063dSJacob Faibussowitsch PetscCall(PetscKernel_A_gets_inverse_A_6(w, shift, allowzeropivot, &zeropivotdetected)); 3507b6c816cSBarry Smith if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 35183287d42SBarry Smith } 35283287d42SBarry Smith 3539566063dSJacob Faibussowitsch PetscCall(PetscFree(rtmp)); 3549566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isicol, &ic)); 3559566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrow, &r)); 35626fbe8dcSKarl Rupp 35706e38f1dSHong Zhang C->ops->solve = MatSolve_SeqBAIJ_6_inplace; 35806e38f1dSHong Zhang C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_6_inplace; 35983287d42SBarry Smith C->assembled = PETSC_TRUE; 36026fbe8dcSKarl Rupp 3619566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(1.333333333333 * 6 * 6 * 6 * b->mbs)); /* from inverting diagonal blocks */ 3623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 36383287d42SBarry Smith } 364bef36659SShri Abhyankar 365d71ae5a4SJacob Faibussowitsch PetscErrorCode MatLUFactorNumeric_SeqBAIJ_6(Mat B, Mat A, const MatFactorInfo *info) 366d71ae5a4SJacob Faibussowitsch { 367bef36659SShri Abhyankar Mat C = B; 368bef36659SShri Abhyankar Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data; 369bef36659SShri Abhyankar IS isrow = b->row, isicol = b->icol; 3705a586d82SBarry Smith const PetscInt *r, *ic; 371bbd65245SShri Abhyankar PetscInt i, j, k, nz, nzL, row; 372bbd65245SShri Abhyankar const PetscInt n = a->mbs, *ai = a->i, *aj = a->j, *bi = b->i, *bj = b->j; 373bbd65245SShri Abhyankar const PetscInt *ajtmp, *bjtmp, *bdiag = b->diag, *pj, bs2 = a->bs2; 374bef36659SShri Abhyankar MatScalar *rtmp, *pc, *mwork, *v, *pv, *aa = a->a; 375bbd65245SShri Abhyankar PetscInt flg; 376182b8fbaSHong Zhang PetscReal shift = info->shiftamount; 377a455e926SHong Zhang PetscBool allowzeropivot, zeropivotdetected; 378bef36659SShri Abhyankar 379bef36659SShri Abhyankar PetscFunctionBegin; 3800164db54SHong Zhang allowzeropivot = PetscNot(A->erroriffailure); 3819566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isrow, &r)); 3829566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isicol, &ic)); 383bef36659SShri Abhyankar 384bef36659SShri Abhyankar /* generate work space needed by the factorization */ 3859566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(bs2 * n, &rtmp, bs2, &mwork)); 3869566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(rtmp, bs2 * n)); 387bef36659SShri Abhyankar 388bef36659SShri Abhyankar for (i = 0; i < n; i++) { 389bef36659SShri Abhyankar /* zero rtmp */ 390bef36659SShri Abhyankar /* L part */ 391bef36659SShri Abhyankar nz = bi[i + 1] - bi[i]; 392bef36659SShri Abhyankar bjtmp = bj + bi[i]; 39348a46eb9SPierre Jolivet for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2)); 394bef36659SShri Abhyankar 395bef36659SShri Abhyankar /* U part */ 3966506fda5SShri Abhyankar nz = bdiag[i] - bdiag[i + 1]; 3976506fda5SShri Abhyankar bjtmp = bj + bdiag[i + 1] + 1; 39848a46eb9SPierre Jolivet for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2)); 3996506fda5SShri Abhyankar 4006506fda5SShri Abhyankar /* load in initial (unfactored row) */ 4016506fda5SShri Abhyankar nz = ai[r[i] + 1] - ai[r[i]]; 4026506fda5SShri Abhyankar ajtmp = aj + ai[r[i]]; 4036506fda5SShri Abhyankar v = aa + bs2 * ai[r[i]]; 40448a46eb9SPierre Jolivet for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(rtmp + bs2 * ic[ajtmp[j]], v + bs2 * j, bs2)); 4056506fda5SShri Abhyankar 4066506fda5SShri Abhyankar /* elimination */ 4076506fda5SShri Abhyankar bjtmp = bj + bi[i]; 4086506fda5SShri Abhyankar nzL = bi[i + 1] - bi[i]; 4096506fda5SShri Abhyankar for (k = 0; k < nzL; k++) { 4106506fda5SShri Abhyankar row = bjtmp[k]; 4116506fda5SShri Abhyankar pc = rtmp + bs2 * row; 412c35f09e5SBarry Smith for (flg = 0, j = 0; j < bs2; j++) { 413c35f09e5SBarry Smith if (pc[j] != 0.0) { 414c35f09e5SBarry Smith flg = 1; 415c35f09e5SBarry Smith break; 416c35f09e5SBarry Smith } 417c35f09e5SBarry Smith } 4186506fda5SShri Abhyankar if (flg) { 4196506fda5SShri Abhyankar pv = b->a + bs2 * bdiag[row]; 42096b95a6bSBarry Smith /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */ 4219566063dSJacob Faibussowitsch PetscCall(PetscKernel_A_gets_A_times_B_6(pc, pv, mwork)); 4226506fda5SShri Abhyankar 423a5b23f4aSJose E. Roman pj = b->j + bdiag[row + 1] + 1; /* beginning of U(row,:) */ 4246506fda5SShri Abhyankar pv = b->a + bs2 * (bdiag[row + 1] + 1); 4256506fda5SShri Abhyankar nz = bdiag[row] - bdiag[row + 1] - 1; /* num of entries inU(row,:), excluding diag */ 4266506fda5SShri Abhyankar for (j = 0; j < nz; j++) { 42796b95a6bSBarry Smith /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */ 4286506fda5SShri Abhyankar /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */ 4296506fda5SShri Abhyankar v = rtmp + bs2 * pj[j]; 4309566063dSJacob Faibussowitsch PetscCall(PetscKernel_A_gets_A_minus_B_times_C_6(v, pc, pv)); 4316506fda5SShri Abhyankar pv += bs2; 4326506fda5SShri Abhyankar } 4339566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(432.0 * nz + 396)); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */ 4346506fda5SShri Abhyankar } 4356506fda5SShri Abhyankar } 4366506fda5SShri Abhyankar 4376506fda5SShri Abhyankar /* finished row so stick it into b->a */ 4386506fda5SShri Abhyankar /* L part */ 4396506fda5SShri Abhyankar pv = b->a + bs2 * bi[i]; 4406506fda5SShri Abhyankar pj = b->j + bi[i]; 4416506fda5SShri Abhyankar nz = bi[i + 1] - bi[i]; 44248a46eb9SPierre Jolivet for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2)); 4436506fda5SShri Abhyankar 444a5b23f4aSJose E. Roman /* Mark diagonal and invert diagonal for simpler triangular solves */ 4456506fda5SShri Abhyankar pv = b->a + bs2 * bdiag[i]; 4466506fda5SShri Abhyankar pj = b->j + bdiag[i]; 4479566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(pv, rtmp + bs2 * pj[0], bs2)); 4489566063dSJacob Faibussowitsch PetscCall(PetscKernel_A_gets_inverse_A_6(pv, shift, allowzeropivot, &zeropivotdetected)); 4497b6c816cSBarry Smith if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 4506506fda5SShri Abhyankar 4516506fda5SShri Abhyankar /* U part */ 4526506fda5SShri Abhyankar pv = b->a + bs2 * (bdiag[i + 1] + 1); 4536506fda5SShri Abhyankar pj = b->j + bdiag[i + 1] + 1; 4546506fda5SShri Abhyankar nz = bdiag[i] - bdiag[i + 1] - 1; 45548a46eb9SPierre Jolivet for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2)); 4566506fda5SShri Abhyankar } 4576506fda5SShri Abhyankar 4589566063dSJacob Faibussowitsch PetscCall(PetscFree2(rtmp, mwork)); 4599566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isicol, &ic)); 4609566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrow, &r)); 46126fbe8dcSKarl Rupp 4624dd39f65SShri Abhyankar C->ops->solve = MatSolve_SeqBAIJ_6; 4634dd39f65SShri Abhyankar C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_6; 4646506fda5SShri Abhyankar C->assembled = PETSC_TRUE; 46526fbe8dcSKarl Rupp 4669566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(1.333333333333 * 6 * 6 * 6 * n)); /* from inverting diagonal blocks */ 4673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4686506fda5SShri Abhyankar } 4696506fda5SShri Abhyankar 470*421480d9SBarry Smith PetscErrorCode MatILUFactorNumeric_SeqBAIJ_6_NaturalOrdering_inplace(Mat C, Mat A, const MatFactorInfo *info) 471d71ae5a4SJacob Faibussowitsch { 472bef36659SShri Abhyankar Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data; 473bef36659SShri Abhyankar PetscInt i, j, n = a->mbs, *bi = b->i, *bj = b->j; 474bef36659SShri Abhyankar PetscInt *ajtmpold, *ajtmp, nz, row; 475*421480d9SBarry Smith PetscInt *ai = a->i, *aj = a->j, *pj; 476*421480d9SBarry Smith const PetscInt *diag_offset; 477bef36659SShri Abhyankar MatScalar *pv, *v, *rtmp, *pc, *w, *x; 478bef36659SShri Abhyankar MatScalar x1, x2, x3, x4, x5, x6, x7, x8, x9, x10, x11, x12, x13, x14, x15; 479bef36659SShri Abhyankar MatScalar x16, x17, x18, x19, x20, x21, x22, x23, x24, x25; 480bef36659SShri Abhyankar MatScalar p1, p2, p3, p4, p5, p6, p7, p8, p9, p10, p11, p12, p13, p14, p15; 481bef36659SShri Abhyankar MatScalar p16, p17, p18, p19, p20, p21, p22, p23, p24, p25; 482bef36659SShri Abhyankar MatScalar m1, m2, m3, m4, m5, m6, m7, m8, m9, m10, m11, m12, m13, m14, m15; 483bef36659SShri Abhyankar MatScalar m16, m17, m18, m19, m20, m21, m22, m23, m24, m25; 484bef36659SShri Abhyankar MatScalar p26, p27, p28, p29, p30, p31, p32, p33, p34, p35, p36; 485bef36659SShri Abhyankar MatScalar x26, x27, x28, x29, x30, x31, x32, x33, x34, x35, x36; 486bef36659SShri Abhyankar MatScalar m26, m27, m28, m29, m30, m31, m32, m33, m34, m35, m36; 487bef36659SShri Abhyankar MatScalar *ba = b->a, *aa = a->a; 488182b8fbaSHong Zhang PetscReal shift = info->shiftamount; 489a455e926SHong Zhang PetscBool allowzeropivot, zeropivotdetected; 490bef36659SShri Abhyankar 491bef36659SShri Abhyankar PetscFunctionBegin; 492*421480d9SBarry Smith /* Since A is C and C is labeled as a factored matrix we need to lie to MatGetDiagonalMarkers_SeqBAIJ() to get it to compute the diagonals */ 493*421480d9SBarry Smith A->factortype = MAT_FACTOR_NONE; 494*421480d9SBarry Smith PetscCall(MatGetDiagonalMarkers_SeqBAIJ(A, &diag_offset, NULL)); 495*421480d9SBarry Smith A->factortype = MAT_FACTOR_ILU; 4960164db54SHong Zhang allowzeropivot = PetscNot(A->erroriffailure); 4979566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(36 * (n + 1), &rtmp)); 498bef36659SShri Abhyankar for (i = 0; i < n; i++) { 499bef36659SShri Abhyankar nz = bi[i + 1] - bi[i]; 500bef36659SShri Abhyankar ajtmp = bj + bi[i]; 501bef36659SShri Abhyankar for (j = 0; j < nz; j++) { 502bef36659SShri Abhyankar x = rtmp + 36 * ajtmp[j]; 503bef36659SShri Abhyankar x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = x[9] = 0.0; 504bef36659SShri Abhyankar x[10] = x[11] = x[12] = x[13] = x[14] = x[15] = x[16] = x[17] = 0.0; 505bef36659SShri Abhyankar x[18] = x[19] = x[20] = x[21] = x[22] = x[23] = x[24] = x[25] = 0.0; 506bef36659SShri Abhyankar x[26] = x[27] = x[28] = x[29] = x[30] = x[31] = x[32] = x[33] = 0.0; 507bef36659SShri Abhyankar x[34] = x[35] = 0.0; 508bef36659SShri Abhyankar } 509bef36659SShri Abhyankar /* load in initial (unfactored row) */ 510bef36659SShri Abhyankar nz = ai[i + 1] - ai[i]; 511bef36659SShri Abhyankar ajtmpold = aj + ai[i]; 512bef36659SShri Abhyankar v = aa + 36 * ai[i]; 513bef36659SShri Abhyankar for (j = 0; j < nz; j++) { 514bef36659SShri Abhyankar x = rtmp + 36 * ajtmpold[j]; 5159371c9d4SSatish Balay x[0] = v[0]; 5169371c9d4SSatish Balay x[1] = v[1]; 5179371c9d4SSatish Balay x[2] = v[2]; 5189371c9d4SSatish Balay x[3] = v[3]; 5199371c9d4SSatish Balay x[4] = v[4]; 5209371c9d4SSatish Balay x[5] = v[5]; 5219371c9d4SSatish Balay x[6] = v[6]; 5229371c9d4SSatish Balay x[7] = v[7]; 5239371c9d4SSatish Balay x[8] = v[8]; 5249371c9d4SSatish Balay x[9] = v[9]; 5259371c9d4SSatish Balay x[10] = v[10]; 5269371c9d4SSatish Balay x[11] = v[11]; 5279371c9d4SSatish Balay x[12] = v[12]; 5289371c9d4SSatish Balay x[13] = v[13]; 5299371c9d4SSatish Balay x[14] = v[14]; 5309371c9d4SSatish Balay x[15] = v[15]; 5319371c9d4SSatish Balay x[16] = v[16]; 5329371c9d4SSatish Balay x[17] = v[17]; 5339371c9d4SSatish Balay x[18] = v[18]; 5349371c9d4SSatish Balay x[19] = v[19]; 5359371c9d4SSatish Balay x[20] = v[20]; 5369371c9d4SSatish Balay x[21] = v[21]; 5379371c9d4SSatish Balay x[22] = v[22]; 5389371c9d4SSatish Balay x[23] = v[23]; 5399371c9d4SSatish Balay x[24] = v[24]; 5409371c9d4SSatish Balay x[25] = v[25]; 5419371c9d4SSatish Balay x[26] = v[26]; 5429371c9d4SSatish Balay x[27] = v[27]; 5439371c9d4SSatish Balay x[28] = v[28]; 5449371c9d4SSatish Balay x[29] = v[29]; 5459371c9d4SSatish Balay x[30] = v[30]; 5469371c9d4SSatish Balay x[31] = v[31]; 5479371c9d4SSatish Balay x[32] = v[32]; 5489371c9d4SSatish Balay x[33] = v[33]; 5499371c9d4SSatish Balay x[34] = v[34]; 5509371c9d4SSatish Balay x[35] = v[35]; 551bef36659SShri Abhyankar v += 36; 552bef36659SShri Abhyankar } 553bef36659SShri Abhyankar row = *ajtmp++; 554bef36659SShri Abhyankar while (row < i) { 555bef36659SShri Abhyankar pc = rtmp + 36 * row; 5569371c9d4SSatish Balay p1 = pc[0]; 5579371c9d4SSatish Balay p2 = pc[1]; 5589371c9d4SSatish Balay p3 = pc[2]; 5599371c9d4SSatish Balay p4 = pc[3]; 5609371c9d4SSatish Balay p5 = pc[4]; 5619371c9d4SSatish Balay p6 = pc[5]; 5629371c9d4SSatish Balay p7 = pc[6]; 5639371c9d4SSatish Balay p8 = pc[7]; 5649371c9d4SSatish Balay p9 = pc[8]; 5659371c9d4SSatish Balay p10 = pc[9]; 5669371c9d4SSatish Balay p11 = pc[10]; 5679371c9d4SSatish Balay p12 = pc[11]; 5689371c9d4SSatish Balay p13 = pc[12]; 5699371c9d4SSatish Balay p14 = pc[13]; 5709371c9d4SSatish Balay p15 = pc[14]; 5719371c9d4SSatish Balay p16 = pc[15]; 5729371c9d4SSatish Balay p17 = pc[16]; 5739371c9d4SSatish Balay p18 = pc[17]; 5749371c9d4SSatish Balay p19 = pc[18]; 5759371c9d4SSatish Balay p20 = pc[19]; 5769371c9d4SSatish Balay p21 = pc[20]; 5779371c9d4SSatish Balay p22 = pc[21]; 5789371c9d4SSatish Balay p23 = pc[22]; 5799371c9d4SSatish Balay p24 = pc[23]; 5809371c9d4SSatish Balay p25 = pc[24]; 5819371c9d4SSatish Balay p26 = pc[25]; 5829371c9d4SSatish Balay p27 = pc[26]; 5839371c9d4SSatish Balay p28 = pc[27]; 5849371c9d4SSatish Balay p29 = pc[28]; 5859371c9d4SSatish Balay p30 = pc[29]; 5869371c9d4SSatish Balay p31 = pc[30]; 5879371c9d4SSatish Balay p32 = pc[31]; 5889371c9d4SSatish Balay p33 = pc[32]; 5899371c9d4SSatish Balay p34 = pc[33]; 5909371c9d4SSatish Balay p35 = pc[34]; 5919371c9d4SSatish Balay p36 = pc[35]; 5929371c9d4SSatish Balay if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 || p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0 || p10 != 0.0 || p11 != 0.0 || p12 != 0.0 || p13 != 0.0 || p14 != 0.0 || p15 != 0.0 || p16 != 0.0 || p17 != 0.0 || p18 != 0.0 || p19 != 0.0 || p20 != 0.0 || p21 != 0.0 || p22 != 0.0 || p23 != 0.0 || p24 != 0.0 || p25 != 0.0 || p26 != 0.0 || p27 != 0.0 || p28 != 0.0 || p29 != 0.0 || p30 != 0.0 || p31 != 0.0 || p32 != 0.0 || p33 != 0.0 || p34 != 0.0 || p35 != 0.0 || p36 != 0.0) { 593bef36659SShri Abhyankar pv = ba + 36 * diag_offset[row]; 594bef36659SShri Abhyankar pj = bj + diag_offset[row] + 1; 5959371c9d4SSatish Balay x1 = pv[0]; 5969371c9d4SSatish Balay x2 = pv[1]; 5979371c9d4SSatish Balay x3 = pv[2]; 5989371c9d4SSatish Balay x4 = pv[3]; 5999371c9d4SSatish Balay x5 = pv[4]; 6009371c9d4SSatish Balay x6 = pv[5]; 6019371c9d4SSatish Balay x7 = pv[6]; 6029371c9d4SSatish Balay x8 = pv[7]; 6039371c9d4SSatish Balay x9 = pv[8]; 6049371c9d4SSatish Balay x10 = pv[9]; 6059371c9d4SSatish Balay x11 = pv[10]; 6069371c9d4SSatish Balay x12 = pv[11]; 6079371c9d4SSatish Balay x13 = pv[12]; 6089371c9d4SSatish Balay x14 = pv[13]; 6099371c9d4SSatish Balay x15 = pv[14]; 6109371c9d4SSatish Balay x16 = pv[15]; 6119371c9d4SSatish Balay x17 = pv[16]; 6129371c9d4SSatish Balay x18 = pv[17]; 6139371c9d4SSatish Balay x19 = pv[18]; 6149371c9d4SSatish Balay x20 = pv[19]; 6159371c9d4SSatish Balay x21 = pv[20]; 6169371c9d4SSatish Balay x22 = pv[21]; 6179371c9d4SSatish Balay x23 = pv[22]; 6189371c9d4SSatish Balay x24 = pv[23]; 6199371c9d4SSatish Balay x25 = pv[24]; 6209371c9d4SSatish Balay x26 = pv[25]; 6219371c9d4SSatish Balay x27 = pv[26]; 6229371c9d4SSatish Balay x28 = pv[27]; 6239371c9d4SSatish Balay x29 = pv[28]; 6249371c9d4SSatish Balay x30 = pv[29]; 6259371c9d4SSatish Balay x31 = pv[30]; 6269371c9d4SSatish Balay x32 = pv[31]; 6279371c9d4SSatish Balay x33 = pv[32]; 6289371c9d4SSatish Balay x34 = pv[33]; 6299371c9d4SSatish Balay x35 = pv[34]; 6309371c9d4SSatish Balay x36 = pv[35]; 631bef36659SShri Abhyankar pc[0] = m1 = p1 * x1 + p7 * x2 + p13 * x3 + p19 * x4 + p25 * x5 + p31 * x6; 632bef36659SShri Abhyankar pc[1] = m2 = p2 * x1 + p8 * x2 + p14 * x3 + p20 * x4 + p26 * x5 + p32 * x6; 633bef36659SShri Abhyankar pc[2] = m3 = p3 * x1 + p9 * x2 + p15 * x3 + p21 * x4 + p27 * x5 + p33 * x6; 634bef36659SShri Abhyankar pc[3] = m4 = p4 * x1 + p10 * x2 + p16 * x3 + p22 * x4 + p28 * x5 + p34 * x6; 635bef36659SShri Abhyankar pc[4] = m5 = p5 * x1 + p11 * x2 + p17 * x3 + p23 * x4 + p29 * x5 + p35 * x6; 636bef36659SShri Abhyankar pc[5] = m6 = p6 * x1 + p12 * x2 + p18 * x3 + p24 * x4 + p30 * x5 + p36 * x6; 637bef36659SShri Abhyankar 638bef36659SShri Abhyankar pc[6] = m7 = p1 * x7 + p7 * x8 + p13 * x9 + p19 * x10 + p25 * x11 + p31 * x12; 639bef36659SShri Abhyankar pc[7] = m8 = p2 * x7 + p8 * x8 + p14 * x9 + p20 * x10 + p26 * x11 + p32 * x12; 640bef36659SShri Abhyankar pc[8] = m9 = p3 * x7 + p9 * x8 + p15 * x9 + p21 * x10 + p27 * x11 + p33 * x12; 641bef36659SShri Abhyankar pc[9] = m10 = p4 * x7 + p10 * x8 + p16 * x9 + p22 * x10 + p28 * x11 + p34 * x12; 642bef36659SShri Abhyankar pc[10] = m11 = p5 * x7 + p11 * x8 + p17 * x9 + p23 * x10 + p29 * x11 + p35 * x12; 643bef36659SShri Abhyankar pc[11] = m12 = p6 * x7 + p12 * x8 + p18 * x9 + p24 * x10 + p30 * x11 + p36 * x12; 644bef36659SShri Abhyankar 645bef36659SShri Abhyankar pc[12] = m13 = p1 * x13 + p7 * x14 + p13 * x15 + p19 * x16 + p25 * x17 + p31 * x18; 646bef36659SShri Abhyankar pc[13] = m14 = p2 * x13 + p8 * x14 + p14 * x15 + p20 * x16 + p26 * x17 + p32 * x18; 647bef36659SShri Abhyankar pc[14] = m15 = p3 * x13 + p9 * x14 + p15 * x15 + p21 * x16 + p27 * x17 + p33 * x18; 648bef36659SShri Abhyankar pc[15] = m16 = p4 * x13 + p10 * x14 + p16 * x15 + p22 * x16 + p28 * x17 + p34 * x18; 649bef36659SShri Abhyankar pc[16] = m17 = p5 * x13 + p11 * x14 + p17 * x15 + p23 * x16 + p29 * x17 + p35 * x18; 650bef36659SShri Abhyankar pc[17] = m18 = p6 * x13 + p12 * x14 + p18 * x15 + p24 * x16 + p30 * x17 + p36 * x18; 651bef36659SShri Abhyankar 652bef36659SShri Abhyankar pc[18] = m19 = p1 * x19 + p7 * x20 + p13 * x21 + p19 * x22 + p25 * x23 + p31 * x24; 653bef36659SShri Abhyankar pc[19] = m20 = p2 * x19 + p8 * x20 + p14 * x21 + p20 * x22 + p26 * x23 + p32 * x24; 654bef36659SShri Abhyankar pc[20] = m21 = p3 * x19 + p9 * x20 + p15 * x21 + p21 * x22 + p27 * x23 + p33 * x24; 655bef36659SShri Abhyankar pc[21] = m22 = p4 * x19 + p10 * x20 + p16 * x21 + p22 * x22 + p28 * x23 + p34 * x24; 656bef36659SShri Abhyankar pc[22] = m23 = p5 * x19 + p11 * x20 + p17 * x21 + p23 * x22 + p29 * x23 + p35 * x24; 657bef36659SShri Abhyankar pc[23] = m24 = p6 * x19 + p12 * x20 + p18 * x21 + p24 * x22 + p30 * x23 + p36 * x24; 658bef36659SShri Abhyankar 659bef36659SShri Abhyankar pc[24] = m25 = p1 * x25 + p7 * x26 + p13 * x27 + p19 * x28 + p25 * x29 + p31 * x30; 660bef36659SShri Abhyankar pc[25] = m26 = p2 * x25 + p8 * x26 + p14 * x27 + p20 * x28 + p26 * x29 + p32 * x30; 661bef36659SShri Abhyankar pc[26] = m27 = p3 * x25 + p9 * x26 + p15 * x27 + p21 * x28 + p27 * x29 + p33 * x30; 662bef36659SShri Abhyankar pc[27] = m28 = p4 * x25 + p10 * x26 + p16 * x27 + p22 * x28 + p28 * x29 + p34 * x30; 663bef36659SShri Abhyankar pc[28] = m29 = p5 * x25 + p11 * x26 + p17 * x27 + p23 * x28 + p29 * x29 + p35 * x30; 664bef36659SShri Abhyankar pc[29] = m30 = p6 * x25 + p12 * x26 + p18 * x27 + p24 * x28 + p30 * x29 + p36 * x30; 665bef36659SShri Abhyankar 666bef36659SShri Abhyankar pc[30] = m31 = p1 * x31 + p7 * x32 + p13 * x33 + p19 * x34 + p25 * x35 + p31 * x36; 667bef36659SShri Abhyankar pc[31] = m32 = p2 * x31 + p8 * x32 + p14 * x33 + p20 * x34 + p26 * x35 + p32 * x36; 668bef36659SShri Abhyankar pc[32] = m33 = p3 * x31 + p9 * x32 + p15 * x33 + p21 * x34 + p27 * x35 + p33 * x36; 669bef36659SShri Abhyankar pc[33] = m34 = p4 * x31 + p10 * x32 + p16 * x33 + p22 * x34 + p28 * x35 + p34 * x36; 670bef36659SShri Abhyankar pc[34] = m35 = p5 * x31 + p11 * x32 + p17 * x33 + p23 * x34 + p29 * x35 + p35 * x36; 671bef36659SShri Abhyankar pc[35] = m36 = p6 * x31 + p12 * x32 + p18 * x33 + p24 * x34 + p30 * x35 + p36 * x36; 672bef36659SShri Abhyankar 673bef36659SShri Abhyankar nz = bi[row + 1] - diag_offset[row] - 1; 674bef36659SShri Abhyankar pv += 36; 675bef36659SShri Abhyankar for (j = 0; j < nz; j++) { 6769371c9d4SSatish Balay x1 = pv[0]; 6779371c9d4SSatish Balay x2 = pv[1]; 6789371c9d4SSatish Balay x3 = pv[2]; 6799371c9d4SSatish Balay x4 = pv[3]; 6809371c9d4SSatish Balay x5 = pv[4]; 6819371c9d4SSatish Balay x6 = pv[5]; 6829371c9d4SSatish Balay x7 = pv[6]; 6839371c9d4SSatish Balay x8 = pv[7]; 6849371c9d4SSatish Balay x9 = pv[8]; 6859371c9d4SSatish Balay x10 = pv[9]; 6869371c9d4SSatish Balay x11 = pv[10]; 6879371c9d4SSatish Balay x12 = pv[11]; 6889371c9d4SSatish Balay x13 = pv[12]; 6899371c9d4SSatish Balay x14 = pv[13]; 6909371c9d4SSatish Balay x15 = pv[14]; 6919371c9d4SSatish Balay x16 = pv[15]; 6929371c9d4SSatish Balay x17 = pv[16]; 6939371c9d4SSatish Balay x18 = pv[17]; 6949371c9d4SSatish Balay x19 = pv[18]; 6959371c9d4SSatish Balay x20 = pv[19]; 6969371c9d4SSatish Balay x21 = pv[20]; 6979371c9d4SSatish Balay x22 = pv[21]; 6989371c9d4SSatish Balay x23 = pv[22]; 6999371c9d4SSatish Balay x24 = pv[23]; 7009371c9d4SSatish Balay x25 = pv[24]; 7019371c9d4SSatish Balay x26 = pv[25]; 7029371c9d4SSatish Balay x27 = pv[26]; 7039371c9d4SSatish Balay x28 = pv[27]; 7049371c9d4SSatish Balay x29 = pv[28]; 7059371c9d4SSatish Balay x30 = pv[29]; 7069371c9d4SSatish Balay x31 = pv[30]; 7079371c9d4SSatish Balay x32 = pv[31]; 7089371c9d4SSatish Balay x33 = pv[32]; 7099371c9d4SSatish Balay x34 = pv[33]; 7109371c9d4SSatish Balay x35 = pv[34]; 7119371c9d4SSatish Balay x36 = pv[35]; 712bef36659SShri Abhyankar x = rtmp + 36 * pj[j]; 713bef36659SShri Abhyankar x[0] -= m1 * x1 + m7 * x2 + m13 * x3 + m19 * x4 + m25 * x5 + m31 * x6; 714bef36659SShri Abhyankar x[1] -= m2 * x1 + m8 * x2 + m14 * x3 + m20 * x4 + m26 * x5 + m32 * x6; 715bef36659SShri Abhyankar x[2] -= m3 * x1 + m9 * x2 + m15 * x3 + m21 * x4 + m27 * x5 + m33 * x6; 716bef36659SShri Abhyankar x[3] -= m4 * x1 + m10 * x2 + m16 * x3 + m22 * x4 + m28 * x5 + m34 * x6; 717bef36659SShri Abhyankar x[4] -= m5 * x1 + m11 * x2 + m17 * x3 + m23 * x4 + m29 * x5 + m35 * x6; 718bef36659SShri Abhyankar x[5] -= m6 * x1 + m12 * x2 + m18 * x3 + m24 * x4 + m30 * x5 + m36 * x6; 719bef36659SShri Abhyankar 720bef36659SShri Abhyankar x[6] -= m1 * x7 + m7 * x8 + m13 * x9 + m19 * x10 + m25 * x11 + m31 * x12; 721bef36659SShri Abhyankar x[7] -= m2 * x7 + m8 * x8 + m14 * x9 + m20 * x10 + m26 * x11 + m32 * x12; 722bef36659SShri Abhyankar x[8] -= m3 * x7 + m9 * x8 + m15 * x9 + m21 * x10 + m27 * x11 + m33 * x12; 723bef36659SShri Abhyankar x[9] -= m4 * x7 + m10 * x8 + m16 * x9 + m22 * x10 + m28 * x11 + m34 * x12; 724bef36659SShri Abhyankar x[10] -= m5 * x7 + m11 * x8 + m17 * x9 + m23 * x10 + m29 * x11 + m35 * x12; 725bef36659SShri Abhyankar x[11] -= m6 * x7 + m12 * x8 + m18 * x9 + m24 * x10 + m30 * x11 + m36 * x12; 726bef36659SShri Abhyankar 727bef36659SShri Abhyankar x[12] -= m1 * x13 + m7 * x14 + m13 * x15 + m19 * x16 + m25 * x17 + m31 * x18; 728bef36659SShri Abhyankar x[13] -= m2 * x13 + m8 * x14 + m14 * x15 + m20 * x16 + m26 * x17 + m32 * x18; 729bef36659SShri Abhyankar x[14] -= m3 * x13 + m9 * x14 + m15 * x15 + m21 * x16 + m27 * x17 + m33 * x18; 730bef36659SShri Abhyankar x[15] -= m4 * x13 + m10 * x14 + m16 * x15 + m22 * x16 + m28 * x17 + m34 * x18; 731bef36659SShri Abhyankar x[16] -= m5 * x13 + m11 * x14 + m17 * x15 + m23 * x16 + m29 * x17 + m35 * x18; 732bef36659SShri Abhyankar x[17] -= m6 * x13 + m12 * x14 + m18 * x15 + m24 * x16 + m30 * x17 + m36 * x18; 733bef36659SShri Abhyankar 734bef36659SShri Abhyankar x[18] -= m1 * x19 + m7 * x20 + m13 * x21 + m19 * x22 + m25 * x23 + m31 * x24; 735bef36659SShri Abhyankar x[19] -= m2 * x19 + m8 * x20 + m14 * x21 + m20 * x22 + m26 * x23 + m32 * x24; 736bef36659SShri Abhyankar x[20] -= m3 * x19 + m9 * x20 + m15 * x21 + m21 * x22 + m27 * x23 + m33 * x24; 737bef36659SShri Abhyankar x[21] -= m4 * x19 + m10 * x20 + m16 * x21 + m22 * x22 + m28 * x23 + m34 * x24; 738bef36659SShri Abhyankar x[22] -= m5 * x19 + m11 * x20 + m17 * x21 + m23 * x22 + m29 * x23 + m35 * x24; 739bef36659SShri Abhyankar x[23] -= m6 * x19 + m12 * x20 + m18 * x21 + m24 * x22 + m30 * x23 + m36 * x24; 740bef36659SShri Abhyankar 741bef36659SShri Abhyankar x[24] -= m1 * x25 + m7 * x26 + m13 * x27 + m19 * x28 + m25 * x29 + m31 * x30; 742bef36659SShri Abhyankar x[25] -= m2 * x25 + m8 * x26 + m14 * x27 + m20 * x28 + m26 * x29 + m32 * x30; 743bef36659SShri Abhyankar x[26] -= m3 * x25 + m9 * x26 + m15 * x27 + m21 * x28 + m27 * x29 + m33 * x30; 744bef36659SShri Abhyankar x[27] -= m4 * x25 + m10 * x26 + m16 * x27 + m22 * x28 + m28 * x29 + m34 * x30; 745bef36659SShri Abhyankar x[28] -= m5 * x25 + m11 * x26 + m17 * x27 + m23 * x28 + m29 * x29 + m35 * x30; 746bef36659SShri Abhyankar x[29] -= m6 * x25 + m12 * x26 + m18 * x27 + m24 * x28 + m30 * x29 + m36 * x30; 747bef36659SShri Abhyankar 748bef36659SShri Abhyankar x[30] -= m1 * x31 + m7 * x32 + m13 * x33 + m19 * x34 + m25 * x35 + m31 * x36; 749bef36659SShri Abhyankar x[31] -= m2 * x31 + m8 * x32 + m14 * x33 + m20 * x34 + m26 * x35 + m32 * x36; 750bef36659SShri Abhyankar x[32] -= m3 * x31 + m9 * x32 + m15 * x33 + m21 * x34 + m27 * x35 + m33 * x36; 751bef36659SShri Abhyankar x[33] -= m4 * x31 + m10 * x32 + m16 * x33 + m22 * x34 + m28 * x35 + m34 * x36; 752bef36659SShri Abhyankar x[34] -= m5 * x31 + m11 * x32 + m17 * x33 + m23 * x34 + m29 * x35 + m35 * x36; 753bef36659SShri Abhyankar x[35] -= m6 * x31 + m12 * x32 + m18 * x33 + m24 * x34 + m30 * x35 + m36 * x36; 754bef36659SShri Abhyankar 755bef36659SShri Abhyankar pv += 36; 756bef36659SShri Abhyankar } 7579566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(432.0 * nz + 396.0)); 758bef36659SShri Abhyankar } 759bef36659SShri Abhyankar row = *ajtmp++; 760bef36659SShri Abhyankar } 761bef36659SShri Abhyankar /* finished row so stick it into b->a */ 762bef36659SShri Abhyankar pv = ba + 36 * bi[i]; 763bef36659SShri Abhyankar pj = bj + bi[i]; 764bef36659SShri Abhyankar nz = bi[i + 1] - bi[i]; 765bef36659SShri Abhyankar for (j = 0; j < nz; j++) { 766bef36659SShri Abhyankar x = rtmp + 36 * pj[j]; 7679371c9d4SSatish Balay pv[0] = x[0]; 7689371c9d4SSatish Balay pv[1] = x[1]; 7699371c9d4SSatish Balay pv[2] = x[2]; 7709371c9d4SSatish Balay pv[3] = x[3]; 7719371c9d4SSatish Balay pv[4] = x[4]; 7729371c9d4SSatish Balay pv[5] = x[5]; 7739371c9d4SSatish Balay pv[6] = x[6]; 7749371c9d4SSatish Balay pv[7] = x[7]; 7759371c9d4SSatish Balay pv[8] = x[8]; 7769371c9d4SSatish Balay pv[9] = x[9]; 7779371c9d4SSatish Balay pv[10] = x[10]; 7789371c9d4SSatish Balay pv[11] = x[11]; 7799371c9d4SSatish Balay pv[12] = x[12]; 7809371c9d4SSatish Balay pv[13] = x[13]; 7819371c9d4SSatish Balay pv[14] = x[14]; 7829371c9d4SSatish Balay pv[15] = x[15]; 7839371c9d4SSatish Balay pv[16] = x[16]; 7849371c9d4SSatish Balay pv[17] = x[17]; 7859371c9d4SSatish Balay pv[18] = x[18]; 7869371c9d4SSatish Balay pv[19] = x[19]; 7879371c9d4SSatish Balay pv[20] = x[20]; 7889371c9d4SSatish Balay pv[21] = x[21]; 7899371c9d4SSatish Balay pv[22] = x[22]; 7909371c9d4SSatish Balay pv[23] = x[23]; 7919371c9d4SSatish Balay pv[24] = x[24]; 7929371c9d4SSatish Balay pv[25] = x[25]; 7939371c9d4SSatish Balay pv[26] = x[26]; 7949371c9d4SSatish Balay pv[27] = x[27]; 7959371c9d4SSatish Balay pv[28] = x[28]; 7969371c9d4SSatish Balay pv[29] = x[29]; 7979371c9d4SSatish Balay pv[30] = x[30]; 7989371c9d4SSatish Balay pv[31] = x[31]; 7999371c9d4SSatish Balay pv[32] = x[32]; 8009371c9d4SSatish Balay pv[33] = x[33]; 8019371c9d4SSatish Balay pv[34] = x[34]; 8029371c9d4SSatish Balay pv[35] = x[35]; 803bef36659SShri Abhyankar pv += 36; 804bef36659SShri Abhyankar } 805bef36659SShri Abhyankar /* invert diagonal block */ 806bef36659SShri Abhyankar w = ba + 36 * diag_offset[i]; 8079566063dSJacob Faibussowitsch PetscCall(PetscKernel_A_gets_inverse_A_6(w, shift, allowzeropivot, &zeropivotdetected)); 8087b6c816cSBarry Smith if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 809bef36659SShri Abhyankar } 810bef36659SShri Abhyankar 8119566063dSJacob Faibussowitsch PetscCall(PetscFree(rtmp)); 81226fbe8dcSKarl Rupp 81306e38f1dSHong Zhang C->ops->solve = MatSolve_SeqBAIJ_6_NaturalOrdering_inplace; 81406e38f1dSHong Zhang C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_6_NaturalOrdering_inplace; 815bef36659SShri Abhyankar C->assembled = PETSC_TRUE; 81626fbe8dcSKarl Rupp 8179566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(1.333333333333 * 6 * 6 * 6 * b->mbs)); /* from inverting diagonal blocks */ 8183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 819bef36659SShri Abhyankar } 820bef36659SShri Abhyankar 821d71ae5a4SJacob Faibussowitsch PetscErrorCode MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering(Mat B, Mat A, const MatFactorInfo *info) 822d71ae5a4SJacob Faibussowitsch { 823bef36659SShri Abhyankar Mat C = B; 824bef36659SShri Abhyankar Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data; 825bbd65245SShri Abhyankar PetscInt i, j, k, nz, nzL, row; 826bbd65245SShri Abhyankar const PetscInt n = a->mbs, *ai = a->i, *aj = a->j, *bi = b->i, *bj = b->j; 827bbd65245SShri Abhyankar const PetscInt *ajtmp, *bjtmp, *bdiag = b->diag, *pj, bs2 = a->bs2; 828bef36659SShri Abhyankar MatScalar *rtmp, *pc, *mwork, *v, *pv, *aa = a->a; 829bbd65245SShri Abhyankar PetscInt flg; 830182b8fbaSHong Zhang PetscReal shift = info->shiftamount; 831a455e926SHong Zhang PetscBool allowzeropivot, zeropivotdetected; 832bef36659SShri Abhyankar 833bef36659SShri Abhyankar PetscFunctionBegin; 8340164db54SHong Zhang allowzeropivot = PetscNot(A->erroriffailure); 8350164db54SHong Zhang 836bef36659SShri Abhyankar /* generate work space needed by the factorization */ 8379566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(bs2 * n, &rtmp, bs2, &mwork)); 8389566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(rtmp, bs2 * n)); 839bef36659SShri Abhyankar 840bef36659SShri Abhyankar for (i = 0; i < n; i++) { 841bef36659SShri Abhyankar /* zero rtmp */ 842bef36659SShri Abhyankar /* L part */ 843bef36659SShri Abhyankar nz = bi[i + 1] - bi[i]; 844bef36659SShri Abhyankar bjtmp = bj + bi[i]; 84548a46eb9SPierre Jolivet for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2)); 846bef36659SShri Abhyankar 847bef36659SShri Abhyankar /* U part */ 84853cca76cSShri Abhyankar nz = bdiag[i] - bdiag[i + 1]; 84953cca76cSShri Abhyankar bjtmp = bj + bdiag[i + 1] + 1; 85048a46eb9SPierre Jolivet for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2)); 85153cca76cSShri Abhyankar 85253cca76cSShri Abhyankar /* load in initial (unfactored row) */ 85353cca76cSShri Abhyankar nz = ai[i + 1] - ai[i]; 85453cca76cSShri Abhyankar ajtmp = aj + ai[i]; 85553cca76cSShri Abhyankar v = aa + bs2 * ai[i]; 85648a46eb9SPierre Jolivet for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(rtmp + bs2 * ajtmp[j], v + bs2 * j, bs2)); 85753cca76cSShri Abhyankar 85853cca76cSShri Abhyankar /* elimination */ 85953cca76cSShri Abhyankar bjtmp = bj + bi[i]; 86053cca76cSShri Abhyankar nzL = bi[i + 1] - bi[i]; 86153cca76cSShri Abhyankar for (k = 0; k < nzL; k++) { 86253cca76cSShri Abhyankar row = bjtmp[k]; 86353cca76cSShri Abhyankar pc = rtmp + bs2 * row; 864c35f09e5SBarry Smith for (flg = 0, j = 0; j < bs2; j++) { 865c35f09e5SBarry Smith if (pc[j] != 0.0) { 866c35f09e5SBarry Smith flg = 1; 867c35f09e5SBarry Smith break; 868c35f09e5SBarry Smith } 869c35f09e5SBarry Smith } 87053cca76cSShri Abhyankar if (flg) { 87153cca76cSShri Abhyankar pv = b->a + bs2 * bdiag[row]; 87296b95a6bSBarry Smith /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */ 8739566063dSJacob Faibussowitsch PetscCall(PetscKernel_A_gets_A_times_B_6(pc, pv, mwork)); 87453cca76cSShri Abhyankar 875a5b23f4aSJose E. Roman pj = b->j + bdiag[row + 1] + 1; /* beginning of U(row,:) */ 87653cca76cSShri Abhyankar pv = b->a + bs2 * (bdiag[row + 1] + 1); 87753cca76cSShri Abhyankar nz = bdiag[row] - bdiag[row + 1] - 1; /* num of entries inU(row,:), excluding diag */ 87853cca76cSShri Abhyankar for (j = 0; j < nz; j++) { 87996b95a6bSBarry Smith /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */ 88053cca76cSShri Abhyankar /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */ 88153cca76cSShri Abhyankar v = rtmp + bs2 * pj[j]; 8829566063dSJacob Faibussowitsch PetscCall(PetscKernel_A_gets_A_minus_B_times_C_6(v, pc, pv)); 88353cca76cSShri Abhyankar pv += bs2; 88453cca76cSShri Abhyankar } 8859566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(432.0 * nz + 396)); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */ 88653cca76cSShri Abhyankar } 88753cca76cSShri Abhyankar } 88853cca76cSShri Abhyankar 88953cca76cSShri Abhyankar /* finished row so stick it into b->a */ 89053cca76cSShri Abhyankar /* L part */ 89153cca76cSShri Abhyankar pv = b->a + bs2 * bi[i]; 89253cca76cSShri Abhyankar pj = b->j + bi[i]; 89353cca76cSShri Abhyankar nz = bi[i + 1] - bi[i]; 89448a46eb9SPierre Jolivet for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2)); 89553cca76cSShri Abhyankar 896a5b23f4aSJose E. Roman /* Mark diagonal and invert diagonal for simpler triangular solves */ 89753cca76cSShri Abhyankar pv = b->a + bs2 * bdiag[i]; 89853cca76cSShri Abhyankar pj = b->j + bdiag[i]; 8999566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(pv, rtmp + bs2 * pj[0], bs2)); 9009566063dSJacob Faibussowitsch PetscCall(PetscKernel_A_gets_inverse_A_6(pv, shift, allowzeropivot, &zeropivotdetected)); 9017b6c816cSBarry Smith if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 90253cca76cSShri Abhyankar 90353cca76cSShri Abhyankar /* U part */ 90453cca76cSShri Abhyankar pv = b->a + bs2 * (bdiag[i + 1] + 1); 90553cca76cSShri Abhyankar pj = b->j + bdiag[i + 1] + 1; 90653cca76cSShri Abhyankar nz = bdiag[i] - bdiag[i + 1] - 1; 90748a46eb9SPierre Jolivet for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2)); 90853cca76cSShri Abhyankar } 9099566063dSJacob Faibussowitsch PetscCall(PetscFree2(rtmp, mwork)); 91026fbe8dcSKarl Rupp 9114dd39f65SShri Abhyankar C->ops->solve = MatSolve_SeqBAIJ_6_NaturalOrdering; 9124dd39f65SShri Abhyankar C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_6_NaturalOrdering; 91353cca76cSShri Abhyankar C->assembled = PETSC_TRUE; 91426fbe8dcSKarl Rupp 9159566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(1.333333333333 * 6 * 6 * 6 * n)); /* from inverting diagonal blocks */ 9163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 91753cca76cSShri Abhyankar } 918