xref: /petsc/src/mat/impls/baij/seq/baijfact7.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1be1d678aSKris Buschelman 
283287d42SBarry Smith /*
383287d42SBarry Smith     Factorization code for BAIJ format.
483287d42SBarry Smith */
5c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h>
6af0996ceSBarry Smith #include <petsc/private/kernels/blockinvert.h>
783287d42SBarry Smith 
883287d42SBarry Smith /* ------------------------------------------------------------*/
983287d42SBarry Smith /*
1083287d42SBarry Smith       Version for when blocks are 6 by 6
1183287d42SBarry Smith */
12*d71ae5a4SJacob Faibussowitsch PetscErrorCode MatLUFactorNumeric_SeqBAIJ_6_inplace(Mat C, Mat A, const MatFactorInfo *info)
13*d71ae5a4SJacob Faibussowitsch {
1483287d42SBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data;
1583287d42SBarry Smith   IS              isrow = b->row, isicol = b->icol;
165d0c19d7SBarry Smith   const PetscInt *ajtmpold, *ajtmp, *diag_offset = b->diag, *r, *ic, *bi = b->i, *bj = b->j, *ai = a->i, *aj = a->j, *pj;
175d0c19d7SBarry Smith   PetscInt        nz, row, i, j, n = a->mbs, idx;
1883287d42SBarry Smith   MatScalar      *pv, *v, *rtmp, *pc, *w, *x;
1983287d42SBarry Smith   MatScalar       p1, p2, p3, p4, m1, m2, m3, m4, m5, m6, m7, m8, m9, x1, x2, x3, x4;
2083287d42SBarry Smith   MatScalar       p5, p6, p7, p8, p9, x5, x6, x7, x8, x9, x10, x11, x12, x13, x14, x15, x16;
2183287d42SBarry Smith   MatScalar       x17, x18, x19, x20, x21, x22, x23, x24, x25, p10, p11, p12, p13, p14;
2283287d42SBarry Smith   MatScalar       p15, p16, p17, p18, p19, p20, p21, p22, p23, p24, p25, m10, m11, m12;
2383287d42SBarry Smith   MatScalar       m13, m14, m15, m16, m17, m18, m19, m20, m21, m22, m23, m24, m25;
2483287d42SBarry Smith   MatScalar       p26, p27, p28, p29, p30, p31, p32, p33, p34, p35, p36;
2583287d42SBarry Smith   MatScalar       x26, x27, x28, x29, x30, x31, x32, x33, x34, x35, x36;
2683287d42SBarry Smith   MatScalar       m26, m27, m28, m29, m30, m31, m32, m33, m34, m35, m36;
2783287d42SBarry Smith   MatScalar      *ba = b->a, *aa = a->a;
28182b8fbaSHong Zhang   PetscReal       shift = info->shiftamount;
29a455e926SHong Zhang   PetscBool       allowzeropivot, zeropivotdetected;
3083287d42SBarry Smith 
3183287d42SBarry Smith   PetscFunctionBegin;
320164db54SHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
339566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isrow, &r));
349566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isicol, &ic));
359566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(36 * (n + 1), &rtmp));
3683287d42SBarry Smith 
3783287d42SBarry Smith   for (i = 0; i < n; i++) {
3883287d42SBarry Smith     nz    = bi[i + 1] - bi[i];
3983287d42SBarry Smith     ajtmp = bj + bi[i];
4083287d42SBarry Smith     for (j = 0; j < nz; j++) {
4183287d42SBarry Smith       x    = rtmp + 36 * ajtmp[j];
4283287d42SBarry Smith       x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = x[9] = 0.0;
4383287d42SBarry Smith       x[10] = x[11] = x[12] = x[13] = x[14] = x[15] = x[16] = x[17] = 0.0;
4483287d42SBarry Smith       x[18] = x[19] = x[20] = x[21] = x[22] = x[23] = x[24] = x[25] = 0.0;
4583287d42SBarry Smith       x[26] = x[27] = x[28] = x[29] = x[30] = x[31] = x[32] = x[33] = 0.0;
4683287d42SBarry Smith       x[34] = x[35] = 0.0;
4783287d42SBarry Smith     }
4883287d42SBarry Smith     /* load in initial (unfactored row) */
4983287d42SBarry Smith     idx      = r[i];
5083287d42SBarry Smith     nz       = ai[idx + 1] - ai[idx];
5183287d42SBarry Smith     ajtmpold = aj + ai[idx];
5283287d42SBarry Smith     v        = aa + 36 * ai[idx];
5383287d42SBarry Smith     for (j = 0; j < nz; j++) {
5483287d42SBarry Smith       x     = rtmp + 36 * ic[ajtmpold[j]];
559371c9d4SSatish Balay       x[0]  = v[0];
569371c9d4SSatish Balay       x[1]  = v[1];
579371c9d4SSatish Balay       x[2]  = v[2];
589371c9d4SSatish Balay       x[3]  = v[3];
599371c9d4SSatish Balay       x[4]  = v[4];
609371c9d4SSatish Balay       x[5]  = v[5];
619371c9d4SSatish Balay       x[6]  = v[6];
629371c9d4SSatish Balay       x[7]  = v[7];
639371c9d4SSatish Balay       x[8]  = v[8];
649371c9d4SSatish Balay       x[9]  = v[9];
659371c9d4SSatish Balay       x[10] = v[10];
669371c9d4SSatish Balay       x[11] = v[11];
679371c9d4SSatish Balay       x[12] = v[12];
689371c9d4SSatish Balay       x[13] = v[13];
699371c9d4SSatish Balay       x[14] = v[14];
709371c9d4SSatish Balay       x[15] = v[15];
719371c9d4SSatish Balay       x[16] = v[16];
729371c9d4SSatish Balay       x[17] = v[17];
739371c9d4SSatish Balay       x[18] = v[18];
749371c9d4SSatish Balay       x[19] = v[19];
759371c9d4SSatish Balay       x[20] = v[20];
769371c9d4SSatish Balay       x[21] = v[21];
779371c9d4SSatish Balay       x[22] = v[22];
789371c9d4SSatish Balay       x[23] = v[23];
799371c9d4SSatish Balay       x[24] = v[24];
809371c9d4SSatish Balay       x[25] = v[25];
819371c9d4SSatish Balay       x[26] = v[26];
829371c9d4SSatish Balay       x[27] = v[27];
839371c9d4SSatish Balay       x[28] = v[28];
849371c9d4SSatish Balay       x[29] = v[29];
859371c9d4SSatish Balay       x[30] = v[30];
869371c9d4SSatish Balay       x[31] = v[31];
879371c9d4SSatish Balay       x[32] = v[32];
889371c9d4SSatish Balay       x[33] = v[33];
899371c9d4SSatish Balay       x[34] = v[34];
909371c9d4SSatish Balay       x[35] = v[35];
9183287d42SBarry Smith       v += 36;
9283287d42SBarry Smith     }
9383287d42SBarry Smith     row = *ajtmp++;
9483287d42SBarry Smith     while (row < i) {
9583287d42SBarry Smith       pc  = rtmp + 36 * row;
969371c9d4SSatish Balay       p1  = pc[0];
979371c9d4SSatish Balay       p2  = pc[1];
989371c9d4SSatish Balay       p3  = pc[2];
999371c9d4SSatish Balay       p4  = pc[3];
1009371c9d4SSatish Balay       p5  = pc[4];
1019371c9d4SSatish Balay       p6  = pc[5];
1029371c9d4SSatish Balay       p7  = pc[6];
1039371c9d4SSatish Balay       p8  = pc[7];
1049371c9d4SSatish Balay       p9  = pc[8];
1059371c9d4SSatish Balay       p10 = pc[9];
1069371c9d4SSatish Balay       p11 = pc[10];
1079371c9d4SSatish Balay       p12 = pc[11];
1089371c9d4SSatish Balay       p13 = pc[12];
1099371c9d4SSatish Balay       p14 = pc[13];
1109371c9d4SSatish Balay       p15 = pc[14];
1119371c9d4SSatish Balay       p16 = pc[15];
1129371c9d4SSatish Balay       p17 = pc[16];
1139371c9d4SSatish Balay       p18 = pc[17];
1149371c9d4SSatish Balay       p19 = pc[18];
1159371c9d4SSatish Balay       p20 = pc[19];
1169371c9d4SSatish Balay       p21 = pc[20];
1179371c9d4SSatish Balay       p22 = pc[21];
1189371c9d4SSatish Balay       p23 = pc[22];
1199371c9d4SSatish Balay       p24 = pc[23];
1209371c9d4SSatish Balay       p25 = pc[24];
1219371c9d4SSatish Balay       p26 = pc[25];
1229371c9d4SSatish Balay       p27 = pc[26];
1239371c9d4SSatish Balay       p28 = pc[27];
1249371c9d4SSatish Balay       p29 = pc[28];
1259371c9d4SSatish Balay       p30 = pc[29];
1269371c9d4SSatish Balay       p31 = pc[30];
1279371c9d4SSatish Balay       p32 = pc[31];
1289371c9d4SSatish Balay       p33 = pc[32];
1299371c9d4SSatish Balay       p34 = pc[33];
1309371c9d4SSatish Balay       p35 = pc[34];
1319371c9d4SSatish Balay       p36 = pc[35];
1329371c9d4SSatish 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) {
13383287d42SBarry Smith         pv    = ba + 36 * diag_offset[row];
13483287d42SBarry Smith         pj    = bj + diag_offset[row] + 1;
1359371c9d4SSatish Balay         x1    = pv[0];
1369371c9d4SSatish Balay         x2    = pv[1];
1379371c9d4SSatish Balay         x3    = pv[2];
1389371c9d4SSatish Balay         x4    = pv[3];
1399371c9d4SSatish Balay         x5    = pv[4];
1409371c9d4SSatish Balay         x6    = pv[5];
1419371c9d4SSatish Balay         x7    = pv[6];
1429371c9d4SSatish Balay         x8    = pv[7];
1439371c9d4SSatish Balay         x9    = pv[8];
1449371c9d4SSatish Balay         x10   = pv[9];
1459371c9d4SSatish Balay         x11   = pv[10];
1469371c9d4SSatish Balay         x12   = pv[11];
1479371c9d4SSatish Balay         x13   = pv[12];
1489371c9d4SSatish Balay         x14   = pv[13];
1499371c9d4SSatish Balay         x15   = pv[14];
1509371c9d4SSatish Balay         x16   = pv[15];
1519371c9d4SSatish Balay         x17   = pv[16];
1529371c9d4SSatish Balay         x18   = pv[17];
1539371c9d4SSatish Balay         x19   = pv[18];
1549371c9d4SSatish Balay         x20   = pv[19];
1559371c9d4SSatish Balay         x21   = pv[20];
1569371c9d4SSatish Balay         x22   = pv[21];
1579371c9d4SSatish Balay         x23   = pv[22];
1589371c9d4SSatish Balay         x24   = pv[23];
1599371c9d4SSatish Balay         x25   = pv[24];
1609371c9d4SSatish Balay         x26   = pv[25];
1619371c9d4SSatish Balay         x27   = pv[26];
1629371c9d4SSatish Balay         x28   = pv[27];
1639371c9d4SSatish Balay         x29   = pv[28];
1649371c9d4SSatish Balay         x30   = pv[29];
1659371c9d4SSatish Balay         x31   = pv[30];
1669371c9d4SSatish Balay         x32   = pv[31];
1679371c9d4SSatish Balay         x33   = pv[32];
1689371c9d4SSatish Balay         x34   = pv[33];
1699371c9d4SSatish Balay         x35   = pv[34];
1709371c9d4SSatish Balay         x36   = pv[35];
17183287d42SBarry Smith         pc[0] = m1 = p1 * x1 + p7 * x2 + p13 * x3 + p19 * x4 + p25 * x5 + p31 * x6;
17283287d42SBarry Smith         pc[1] = m2 = p2 * x1 + p8 * x2 + p14 * x3 + p20 * x4 + p26 * x5 + p32 * x6;
17383287d42SBarry Smith         pc[2] = m3 = p3 * x1 + p9 * x2 + p15 * x3 + p21 * x4 + p27 * x5 + p33 * x6;
17483287d42SBarry Smith         pc[3] = m4 = p4 * x1 + p10 * x2 + p16 * x3 + p22 * x4 + p28 * x5 + p34 * x6;
17583287d42SBarry Smith         pc[4] = m5 = p5 * x1 + p11 * x2 + p17 * x3 + p23 * x4 + p29 * x5 + p35 * x6;
17683287d42SBarry Smith         pc[5] = m6 = p6 * x1 + p12 * x2 + p18 * x3 + p24 * x4 + p30 * x5 + p36 * x6;
17783287d42SBarry Smith 
17883287d42SBarry Smith         pc[6] = m7 = p1 * x7 + p7 * x8 + p13 * x9 + p19 * x10 + p25 * x11 + p31 * x12;
17983287d42SBarry Smith         pc[7] = m8 = p2 * x7 + p8 * x8 + p14 * x9 + p20 * x10 + p26 * x11 + p32 * x12;
18083287d42SBarry Smith         pc[8] = m9 = p3 * x7 + p9 * x8 + p15 * x9 + p21 * x10 + p27 * x11 + p33 * x12;
18183287d42SBarry Smith         pc[9] = m10 = p4 * x7 + p10 * x8 + p16 * x9 + p22 * x10 + p28 * x11 + p34 * x12;
18283287d42SBarry Smith         pc[10] = m11 = p5 * x7 + p11 * x8 + p17 * x9 + p23 * x10 + p29 * x11 + p35 * x12;
18383287d42SBarry Smith         pc[11] = m12 = p6 * x7 + p12 * x8 + p18 * x9 + p24 * x10 + p30 * x11 + p36 * x12;
18483287d42SBarry Smith 
18583287d42SBarry Smith         pc[12] = m13 = p1 * x13 + p7 * x14 + p13 * x15 + p19 * x16 + p25 * x17 + p31 * x18;
18683287d42SBarry Smith         pc[13] = m14 = p2 * x13 + p8 * x14 + p14 * x15 + p20 * x16 + p26 * x17 + p32 * x18;
18783287d42SBarry Smith         pc[14] = m15 = p3 * x13 + p9 * x14 + p15 * x15 + p21 * x16 + p27 * x17 + p33 * x18;
18883287d42SBarry Smith         pc[15] = m16 = p4 * x13 + p10 * x14 + p16 * x15 + p22 * x16 + p28 * x17 + p34 * x18;
18983287d42SBarry Smith         pc[16] = m17 = p5 * x13 + p11 * x14 + p17 * x15 + p23 * x16 + p29 * x17 + p35 * x18;
19083287d42SBarry Smith         pc[17] = m18 = p6 * x13 + p12 * x14 + p18 * x15 + p24 * x16 + p30 * x17 + p36 * x18;
19183287d42SBarry Smith 
19283287d42SBarry Smith         pc[18] = m19 = p1 * x19 + p7 * x20 + p13 * x21 + p19 * x22 + p25 * x23 + p31 * x24;
19383287d42SBarry Smith         pc[19] = m20 = p2 * x19 + p8 * x20 + p14 * x21 + p20 * x22 + p26 * x23 + p32 * x24;
19483287d42SBarry Smith         pc[20] = m21 = p3 * x19 + p9 * x20 + p15 * x21 + p21 * x22 + p27 * x23 + p33 * x24;
19583287d42SBarry Smith         pc[21] = m22 = p4 * x19 + p10 * x20 + p16 * x21 + p22 * x22 + p28 * x23 + p34 * x24;
19683287d42SBarry Smith         pc[22] = m23 = p5 * x19 + p11 * x20 + p17 * x21 + p23 * x22 + p29 * x23 + p35 * x24;
19783287d42SBarry Smith         pc[23] = m24 = p6 * x19 + p12 * x20 + p18 * x21 + p24 * x22 + p30 * x23 + p36 * x24;
19883287d42SBarry Smith 
19983287d42SBarry Smith         pc[24] = m25 = p1 * x25 + p7 * x26 + p13 * x27 + p19 * x28 + p25 * x29 + p31 * x30;
20083287d42SBarry Smith         pc[25] = m26 = p2 * x25 + p8 * x26 + p14 * x27 + p20 * x28 + p26 * x29 + p32 * x30;
20183287d42SBarry Smith         pc[26] = m27 = p3 * x25 + p9 * x26 + p15 * x27 + p21 * x28 + p27 * x29 + p33 * x30;
20283287d42SBarry Smith         pc[27] = m28 = p4 * x25 + p10 * x26 + p16 * x27 + p22 * x28 + p28 * x29 + p34 * x30;
20383287d42SBarry Smith         pc[28] = m29 = p5 * x25 + p11 * x26 + p17 * x27 + p23 * x28 + p29 * x29 + p35 * x30;
20483287d42SBarry Smith         pc[29] = m30 = p6 * x25 + p12 * x26 + p18 * x27 + p24 * x28 + p30 * x29 + p36 * x30;
20583287d42SBarry Smith 
20683287d42SBarry Smith         pc[30] = m31 = p1 * x31 + p7 * x32 + p13 * x33 + p19 * x34 + p25 * x35 + p31 * x36;
20783287d42SBarry Smith         pc[31] = m32 = p2 * x31 + p8 * x32 + p14 * x33 + p20 * x34 + p26 * x35 + p32 * x36;
20883287d42SBarry Smith         pc[32] = m33 = p3 * x31 + p9 * x32 + p15 * x33 + p21 * x34 + p27 * x35 + p33 * x36;
20983287d42SBarry Smith         pc[33] = m34 = p4 * x31 + p10 * x32 + p16 * x33 + p22 * x34 + p28 * x35 + p34 * x36;
21083287d42SBarry Smith         pc[34] = m35 = p5 * x31 + p11 * x32 + p17 * x33 + p23 * x34 + p29 * x35 + p35 * x36;
21183287d42SBarry Smith         pc[35] = m36 = p6 * x31 + p12 * x32 + p18 * x33 + p24 * x34 + p30 * x35 + p36 * x36;
21283287d42SBarry Smith 
21383287d42SBarry Smith         nz = bi[row + 1] - diag_offset[row] - 1;
21483287d42SBarry Smith         pv += 36;
21583287d42SBarry Smith         for (j = 0; j < nz; j++) {
2169371c9d4SSatish Balay           x1  = pv[0];
2179371c9d4SSatish Balay           x2  = pv[1];
2189371c9d4SSatish Balay           x3  = pv[2];
2199371c9d4SSatish Balay           x4  = pv[3];
2209371c9d4SSatish Balay           x5  = pv[4];
2219371c9d4SSatish Balay           x6  = pv[5];
2229371c9d4SSatish Balay           x7  = pv[6];
2239371c9d4SSatish Balay           x8  = pv[7];
2249371c9d4SSatish Balay           x9  = pv[8];
2259371c9d4SSatish Balay           x10 = pv[9];
2269371c9d4SSatish Balay           x11 = pv[10];
2279371c9d4SSatish Balay           x12 = pv[11];
2289371c9d4SSatish Balay           x13 = pv[12];
2299371c9d4SSatish Balay           x14 = pv[13];
2309371c9d4SSatish Balay           x15 = pv[14];
2319371c9d4SSatish Balay           x16 = pv[15];
2329371c9d4SSatish Balay           x17 = pv[16];
2339371c9d4SSatish Balay           x18 = pv[17];
2349371c9d4SSatish Balay           x19 = pv[18];
2359371c9d4SSatish Balay           x20 = pv[19];
2369371c9d4SSatish Balay           x21 = pv[20];
2379371c9d4SSatish Balay           x22 = pv[21];
2389371c9d4SSatish Balay           x23 = pv[22];
2399371c9d4SSatish Balay           x24 = pv[23];
2409371c9d4SSatish Balay           x25 = pv[24];
2419371c9d4SSatish Balay           x26 = pv[25];
2429371c9d4SSatish Balay           x27 = pv[26];
2439371c9d4SSatish Balay           x28 = pv[27];
2449371c9d4SSatish Balay           x29 = pv[28];
2459371c9d4SSatish Balay           x30 = pv[29];
2469371c9d4SSatish Balay           x31 = pv[30];
2479371c9d4SSatish Balay           x32 = pv[31];
2489371c9d4SSatish Balay           x33 = pv[32];
2499371c9d4SSatish Balay           x34 = pv[33];
2509371c9d4SSatish Balay           x35 = pv[34];
2519371c9d4SSatish Balay           x36 = pv[35];
25283287d42SBarry Smith           x   = rtmp + 36 * pj[j];
25383287d42SBarry Smith           x[0] -= m1 * x1 + m7 * x2 + m13 * x3 + m19 * x4 + m25 * x5 + m31 * x6;
25483287d42SBarry Smith           x[1] -= m2 * x1 + m8 * x2 + m14 * x3 + m20 * x4 + m26 * x5 + m32 * x6;
25583287d42SBarry Smith           x[2] -= m3 * x1 + m9 * x2 + m15 * x3 + m21 * x4 + m27 * x5 + m33 * x6;
25683287d42SBarry Smith           x[3] -= m4 * x1 + m10 * x2 + m16 * x3 + m22 * x4 + m28 * x5 + m34 * x6;
25783287d42SBarry Smith           x[4] -= m5 * x1 + m11 * x2 + m17 * x3 + m23 * x4 + m29 * x5 + m35 * x6;
25883287d42SBarry Smith           x[5] -= m6 * x1 + m12 * x2 + m18 * x3 + m24 * x4 + m30 * x5 + m36 * x6;
25983287d42SBarry Smith 
26083287d42SBarry Smith           x[6] -= m1 * x7 + m7 * x8 + m13 * x9 + m19 * x10 + m25 * x11 + m31 * x12;
26183287d42SBarry Smith           x[7] -= m2 * x7 + m8 * x8 + m14 * x9 + m20 * x10 + m26 * x11 + m32 * x12;
26283287d42SBarry Smith           x[8] -= m3 * x7 + m9 * x8 + m15 * x9 + m21 * x10 + m27 * x11 + m33 * x12;
26383287d42SBarry Smith           x[9] -= m4 * x7 + m10 * x8 + m16 * x9 + m22 * x10 + m28 * x11 + m34 * x12;
26483287d42SBarry Smith           x[10] -= m5 * x7 + m11 * x8 + m17 * x9 + m23 * x10 + m29 * x11 + m35 * x12;
26583287d42SBarry Smith           x[11] -= m6 * x7 + m12 * x8 + m18 * x9 + m24 * x10 + m30 * x11 + m36 * x12;
26683287d42SBarry Smith 
26783287d42SBarry Smith           x[12] -= m1 * x13 + m7 * x14 + m13 * x15 + m19 * x16 + m25 * x17 + m31 * x18;
26883287d42SBarry Smith           x[13] -= m2 * x13 + m8 * x14 + m14 * x15 + m20 * x16 + m26 * x17 + m32 * x18;
26983287d42SBarry Smith           x[14] -= m3 * x13 + m9 * x14 + m15 * x15 + m21 * x16 + m27 * x17 + m33 * x18;
27083287d42SBarry Smith           x[15] -= m4 * x13 + m10 * x14 + m16 * x15 + m22 * x16 + m28 * x17 + m34 * x18;
27183287d42SBarry Smith           x[16] -= m5 * x13 + m11 * x14 + m17 * x15 + m23 * x16 + m29 * x17 + m35 * x18;
27283287d42SBarry Smith           x[17] -= m6 * x13 + m12 * x14 + m18 * x15 + m24 * x16 + m30 * x17 + m36 * x18;
27383287d42SBarry Smith 
27483287d42SBarry Smith           x[18] -= m1 * x19 + m7 * x20 + m13 * x21 + m19 * x22 + m25 * x23 + m31 * x24;
27583287d42SBarry Smith           x[19] -= m2 * x19 + m8 * x20 + m14 * x21 + m20 * x22 + m26 * x23 + m32 * x24;
27683287d42SBarry Smith           x[20] -= m3 * x19 + m9 * x20 + m15 * x21 + m21 * x22 + m27 * x23 + m33 * x24;
27783287d42SBarry Smith           x[21] -= m4 * x19 + m10 * x20 + m16 * x21 + m22 * x22 + m28 * x23 + m34 * x24;
27883287d42SBarry Smith           x[22] -= m5 * x19 + m11 * x20 + m17 * x21 + m23 * x22 + m29 * x23 + m35 * x24;
27983287d42SBarry Smith           x[23] -= m6 * x19 + m12 * x20 + m18 * x21 + m24 * x22 + m30 * x23 + m36 * x24;
28083287d42SBarry Smith 
28183287d42SBarry Smith           x[24] -= m1 * x25 + m7 * x26 + m13 * x27 + m19 * x28 + m25 * x29 + m31 * x30;
28283287d42SBarry Smith           x[25] -= m2 * x25 + m8 * x26 + m14 * x27 + m20 * x28 + m26 * x29 + m32 * x30;
28383287d42SBarry Smith           x[26] -= m3 * x25 + m9 * x26 + m15 * x27 + m21 * x28 + m27 * x29 + m33 * x30;
28483287d42SBarry Smith           x[27] -= m4 * x25 + m10 * x26 + m16 * x27 + m22 * x28 + m28 * x29 + m34 * x30;
28583287d42SBarry Smith           x[28] -= m5 * x25 + m11 * x26 + m17 * x27 + m23 * x28 + m29 * x29 + m35 * x30;
28683287d42SBarry Smith           x[29] -= m6 * x25 + m12 * x26 + m18 * x27 + m24 * x28 + m30 * x29 + m36 * x30;
28783287d42SBarry Smith 
28883287d42SBarry Smith           x[30] -= m1 * x31 + m7 * x32 + m13 * x33 + m19 * x34 + m25 * x35 + m31 * x36;
28983287d42SBarry Smith           x[31] -= m2 * x31 + m8 * x32 + m14 * x33 + m20 * x34 + m26 * x35 + m32 * x36;
29083287d42SBarry Smith           x[32] -= m3 * x31 + m9 * x32 + m15 * x33 + m21 * x34 + m27 * x35 + m33 * x36;
29183287d42SBarry Smith           x[33] -= m4 * x31 + m10 * x32 + m16 * x33 + m22 * x34 + m28 * x35 + m34 * x36;
29283287d42SBarry Smith           x[34] -= m5 * x31 + m11 * x32 + m17 * x33 + m23 * x34 + m29 * x35 + m35 * x36;
29383287d42SBarry Smith           x[35] -= m6 * x31 + m12 * x32 + m18 * x33 + m24 * x34 + m30 * x35 + m36 * x36;
29483287d42SBarry Smith 
29583287d42SBarry Smith           pv += 36;
29683287d42SBarry Smith         }
2979566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(432.0 * nz + 396.0));
29883287d42SBarry Smith       }
29983287d42SBarry Smith       row = *ajtmp++;
30083287d42SBarry Smith     }
30183287d42SBarry Smith     /* finished row so stick it into b->a */
30283287d42SBarry Smith     pv = ba + 36 * bi[i];
30383287d42SBarry Smith     pj = bj + bi[i];
30483287d42SBarry Smith     nz = bi[i + 1] - bi[i];
30583287d42SBarry Smith     for (j = 0; j < nz; j++) {
30683287d42SBarry Smith       x      = rtmp + 36 * pj[j];
3079371c9d4SSatish Balay       pv[0]  = x[0];
3089371c9d4SSatish Balay       pv[1]  = x[1];
3099371c9d4SSatish Balay       pv[2]  = x[2];
3109371c9d4SSatish Balay       pv[3]  = x[3];
3119371c9d4SSatish Balay       pv[4]  = x[4];
3129371c9d4SSatish Balay       pv[5]  = x[5];
3139371c9d4SSatish Balay       pv[6]  = x[6];
3149371c9d4SSatish Balay       pv[7]  = x[7];
3159371c9d4SSatish Balay       pv[8]  = x[8];
3169371c9d4SSatish Balay       pv[9]  = x[9];
3179371c9d4SSatish Balay       pv[10] = x[10];
3189371c9d4SSatish Balay       pv[11] = x[11];
3199371c9d4SSatish Balay       pv[12] = x[12];
3209371c9d4SSatish Balay       pv[13] = x[13];
3219371c9d4SSatish Balay       pv[14] = x[14];
3229371c9d4SSatish Balay       pv[15] = x[15];
3239371c9d4SSatish Balay       pv[16] = x[16];
3249371c9d4SSatish Balay       pv[17] = x[17];
3259371c9d4SSatish Balay       pv[18] = x[18];
3269371c9d4SSatish Balay       pv[19] = x[19];
3279371c9d4SSatish Balay       pv[20] = x[20];
3289371c9d4SSatish Balay       pv[21] = x[21];
3299371c9d4SSatish Balay       pv[22] = x[22];
3309371c9d4SSatish Balay       pv[23] = x[23];
3319371c9d4SSatish Balay       pv[24] = x[24];
3329371c9d4SSatish Balay       pv[25] = x[25];
3339371c9d4SSatish Balay       pv[26] = x[26];
3349371c9d4SSatish Balay       pv[27] = x[27];
3359371c9d4SSatish Balay       pv[28] = x[28];
3369371c9d4SSatish Balay       pv[29] = x[29];
3379371c9d4SSatish Balay       pv[30] = x[30];
3389371c9d4SSatish Balay       pv[31] = x[31];
3399371c9d4SSatish Balay       pv[32] = x[32];
3409371c9d4SSatish Balay       pv[33] = x[33];
3419371c9d4SSatish Balay       pv[34] = x[34];
3429371c9d4SSatish Balay       pv[35] = x[35];
34383287d42SBarry Smith       pv += 36;
34483287d42SBarry Smith     }
34583287d42SBarry Smith     /* invert diagonal block */
34683287d42SBarry Smith     w = ba + 36 * diag_offset[i];
3479566063dSJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_6(w, shift, allowzeropivot, &zeropivotdetected));
3487b6c816cSBarry Smith     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
34983287d42SBarry Smith   }
35083287d42SBarry Smith 
3519566063dSJacob Faibussowitsch   PetscCall(PetscFree(rtmp));
3529566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isicol, &ic));
3539566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow, &r));
35426fbe8dcSKarl Rupp 
35506e38f1dSHong Zhang   C->ops->solve          = MatSolve_SeqBAIJ_6_inplace;
35606e38f1dSHong Zhang   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_6_inplace;
35783287d42SBarry Smith   C->assembled           = PETSC_TRUE;
35826fbe8dcSKarl Rupp 
3599566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(1.333333333333 * 6 * 6 * 6 * b->mbs)); /* from inverting diagonal blocks */
36083287d42SBarry Smith   PetscFunctionReturn(0);
36183287d42SBarry Smith }
362bef36659SShri Abhyankar 
363*d71ae5a4SJacob Faibussowitsch PetscErrorCode MatLUFactorNumeric_SeqBAIJ_6(Mat B, Mat A, const MatFactorInfo *info)
364*d71ae5a4SJacob Faibussowitsch {
365bef36659SShri Abhyankar   Mat             C = B;
366bef36659SShri Abhyankar   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data;
367bef36659SShri Abhyankar   IS              isrow = b->row, isicol = b->icol;
3685a586d82SBarry Smith   const PetscInt *r, *ic;
369bbd65245SShri Abhyankar   PetscInt        i, j, k, nz, nzL, row;
370bbd65245SShri Abhyankar   const PetscInt  n = a->mbs, *ai = a->i, *aj = a->j, *bi = b->i, *bj = b->j;
371bbd65245SShri Abhyankar   const PetscInt *ajtmp, *bjtmp, *bdiag = b->diag, *pj, bs2 = a->bs2;
372bef36659SShri Abhyankar   MatScalar      *rtmp, *pc, *mwork, *v, *pv, *aa = a->a;
373bbd65245SShri Abhyankar   PetscInt        flg;
374182b8fbaSHong Zhang   PetscReal       shift = info->shiftamount;
375a455e926SHong Zhang   PetscBool       allowzeropivot, zeropivotdetected;
376bef36659SShri Abhyankar 
377bef36659SShri Abhyankar   PetscFunctionBegin;
3780164db54SHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
3799566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isrow, &r));
3809566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isicol, &ic));
381bef36659SShri Abhyankar 
382bef36659SShri Abhyankar   /* generate work space needed by the factorization */
3839566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(bs2 * n, &rtmp, bs2, &mwork));
3849566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(rtmp, bs2 * n));
385bef36659SShri Abhyankar 
386bef36659SShri Abhyankar   for (i = 0; i < n; i++) {
387bef36659SShri Abhyankar     /* zero rtmp */
388bef36659SShri Abhyankar     /* L part */
389bef36659SShri Abhyankar     nz    = bi[i + 1] - bi[i];
390bef36659SShri Abhyankar     bjtmp = bj + bi[i];
39148a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2));
392bef36659SShri Abhyankar 
393bef36659SShri Abhyankar     /* U part */
3946506fda5SShri Abhyankar     nz    = bdiag[i] - bdiag[i + 1];
3956506fda5SShri Abhyankar     bjtmp = bj + bdiag[i + 1] + 1;
39648a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2));
3976506fda5SShri Abhyankar 
3986506fda5SShri Abhyankar     /* load in initial (unfactored row) */
3996506fda5SShri Abhyankar     nz    = ai[r[i] + 1] - ai[r[i]];
4006506fda5SShri Abhyankar     ajtmp = aj + ai[r[i]];
4016506fda5SShri Abhyankar     v     = aa + bs2 * ai[r[i]];
40248a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(rtmp + bs2 * ic[ajtmp[j]], v + bs2 * j, bs2));
4036506fda5SShri Abhyankar 
4046506fda5SShri Abhyankar     /* elimination */
4056506fda5SShri Abhyankar     bjtmp = bj + bi[i];
4066506fda5SShri Abhyankar     nzL   = bi[i + 1] - bi[i];
4076506fda5SShri Abhyankar     for (k = 0; k < nzL; k++) {
4086506fda5SShri Abhyankar       row = bjtmp[k];
4096506fda5SShri Abhyankar       pc  = rtmp + bs2 * row;
410c35f09e5SBarry Smith       for (flg = 0, j = 0; j < bs2; j++) {
411c35f09e5SBarry Smith         if (pc[j] != 0.0) {
412c35f09e5SBarry Smith           flg = 1;
413c35f09e5SBarry Smith           break;
414c35f09e5SBarry Smith         }
415c35f09e5SBarry Smith       }
4166506fda5SShri Abhyankar       if (flg) {
4176506fda5SShri Abhyankar         pv = b->a + bs2 * bdiag[row];
41896b95a6bSBarry Smith         /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
4199566063dSJacob Faibussowitsch         PetscCall(PetscKernel_A_gets_A_times_B_6(pc, pv, mwork));
4206506fda5SShri Abhyankar 
421a5b23f4aSJose E. Roman         pj = b->j + bdiag[row + 1] + 1; /* beginning of U(row,:) */
4226506fda5SShri Abhyankar         pv = b->a + bs2 * (bdiag[row + 1] + 1);
4236506fda5SShri Abhyankar         nz = bdiag[row] - bdiag[row + 1] - 1; /* num of entries inU(row,:), excluding diag */
4246506fda5SShri Abhyankar         for (j = 0; j < nz; j++) {
42596b95a6bSBarry Smith           /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
4266506fda5SShri Abhyankar           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
4276506fda5SShri Abhyankar           v = rtmp + bs2 * pj[j];
4289566063dSJacob Faibussowitsch           PetscCall(PetscKernel_A_gets_A_minus_B_times_C_6(v, pc, pv));
4296506fda5SShri Abhyankar           pv += bs2;
4306506fda5SShri Abhyankar         }
4319566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(432.0 * nz + 396)); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
4326506fda5SShri Abhyankar       }
4336506fda5SShri Abhyankar     }
4346506fda5SShri Abhyankar 
4356506fda5SShri Abhyankar     /* finished row so stick it into b->a */
4366506fda5SShri Abhyankar     /* L part */
4376506fda5SShri Abhyankar     pv = b->a + bs2 * bi[i];
4386506fda5SShri Abhyankar     pj = b->j + bi[i];
4396506fda5SShri Abhyankar     nz = bi[i + 1] - bi[i];
44048a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2));
4416506fda5SShri Abhyankar 
442a5b23f4aSJose E. Roman     /* Mark diagonal and invert diagonal for simpler triangular solves */
4436506fda5SShri Abhyankar     pv = b->a + bs2 * bdiag[i];
4446506fda5SShri Abhyankar     pj = b->j + bdiag[i];
4459566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(pv, rtmp + bs2 * pj[0], bs2));
4469566063dSJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_6(pv, shift, allowzeropivot, &zeropivotdetected));
4477b6c816cSBarry Smith     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
4486506fda5SShri Abhyankar 
4496506fda5SShri Abhyankar     /* U part */
4506506fda5SShri Abhyankar     pv = b->a + bs2 * (bdiag[i + 1] + 1);
4516506fda5SShri Abhyankar     pj = b->j + bdiag[i + 1] + 1;
4526506fda5SShri Abhyankar     nz = bdiag[i] - bdiag[i + 1] - 1;
45348a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2));
4546506fda5SShri Abhyankar   }
4556506fda5SShri Abhyankar 
4569566063dSJacob Faibussowitsch   PetscCall(PetscFree2(rtmp, mwork));
4579566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isicol, &ic));
4589566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow, &r));
45926fbe8dcSKarl Rupp 
4604dd39f65SShri Abhyankar   C->ops->solve          = MatSolve_SeqBAIJ_6;
4614dd39f65SShri Abhyankar   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_6;
4626506fda5SShri Abhyankar   C->assembled           = PETSC_TRUE;
46326fbe8dcSKarl Rupp 
4649566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(1.333333333333 * 6 * 6 * 6 * n)); /* from inverting diagonal blocks */
4656506fda5SShri Abhyankar   PetscFunctionReturn(0);
4666506fda5SShri Abhyankar }
4676506fda5SShri Abhyankar 
468*d71ae5a4SJacob Faibussowitsch PetscErrorCode MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering_inplace(Mat C, Mat A, const MatFactorInfo *info)
469*d71ae5a4SJacob Faibussowitsch {
470bef36659SShri Abhyankar   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data;
471bef36659SShri Abhyankar   PetscInt     i, j, n = a->mbs, *bi = b->i, *bj = b->j;
472bef36659SShri Abhyankar   PetscInt    *ajtmpold, *ajtmp, nz, row;
473bef36659SShri Abhyankar   PetscInt    *diag_offset = b->diag, *ai = a->i, *aj = a->j, *pj;
474bef36659SShri Abhyankar   MatScalar   *pv, *v, *rtmp, *pc, *w, *x;
475bef36659SShri Abhyankar   MatScalar    x1, x2, x3, x4, x5, x6, x7, x8, x9, x10, x11, x12, x13, x14, x15;
476bef36659SShri Abhyankar   MatScalar    x16, x17, x18, x19, x20, x21, x22, x23, x24, x25;
477bef36659SShri Abhyankar   MatScalar    p1, p2, p3, p4, p5, p6, p7, p8, p9, p10, p11, p12, p13, p14, p15;
478bef36659SShri Abhyankar   MatScalar    p16, p17, p18, p19, p20, p21, p22, p23, p24, p25;
479bef36659SShri Abhyankar   MatScalar    m1, m2, m3, m4, m5, m6, m7, m8, m9, m10, m11, m12, m13, m14, m15;
480bef36659SShri Abhyankar   MatScalar    m16, m17, m18, m19, m20, m21, m22, m23, m24, m25;
481bef36659SShri Abhyankar   MatScalar    p26, p27, p28, p29, p30, p31, p32, p33, p34, p35, p36;
482bef36659SShri Abhyankar   MatScalar    x26, x27, x28, x29, x30, x31, x32, x33, x34, x35, x36;
483bef36659SShri Abhyankar   MatScalar    m26, m27, m28, m29, m30, m31, m32, m33, m34, m35, m36;
484bef36659SShri Abhyankar   MatScalar   *ba = b->a, *aa = a->a;
485182b8fbaSHong Zhang   PetscReal    shift = info->shiftamount;
486a455e926SHong Zhang   PetscBool    allowzeropivot, zeropivotdetected;
487bef36659SShri Abhyankar 
488bef36659SShri Abhyankar   PetscFunctionBegin;
4890164db54SHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
4909566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(36 * (n + 1), &rtmp));
491bef36659SShri Abhyankar   for (i = 0; i < n; i++) {
492bef36659SShri Abhyankar     nz    = bi[i + 1] - bi[i];
493bef36659SShri Abhyankar     ajtmp = bj + bi[i];
494bef36659SShri Abhyankar     for (j = 0; j < nz; j++) {
495bef36659SShri Abhyankar       x    = rtmp + 36 * ajtmp[j];
496bef36659SShri Abhyankar       x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = x[9] = 0.0;
497bef36659SShri Abhyankar       x[10] = x[11] = x[12] = x[13] = x[14] = x[15] = x[16] = x[17] = 0.0;
498bef36659SShri Abhyankar       x[18] = x[19] = x[20] = x[21] = x[22] = x[23] = x[24] = x[25] = 0.0;
499bef36659SShri Abhyankar       x[26] = x[27] = x[28] = x[29] = x[30] = x[31] = x[32] = x[33] = 0.0;
500bef36659SShri Abhyankar       x[34] = x[35] = 0.0;
501bef36659SShri Abhyankar     }
502bef36659SShri Abhyankar     /* load in initial (unfactored row) */
503bef36659SShri Abhyankar     nz       = ai[i + 1] - ai[i];
504bef36659SShri Abhyankar     ajtmpold = aj + ai[i];
505bef36659SShri Abhyankar     v        = aa + 36 * ai[i];
506bef36659SShri Abhyankar     for (j = 0; j < nz; j++) {
507bef36659SShri Abhyankar       x     = rtmp + 36 * ajtmpold[j];
5089371c9d4SSatish Balay       x[0]  = v[0];
5099371c9d4SSatish Balay       x[1]  = v[1];
5109371c9d4SSatish Balay       x[2]  = v[2];
5119371c9d4SSatish Balay       x[3]  = v[3];
5129371c9d4SSatish Balay       x[4]  = v[4];
5139371c9d4SSatish Balay       x[5]  = v[5];
5149371c9d4SSatish Balay       x[6]  = v[6];
5159371c9d4SSatish Balay       x[7]  = v[7];
5169371c9d4SSatish Balay       x[8]  = v[8];
5179371c9d4SSatish Balay       x[9]  = v[9];
5189371c9d4SSatish Balay       x[10] = v[10];
5199371c9d4SSatish Balay       x[11] = v[11];
5209371c9d4SSatish Balay       x[12] = v[12];
5219371c9d4SSatish Balay       x[13] = v[13];
5229371c9d4SSatish Balay       x[14] = v[14];
5239371c9d4SSatish Balay       x[15] = v[15];
5249371c9d4SSatish Balay       x[16] = v[16];
5259371c9d4SSatish Balay       x[17] = v[17];
5269371c9d4SSatish Balay       x[18] = v[18];
5279371c9d4SSatish Balay       x[19] = v[19];
5289371c9d4SSatish Balay       x[20] = v[20];
5299371c9d4SSatish Balay       x[21] = v[21];
5309371c9d4SSatish Balay       x[22] = v[22];
5319371c9d4SSatish Balay       x[23] = v[23];
5329371c9d4SSatish Balay       x[24] = v[24];
5339371c9d4SSatish Balay       x[25] = v[25];
5349371c9d4SSatish Balay       x[26] = v[26];
5359371c9d4SSatish Balay       x[27] = v[27];
5369371c9d4SSatish Balay       x[28] = v[28];
5379371c9d4SSatish Balay       x[29] = v[29];
5389371c9d4SSatish Balay       x[30] = v[30];
5399371c9d4SSatish Balay       x[31] = v[31];
5409371c9d4SSatish Balay       x[32] = v[32];
5419371c9d4SSatish Balay       x[33] = v[33];
5429371c9d4SSatish Balay       x[34] = v[34];
5439371c9d4SSatish Balay       x[35] = v[35];
544bef36659SShri Abhyankar       v += 36;
545bef36659SShri Abhyankar     }
546bef36659SShri Abhyankar     row = *ajtmp++;
547bef36659SShri Abhyankar     while (row < i) {
548bef36659SShri Abhyankar       pc  = rtmp + 36 * row;
5499371c9d4SSatish Balay       p1  = pc[0];
5509371c9d4SSatish Balay       p2  = pc[1];
5519371c9d4SSatish Balay       p3  = pc[2];
5529371c9d4SSatish Balay       p4  = pc[3];
5539371c9d4SSatish Balay       p5  = pc[4];
5549371c9d4SSatish Balay       p6  = pc[5];
5559371c9d4SSatish Balay       p7  = pc[6];
5569371c9d4SSatish Balay       p8  = pc[7];
5579371c9d4SSatish Balay       p9  = pc[8];
5589371c9d4SSatish Balay       p10 = pc[9];
5599371c9d4SSatish Balay       p11 = pc[10];
5609371c9d4SSatish Balay       p12 = pc[11];
5619371c9d4SSatish Balay       p13 = pc[12];
5629371c9d4SSatish Balay       p14 = pc[13];
5639371c9d4SSatish Balay       p15 = pc[14];
5649371c9d4SSatish Balay       p16 = pc[15];
5659371c9d4SSatish Balay       p17 = pc[16];
5669371c9d4SSatish Balay       p18 = pc[17];
5679371c9d4SSatish Balay       p19 = pc[18];
5689371c9d4SSatish Balay       p20 = pc[19];
5699371c9d4SSatish Balay       p21 = pc[20];
5709371c9d4SSatish Balay       p22 = pc[21];
5719371c9d4SSatish Balay       p23 = pc[22];
5729371c9d4SSatish Balay       p24 = pc[23];
5739371c9d4SSatish Balay       p25 = pc[24];
5749371c9d4SSatish Balay       p26 = pc[25];
5759371c9d4SSatish Balay       p27 = pc[26];
5769371c9d4SSatish Balay       p28 = pc[27];
5779371c9d4SSatish Balay       p29 = pc[28];
5789371c9d4SSatish Balay       p30 = pc[29];
5799371c9d4SSatish Balay       p31 = pc[30];
5809371c9d4SSatish Balay       p32 = pc[31];
5819371c9d4SSatish Balay       p33 = pc[32];
5829371c9d4SSatish Balay       p34 = pc[33];
5839371c9d4SSatish Balay       p35 = pc[34];
5849371c9d4SSatish Balay       p36 = pc[35];
5859371c9d4SSatish 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) {
586bef36659SShri Abhyankar         pv    = ba + 36 * diag_offset[row];
587bef36659SShri Abhyankar         pj    = bj + diag_offset[row] + 1;
5889371c9d4SSatish Balay         x1    = pv[0];
5899371c9d4SSatish Balay         x2    = pv[1];
5909371c9d4SSatish Balay         x3    = pv[2];
5919371c9d4SSatish Balay         x4    = pv[3];
5929371c9d4SSatish Balay         x5    = pv[4];
5939371c9d4SSatish Balay         x6    = pv[5];
5949371c9d4SSatish Balay         x7    = pv[6];
5959371c9d4SSatish Balay         x8    = pv[7];
5969371c9d4SSatish Balay         x9    = pv[8];
5979371c9d4SSatish Balay         x10   = pv[9];
5989371c9d4SSatish Balay         x11   = pv[10];
5999371c9d4SSatish Balay         x12   = pv[11];
6009371c9d4SSatish Balay         x13   = pv[12];
6019371c9d4SSatish Balay         x14   = pv[13];
6029371c9d4SSatish Balay         x15   = pv[14];
6039371c9d4SSatish Balay         x16   = pv[15];
6049371c9d4SSatish Balay         x17   = pv[16];
6059371c9d4SSatish Balay         x18   = pv[17];
6069371c9d4SSatish Balay         x19   = pv[18];
6079371c9d4SSatish Balay         x20   = pv[19];
6089371c9d4SSatish Balay         x21   = pv[20];
6099371c9d4SSatish Balay         x22   = pv[21];
6109371c9d4SSatish Balay         x23   = pv[22];
6119371c9d4SSatish Balay         x24   = pv[23];
6129371c9d4SSatish Balay         x25   = pv[24];
6139371c9d4SSatish Balay         x26   = pv[25];
6149371c9d4SSatish Balay         x27   = pv[26];
6159371c9d4SSatish Balay         x28   = pv[27];
6169371c9d4SSatish Balay         x29   = pv[28];
6179371c9d4SSatish Balay         x30   = pv[29];
6189371c9d4SSatish Balay         x31   = pv[30];
6199371c9d4SSatish Balay         x32   = pv[31];
6209371c9d4SSatish Balay         x33   = pv[32];
6219371c9d4SSatish Balay         x34   = pv[33];
6229371c9d4SSatish Balay         x35   = pv[34];
6239371c9d4SSatish Balay         x36   = pv[35];
624bef36659SShri Abhyankar         pc[0] = m1 = p1 * x1 + p7 * x2 + p13 * x3 + p19 * x4 + p25 * x5 + p31 * x6;
625bef36659SShri Abhyankar         pc[1] = m2 = p2 * x1 + p8 * x2 + p14 * x3 + p20 * x4 + p26 * x5 + p32 * x6;
626bef36659SShri Abhyankar         pc[2] = m3 = p3 * x1 + p9 * x2 + p15 * x3 + p21 * x4 + p27 * x5 + p33 * x6;
627bef36659SShri Abhyankar         pc[3] = m4 = p4 * x1 + p10 * x2 + p16 * x3 + p22 * x4 + p28 * x5 + p34 * x6;
628bef36659SShri Abhyankar         pc[4] = m5 = p5 * x1 + p11 * x2 + p17 * x3 + p23 * x4 + p29 * x5 + p35 * x6;
629bef36659SShri Abhyankar         pc[5] = m6 = p6 * x1 + p12 * x2 + p18 * x3 + p24 * x4 + p30 * x5 + p36 * x6;
630bef36659SShri Abhyankar 
631bef36659SShri Abhyankar         pc[6] = m7 = p1 * x7 + p7 * x8 + p13 * x9 + p19 * x10 + p25 * x11 + p31 * x12;
632bef36659SShri Abhyankar         pc[7] = m8 = p2 * x7 + p8 * x8 + p14 * x9 + p20 * x10 + p26 * x11 + p32 * x12;
633bef36659SShri Abhyankar         pc[8] = m9 = p3 * x7 + p9 * x8 + p15 * x9 + p21 * x10 + p27 * x11 + p33 * x12;
634bef36659SShri Abhyankar         pc[9] = m10 = p4 * x7 + p10 * x8 + p16 * x9 + p22 * x10 + p28 * x11 + p34 * x12;
635bef36659SShri Abhyankar         pc[10] = m11 = p5 * x7 + p11 * x8 + p17 * x9 + p23 * x10 + p29 * x11 + p35 * x12;
636bef36659SShri Abhyankar         pc[11] = m12 = p6 * x7 + p12 * x8 + p18 * x9 + p24 * x10 + p30 * x11 + p36 * x12;
637bef36659SShri Abhyankar 
638bef36659SShri Abhyankar         pc[12] = m13 = p1 * x13 + p7 * x14 + p13 * x15 + p19 * x16 + p25 * x17 + p31 * x18;
639bef36659SShri Abhyankar         pc[13] = m14 = p2 * x13 + p8 * x14 + p14 * x15 + p20 * x16 + p26 * x17 + p32 * x18;
640bef36659SShri Abhyankar         pc[14] = m15 = p3 * x13 + p9 * x14 + p15 * x15 + p21 * x16 + p27 * x17 + p33 * x18;
641bef36659SShri Abhyankar         pc[15] = m16 = p4 * x13 + p10 * x14 + p16 * x15 + p22 * x16 + p28 * x17 + p34 * x18;
642bef36659SShri Abhyankar         pc[16] = m17 = p5 * x13 + p11 * x14 + p17 * x15 + p23 * x16 + p29 * x17 + p35 * x18;
643bef36659SShri Abhyankar         pc[17] = m18 = p6 * x13 + p12 * x14 + p18 * x15 + p24 * x16 + p30 * x17 + p36 * x18;
644bef36659SShri Abhyankar 
645bef36659SShri Abhyankar         pc[18] = m19 = p1 * x19 + p7 * x20 + p13 * x21 + p19 * x22 + p25 * x23 + p31 * x24;
646bef36659SShri Abhyankar         pc[19] = m20 = p2 * x19 + p8 * x20 + p14 * x21 + p20 * x22 + p26 * x23 + p32 * x24;
647bef36659SShri Abhyankar         pc[20] = m21 = p3 * x19 + p9 * x20 + p15 * x21 + p21 * x22 + p27 * x23 + p33 * x24;
648bef36659SShri Abhyankar         pc[21] = m22 = p4 * x19 + p10 * x20 + p16 * x21 + p22 * x22 + p28 * x23 + p34 * x24;
649bef36659SShri Abhyankar         pc[22] = m23 = p5 * x19 + p11 * x20 + p17 * x21 + p23 * x22 + p29 * x23 + p35 * x24;
650bef36659SShri Abhyankar         pc[23] = m24 = p6 * x19 + p12 * x20 + p18 * x21 + p24 * x22 + p30 * x23 + p36 * x24;
651bef36659SShri Abhyankar 
652bef36659SShri Abhyankar         pc[24] = m25 = p1 * x25 + p7 * x26 + p13 * x27 + p19 * x28 + p25 * x29 + p31 * x30;
653bef36659SShri Abhyankar         pc[25] = m26 = p2 * x25 + p8 * x26 + p14 * x27 + p20 * x28 + p26 * x29 + p32 * x30;
654bef36659SShri Abhyankar         pc[26] = m27 = p3 * x25 + p9 * x26 + p15 * x27 + p21 * x28 + p27 * x29 + p33 * x30;
655bef36659SShri Abhyankar         pc[27] = m28 = p4 * x25 + p10 * x26 + p16 * x27 + p22 * x28 + p28 * x29 + p34 * x30;
656bef36659SShri Abhyankar         pc[28] = m29 = p5 * x25 + p11 * x26 + p17 * x27 + p23 * x28 + p29 * x29 + p35 * x30;
657bef36659SShri Abhyankar         pc[29] = m30 = p6 * x25 + p12 * x26 + p18 * x27 + p24 * x28 + p30 * x29 + p36 * x30;
658bef36659SShri Abhyankar 
659bef36659SShri Abhyankar         pc[30] = m31 = p1 * x31 + p7 * x32 + p13 * x33 + p19 * x34 + p25 * x35 + p31 * x36;
660bef36659SShri Abhyankar         pc[31] = m32 = p2 * x31 + p8 * x32 + p14 * x33 + p20 * x34 + p26 * x35 + p32 * x36;
661bef36659SShri Abhyankar         pc[32] = m33 = p3 * x31 + p9 * x32 + p15 * x33 + p21 * x34 + p27 * x35 + p33 * x36;
662bef36659SShri Abhyankar         pc[33] = m34 = p4 * x31 + p10 * x32 + p16 * x33 + p22 * x34 + p28 * x35 + p34 * x36;
663bef36659SShri Abhyankar         pc[34] = m35 = p5 * x31 + p11 * x32 + p17 * x33 + p23 * x34 + p29 * x35 + p35 * x36;
664bef36659SShri Abhyankar         pc[35] = m36 = p6 * x31 + p12 * x32 + p18 * x33 + p24 * x34 + p30 * x35 + p36 * x36;
665bef36659SShri Abhyankar 
666bef36659SShri Abhyankar         nz = bi[row + 1] - diag_offset[row] - 1;
667bef36659SShri Abhyankar         pv += 36;
668bef36659SShri Abhyankar         for (j = 0; j < nz; j++) {
6699371c9d4SSatish Balay           x1  = pv[0];
6709371c9d4SSatish Balay           x2  = pv[1];
6719371c9d4SSatish Balay           x3  = pv[2];
6729371c9d4SSatish Balay           x4  = pv[3];
6739371c9d4SSatish Balay           x5  = pv[4];
6749371c9d4SSatish Balay           x6  = pv[5];
6759371c9d4SSatish Balay           x7  = pv[6];
6769371c9d4SSatish Balay           x8  = pv[7];
6779371c9d4SSatish Balay           x9  = pv[8];
6789371c9d4SSatish Balay           x10 = pv[9];
6799371c9d4SSatish Balay           x11 = pv[10];
6809371c9d4SSatish Balay           x12 = pv[11];
6819371c9d4SSatish Balay           x13 = pv[12];
6829371c9d4SSatish Balay           x14 = pv[13];
6839371c9d4SSatish Balay           x15 = pv[14];
6849371c9d4SSatish Balay           x16 = pv[15];
6859371c9d4SSatish Balay           x17 = pv[16];
6869371c9d4SSatish Balay           x18 = pv[17];
6879371c9d4SSatish Balay           x19 = pv[18];
6889371c9d4SSatish Balay           x20 = pv[19];
6899371c9d4SSatish Balay           x21 = pv[20];
6909371c9d4SSatish Balay           x22 = pv[21];
6919371c9d4SSatish Balay           x23 = pv[22];
6929371c9d4SSatish Balay           x24 = pv[23];
6939371c9d4SSatish Balay           x25 = pv[24];
6949371c9d4SSatish Balay           x26 = pv[25];
6959371c9d4SSatish Balay           x27 = pv[26];
6969371c9d4SSatish Balay           x28 = pv[27];
6979371c9d4SSatish Balay           x29 = pv[28];
6989371c9d4SSatish Balay           x30 = pv[29];
6999371c9d4SSatish Balay           x31 = pv[30];
7009371c9d4SSatish Balay           x32 = pv[31];
7019371c9d4SSatish Balay           x33 = pv[32];
7029371c9d4SSatish Balay           x34 = pv[33];
7039371c9d4SSatish Balay           x35 = pv[34];
7049371c9d4SSatish Balay           x36 = pv[35];
705bef36659SShri Abhyankar           x   = rtmp + 36 * pj[j];
706bef36659SShri Abhyankar           x[0] -= m1 * x1 + m7 * x2 + m13 * x3 + m19 * x4 + m25 * x5 + m31 * x6;
707bef36659SShri Abhyankar           x[1] -= m2 * x1 + m8 * x2 + m14 * x3 + m20 * x4 + m26 * x5 + m32 * x6;
708bef36659SShri Abhyankar           x[2] -= m3 * x1 + m9 * x2 + m15 * x3 + m21 * x4 + m27 * x5 + m33 * x6;
709bef36659SShri Abhyankar           x[3] -= m4 * x1 + m10 * x2 + m16 * x3 + m22 * x4 + m28 * x5 + m34 * x6;
710bef36659SShri Abhyankar           x[4] -= m5 * x1 + m11 * x2 + m17 * x3 + m23 * x4 + m29 * x5 + m35 * x6;
711bef36659SShri Abhyankar           x[5] -= m6 * x1 + m12 * x2 + m18 * x3 + m24 * x4 + m30 * x5 + m36 * x6;
712bef36659SShri Abhyankar 
713bef36659SShri Abhyankar           x[6] -= m1 * x7 + m7 * x8 + m13 * x9 + m19 * x10 + m25 * x11 + m31 * x12;
714bef36659SShri Abhyankar           x[7] -= m2 * x7 + m8 * x8 + m14 * x9 + m20 * x10 + m26 * x11 + m32 * x12;
715bef36659SShri Abhyankar           x[8] -= m3 * x7 + m9 * x8 + m15 * x9 + m21 * x10 + m27 * x11 + m33 * x12;
716bef36659SShri Abhyankar           x[9] -= m4 * x7 + m10 * x8 + m16 * x9 + m22 * x10 + m28 * x11 + m34 * x12;
717bef36659SShri Abhyankar           x[10] -= m5 * x7 + m11 * x8 + m17 * x9 + m23 * x10 + m29 * x11 + m35 * x12;
718bef36659SShri Abhyankar           x[11] -= m6 * x7 + m12 * x8 + m18 * x9 + m24 * x10 + m30 * x11 + m36 * x12;
719bef36659SShri Abhyankar 
720bef36659SShri Abhyankar           x[12] -= m1 * x13 + m7 * x14 + m13 * x15 + m19 * x16 + m25 * x17 + m31 * x18;
721bef36659SShri Abhyankar           x[13] -= m2 * x13 + m8 * x14 + m14 * x15 + m20 * x16 + m26 * x17 + m32 * x18;
722bef36659SShri Abhyankar           x[14] -= m3 * x13 + m9 * x14 + m15 * x15 + m21 * x16 + m27 * x17 + m33 * x18;
723bef36659SShri Abhyankar           x[15] -= m4 * x13 + m10 * x14 + m16 * x15 + m22 * x16 + m28 * x17 + m34 * x18;
724bef36659SShri Abhyankar           x[16] -= m5 * x13 + m11 * x14 + m17 * x15 + m23 * x16 + m29 * x17 + m35 * x18;
725bef36659SShri Abhyankar           x[17] -= m6 * x13 + m12 * x14 + m18 * x15 + m24 * x16 + m30 * x17 + m36 * x18;
726bef36659SShri Abhyankar 
727bef36659SShri Abhyankar           x[18] -= m1 * x19 + m7 * x20 + m13 * x21 + m19 * x22 + m25 * x23 + m31 * x24;
728bef36659SShri Abhyankar           x[19] -= m2 * x19 + m8 * x20 + m14 * x21 + m20 * x22 + m26 * x23 + m32 * x24;
729bef36659SShri Abhyankar           x[20] -= m3 * x19 + m9 * x20 + m15 * x21 + m21 * x22 + m27 * x23 + m33 * x24;
730bef36659SShri Abhyankar           x[21] -= m4 * x19 + m10 * x20 + m16 * x21 + m22 * x22 + m28 * x23 + m34 * x24;
731bef36659SShri Abhyankar           x[22] -= m5 * x19 + m11 * x20 + m17 * x21 + m23 * x22 + m29 * x23 + m35 * x24;
732bef36659SShri Abhyankar           x[23] -= m6 * x19 + m12 * x20 + m18 * x21 + m24 * x22 + m30 * x23 + m36 * x24;
733bef36659SShri Abhyankar 
734bef36659SShri Abhyankar           x[24] -= m1 * x25 + m7 * x26 + m13 * x27 + m19 * x28 + m25 * x29 + m31 * x30;
735bef36659SShri Abhyankar           x[25] -= m2 * x25 + m8 * x26 + m14 * x27 + m20 * x28 + m26 * x29 + m32 * x30;
736bef36659SShri Abhyankar           x[26] -= m3 * x25 + m9 * x26 + m15 * x27 + m21 * x28 + m27 * x29 + m33 * x30;
737bef36659SShri Abhyankar           x[27] -= m4 * x25 + m10 * x26 + m16 * x27 + m22 * x28 + m28 * x29 + m34 * x30;
738bef36659SShri Abhyankar           x[28] -= m5 * x25 + m11 * x26 + m17 * x27 + m23 * x28 + m29 * x29 + m35 * x30;
739bef36659SShri Abhyankar           x[29] -= m6 * x25 + m12 * x26 + m18 * x27 + m24 * x28 + m30 * x29 + m36 * x30;
740bef36659SShri Abhyankar 
741bef36659SShri Abhyankar           x[30] -= m1 * x31 + m7 * x32 + m13 * x33 + m19 * x34 + m25 * x35 + m31 * x36;
742bef36659SShri Abhyankar           x[31] -= m2 * x31 + m8 * x32 + m14 * x33 + m20 * x34 + m26 * x35 + m32 * x36;
743bef36659SShri Abhyankar           x[32] -= m3 * x31 + m9 * x32 + m15 * x33 + m21 * x34 + m27 * x35 + m33 * x36;
744bef36659SShri Abhyankar           x[33] -= m4 * x31 + m10 * x32 + m16 * x33 + m22 * x34 + m28 * x35 + m34 * x36;
745bef36659SShri Abhyankar           x[34] -= m5 * x31 + m11 * x32 + m17 * x33 + m23 * x34 + m29 * x35 + m35 * x36;
746bef36659SShri Abhyankar           x[35] -= m6 * x31 + m12 * x32 + m18 * x33 + m24 * x34 + m30 * x35 + m36 * x36;
747bef36659SShri Abhyankar 
748bef36659SShri Abhyankar           pv += 36;
749bef36659SShri Abhyankar         }
7509566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(432.0 * nz + 396.0));
751bef36659SShri Abhyankar       }
752bef36659SShri Abhyankar       row = *ajtmp++;
753bef36659SShri Abhyankar     }
754bef36659SShri Abhyankar     /* finished row so stick it into b->a */
755bef36659SShri Abhyankar     pv = ba + 36 * bi[i];
756bef36659SShri Abhyankar     pj = bj + bi[i];
757bef36659SShri Abhyankar     nz = bi[i + 1] - bi[i];
758bef36659SShri Abhyankar     for (j = 0; j < nz; j++) {
759bef36659SShri Abhyankar       x      = rtmp + 36 * pj[j];
7609371c9d4SSatish Balay       pv[0]  = x[0];
7619371c9d4SSatish Balay       pv[1]  = x[1];
7629371c9d4SSatish Balay       pv[2]  = x[2];
7639371c9d4SSatish Balay       pv[3]  = x[3];
7649371c9d4SSatish Balay       pv[4]  = x[4];
7659371c9d4SSatish Balay       pv[5]  = x[5];
7669371c9d4SSatish Balay       pv[6]  = x[6];
7679371c9d4SSatish Balay       pv[7]  = x[7];
7689371c9d4SSatish Balay       pv[8]  = x[8];
7699371c9d4SSatish Balay       pv[9]  = x[9];
7709371c9d4SSatish Balay       pv[10] = x[10];
7719371c9d4SSatish Balay       pv[11] = x[11];
7729371c9d4SSatish Balay       pv[12] = x[12];
7739371c9d4SSatish Balay       pv[13] = x[13];
7749371c9d4SSatish Balay       pv[14] = x[14];
7759371c9d4SSatish Balay       pv[15] = x[15];
7769371c9d4SSatish Balay       pv[16] = x[16];
7779371c9d4SSatish Balay       pv[17] = x[17];
7789371c9d4SSatish Balay       pv[18] = x[18];
7799371c9d4SSatish Balay       pv[19] = x[19];
7809371c9d4SSatish Balay       pv[20] = x[20];
7819371c9d4SSatish Balay       pv[21] = x[21];
7829371c9d4SSatish Balay       pv[22] = x[22];
7839371c9d4SSatish Balay       pv[23] = x[23];
7849371c9d4SSatish Balay       pv[24] = x[24];
7859371c9d4SSatish Balay       pv[25] = x[25];
7869371c9d4SSatish Balay       pv[26] = x[26];
7879371c9d4SSatish Balay       pv[27] = x[27];
7889371c9d4SSatish Balay       pv[28] = x[28];
7899371c9d4SSatish Balay       pv[29] = x[29];
7909371c9d4SSatish Balay       pv[30] = x[30];
7919371c9d4SSatish Balay       pv[31] = x[31];
7929371c9d4SSatish Balay       pv[32] = x[32];
7939371c9d4SSatish Balay       pv[33] = x[33];
7949371c9d4SSatish Balay       pv[34] = x[34];
7959371c9d4SSatish Balay       pv[35] = x[35];
796bef36659SShri Abhyankar       pv += 36;
797bef36659SShri Abhyankar     }
798bef36659SShri Abhyankar     /* invert diagonal block */
799bef36659SShri Abhyankar     w = ba + 36 * diag_offset[i];
8009566063dSJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_6(w, shift, allowzeropivot, &zeropivotdetected));
8017b6c816cSBarry Smith     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
802bef36659SShri Abhyankar   }
803bef36659SShri Abhyankar 
8049566063dSJacob Faibussowitsch   PetscCall(PetscFree(rtmp));
80526fbe8dcSKarl Rupp 
80606e38f1dSHong Zhang   C->ops->solve          = MatSolve_SeqBAIJ_6_NaturalOrdering_inplace;
80706e38f1dSHong Zhang   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_6_NaturalOrdering_inplace;
808bef36659SShri Abhyankar   C->assembled           = PETSC_TRUE;
80926fbe8dcSKarl Rupp 
8109566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(1.333333333333 * 6 * 6 * 6 * b->mbs)); /* from inverting diagonal blocks */
811bef36659SShri Abhyankar   PetscFunctionReturn(0);
812bef36659SShri Abhyankar }
813bef36659SShri Abhyankar 
814*d71ae5a4SJacob Faibussowitsch PetscErrorCode MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering(Mat B, Mat A, const MatFactorInfo *info)
815*d71ae5a4SJacob Faibussowitsch {
816bef36659SShri Abhyankar   Mat             C = B;
817bef36659SShri Abhyankar   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data;
818bbd65245SShri Abhyankar   PetscInt        i, j, k, nz, nzL, row;
819bbd65245SShri Abhyankar   const PetscInt  n = a->mbs, *ai = a->i, *aj = a->j, *bi = b->i, *bj = b->j;
820bbd65245SShri Abhyankar   const PetscInt *ajtmp, *bjtmp, *bdiag = b->diag, *pj, bs2 = a->bs2;
821bef36659SShri Abhyankar   MatScalar      *rtmp, *pc, *mwork, *v, *pv, *aa = a->a;
822bbd65245SShri Abhyankar   PetscInt        flg;
823182b8fbaSHong Zhang   PetscReal       shift = info->shiftamount;
824a455e926SHong Zhang   PetscBool       allowzeropivot, zeropivotdetected;
825bef36659SShri Abhyankar 
826bef36659SShri Abhyankar   PetscFunctionBegin;
8270164db54SHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
8280164db54SHong Zhang 
829bef36659SShri Abhyankar   /* generate work space needed by the factorization */
8309566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(bs2 * n, &rtmp, bs2, &mwork));
8319566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(rtmp, bs2 * n));
832bef36659SShri Abhyankar 
833bef36659SShri Abhyankar   for (i = 0; i < n; i++) {
834bef36659SShri Abhyankar     /* zero rtmp */
835bef36659SShri Abhyankar     /* L part */
836bef36659SShri Abhyankar     nz    = bi[i + 1] - bi[i];
837bef36659SShri Abhyankar     bjtmp = bj + bi[i];
83848a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2));
839bef36659SShri Abhyankar 
840bef36659SShri Abhyankar     /* U part */
84153cca76cSShri Abhyankar     nz    = bdiag[i] - bdiag[i + 1];
84253cca76cSShri Abhyankar     bjtmp = bj + bdiag[i + 1] + 1;
84348a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2));
84453cca76cSShri Abhyankar 
84553cca76cSShri Abhyankar     /* load in initial (unfactored row) */
84653cca76cSShri Abhyankar     nz    = ai[i + 1] - ai[i];
84753cca76cSShri Abhyankar     ajtmp = aj + ai[i];
84853cca76cSShri Abhyankar     v     = aa + bs2 * ai[i];
84948a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(rtmp + bs2 * ajtmp[j], v + bs2 * j, bs2));
85053cca76cSShri Abhyankar 
85153cca76cSShri Abhyankar     /* elimination */
85253cca76cSShri Abhyankar     bjtmp = bj + bi[i];
85353cca76cSShri Abhyankar     nzL   = bi[i + 1] - bi[i];
85453cca76cSShri Abhyankar     for (k = 0; k < nzL; k++) {
85553cca76cSShri Abhyankar       row = bjtmp[k];
85653cca76cSShri Abhyankar       pc  = rtmp + bs2 * row;
857c35f09e5SBarry Smith       for (flg = 0, j = 0; j < bs2; j++) {
858c35f09e5SBarry Smith         if (pc[j] != 0.0) {
859c35f09e5SBarry Smith           flg = 1;
860c35f09e5SBarry Smith           break;
861c35f09e5SBarry Smith         }
862c35f09e5SBarry Smith       }
86353cca76cSShri Abhyankar       if (flg) {
86453cca76cSShri Abhyankar         pv = b->a + bs2 * bdiag[row];
86596b95a6bSBarry Smith         /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
8669566063dSJacob Faibussowitsch         PetscCall(PetscKernel_A_gets_A_times_B_6(pc, pv, mwork));
86753cca76cSShri Abhyankar 
868a5b23f4aSJose E. Roman         pj = b->j + bdiag[row + 1] + 1; /* beginning of U(row,:) */
86953cca76cSShri Abhyankar         pv = b->a + bs2 * (bdiag[row + 1] + 1);
87053cca76cSShri Abhyankar         nz = bdiag[row] - bdiag[row + 1] - 1; /* num of entries inU(row,:), excluding diag */
87153cca76cSShri Abhyankar         for (j = 0; j < nz; j++) {
87296b95a6bSBarry Smith           /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
87353cca76cSShri Abhyankar           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
87453cca76cSShri Abhyankar           v = rtmp + bs2 * pj[j];
8759566063dSJacob Faibussowitsch           PetscCall(PetscKernel_A_gets_A_minus_B_times_C_6(v, pc, pv));
87653cca76cSShri Abhyankar           pv += bs2;
87753cca76cSShri Abhyankar         }
8789566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(432.0 * nz + 396)); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
87953cca76cSShri Abhyankar       }
88053cca76cSShri Abhyankar     }
88153cca76cSShri Abhyankar 
88253cca76cSShri Abhyankar     /* finished row so stick it into b->a */
88353cca76cSShri Abhyankar     /* L part */
88453cca76cSShri Abhyankar     pv = b->a + bs2 * bi[i];
88553cca76cSShri Abhyankar     pj = b->j + bi[i];
88653cca76cSShri Abhyankar     nz = bi[i + 1] - bi[i];
88748a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2));
88853cca76cSShri Abhyankar 
889a5b23f4aSJose E. Roman     /* Mark diagonal and invert diagonal for simpler triangular solves */
89053cca76cSShri Abhyankar     pv = b->a + bs2 * bdiag[i];
89153cca76cSShri Abhyankar     pj = b->j + bdiag[i];
8929566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(pv, rtmp + bs2 * pj[0], bs2));
8939566063dSJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_6(pv, shift, allowzeropivot, &zeropivotdetected));
8947b6c816cSBarry Smith     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
89553cca76cSShri Abhyankar 
89653cca76cSShri Abhyankar     /* U part */
89753cca76cSShri Abhyankar     pv = b->a + bs2 * (bdiag[i + 1] + 1);
89853cca76cSShri Abhyankar     pj = b->j + bdiag[i + 1] + 1;
89953cca76cSShri Abhyankar     nz = bdiag[i] - bdiag[i + 1] - 1;
90048a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2));
90153cca76cSShri Abhyankar   }
9029566063dSJacob Faibussowitsch   PetscCall(PetscFree2(rtmp, mwork));
90326fbe8dcSKarl Rupp 
9044dd39f65SShri Abhyankar   C->ops->solve          = MatSolve_SeqBAIJ_6_NaturalOrdering;
9054dd39f65SShri Abhyankar   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_6_NaturalOrdering;
90653cca76cSShri Abhyankar   C->assembled           = PETSC_TRUE;
90726fbe8dcSKarl Rupp 
9089566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(1.333333333333 * 6 * 6 * 6 * n)); /* from inverting diagonal blocks */
90953cca76cSShri Abhyankar   PetscFunctionReturn(0);
91053cca76cSShri Abhyankar }
911