xref: /petsc/src/mat/impls/baij/seq/baijfact9.c (revision e8c0849ab8fe171bed529bea27238c9b402db591)
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 5 by 5
983287d42SBarry Smith */
MatILUFactorNumeric_SeqBAIJ_5_inplace(Mat C,Mat A,const MatFactorInfo * info)10*421480d9SBarry Smith PetscErrorCode MatILUFactorNumeric_SeqBAIJ_5_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;
14*421480d9SBarry Smith   const PetscInt  *r, *ic;
15*421480d9SBarry Smith   PetscInt        *bi = b->i, *bj = b->j, *ajtmpold, *ajtmp;
16766f9fbaSBarry Smith   PetscInt         i, j, n = a->mbs, nz, row, idx, ipvt[5];
17*421480d9SBarry Smith   const PetscInt  *diag_offset;
18*421480d9SBarry Smith   PetscInt        *ai = a->i, *aj = a->j, *pj;
198397fe1aSBarry Smith   MatScalar       *w, *pv, *rtmp, *x, *pc;
208397fe1aSBarry Smith   const MatScalar *v, *aa = a->a;
2183287d42SBarry Smith   MatScalar        p1, p2, p3, p4, m1, m2, m3, m4, m5, m6, m7, m8, m9, x1, x2, x3, x4;
2283287d42SBarry Smith   MatScalar        p5, p6, p7, p8, p9, x5, x6, x7, x8, x9, x10, x11, x12, x13, x14, x15, x16;
2383287d42SBarry Smith   MatScalar        x17, x18, x19, x20, x21, x22, x23, x24, x25, p10, p11, p12, p13, p14;
2483287d42SBarry Smith   MatScalar        p15, p16, p17, p18, p19, p20, p21, p22, p23, p24, p25, m10, m11, m12;
2583287d42SBarry Smith   MatScalar        m13, m14, m15, m16, m17, m18, m19, m20, m21, m22, m23, m24, m25;
26766f9fbaSBarry Smith   MatScalar       *ba    = b->a, work[25];
27182b8fbaSHong Zhang   PetscReal        shift = info->shiftamount;
28a455e926SHong Zhang   PetscBool        allowzeropivot, zeropivotdetected;
2983287d42SBarry Smith 
3083287d42SBarry Smith   PetscFunctionBegin;
31*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 */
32*421480d9SBarry Smith   A->factortype = MAT_FACTOR_NONE;
33*421480d9SBarry Smith   PetscCall(MatGetDiagonalMarkers_SeqBAIJ(A, &diag_offset, NULL));
34*421480d9SBarry Smith   A->factortype  = MAT_FACTOR_ILU;
350164db54SHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
369566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isrow, &r));
379566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isicol, &ic));
389566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(25 * (n + 1), &rtmp));
3983287d42SBarry Smith 
408397fe1aSBarry Smith #define PETSC_USE_MEMZERO 1
418397fe1aSBarry Smith #define PETSC_USE_MEMCPY  1
428397fe1aSBarry Smith 
4383287d42SBarry Smith   for (i = 0; i < n; i++) {
4483287d42SBarry Smith     nz    = bi[i + 1] - bi[i];
4583287d42SBarry Smith     ajtmp = bj + bi[i];
4683287d42SBarry Smith     for (j = 0; j < nz; j++) {
478397fe1aSBarry Smith #if defined(PETSC_USE_MEMZERO)
489566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(rtmp + 25 * ajtmp[j], 25));
498397fe1aSBarry Smith #else
5083287d42SBarry Smith       x    = rtmp + 25 * ajtmp[j];
5183287d42SBarry Smith       x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = x[9] = 0.0;
5283287d42SBarry Smith       x[10] = x[11] = x[12] = x[13] = x[14] = x[15] = x[16] = x[17] = 0.0;
5383287d42SBarry Smith       x[18] = x[19] = x[20] = x[21] = x[22] = x[23] = x[24] = 0.0;
548397fe1aSBarry Smith #endif
5583287d42SBarry Smith     }
5683287d42SBarry Smith     /* load in initial (unfactored row) */
5783287d42SBarry Smith     idx      = r[i];
5883287d42SBarry Smith     nz       = ai[idx + 1] - ai[idx];
5983287d42SBarry Smith     ajtmpold = aj + ai[idx];
6083287d42SBarry Smith     v        = aa + 25 * ai[idx];
6183287d42SBarry Smith     for (j = 0; j < nz; j++) {
628397fe1aSBarry Smith #if defined(PETSC_USE_MEMCPY)
639566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(rtmp + 25 * ic[ajtmpold[j]], v, 25));
648397fe1aSBarry Smith #else
6583287d42SBarry Smith       x     = rtmp + 25 * ic[ajtmpold[j]];
669371c9d4SSatish Balay       x[0]  = v[0];
679371c9d4SSatish Balay       x[1]  = v[1];
689371c9d4SSatish Balay       x[2]  = v[2];
699371c9d4SSatish Balay       x[3]  = v[3];
709371c9d4SSatish Balay       x[4]  = v[4];
719371c9d4SSatish Balay       x[5]  = v[5];
729371c9d4SSatish Balay       x[6]  = v[6];
739371c9d4SSatish Balay       x[7]  = v[7];
749371c9d4SSatish Balay       x[8]  = v[8];
759371c9d4SSatish Balay       x[9]  = v[9];
769371c9d4SSatish Balay       x[10] = v[10];
779371c9d4SSatish Balay       x[11] = v[11];
789371c9d4SSatish Balay       x[12] = v[12];
799371c9d4SSatish Balay       x[13] = v[13];
809371c9d4SSatish Balay       x[14] = v[14];
819371c9d4SSatish Balay       x[15] = v[15];
829371c9d4SSatish Balay       x[16] = v[16];
839371c9d4SSatish Balay       x[17] = v[17];
849371c9d4SSatish Balay       x[18] = v[18];
859371c9d4SSatish Balay       x[19] = v[19];
869371c9d4SSatish Balay       x[20] = v[20];
879371c9d4SSatish Balay       x[21] = v[21];
889371c9d4SSatish Balay       x[22] = v[22];
899371c9d4SSatish Balay       x[23] = v[23];
909371c9d4SSatish Balay       x[24] = v[24];
918397fe1aSBarry Smith #endif
9283287d42SBarry Smith       v += 25;
9383287d42SBarry Smith     }
9483287d42SBarry Smith     row = *ajtmp++;
9583287d42SBarry Smith     while (row < i) {
9683287d42SBarry Smith       pc  = rtmp + 25 * row;
979371c9d4SSatish Balay       p1  = pc[0];
989371c9d4SSatish Balay       p2  = pc[1];
999371c9d4SSatish Balay       p3  = pc[2];
1009371c9d4SSatish Balay       p4  = pc[3];
1019371c9d4SSatish Balay       p5  = pc[4];
1029371c9d4SSatish Balay       p6  = pc[5];
1039371c9d4SSatish Balay       p7  = pc[6];
1049371c9d4SSatish Balay       p8  = pc[7];
1059371c9d4SSatish Balay       p9  = pc[8];
1069371c9d4SSatish Balay       p10 = pc[9];
1079371c9d4SSatish Balay       p11 = pc[10];
1089371c9d4SSatish Balay       p12 = pc[11];
1099371c9d4SSatish Balay       p13 = pc[12];
1109371c9d4SSatish Balay       p14 = pc[13];
1119371c9d4SSatish Balay       p15 = pc[14];
1129371c9d4SSatish Balay       p16 = pc[15];
1139371c9d4SSatish Balay       p17 = pc[16];
1149371c9d4SSatish Balay       p18 = pc[17];
1159371c9d4SSatish Balay       p19 = pc[18];
1169371c9d4SSatish Balay       p20 = pc[19];
1179371c9d4SSatish Balay       p21 = pc[20];
1189371c9d4SSatish Balay       p22 = pc[21];
1199371c9d4SSatish Balay       p23 = pc[22];
1209371c9d4SSatish Balay       p24 = pc[23];
12183287d42SBarry Smith       p25 = pc[24];
1229371c9d4SSatish 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) {
12383287d42SBarry Smith         pv    = ba + 25 * diag_offset[row];
12483287d42SBarry Smith         pj    = bj + diag_offset[row] + 1;
1259371c9d4SSatish Balay         x1    = pv[0];
1269371c9d4SSatish Balay         x2    = pv[1];
1279371c9d4SSatish Balay         x3    = pv[2];
1289371c9d4SSatish Balay         x4    = pv[3];
1299371c9d4SSatish Balay         x5    = pv[4];
1309371c9d4SSatish Balay         x6    = pv[5];
1319371c9d4SSatish Balay         x7    = pv[6];
1329371c9d4SSatish Balay         x8    = pv[7];
1339371c9d4SSatish Balay         x9    = pv[8];
1349371c9d4SSatish Balay         x10   = pv[9];
1359371c9d4SSatish Balay         x11   = pv[10];
1369371c9d4SSatish Balay         x12   = pv[11];
1379371c9d4SSatish Balay         x13   = pv[12];
1389371c9d4SSatish Balay         x14   = pv[13];
1399371c9d4SSatish Balay         x15   = pv[14];
1409371c9d4SSatish Balay         x16   = pv[15];
1419371c9d4SSatish Balay         x17   = pv[16];
1429371c9d4SSatish Balay         x18   = pv[17];
1439371c9d4SSatish Balay         x19   = pv[18];
1449371c9d4SSatish Balay         x20   = pv[19];
1459371c9d4SSatish Balay         x21   = pv[20];
1469371c9d4SSatish Balay         x22   = pv[21];
1479371c9d4SSatish Balay         x23   = pv[22];
1489371c9d4SSatish Balay         x24   = pv[23];
1499371c9d4SSatish Balay         x25   = pv[24];
15083287d42SBarry Smith         pc[0] = m1 = p1 * x1 + p6 * x2 + p11 * x3 + p16 * x4 + p21 * x5;
15183287d42SBarry Smith         pc[1] = m2 = p2 * x1 + p7 * x2 + p12 * x3 + p17 * x4 + p22 * x5;
15283287d42SBarry Smith         pc[2] = m3 = p3 * x1 + p8 * x2 + p13 * x3 + p18 * x4 + p23 * x5;
15383287d42SBarry Smith         pc[3] = m4 = p4 * x1 + p9 * x2 + p14 * x3 + p19 * x4 + p24 * x5;
15483287d42SBarry Smith         pc[4] = m5 = p5 * x1 + p10 * x2 + p15 * x3 + p20 * x4 + p25 * x5;
15583287d42SBarry Smith 
15683287d42SBarry Smith         pc[5] = m6 = p1 * x6 + p6 * x7 + p11 * x8 + p16 * x9 + p21 * x10;
15783287d42SBarry Smith         pc[6] = m7 = p2 * x6 + p7 * x7 + p12 * x8 + p17 * x9 + p22 * x10;
15883287d42SBarry Smith         pc[7] = m8 = p3 * x6 + p8 * x7 + p13 * x8 + p18 * x9 + p23 * x10;
15983287d42SBarry Smith         pc[8] = m9 = p4 * x6 + p9 * x7 + p14 * x8 + p19 * x9 + p24 * x10;
16083287d42SBarry Smith         pc[9] = m10 = p5 * x6 + p10 * x7 + p15 * x8 + p20 * x9 + p25 * x10;
16183287d42SBarry Smith 
16283287d42SBarry Smith         pc[10] = m11 = p1 * x11 + p6 * x12 + p11 * x13 + p16 * x14 + p21 * x15;
16383287d42SBarry Smith         pc[11] = m12 = p2 * x11 + p7 * x12 + p12 * x13 + p17 * x14 + p22 * x15;
16483287d42SBarry Smith         pc[12] = m13 = p3 * x11 + p8 * x12 + p13 * x13 + p18 * x14 + p23 * x15;
16583287d42SBarry Smith         pc[13] = m14 = p4 * x11 + p9 * x12 + p14 * x13 + p19 * x14 + p24 * x15;
16683287d42SBarry Smith         pc[14] = m15 = p5 * x11 + p10 * x12 + p15 * x13 + p20 * x14 + p25 * x15;
16783287d42SBarry Smith 
16883287d42SBarry Smith         pc[15] = m16 = p1 * x16 + p6 * x17 + p11 * x18 + p16 * x19 + p21 * x20;
16983287d42SBarry Smith         pc[16] = m17 = p2 * x16 + p7 * x17 + p12 * x18 + p17 * x19 + p22 * x20;
17083287d42SBarry Smith         pc[17] = m18 = p3 * x16 + p8 * x17 + p13 * x18 + p18 * x19 + p23 * x20;
17183287d42SBarry Smith         pc[18] = m19 = p4 * x16 + p9 * x17 + p14 * x18 + p19 * x19 + p24 * x20;
17283287d42SBarry Smith         pc[19] = m20 = p5 * x16 + p10 * x17 + p15 * x18 + p20 * x19 + p25 * x20;
17383287d42SBarry Smith 
17483287d42SBarry Smith         pc[20] = m21 = p1 * x21 + p6 * x22 + p11 * x23 + p16 * x24 + p21 * x25;
17583287d42SBarry Smith         pc[21] = m22 = p2 * x21 + p7 * x22 + p12 * x23 + p17 * x24 + p22 * x25;
17683287d42SBarry Smith         pc[22] = m23 = p3 * x21 + p8 * x22 + p13 * x23 + p18 * x24 + p23 * x25;
17783287d42SBarry Smith         pc[23] = m24 = p4 * x21 + p9 * x22 + p14 * x23 + p19 * x24 + p24 * x25;
17883287d42SBarry Smith         pc[24] = m25 = p5 * x21 + p10 * x22 + p15 * x23 + p20 * x24 + p25 * x25;
17983287d42SBarry Smith 
18083287d42SBarry Smith         nz = bi[row + 1] - diag_offset[row] - 1;
18183287d42SBarry Smith         pv += 25;
18283287d42SBarry Smith         for (j = 0; j < nz; j++) {
1839371c9d4SSatish Balay           x1  = pv[0];
1849371c9d4SSatish Balay           x2  = pv[1];
1859371c9d4SSatish Balay           x3  = pv[2];
1869371c9d4SSatish Balay           x4  = pv[3];
1879371c9d4SSatish Balay           x5  = pv[4];
1889371c9d4SSatish Balay           x6  = pv[5];
1899371c9d4SSatish Balay           x7  = pv[6];
1909371c9d4SSatish Balay           x8  = pv[7];
1919371c9d4SSatish Balay           x9  = pv[8];
1929371c9d4SSatish Balay           x10 = pv[9];
1939371c9d4SSatish Balay           x11 = pv[10];
1949371c9d4SSatish Balay           x12 = pv[11];
1959371c9d4SSatish Balay           x13 = pv[12];
1969371c9d4SSatish Balay           x14 = pv[13];
1979371c9d4SSatish Balay           x15 = pv[14];
1989371c9d4SSatish Balay           x16 = pv[15];
1999371c9d4SSatish Balay           x17 = pv[16];
2009371c9d4SSatish Balay           x18 = pv[17];
2019371c9d4SSatish Balay           x19 = pv[18];
2029371c9d4SSatish Balay           x20 = pv[19];
2039371c9d4SSatish Balay           x21 = pv[20];
2049371c9d4SSatish Balay           x22 = pv[21];
2059371c9d4SSatish Balay           x23 = pv[22];
2069371c9d4SSatish Balay           x24 = pv[23];
2079371c9d4SSatish Balay           x25 = pv[24];
20883287d42SBarry Smith           x   = rtmp + 25 * pj[j];
20983287d42SBarry Smith           x[0] -= m1 * x1 + m6 * x2 + m11 * x3 + m16 * x4 + m21 * x5;
21083287d42SBarry Smith           x[1] -= m2 * x1 + m7 * x2 + m12 * x3 + m17 * x4 + m22 * x5;
21183287d42SBarry Smith           x[2] -= m3 * x1 + m8 * x2 + m13 * x3 + m18 * x4 + m23 * x5;
21283287d42SBarry Smith           x[3] -= m4 * x1 + m9 * x2 + m14 * x3 + m19 * x4 + m24 * x5;
21383287d42SBarry Smith           x[4] -= m5 * x1 + m10 * x2 + m15 * x3 + m20 * x4 + m25 * x5;
21483287d42SBarry Smith 
21583287d42SBarry Smith           x[5] -= m1 * x6 + m6 * x7 + m11 * x8 + m16 * x9 + m21 * x10;
21683287d42SBarry Smith           x[6] -= m2 * x6 + m7 * x7 + m12 * x8 + m17 * x9 + m22 * x10;
21783287d42SBarry Smith           x[7] -= m3 * x6 + m8 * x7 + m13 * x8 + m18 * x9 + m23 * x10;
21883287d42SBarry Smith           x[8] -= m4 * x6 + m9 * x7 + m14 * x8 + m19 * x9 + m24 * x10;
21983287d42SBarry Smith           x[9] -= m5 * x6 + m10 * x7 + m15 * x8 + m20 * x9 + m25 * x10;
22083287d42SBarry Smith 
22183287d42SBarry Smith           x[10] -= m1 * x11 + m6 * x12 + m11 * x13 + m16 * x14 + m21 * x15;
22283287d42SBarry Smith           x[11] -= m2 * x11 + m7 * x12 + m12 * x13 + m17 * x14 + m22 * x15;
22383287d42SBarry Smith           x[12] -= m3 * x11 + m8 * x12 + m13 * x13 + m18 * x14 + m23 * x15;
22483287d42SBarry Smith           x[13] -= m4 * x11 + m9 * x12 + m14 * x13 + m19 * x14 + m24 * x15;
22583287d42SBarry Smith           x[14] -= m5 * x11 + m10 * x12 + m15 * x13 + m20 * x14 + m25 * x15;
22683287d42SBarry Smith 
22783287d42SBarry Smith           x[15] -= m1 * x16 + m6 * x17 + m11 * x18 + m16 * x19 + m21 * x20;
22883287d42SBarry Smith           x[16] -= m2 * x16 + m7 * x17 + m12 * x18 + m17 * x19 + m22 * x20;
22983287d42SBarry Smith           x[17] -= m3 * x16 + m8 * x17 + m13 * x18 + m18 * x19 + m23 * x20;
23083287d42SBarry Smith           x[18] -= m4 * x16 + m9 * x17 + m14 * x18 + m19 * x19 + m24 * x20;
23183287d42SBarry Smith           x[19] -= m5 * x16 + m10 * x17 + m15 * x18 + m20 * x19 + m25 * x20;
23283287d42SBarry Smith 
23383287d42SBarry Smith           x[20] -= m1 * x21 + m6 * x22 + m11 * x23 + m16 * x24 + m21 * x25;
23483287d42SBarry Smith           x[21] -= m2 * x21 + m7 * x22 + m12 * x23 + m17 * x24 + m22 * x25;
23583287d42SBarry Smith           x[22] -= m3 * x21 + m8 * x22 + m13 * x23 + m18 * x24 + m23 * x25;
23683287d42SBarry Smith           x[23] -= m4 * x21 + m9 * x22 + m14 * x23 + m19 * x24 + m24 * x25;
23783287d42SBarry Smith           x[24] -= m5 * x21 + m10 * x22 + m15 * x23 + m20 * x24 + m25 * x25;
23883287d42SBarry Smith 
23983287d42SBarry Smith           pv += 25;
24083287d42SBarry Smith         }
2419566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(250.0 * nz + 225.0));
24283287d42SBarry Smith       }
24383287d42SBarry Smith       row = *ajtmp++;
24483287d42SBarry Smith     }
24583287d42SBarry Smith     /* finished row so stick it into b->a */
24683287d42SBarry Smith     pv = ba + 25 * bi[i];
24783287d42SBarry Smith     pj = bj + bi[i];
24883287d42SBarry Smith     nz = bi[i + 1] - bi[i];
24983287d42SBarry Smith     for (j = 0; j < nz; j++) {
2508397fe1aSBarry Smith #if defined(PETSC_USE_MEMCPY)
2519566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(pv, rtmp + 25 * pj[j], 25));
2528397fe1aSBarry Smith #else
25383287d42SBarry Smith       x      = rtmp + 25 * pj[j];
2549371c9d4SSatish Balay       pv[0]  = x[0];
2559371c9d4SSatish Balay       pv[1]  = x[1];
2569371c9d4SSatish Balay       pv[2]  = x[2];
2579371c9d4SSatish Balay       pv[3]  = x[3];
2589371c9d4SSatish Balay       pv[4]  = x[4];
2599371c9d4SSatish Balay       pv[5]  = x[5];
2609371c9d4SSatish Balay       pv[6]  = x[6];
2619371c9d4SSatish Balay       pv[7]  = x[7];
2629371c9d4SSatish Balay       pv[8]  = x[8];
2639371c9d4SSatish Balay       pv[9]  = x[9];
2649371c9d4SSatish Balay       pv[10] = x[10];
2659371c9d4SSatish Balay       pv[11] = x[11];
2669371c9d4SSatish Balay       pv[12] = x[12];
2679371c9d4SSatish Balay       pv[13] = x[13];
2689371c9d4SSatish Balay       pv[14] = x[14];
2699371c9d4SSatish Balay       pv[15] = x[15];
2709371c9d4SSatish Balay       pv[16] = x[16];
2719371c9d4SSatish Balay       pv[17] = x[17];
2729371c9d4SSatish Balay       pv[18] = x[18];
2739371c9d4SSatish Balay       pv[19] = x[19];
2749371c9d4SSatish Balay       pv[20] = x[20];
2759371c9d4SSatish Balay       pv[21] = x[21];
2769371c9d4SSatish Balay       pv[22] = x[22];
2779371c9d4SSatish Balay       pv[23] = x[23];
2789371c9d4SSatish Balay       pv[24] = x[24];
2798397fe1aSBarry Smith #endif
28083287d42SBarry Smith       pv += 25;
28183287d42SBarry Smith     }
28283287d42SBarry Smith     /* invert diagonal block */
28383287d42SBarry Smith     w = ba + 25 * diag_offset[i];
2849566063dSJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_5(w, ipvt, work, shift, allowzeropivot, &zeropivotdetected));
2857b6c816cSBarry Smith     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
28683287d42SBarry Smith   }
28783287d42SBarry Smith 
2889566063dSJacob Faibussowitsch   PetscCall(PetscFree(rtmp));
2899566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isicol, &ic));
2909566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow, &r));
29126fbe8dcSKarl Rupp 
29206e38f1dSHong Zhang   C->ops->solve          = MatSolve_SeqBAIJ_5_inplace;
29306e38f1dSHong Zhang   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_5_inplace;
29483287d42SBarry Smith   C->assembled           = PETSC_TRUE;
29526fbe8dcSKarl Rupp 
2969566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(1.333333333333 * 5 * 5 * 5 * b->mbs)); /* from inverting diagonal blocks */
2973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
29883287d42SBarry Smith }
2990deeaf61SShri Abhyankar 
3004dd39f65SShri Abhyankar /* MatLUFactorNumeric_SeqBAIJ_5 -
3014dd39f65SShri Abhyankar      copied from MatLUFactorNumeric_SeqBAIJ_N_inplace() and manually re-implemented
30296b95a6bSBarry Smith        PetscKernel_A_gets_A_times_B()
30396b95a6bSBarry Smith        PetscKernel_A_gets_A_minus_B_times_C()
30496b95a6bSBarry Smith        PetscKernel_A_gets_inverse_A()
3050deeaf61SShri Abhyankar */
306c0c7eb62SShri Abhyankar 
MatLUFactorNumeric_SeqBAIJ_5(Mat B,Mat A,const MatFactorInfo * info)307d71ae5a4SJacob Faibussowitsch PetscErrorCode MatLUFactorNumeric_SeqBAIJ_5(Mat B, Mat A, const MatFactorInfo *info)
308d71ae5a4SJacob Faibussowitsch {
3090deeaf61SShri Abhyankar   Mat             C = B;
3100deeaf61SShri Abhyankar   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data;
3110deeaf61SShri Abhyankar   IS              isrow = b->row, isicol = b->icol;
3125a586d82SBarry Smith   const PetscInt *r, *ic;
313bbd65245SShri Abhyankar   PetscInt        i, j, k, nz, nzL, row;
314bbd65245SShri Abhyankar   const PetscInt  n = a->mbs, *ai = a->i, *aj = a->j, *bi = b->i, *bj = b->j;
315bbd65245SShri Abhyankar   const PetscInt *ajtmp, *bjtmp, *bdiag = b->diag, *pj, bs2 = a->bs2;
316766f9fbaSBarry Smith   MatScalar      *rtmp, *pc, *mwork, *v, *pv, *aa = a->a, work[25];
317bbd65245SShri Abhyankar   PetscInt        flg, ipvt[5];
318182b8fbaSHong Zhang   PetscReal       shift = info->shiftamount;
319a455e926SHong Zhang   PetscBool       allowzeropivot, zeropivotdetected;
3200deeaf61SShri Abhyankar 
3210deeaf61SShri Abhyankar   PetscFunctionBegin;
3220164db54SHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
3239566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isrow, &r));
3249566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isicol, &ic));
3250deeaf61SShri Abhyankar 
3260deeaf61SShri Abhyankar   /* generate work space needed by the factorization */
3279566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(bs2 * n, &rtmp, bs2, &mwork));
3289566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(rtmp, bs2 * n));
3290deeaf61SShri Abhyankar 
3300deeaf61SShri Abhyankar   for (i = 0; i < n; i++) {
3310deeaf61SShri Abhyankar     /* zero rtmp */
3320deeaf61SShri Abhyankar     /* L part */
3330deeaf61SShri Abhyankar     nz    = bi[i + 1] - bi[i];
3340deeaf61SShri Abhyankar     bjtmp = bj + bi[i];
33548a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2));
3360deeaf61SShri Abhyankar 
3370deeaf61SShri Abhyankar     /* U part */
33878bb4007SShri Abhyankar     nz    = bdiag[i] - bdiag[i + 1];
33978bb4007SShri Abhyankar     bjtmp = bj + bdiag[i + 1] + 1;
34048a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2));
34178bb4007SShri Abhyankar 
34278bb4007SShri Abhyankar     /* load in initial (unfactored row) */
34378bb4007SShri Abhyankar     nz    = ai[r[i] + 1] - ai[r[i]];
34478bb4007SShri Abhyankar     ajtmp = aj + ai[r[i]];
34578bb4007SShri Abhyankar     v     = aa + bs2 * ai[r[i]];
34648a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(rtmp + bs2 * ic[ajtmp[j]], v + bs2 * j, bs2));
34778bb4007SShri Abhyankar 
34878bb4007SShri Abhyankar     /* elimination */
34978bb4007SShri Abhyankar     bjtmp = bj + bi[i];
35078bb4007SShri Abhyankar     nzL   = bi[i + 1] - bi[i];
35178bb4007SShri Abhyankar     for (k = 0; k < nzL; k++) {
35278bb4007SShri Abhyankar       row = bjtmp[k];
35378bb4007SShri Abhyankar       pc  = rtmp + bs2 * row;
354c35f09e5SBarry Smith       for (flg = 0, j = 0; j < bs2; j++) {
355c35f09e5SBarry Smith         if (pc[j] != 0.0) {
356c35f09e5SBarry Smith           flg = 1;
357c35f09e5SBarry Smith           break;
358c35f09e5SBarry Smith         }
359c35f09e5SBarry Smith       }
36078bb4007SShri Abhyankar       if (flg) {
36178bb4007SShri Abhyankar         pv = b->a + bs2 * bdiag[row];
36296b95a6bSBarry Smith         /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
3639566063dSJacob Faibussowitsch         PetscCall(PetscKernel_A_gets_A_times_B_5(pc, pv, mwork));
36478bb4007SShri Abhyankar 
365a5b23f4aSJose E. Roman         pj = b->j + bdiag[row + 1] + 1; /* beginning of U(row,:) */
36678bb4007SShri Abhyankar         pv = b->a + bs2 * (bdiag[row + 1] + 1);
36778bb4007SShri Abhyankar         nz = bdiag[row] - bdiag[row + 1] - 1; /* num of entries inU(row,:), excluding diag */
36878bb4007SShri Abhyankar         for (j = 0; j < nz; j++) {
36996b95a6bSBarry Smith           /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
37078bb4007SShri Abhyankar           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
37178bb4007SShri Abhyankar           v = rtmp + bs2 * pj[j];
3729566063dSJacob Faibussowitsch           PetscCall(PetscKernel_A_gets_A_minus_B_times_C_5(v, pc, pv));
37378bb4007SShri Abhyankar           pv += bs2;
37478bb4007SShri Abhyankar         }
3759566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(250.0 * nz + 225)); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
37678bb4007SShri Abhyankar       }
37778bb4007SShri Abhyankar     }
37878bb4007SShri Abhyankar 
37978bb4007SShri Abhyankar     /* finished row so stick it into b->a */
38078bb4007SShri Abhyankar     /* L part */
38178bb4007SShri Abhyankar     pv = b->a + bs2 * bi[i];
38278bb4007SShri Abhyankar     pj = b->j + bi[i];
38378bb4007SShri Abhyankar     nz = bi[i + 1] - bi[i];
38448a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2));
38578bb4007SShri Abhyankar 
386a5b23f4aSJose E. Roman     /* Mark diagonal and invert diagonal for simpler triangular solves */
38778bb4007SShri Abhyankar     pv = b->a + bs2 * bdiag[i];
38878bb4007SShri Abhyankar     pj = b->j + bdiag[i];
3899566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(pv, rtmp + bs2 * pj[0], bs2));
3909566063dSJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_5(pv, ipvt, work, shift, allowzeropivot, &zeropivotdetected));
3917b6c816cSBarry Smith     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
39278bb4007SShri Abhyankar 
39378bb4007SShri Abhyankar     /* U part */
39478bb4007SShri Abhyankar     pv = b->a + bs2 * (bdiag[i + 1] + 1);
39578bb4007SShri Abhyankar     pj = b->j + bdiag[i + 1] + 1;
39678bb4007SShri Abhyankar     nz = bdiag[i] - bdiag[i + 1] - 1;
39748a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2));
39878bb4007SShri Abhyankar   }
39978bb4007SShri Abhyankar 
4009566063dSJacob Faibussowitsch   PetscCall(PetscFree2(rtmp, mwork));
4019566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isicol, &ic));
4029566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow, &r));
40326fbe8dcSKarl Rupp 
4044dd39f65SShri Abhyankar   C->ops->solve          = MatSolve_SeqBAIJ_5;
4054dd39f65SShri Abhyankar   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_5;
40678bb4007SShri Abhyankar   C->assembled           = PETSC_TRUE;
40726fbe8dcSKarl Rupp 
4089566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(1.333333333333 * 5 * 5 * 5 * n)); /* from inverting diagonal blocks */
4093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
41078bb4007SShri Abhyankar }
41178bb4007SShri Abhyankar 
4120deeaf61SShri Abhyankar /*
4130deeaf61SShri Abhyankar       Version for when blocks are 5 by 5 Using natural ordering
4140deeaf61SShri Abhyankar */
MatILUFactorNumeric_SeqBAIJ_5_NaturalOrdering_inplace(Mat C,Mat A,const MatFactorInfo * info)415*421480d9SBarry Smith PetscErrorCode MatILUFactorNumeric_SeqBAIJ_5_NaturalOrdering_inplace(Mat C, Mat A, const MatFactorInfo *info)
416d71ae5a4SJacob Faibussowitsch {
4170deeaf61SShri Abhyankar   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data;
418766f9fbaSBarry Smith   PetscInt     i, j, n = a->mbs, *bi = b->i, *bj = b->j, ipvt[5];
4190deeaf61SShri Abhyankar   PetscInt    *ajtmpold, *ajtmp, nz, row;
4200deeaf61SShri Abhyankar   PetscInt    *diag_offset = b->diag, *ai = a->i, *aj = a->j, *pj;
4210deeaf61SShri Abhyankar   MatScalar   *pv, *v, *rtmp, *pc, *w, *x;
4220deeaf61SShri Abhyankar   MatScalar    x1, x2, x3, x4, x5, x6, x7, x8, x9, x10, x11, x12, x13, x14, x15;
4230deeaf61SShri Abhyankar   MatScalar    x16, x17, x18, x19, x20, x21, x22, x23, x24, x25;
4240deeaf61SShri Abhyankar   MatScalar    p1, p2, p3, p4, p5, p6, p7, p8, p9, p10, p11, p12, p13, p14, p15;
4250deeaf61SShri Abhyankar   MatScalar    p16, p17, p18, p19, p20, p21, p22, p23, p24, p25;
4260deeaf61SShri Abhyankar   MatScalar    m1, m2, m3, m4, m5, m6, m7, m8, m9, m10, m11, m12, m13, m14, m15;
4270deeaf61SShri Abhyankar   MatScalar    m16, m17, m18, m19, m20, m21, m22, m23, m24, m25;
428766f9fbaSBarry Smith   MatScalar   *ba = b->a, *aa = a->a, work[25];
429182b8fbaSHong Zhang   PetscReal    shift = info->shiftamount;
430a455e926SHong Zhang   PetscBool    allowzeropivot, zeropivotdetected;
4310deeaf61SShri Abhyankar 
4320deeaf61SShri Abhyankar   PetscFunctionBegin;
4330164db54SHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
4349566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(25 * (n + 1), &rtmp));
4350deeaf61SShri Abhyankar   for (i = 0; i < n; i++) {
4360deeaf61SShri Abhyankar     nz    = bi[i + 1] - bi[i];
4370deeaf61SShri Abhyankar     ajtmp = bj + bi[i];
4380deeaf61SShri Abhyankar     for (j = 0; j < nz; j++) {
4390deeaf61SShri Abhyankar       x    = rtmp + 25 * ajtmp[j];
4400deeaf61SShri Abhyankar       x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = x[9] = 0.0;
4410deeaf61SShri Abhyankar       x[10] = x[11] = x[12] = x[13] = x[14] = x[15] = 0.0;
4420deeaf61SShri Abhyankar       x[16] = x[17] = x[18] = x[19] = x[20] = x[21] = x[22] = x[23] = x[24] = 0.0;
4430deeaf61SShri Abhyankar     }
4440deeaf61SShri Abhyankar     /* load in initial (unfactored row) */
4450deeaf61SShri Abhyankar     nz       = ai[i + 1] - ai[i];
4460deeaf61SShri Abhyankar     ajtmpold = aj + ai[i];
4470deeaf61SShri Abhyankar     v        = aa + 25 * ai[i];
4480deeaf61SShri Abhyankar     for (j = 0; j < nz; j++) {
4490deeaf61SShri Abhyankar       x     = rtmp + 25 * ajtmpold[j];
4509371c9d4SSatish Balay       x[0]  = v[0];
4519371c9d4SSatish Balay       x[1]  = v[1];
4529371c9d4SSatish Balay       x[2]  = v[2];
4539371c9d4SSatish Balay       x[3]  = v[3];
4549371c9d4SSatish Balay       x[4]  = v[4];
4559371c9d4SSatish Balay       x[5]  = v[5];
4569371c9d4SSatish Balay       x[6]  = v[6];
4579371c9d4SSatish Balay       x[7]  = v[7];
4589371c9d4SSatish Balay       x[8]  = v[8];
4599371c9d4SSatish Balay       x[9]  = v[9];
4609371c9d4SSatish Balay       x[10] = v[10];
4619371c9d4SSatish Balay       x[11] = v[11];
4629371c9d4SSatish Balay       x[12] = v[12];
4639371c9d4SSatish Balay       x[13] = v[13];
4649371c9d4SSatish Balay       x[14] = v[14];
4659371c9d4SSatish Balay       x[15] = v[15];
4669371c9d4SSatish Balay       x[16] = v[16];
4679371c9d4SSatish Balay       x[17] = v[17];
4689371c9d4SSatish Balay       x[18] = v[18];
4699371c9d4SSatish Balay       x[19] = v[19];
4709371c9d4SSatish Balay       x[20] = v[20];
4719371c9d4SSatish Balay       x[21] = v[21];
4729371c9d4SSatish Balay       x[22] = v[22];
4739371c9d4SSatish Balay       x[23] = v[23];
4740deeaf61SShri Abhyankar       x[24] = v[24];
4750deeaf61SShri Abhyankar       v += 25;
4760deeaf61SShri Abhyankar     }
4770deeaf61SShri Abhyankar     row = *ajtmp++;
4780deeaf61SShri Abhyankar     while (row < i) {
4790deeaf61SShri Abhyankar       pc  = rtmp + 25 * row;
4809371c9d4SSatish Balay       p1  = pc[0];
4819371c9d4SSatish Balay       p2  = pc[1];
4829371c9d4SSatish Balay       p3  = pc[2];
4839371c9d4SSatish Balay       p4  = pc[3];
4849371c9d4SSatish Balay       p5  = pc[4];
4859371c9d4SSatish Balay       p6  = pc[5];
4869371c9d4SSatish Balay       p7  = pc[6];
4879371c9d4SSatish Balay       p8  = pc[7];
4889371c9d4SSatish Balay       p9  = pc[8];
4899371c9d4SSatish Balay       p10 = pc[9];
4909371c9d4SSatish Balay       p11 = pc[10];
4919371c9d4SSatish Balay       p12 = pc[11];
4929371c9d4SSatish Balay       p13 = pc[12];
4939371c9d4SSatish Balay       p14 = pc[13];
4949371c9d4SSatish Balay       p15 = pc[14];
4959371c9d4SSatish Balay       p16 = pc[15];
4969371c9d4SSatish Balay       p17 = pc[16];
4979371c9d4SSatish Balay       p18 = pc[17];
4989371c9d4SSatish Balay       p19 = pc[18];
4999371c9d4SSatish Balay       p20 = pc[19];
5009371c9d4SSatish Balay       p21 = pc[20];
5019371c9d4SSatish Balay       p22 = pc[21];
5029371c9d4SSatish Balay       p23 = pc[22];
5039371c9d4SSatish Balay       p24 = pc[23];
5049371c9d4SSatish Balay       p25 = pc[24];
5059371c9d4SSatish 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) {
5060deeaf61SShri Abhyankar         pv    = ba + 25 * diag_offset[row];
5070deeaf61SShri Abhyankar         pj    = bj + diag_offset[row] + 1;
5089371c9d4SSatish Balay         x1    = pv[0];
5099371c9d4SSatish Balay         x2    = pv[1];
5109371c9d4SSatish Balay         x3    = pv[2];
5119371c9d4SSatish Balay         x4    = pv[3];
5129371c9d4SSatish Balay         x5    = pv[4];
5139371c9d4SSatish Balay         x6    = pv[5];
5149371c9d4SSatish Balay         x7    = pv[6];
5159371c9d4SSatish Balay         x8    = pv[7];
5169371c9d4SSatish Balay         x9    = pv[8];
5179371c9d4SSatish Balay         x10   = pv[9];
5189371c9d4SSatish Balay         x11   = pv[10];
5199371c9d4SSatish Balay         x12   = pv[11];
5209371c9d4SSatish Balay         x13   = pv[12];
5219371c9d4SSatish Balay         x14   = pv[13];
5229371c9d4SSatish Balay         x15   = pv[14];
5239371c9d4SSatish Balay         x16   = pv[15];
5249371c9d4SSatish Balay         x17   = pv[16];
5259371c9d4SSatish Balay         x18   = pv[17];
5269371c9d4SSatish Balay         x19   = pv[18];
5279371c9d4SSatish Balay         x20   = pv[19];
5289371c9d4SSatish Balay         x21   = pv[20];
5299371c9d4SSatish Balay         x22   = pv[21];
5309371c9d4SSatish Balay         x23   = pv[22];
5319371c9d4SSatish Balay         x24   = pv[23];
5320deeaf61SShri Abhyankar         x25   = pv[24];
5330deeaf61SShri Abhyankar         pc[0] = m1 = p1 * x1 + p6 * x2 + p11 * x3 + p16 * x4 + p21 * x5;
5340deeaf61SShri Abhyankar         pc[1] = m2 = p2 * x1 + p7 * x2 + p12 * x3 + p17 * x4 + p22 * x5;
5350deeaf61SShri Abhyankar         pc[2] = m3 = p3 * x1 + p8 * x2 + p13 * x3 + p18 * x4 + p23 * x5;
5360deeaf61SShri Abhyankar         pc[3] = m4 = p4 * x1 + p9 * x2 + p14 * x3 + p19 * x4 + p24 * x5;
5370deeaf61SShri Abhyankar         pc[4] = m5 = p5 * x1 + p10 * x2 + p15 * x3 + p20 * x4 + p25 * x5;
5380deeaf61SShri Abhyankar 
5390deeaf61SShri Abhyankar         pc[5] = m6 = p1 * x6 + p6 * x7 + p11 * x8 + p16 * x9 + p21 * x10;
5400deeaf61SShri Abhyankar         pc[6] = m7 = p2 * x6 + p7 * x7 + p12 * x8 + p17 * x9 + p22 * x10;
5410deeaf61SShri Abhyankar         pc[7] = m8 = p3 * x6 + p8 * x7 + p13 * x8 + p18 * x9 + p23 * x10;
5420deeaf61SShri Abhyankar         pc[8] = m9 = p4 * x6 + p9 * x7 + p14 * x8 + p19 * x9 + p24 * x10;
5430deeaf61SShri Abhyankar         pc[9] = m10 = p5 * x6 + p10 * x7 + p15 * x8 + p20 * x9 + p25 * x10;
5440deeaf61SShri Abhyankar 
5450deeaf61SShri Abhyankar         pc[10] = m11 = p1 * x11 + p6 * x12 + p11 * x13 + p16 * x14 + p21 * x15;
5460deeaf61SShri Abhyankar         pc[11] = m12 = p2 * x11 + p7 * x12 + p12 * x13 + p17 * x14 + p22 * x15;
5470deeaf61SShri Abhyankar         pc[12] = m13 = p3 * x11 + p8 * x12 + p13 * x13 + p18 * x14 + p23 * x15;
5480deeaf61SShri Abhyankar         pc[13] = m14 = p4 * x11 + p9 * x12 + p14 * x13 + p19 * x14 + p24 * x15;
5490deeaf61SShri Abhyankar         pc[14] = m15 = p5 * x11 + p10 * x12 + p15 * x13 + p20 * x14 + p25 * x15;
5500deeaf61SShri Abhyankar 
5510deeaf61SShri Abhyankar         pc[15] = m16 = p1 * x16 + p6 * x17 + p11 * x18 + p16 * x19 + p21 * x20;
5520deeaf61SShri Abhyankar         pc[16] = m17 = p2 * x16 + p7 * x17 + p12 * x18 + p17 * x19 + p22 * x20;
5530deeaf61SShri Abhyankar         pc[17] = m18 = p3 * x16 + p8 * x17 + p13 * x18 + p18 * x19 + p23 * x20;
5540deeaf61SShri Abhyankar         pc[18] = m19 = p4 * x16 + p9 * x17 + p14 * x18 + p19 * x19 + p24 * x20;
5550deeaf61SShri Abhyankar         pc[19] = m20 = p5 * x16 + p10 * x17 + p15 * x18 + p20 * x19 + p25 * x20;
5560deeaf61SShri Abhyankar 
5570deeaf61SShri Abhyankar         pc[20] = m21 = p1 * x21 + p6 * x22 + p11 * x23 + p16 * x24 + p21 * x25;
5580deeaf61SShri Abhyankar         pc[21] = m22 = p2 * x21 + p7 * x22 + p12 * x23 + p17 * x24 + p22 * x25;
5590deeaf61SShri Abhyankar         pc[22] = m23 = p3 * x21 + p8 * x22 + p13 * x23 + p18 * x24 + p23 * x25;
5600deeaf61SShri Abhyankar         pc[23] = m24 = p4 * x21 + p9 * x22 + p14 * x23 + p19 * x24 + p24 * x25;
5610deeaf61SShri Abhyankar         pc[24] = m25 = p5 * x21 + p10 * x22 + p15 * x23 + p20 * x24 + p25 * x25;
5620deeaf61SShri Abhyankar 
5630deeaf61SShri Abhyankar         nz = bi[row + 1] - diag_offset[row] - 1;
5640deeaf61SShri Abhyankar         pv += 25;
5650deeaf61SShri Abhyankar         for (j = 0; j < nz; j++) {
5669371c9d4SSatish Balay           x1  = pv[0];
5679371c9d4SSatish Balay           x2  = pv[1];
5689371c9d4SSatish Balay           x3  = pv[2];
5699371c9d4SSatish Balay           x4  = pv[3];
5709371c9d4SSatish Balay           x5  = pv[4];
5719371c9d4SSatish Balay           x6  = pv[5];
5729371c9d4SSatish Balay           x7  = pv[6];
5739371c9d4SSatish Balay           x8  = pv[7];
5749371c9d4SSatish Balay           x9  = pv[8];
5759371c9d4SSatish Balay           x10 = pv[9];
5769371c9d4SSatish Balay           x11 = pv[10];
5779371c9d4SSatish Balay           x12 = pv[11];
5789371c9d4SSatish Balay           x13 = pv[12];
5799371c9d4SSatish Balay           x14 = pv[13];
5809371c9d4SSatish Balay           x15 = pv[14];
5819371c9d4SSatish Balay           x16 = pv[15];
5829371c9d4SSatish Balay           x17 = pv[16];
5839371c9d4SSatish Balay           x18 = pv[17];
5849371c9d4SSatish Balay           x19 = pv[18];
5859371c9d4SSatish Balay           x20 = pv[19];
5869371c9d4SSatish Balay           x21 = pv[20];
5879371c9d4SSatish Balay           x22 = pv[21];
5889371c9d4SSatish Balay           x23 = pv[22];
5899371c9d4SSatish Balay           x24 = pv[23];
5909371c9d4SSatish Balay           x25 = pv[24];
5910deeaf61SShri Abhyankar           x   = rtmp + 25 * pj[j];
5920deeaf61SShri Abhyankar           x[0] -= m1 * x1 + m6 * x2 + m11 * x3 + m16 * x4 + m21 * x5;
5930deeaf61SShri Abhyankar           x[1] -= m2 * x1 + m7 * x2 + m12 * x3 + m17 * x4 + m22 * x5;
5940deeaf61SShri Abhyankar           x[2] -= m3 * x1 + m8 * x2 + m13 * x3 + m18 * x4 + m23 * x5;
5950deeaf61SShri Abhyankar           x[3] -= m4 * x1 + m9 * x2 + m14 * x3 + m19 * x4 + m24 * x5;
5960deeaf61SShri Abhyankar           x[4] -= m5 * x1 + m10 * x2 + m15 * x3 + m20 * x4 + m25 * x5;
5970deeaf61SShri Abhyankar 
5980deeaf61SShri Abhyankar           x[5] -= m1 * x6 + m6 * x7 + m11 * x8 + m16 * x9 + m21 * x10;
5990deeaf61SShri Abhyankar           x[6] -= m2 * x6 + m7 * x7 + m12 * x8 + m17 * x9 + m22 * x10;
6000deeaf61SShri Abhyankar           x[7] -= m3 * x6 + m8 * x7 + m13 * x8 + m18 * x9 + m23 * x10;
6010deeaf61SShri Abhyankar           x[8] -= m4 * x6 + m9 * x7 + m14 * x8 + m19 * x9 + m24 * x10;
6020deeaf61SShri Abhyankar           x[9] -= m5 * x6 + m10 * x7 + m15 * x8 + m20 * x9 + m25 * x10;
6030deeaf61SShri Abhyankar 
6040deeaf61SShri Abhyankar           x[10] -= m1 * x11 + m6 * x12 + m11 * x13 + m16 * x14 + m21 * x15;
6050deeaf61SShri Abhyankar           x[11] -= m2 * x11 + m7 * x12 + m12 * x13 + m17 * x14 + m22 * x15;
6060deeaf61SShri Abhyankar           x[12] -= m3 * x11 + m8 * x12 + m13 * x13 + m18 * x14 + m23 * x15;
6070deeaf61SShri Abhyankar           x[13] -= m4 * x11 + m9 * x12 + m14 * x13 + m19 * x14 + m24 * x15;
6080deeaf61SShri Abhyankar           x[14] -= m5 * x11 + m10 * x12 + m15 * x13 + m20 * x14 + m25 * x15;
6090deeaf61SShri Abhyankar 
6100deeaf61SShri Abhyankar           x[15] -= m1 * x16 + m6 * x17 + m11 * x18 + m16 * x19 + m21 * x20;
6110deeaf61SShri Abhyankar           x[16] -= m2 * x16 + m7 * x17 + m12 * x18 + m17 * x19 + m22 * x20;
6120deeaf61SShri Abhyankar           x[17] -= m3 * x16 + m8 * x17 + m13 * x18 + m18 * x19 + m23 * x20;
6130deeaf61SShri Abhyankar           x[18] -= m4 * x16 + m9 * x17 + m14 * x18 + m19 * x19 + m24 * x20;
6140deeaf61SShri Abhyankar           x[19] -= m5 * x16 + m10 * x17 + m15 * x18 + m20 * x19 + m25 * x20;
6150deeaf61SShri Abhyankar 
6160deeaf61SShri Abhyankar           x[20] -= m1 * x21 + m6 * x22 + m11 * x23 + m16 * x24 + m21 * x25;
6170deeaf61SShri Abhyankar           x[21] -= m2 * x21 + m7 * x22 + m12 * x23 + m17 * x24 + m22 * x25;
6180deeaf61SShri Abhyankar           x[22] -= m3 * x21 + m8 * x22 + m13 * x23 + m18 * x24 + m23 * x25;
6190deeaf61SShri Abhyankar           x[23] -= m4 * x21 + m9 * x22 + m14 * x23 + m19 * x24 + m24 * x25;
6200deeaf61SShri Abhyankar           x[24] -= m5 * x21 + m10 * x22 + m15 * x23 + m20 * x24 + m25 * x25;
6210deeaf61SShri Abhyankar           pv += 25;
6220deeaf61SShri Abhyankar         }
6239566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(250.0 * nz + 225.0));
6240deeaf61SShri Abhyankar       }
6250deeaf61SShri Abhyankar       row = *ajtmp++;
6260deeaf61SShri Abhyankar     }
6270deeaf61SShri Abhyankar     /* finished row so stick it into b->a */
6280deeaf61SShri Abhyankar     pv = ba + 25 * bi[i];
6290deeaf61SShri Abhyankar     pj = bj + bi[i];
6300deeaf61SShri Abhyankar     nz = bi[i + 1] - bi[i];
6310deeaf61SShri Abhyankar     for (j = 0; j < nz; j++) {
6320deeaf61SShri Abhyankar       x      = rtmp + 25 * pj[j];
6339371c9d4SSatish Balay       pv[0]  = x[0];
6349371c9d4SSatish Balay       pv[1]  = x[1];
6359371c9d4SSatish Balay       pv[2]  = x[2];
6369371c9d4SSatish Balay       pv[3]  = x[3];
6379371c9d4SSatish Balay       pv[4]  = x[4];
6389371c9d4SSatish Balay       pv[5]  = x[5];
6399371c9d4SSatish Balay       pv[6]  = x[6];
6409371c9d4SSatish Balay       pv[7]  = x[7];
6419371c9d4SSatish Balay       pv[8]  = x[8];
6429371c9d4SSatish Balay       pv[9]  = x[9];
6439371c9d4SSatish Balay       pv[10] = x[10];
6449371c9d4SSatish Balay       pv[11] = x[11];
6459371c9d4SSatish Balay       pv[12] = x[12];
6469371c9d4SSatish Balay       pv[13] = x[13];
6479371c9d4SSatish Balay       pv[14] = x[14];
6489371c9d4SSatish Balay       pv[15] = x[15];
6499371c9d4SSatish Balay       pv[16] = x[16];
6509371c9d4SSatish Balay       pv[17] = x[17];
6519371c9d4SSatish Balay       pv[18] = x[18];
6529371c9d4SSatish Balay       pv[19] = x[19];
6539371c9d4SSatish Balay       pv[20] = x[20];
6549371c9d4SSatish Balay       pv[21] = x[21];
6559371c9d4SSatish Balay       pv[22] = x[22];
6569371c9d4SSatish Balay       pv[23] = x[23];
6579371c9d4SSatish Balay       pv[24] = x[24];
6580deeaf61SShri Abhyankar       pv += 25;
6590deeaf61SShri Abhyankar     }
6600deeaf61SShri Abhyankar     /* invert diagonal block */
6610deeaf61SShri Abhyankar     w = ba + 25 * diag_offset[i];
6629566063dSJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_5(w, ipvt, work, shift, allowzeropivot, &zeropivotdetected));
6637b6c816cSBarry Smith     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
6640deeaf61SShri Abhyankar   }
6650deeaf61SShri Abhyankar 
6669566063dSJacob Faibussowitsch   PetscCall(PetscFree(rtmp));
66726fbe8dcSKarl Rupp 
66806e38f1dSHong Zhang   C->ops->solve          = MatSolve_SeqBAIJ_5_NaturalOrdering_inplace;
66906e38f1dSHong Zhang   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_5_NaturalOrdering_inplace;
6700deeaf61SShri Abhyankar   C->assembled           = PETSC_TRUE;
67126fbe8dcSKarl Rupp 
6729566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(1.333333333333 * 5 * 5 * 5 * b->mbs)); /* from inverting diagonal blocks */
6733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6740deeaf61SShri Abhyankar }
6750deeaf61SShri Abhyankar 
MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering(Mat B,Mat A,const MatFactorInfo * info)676d71ae5a4SJacob Faibussowitsch PetscErrorCode MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering(Mat B, Mat A, const MatFactorInfo *info)
677d71ae5a4SJacob Faibussowitsch {
6780deeaf61SShri Abhyankar   Mat             C = B;
6790deeaf61SShri Abhyankar   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data;
680bbd65245SShri Abhyankar   PetscInt        i, j, k, nz, nzL, row;
681766f9fbaSBarry Smith   const PetscInt  n = a->mbs, *ai = a->i, *aj = a->j, *bi = b->i, *bj = b->j;
682bbd65245SShri Abhyankar   const PetscInt *ajtmp, *bjtmp, *bdiag = b->diag, *pj, bs2 = a->bs2;
683bbd65245SShri Abhyankar   MatScalar      *rtmp, *pc, *mwork, *v, *vv, *pv, *aa = a->a, work[25];
684bbd65245SShri Abhyankar   PetscInt        flg, ipvt[5];
685182b8fbaSHong Zhang   PetscReal       shift = info->shiftamount;
686a455e926SHong Zhang   PetscBool       allowzeropivot, zeropivotdetected;
6870deeaf61SShri Abhyankar 
6880deeaf61SShri Abhyankar   PetscFunctionBegin;
6890164db54SHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
6900164db54SHong Zhang 
6910deeaf61SShri Abhyankar   /* generate work space needed by the factorization */
6929566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(bs2 * n, &rtmp, bs2, &mwork));
6939566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(rtmp, bs2 * n));
6940deeaf61SShri Abhyankar 
6950deeaf61SShri Abhyankar   for (i = 0; i < n; i++) {
6960deeaf61SShri Abhyankar     /* zero rtmp */
6970deeaf61SShri Abhyankar     /* L part */
6980deeaf61SShri Abhyankar     nz    = bi[i + 1] - bi[i];
6990deeaf61SShri Abhyankar     bjtmp = bj + bi[i];
70048a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2));
7010deeaf61SShri Abhyankar 
7020deeaf61SShri Abhyankar     /* U part */
70353cca76cSShri Abhyankar     nz    = bdiag[i] - bdiag[i + 1];
70453cca76cSShri Abhyankar     bjtmp = bj + bdiag[i + 1] + 1;
70548a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2));
70653cca76cSShri Abhyankar 
70753cca76cSShri Abhyankar     /* load in initial (unfactored row) */
70853cca76cSShri Abhyankar     nz    = ai[i + 1] - ai[i];
70953cca76cSShri Abhyankar     ajtmp = aj + ai[i];
71053cca76cSShri Abhyankar     v     = aa + bs2 * ai[i];
71148a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(rtmp + bs2 * ajtmp[j], v + bs2 * j, bs2));
71253cca76cSShri Abhyankar 
71353cca76cSShri Abhyankar     /* elimination */
71453cca76cSShri Abhyankar     bjtmp = bj + bi[i];
71553cca76cSShri Abhyankar     nzL   = bi[i + 1] - bi[i];
71653cca76cSShri Abhyankar     for (k = 0; k < nzL; k++) {
71753cca76cSShri Abhyankar       row = bjtmp[k];
71853cca76cSShri Abhyankar       pc  = rtmp + bs2 * row;
719c35f09e5SBarry Smith       for (flg = 0, j = 0; j < bs2; j++) {
720c35f09e5SBarry Smith         if (pc[j] != 0.0) {
721c35f09e5SBarry Smith           flg = 1;
722c35f09e5SBarry Smith           break;
723c35f09e5SBarry Smith         }
724c35f09e5SBarry Smith       }
72553cca76cSShri Abhyankar       if (flg) {
72653cca76cSShri Abhyankar         pv = b->a + bs2 * bdiag[row];
72796b95a6bSBarry Smith         /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
7289566063dSJacob Faibussowitsch         PetscCall(PetscKernel_A_gets_A_times_B_5(pc, pv, mwork));
72953cca76cSShri Abhyankar 
730a5b23f4aSJose E. Roman         pj = b->j + bdiag[row + 1] + 1; /* beginning of U(row,:) */
73153cca76cSShri Abhyankar         pv = b->a + bs2 * (bdiag[row + 1] + 1);
73253cca76cSShri Abhyankar         nz = bdiag[row] - bdiag[row + 1] - 1; /* num of entries inU(row,:), excluding diag */
73353cca76cSShri Abhyankar         for (j = 0; j < nz; j++) {
73496b95a6bSBarry Smith           /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
73553cca76cSShri Abhyankar           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
736766f9fbaSBarry Smith           vv = rtmp + bs2 * pj[j];
7379566063dSJacob Faibussowitsch           PetscCall(PetscKernel_A_gets_A_minus_B_times_C_5(vv, pc, pv));
73853cca76cSShri Abhyankar           pv += bs2;
73953cca76cSShri Abhyankar         }
7409566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(250.0 * nz + 225)); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
74153cca76cSShri Abhyankar       }
74253cca76cSShri Abhyankar     }
74353cca76cSShri Abhyankar 
74453cca76cSShri Abhyankar     /* finished row so stick it into b->a */
74553cca76cSShri Abhyankar     /* L part */
74653cca76cSShri Abhyankar     pv = b->a + bs2 * bi[i];
74753cca76cSShri Abhyankar     pj = b->j + bi[i];
74853cca76cSShri Abhyankar     nz = bi[i + 1] - bi[i];
74948a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2));
75053cca76cSShri Abhyankar 
751a5b23f4aSJose E. Roman     /* Mark diagonal and invert diagonal for simpler triangular solves */
75253cca76cSShri Abhyankar     pv = b->a + bs2 * bdiag[i];
75353cca76cSShri Abhyankar     pj = b->j + bdiag[i];
7549566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(pv, rtmp + bs2 * pj[0], bs2));
7559566063dSJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_5(pv, ipvt, work, shift, allowzeropivot, &zeropivotdetected));
7567b6c816cSBarry Smith     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
75753cca76cSShri Abhyankar 
75853cca76cSShri Abhyankar     /* U part */
75953cca76cSShri Abhyankar     pv = b->a + bs2 * (bdiag[i + 1] + 1);
76053cca76cSShri Abhyankar     pj = b->j + bdiag[i + 1] + 1;
76153cca76cSShri Abhyankar     nz = bdiag[i] - bdiag[i + 1] - 1;
76248a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2));
76353cca76cSShri Abhyankar   }
7649566063dSJacob Faibussowitsch   PetscCall(PetscFree2(rtmp, mwork));
76526fbe8dcSKarl Rupp 
7664dd39f65SShri Abhyankar   C->ops->solve          = MatSolve_SeqBAIJ_5_NaturalOrdering;
7674dd39f65SShri Abhyankar   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_5_NaturalOrdering;
76853cca76cSShri Abhyankar   C->assembled           = PETSC_TRUE;
76926fbe8dcSKarl Rupp 
7709566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(1.333333333333 * 5 * 5 * 5 * n)); /* from inverting diagonal blocks */
7713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
77253cca76cSShri Abhyankar }
773c3034d77SJose E. Roman 
774c3034d77SJose E. Roman /*
775c3034d77SJose E. Roman    Version for when blocks are 9 by 9
776c3034d77SJose E. Roman  */
777c3034d77SJose E. Roman #if defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX2__) && defined(__FMA__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES)
778c3034d77SJose E. Roman   #include <immintrin.h>
MatLUFactorNumeric_SeqBAIJ_9_NaturalOrdering(Mat B,Mat A,const MatFactorInfo * info)779c3034d77SJose E. Roman PetscErrorCode MatLUFactorNumeric_SeqBAIJ_9_NaturalOrdering(Mat B, Mat A, const MatFactorInfo *info)
780c3034d77SJose E. Roman {
781c3034d77SJose E. Roman   Mat             C = B;
782c3034d77SJose E. Roman   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data;
783c3034d77SJose E. Roman   PetscInt        i, j, k, nz, nzL, row;
784c3034d77SJose E. Roman   const PetscInt  n = a->mbs, *ai = a->i, *aj = a->j, *bi = b->i, *bj = b->j;
785c3034d77SJose E. Roman   const PetscInt *ajtmp, *bjtmp, *bdiag = b->diag, *pj, bs2 = a->bs2;
786c3034d77SJose E. Roman   MatScalar      *rtmp, *pc, *mwork, *v, *pv, *aa = a->a;
787c3034d77SJose E. Roman   PetscInt        flg;
788c3034d77SJose E. Roman   PetscReal       shift = info->shiftamount;
789c3034d77SJose E. Roman   PetscBool       allowzeropivot, zeropivotdetected;
790c3034d77SJose E. Roman 
791c3034d77SJose E. Roman   PetscFunctionBegin;
792c3034d77SJose E. Roman   allowzeropivot = PetscNot(A->erroriffailure);
793c3034d77SJose E. Roman 
794c3034d77SJose E. Roman   /* generate work space needed by the factorization */
795c3034d77SJose E. Roman   PetscCall(PetscMalloc2(bs2 * n, &rtmp, bs2, &mwork));
796c3034d77SJose E. Roman   PetscCall(PetscArrayzero(rtmp, bs2 * n));
797c3034d77SJose E. Roman 
798c3034d77SJose E. Roman   for (i = 0; i < n; i++) {
799c3034d77SJose E. Roman     /* zero rtmp */
800c3034d77SJose E. Roman     /* L part */
801c3034d77SJose E. Roman     nz    = bi[i + 1] - bi[i];
802c3034d77SJose E. Roman     bjtmp = bj + bi[i];
803c3034d77SJose E. Roman     for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2));
804c3034d77SJose E. Roman 
805c3034d77SJose E. Roman     /* U part */
806c3034d77SJose E. Roman     nz    = bdiag[i] - bdiag[i + 1];
807c3034d77SJose E. Roman     bjtmp = bj + bdiag[i + 1] + 1;
808c3034d77SJose E. Roman     for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2));
809c3034d77SJose E. Roman 
810c3034d77SJose E. Roman     /* load in initial (unfactored row) */
811c3034d77SJose E. Roman     nz    = ai[i + 1] - ai[i];
812c3034d77SJose E. Roman     ajtmp = aj + ai[i];
813c3034d77SJose E. Roman     v     = aa + bs2 * ai[i];
814c3034d77SJose E. Roman     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(rtmp + bs2 * ajtmp[j], v + bs2 * j, bs2));
815c3034d77SJose E. Roman 
816c3034d77SJose E. Roman     /* elimination */
817c3034d77SJose E. Roman     bjtmp = bj + bi[i];
818c3034d77SJose E. Roman     nzL   = bi[i + 1] - bi[i];
819c3034d77SJose E. Roman     for (k = 0; k < nzL; k++) {
820c3034d77SJose E. Roman       row = bjtmp[k];
821c3034d77SJose E. Roman       pc  = rtmp + bs2 * row;
822c3034d77SJose E. Roman       for (flg = 0, j = 0; j < bs2; j++) {
823c3034d77SJose E. Roman         if (pc[j] != 0.0) {
824c3034d77SJose E. Roman           flg = 1;
825c3034d77SJose E. Roman           break;
826c3034d77SJose E. Roman         }
827c3034d77SJose E. Roman       }
828c3034d77SJose E. Roman       if (flg) {
829c3034d77SJose E. Roman         pv = b->a + bs2 * bdiag[row];
830c3034d77SJose E. Roman         /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
831c3034d77SJose E. Roman         PetscCall(PetscKernel_A_gets_A_times_B_9(pc, pv, mwork));
832c3034d77SJose E. Roman 
833c3034d77SJose E. Roman         pj = b->j + bdiag[row + 1] + 1; /* beginning of U(row,:) */
834c3034d77SJose E. Roman         pv = b->a + bs2 * (bdiag[row + 1] + 1);
835c3034d77SJose E. Roman         nz = bdiag[row] - bdiag[row + 1] - 1; /* num of entries inU(row,:), excluding diag */
836c3034d77SJose E. Roman         for (j = 0; j < nz; j++) {
837c3034d77SJose E. Roman           /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
838c3034d77SJose E. Roman           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
839c3034d77SJose E. Roman           v = rtmp + bs2 * pj[j];
840c3034d77SJose E. Roman           PetscCall(PetscKernel_A_gets_A_minus_B_times_C_9(v, pc, pv + 81 * j));
841c3034d77SJose E. Roman           /* pv incremented in PetscKernel_A_gets_A_minus_B_times_C_9 */
842c3034d77SJose E. Roman         }
843c3034d77SJose E. Roman         PetscCall(PetscLogFlops(1458 * nz + 1377)); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
844c3034d77SJose E. Roman       }
845c3034d77SJose E. Roman     }
846c3034d77SJose E. Roman 
847c3034d77SJose E. Roman     /* finished row so stick it into b->a */
848c3034d77SJose E. Roman     /* L part */
849c3034d77SJose E. Roman     pv = b->a + bs2 * bi[i];
850c3034d77SJose E. Roman     pj = b->j + bi[i];
851c3034d77SJose E. Roman     nz = bi[i + 1] - bi[i];
852c3034d77SJose E. Roman     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2));
853c3034d77SJose E. Roman 
854c3034d77SJose E. Roman     /* Mark diagonal and invert diagonal for simpler triangular solves */
855c3034d77SJose E. Roman     pv = b->a + bs2 * bdiag[i];
856c3034d77SJose E. Roman     pj = b->j + bdiag[i];
857c3034d77SJose E. Roman     PetscCall(PetscArraycpy(pv, rtmp + bs2 * pj[0], bs2));
858c3034d77SJose E. Roman     PetscCall(PetscKernel_A_gets_inverse_A_9(pv, shift, allowzeropivot, &zeropivotdetected));
859c3034d77SJose E. Roman     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
860c3034d77SJose E. Roman 
861c3034d77SJose E. Roman     /* U part */
862c3034d77SJose E. Roman     pv = b->a + bs2 * (bdiag[i + 1] + 1);
863c3034d77SJose E. Roman     pj = b->j + bdiag[i + 1] + 1;
864c3034d77SJose E. Roman     nz = bdiag[i] - bdiag[i + 1] - 1;
865c3034d77SJose E. Roman     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2));
866c3034d77SJose E. Roman   }
867c3034d77SJose E. Roman   PetscCall(PetscFree2(rtmp, mwork));
868c3034d77SJose E. Roman 
869c3034d77SJose E. Roman   C->ops->solve          = MatSolve_SeqBAIJ_9_NaturalOrdering;
870c3034d77SJose E. Roman   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_N;
871c3034d77SJose E. Roman   C->assembled           = PETSC_TRUE;
872c3034d77SJose E. Roman 
873c3034d77SJose E. Roman   PetscCall(PetscLogFlops(1.333333333333 * 9 * 9 * 9 * n)); /* from inverting diagonal blocks */
874c3034d77SJose E. Roman   PetscFunctionReturn(PETSC_SUCCESS);
875c3034d77SJose E. Roman }
876c3034d77SJose E. Roman 
MatSolve_SeqBAIJ_9_NaturalOrdering(Mat A,Vec bb,Vec xx)877c3034d77SJose E. Roman PetscErrorCode MatSolve_SeqBAIJ_9_NaturalOrdering(Mat A, Vec bb, Vec xx)
878c3034d77SJose E. Roman {
879c3034d77SJose E. Roman   Mat_SeqBAIJ       *a  = (Mat_SeqBAIJ *)A->data;
880c3034d77SJose E. Roman   const PetscInt    *ai = a->i, *aj = a->j, *adiag = a->diag, *vi;
881c3034d77SJose E. Roman   PetscInt           i, k, n                       = a->mbs;
882c3034d77SJose E. Roman   PetscInt           nz, bs = A->rmap->bs, bs2 = a->bs2;
883c3034d77SJose E. Roman   const MatScalar   *aa = a->a, *v;
884c3034d77SJose E. Roman   PetscScalar       *x, *s, *t, *ls;
885c3034d77SJose E. Roman   const PetscScalar *b;
886c3034d77SJose E. Roman   __m256d            a0, a1, a2, a3, a4, a5, w0, w1, w2, w3, s0, s1, s2, v0, v1, v2, v3;
887c3034d77SJose E. Roman 
888c3034d77SJose E. Roman   PetscFunctionBegin;
889c3034d77SJose E. Roman   PetscCall(VecGetArrayRead(bb, &b));
890c3034d77SJose E. Roman   PetscCall(VecGetArray(xx, &x));
891c3034d77SJose E. Roman   t = a->solve_work;
892c3034d77SJose E. Roman 
893c3034d77SJose E. Roman   /* forward solve the lower triangular */
894c3034d77SJose E. Roman   PetscCall(PetscArraycpy(t, b, bs)); /* copy 1st block of b to t */
895c3034d77SJose E. Roman 
896c3034d77SJose E. Roman   for (i = 1; i < n; i++) {
897c3034d77SJose E. Roman     v  = aa + bs2 * ai[i];
898c3034d77SJose E. Roman     vi = aj + ai[i];
899c3034d77SJose E. Roman     nz = ai[i + 1] - ai[i];
900c3034d77SJose E. Roman     s  = t + bs * i;
901c3034d77SJose E. Roman     PetscCall(PetscArraycpy(s, b + bs * i, bs)); /* copy i_th block of b to t */
902c3034d77SJose E. Roman 
903c3034d77SJose E. Roman     __m256d s0, s1, s2;
904c3034d77SJose E. Roman     s0 = _mm256_loadu_pd(s + 0);
905c3034d77SJose E. Roman     s1 = _mm256_loadu_pd(s + 4);
906c3034d77SJose E. Roman     s2 = _mm256_maskload_pd(s + 8, _mm256_set_epi64x(0LL, 0LL, 0LL, 1LL << 63));
907c3034d77SJose E. Roman 
908c3034d77SJose E. Roman     for (k = 0; k < nz; k++) {
909c3034d77SJose E. Roman       w0 = _mm256_set1_pd((t + bs * vi[k])[0]);
910c3034d77SJose E. Roman       a0 = _mm256_loadu_pd(&v[0]);
911c3034d77SJose E. Roman       s0 = _mm256_fnmadd_pd(a0, w0, s0);
912c3034d77SJose E. Roman       a1 = _mm256_loadu_pd(&v[4]);
913c3034d77SJose E. Roman       s1 = _mm256_fnmadd_pd(a1, w0, s1);
914c3034d77SJose E. Roman       a2 = _mm256_loadu_pd(&v[8]);
915c3034d77SJose E. Roman       s2 = _mm256_fnmadd_pd(a2, w0, s2);
916c3034d77SJose E. Roman 
917c3034d77SJose E. Roman       w1 = _mm256_set1_pd((t + bs * vi[k])[1]);
918c3034d77SJose E. Roman       a3 = _mm256_loadu_pd(&v[9]);
919c3034d77SJose E. Roman       s0 = _mm256_fnmadd_pd(a3, w1, s0);
920c3034d77SJose E. Roman       a4 = _mm256_loadu_pd(&v[13]);
921c3034d77SJose E. Roman       s1 = _mm256_fnmadd_pd(a4, w1, s1);
922c3034d77SJose E. Roman       a5 = _mm256_loadu_pd(&v[17]);
923c3034d77SJose E. Roman       s2 = _mm256_fnmadd_pd(a5, w1, s2);
924c3034d77SJose E. Roman 
925c3034d77SJose E. Roman       w2 = _mm256_set1_pd((t + bs * vi[k])[2]);
926c3034d77SJose E. Roman       a0 = _mm256_loadu_pd(&v[18]);
927c3034d77SJose E. Roman       s0 = _mm256_fnmadd_pd(a0, w2, s0);
928c3034d77SJose E. Roman       a1 = _mm256_loadu_pd(&v[22]);
929c3034d77SJose E. Roman       s1 = _mm256_fnmadd_pd(a1, w2, s1);
930c3034d77SJose E. Roman       a2 = _mm256_loadu_pd(&v[26]);
931c3034d77SJose E. Roman       s2 = _mm256_fnmadd_pd(a2, w2, s2);
932c3034d77SJose E. Roman 
933c3034d77SJose E. Roman       w3 = _mm256_set1_pd((t + bs * vi[k])[3]);
934c3034d77SJose E. Roman       a3 = _mm256_loadu_pd(&v[27]);
935c3034d77SJose E. Roman       s0 = _mm256_fnmadd_pd(a3, w3, s0);
936c3034d77SJose E. Roman       a4 = _mm256_loadu_pd(&v[31]);
937c3034d77SJose E. Roman       s1 = _mm256_fnmadd_pd(a4, w3, s1);
938c3034d77SJose E. Roman       a5 = _mm256_loadu_pd(&v[35]);
939c3034d77SJose E. Roman       s2 = _mm256_fnmadd_pd(a5, w3, s2);
940c3034d77SJose E. Roman 
941c3034d77SJose E. Roman       w0 = _mm256_set1_pd((t + bs * vi[k])[4]);
942c3034d77SJose E. Roman       a0 = _mm256_loadu_pd(&v[36]);
943c3034d77SJose E. Roman       s0 = _mm256_fnmadd_pd(a0, w0, s0);
944c3034d77SJose E. Roman       a1 = _mm256_loadu_pd(&v[40]);
945c3034d77SJose E. Roman       s1 = _mm256_fnmadd_pd(a1, w0, s1);
946c3034d77SJose E. Roman       a2 = _mm256_loadu_pd(&v[44]);
947c3034d77SJose E. Roman       s2 = _mm256_fnmadd_pd(a2, w0, s2);
948c3034d77SJose E. Roman 
949c3034d77SJose E. Roman       w1 = _mm256_set1_pd((t + bs * vi[k])[5]);
950c3034d77SJose E. Roman       a3 = _mm256_loadu_pd(&v[45]);
951c3034d77SJose E. Roman       s0 = _mm256_fnmadd_pd(a3, w1, s0);
952c3034d77SJose E. Roman       a4 = _mm256_loadu_pd(&v[49]);
953c3034d77SJose E. Roman       s1 = _mm256_fnmadd_pd(a4, w1, s1);
954c3034d77SJose E. Roman       a5 = _mm256_loadu_pd(&v[53]);
955c3034d77SJose E. Roman       s2 = _mm256_fnmadd_pd(a5, w1, s2);
956c3034d77SJose E. Roman 
957c3034d77SJose E. Roman       w2 = _mm256_set1_pd((t + bs * vi[k])[6]);
958c3034d77SJose E. Roman       a0 = _mm256_loadu_pd(&v[54]);
959c3034d77SJose E. Roman       s0 = _mm256_fnmadd_pd(a0, w2, s0);
960c3034d77SJose E. Roman       a1 = _mm256_loadu_pd(&v[58]);
961c3034d77SJose E. Roman       s1 = _mm256_fnmadd_pd(a1, w2, s1);
962c3034d77SJose E. Roman       a2 = _mm256_loadu_pd(&v[62]);
963c3034d77SJose E. Roman       s2 = _mm256_fnmadd_pd(a2, w2, s2);
964c3034d77SJose E. Roman 
965c3034d77SJose E. Roman       w3 = _mm256_set1_pd((t + bs * vi[k])[7]);
966c3034d77SJose E. Roman       a3 = _mm256_loadu_pd(&v[63]);
967c3034d77SJose E. Roman       s0 = _mm256_fnmadd_pd(a3, w3, s0);
968c3034d77SJose E. Roman       a4 = _mm256_loadu_pd(&v[67]);
969c3034d77SJose E. Roman       s1 = _mm256_fnmadd_pd(a4, w3, s1);
970c3034d77SJose E. Roman       a5 = _mm256_loadu_pd(&v[71]);
971c3034d77SJose E. Roman       s2 = _mm256_fnmadd_pd(a5, w3, s2);
972c3034d77SJose E. Roman 
973c3034d77SJose E. Roman       w0 = _mm256_set1_pd((t + bs * vi[k])[8]);
974c3034d77SJose E. Roman       a0 = _mm256_loadu_pd(&v[72]);
975c3034d77SJose E. Roman       s0 = _mm256_fnmadd_pd(a0, w0, s0);
976c3034d77SJose E. Roman       a1 = _mm256_loadu_pd(&v[76]);
977c3034d77SJose E. Roman       s1 = _mm256_fnmadd_pd(a1, w0, s1);
978c3034d77SJose E. Roman       a2 = _mm256_maskload_pd(v + 80, _mm256_set_epi64x(0LL, 0LL, 0LL, 1LL << 63));
979c3034d77SJose E. Roman       s2 = _mm256_fnmadd_pd(a2, w0, s2);
980c3034d77SJose E. Roman       v += bs2;
981c3034d77SJose E. Roman     }
982c3034d77SJose E. Roman     _mm256_storeu_pd(&s[0], s0);
983c3034d77SJose E. Roman     _mm256_storeu_pd(&s[4], s1);
984c3034d77SJose E. Roman     _mm256_maskstore_pd(&s[8], _mm256_set_epi64x(0LL, 0LL, 0LL, 1LL << 63), s2);
985c3034d77SJose E. Roman   }
986c3034d77SJose E. Roman 
987c3034d77SJose E. Roman   /* backward solve the upper triangular */
988c3034d77SJose E. Roman   ls = a->solve_work + A->cmap->n;
989c3034d77SJose E. Roman   for (i = n - 1; i >= 0; i--) {
990c3034d77SJose E. Roman     v  = aa + bs2 * (adiag[i + 1] + 1);
991c3034d77SJose E. Roman     vi = aj + adiag[i + 1] + 1;
992c3034d77SJose E. Roman     nz = adiag[i] - adiag[i + 1] - 1;
993c3034d77SJose E. Roman     PetscCall(PetscArraycpy(ls, t + i * bs, bs));
994c3034d77SJose E. Roman 
995c3034d77SJose E. Roman     s0 = _mm256_loadu_pd(ls + 0);
996c3034d77SJose E. Roman     s1 = _mm256_loadu_pd(ls + 4);
997c3034d77SJose E. Roman     s2 = _mm256_maskload_pd(ls + 8, _mm256_set_epi64x(0LL, 0LL, 0LL, 1LL << 63));
998c3034d77SJose E. Roman 
999c3034d77SJose E. Roman     for (k = 0; k < nz; k++) {
1000c3034d77SJose E. Roman       w0 = _mm256_set1_pd((t + bs * vi[k])[0]);
1001c3034d77SJose E. Roman       a0 = _mm256_loadu_pd(&v[0]);
1002c3034d77SJose E. Roman       s0 = _mm256_fnmadd_pd(a0, w0, s0);
1003c3034d77SJose E. Roman       a1 = _mm256_loadu_pd(&v[4]);
1004c3034d77SJose E. Roman       s1 = _mm256_fnmadd_pd(a1, w0, s1);
1005c3034d77SJose E. Roman       a2 = _mm256_loadu_pd(&v[8]);
1006c3034d77SJose E. Roman       s2 = _mm256_fnmadd_pd(a2, w0, s2);
1007c3034d77SJose E. Roman 
1008c3034d77SJose E. Roman       /* v += 9; */
1009c3034d77SJose E. Roman       w1 = _mm256_set1_pd((t + bs * vi[k])[1]);
1010c3034d77SJose E. Roman       a3 = _mm256_loadu_pd(&v[9]);
1011c3034d77SJose E. Roman       s0 = _mm256_fnmadd_pd(a3, w1, s0);
1012c3034d77SJose E. Roman       a4 = _mm256_loadu_pd(&v[13]);
1013c3034d77SJose E. Roman       s1 = _mm256_fnmadd_pd(a4, w1, s1);
1014c3034d77SJose E. Roman       a5 = _mm256_loadu_pd(&v[17]);
1015c3034d77SJose E. Roman       s2 = _mm256_fnmadd_pd(a5, w1, s2);
1016c3034d77SJose E. Roman 
1017c3034d77SJose E. Roman       /* v += 9; */
1018c3034d77SJose E. Roman       w2 = _mm256_set1_pd((t + bs * vi[k])[2]);
1019c3034d77SJose E. Roman       a0 = _mm256_loadu_pd(&v[18]);
1020c3034d77SJose E. Roman       s0 = _mm256_fnmadd_pd(a0, w2, s0);
1021c3034d77SJose E. Roman       a1 = _mm256_loadu_pd(&v[22]);
1022c3034d77SJose E. Roman       s1 = _mm256_fnmadd_pd(a1, w2, s1);
1023c3034d77SJose E. Roman       a2 = _mm256_loadu_pd(&v[26]);
1024c3034d77SJose E. Roman       s2 = _mm256_fnmadd_pd(a2, w2, s2);
1025c3034d77SJose E. Roman 
1026c3034d77SJose E. Roman       /* v += 9; */
1027c3034d77SJose E. Roman       w3 = _mm256_set1_pd((t + bs * vi[k])[3]);
1028c3034d77SJose E. Roman       a3 = _mm256_loadu_pd(&v[27]);
1029c3034d77SJose E. Roman       s0 = _mm256_fnmadd_pd(a3, w3, s0);
1030c3034d77SJose E. Roman       a4 = _mm256_loadu_pd(&v[31]);
1031c3034d77SJose E. Roman       s1 = _mm256_fnmadd_pd(a4, w3, s1);
1032c3034d77SJose E. Roman       a5 = _mm256_loadu_pd(&v[35]);
1033c3034d77SJose E. Roman       s2 = _mm256_fnmadd_pd(a5, w3, s2);
1034c3034d77SJose E. Roman 
1035c3034d77SJose E. Roman       /* v += 9; */
1036c3034d77SJose E. Roman       w0 = _mm256_set1_pd((t + bs * vi[k])[4]);
1037c3034d77SJose E. Roman       a0 = _mm256_loadu_pd(&v[36]);
1038c3034d77SJose E. Roman       s0 = _mm256_fnmadd_pd(a0, w0, s0);
1039c3034d77SJose E. Roman       a1 = _mm256_loadu_pd(&v[40]);
1040c3034d77SJose E. Roman       s1 = _mm256_fnmadd_pd(a1, w0, s1);
1041c3034d77SJose E. Roman       a2 = _mm256_loadu_pd(&v[44]);
1042c3034d77SJose E. Roman       s2 = _mm256_fnmadd_pd(a2, w0, s2);
1043c3034d77SJose E. Roman 
1044c3034d77SJose E. Roman       /* v += 9; */
1045c3034d77SJose E. Roman       w1 = _mm256_set1_pd((t + bs * vi[k])[5]);
1046c3034d77SJose E. Roman       a3 = _mm256_loadu_pd(&v[45]);
1047c3034d77SJose E. Roman       s0 = _mm256_fnmadd_pd(a3, w1, s0);
1048c3034d77SJose E. Roman       a4 = _mm256_loadu_pd(&v[49]);
1049c3034d77SJose E. Roman       s1 = _mm256_fnmadd_pd(a4, w1, s1);
1050c3034d77SJose E. Roman       a5 = _mm256_loadu_pd(&v[53]);
1051c3034d77SJose E. Roman       s2 = _mm256_fnmadd_pd(a5, w1, s2);
1052c3034d77SJose E. Roman 
1053c3034d77SJose E. Roman       /* v += 9; */
1054c3034d77SJose E. Roman       w2 = _mm256_set1_pd((t + bs * vi[k])[6]);
1055c3034d77SJose E. Roman       a0 = _mm256_loadu_pd(&v[54]);
1056c3034d77SJose E. Roman       s0 = _mm256_fnmadd_pd(a0, w2, s0);
1057c3034d77SJose E. Roman       a1 = _mm256_loadu_pd(&v[58]);
1058c3034d77SJose E. Roman       s1 = _mm256_fnmadd_pd(a1, w2, s1);
1059c3034d77SJose E. Roman       a2 = _mm256_loadu_pd(&v[62]);
1060c3034d77SJose E. Roman       s2 = _mm256_fnmadd_pd(a2, w2, s2);
1061c3034d77SJose E. Roman 
1062c3034d77SJose E. Roman       /* v += 9; */
1063c3034d77SJose E. Roman       w3 = _mm256_set1_pd((t + bs * vi[k])[7]);
1064c3034d77SJose E. Roman       a3 = _mm256_loadu_pd(&v[63]);
1065c3034d77SJose E. Roman       s0 = _mm256_fnmadd_pd(a3, w3, s0);
1066c3034d77SJose E. Roman       a4 = _mm256_loadu_pd(&v[67]);
1067c3034d77SJose E. Roman       s1 = _mm256_fnmadd_pd(a4, w3, s1);
1068c3034d77SJose E. Roman       a5 = _mm256_loadu_pd(&v[71]);
1069c3034d77SJose E. Roman       s2 = _mm256_fnmadd_pd(a5, w3, s2);
1070c3034d77SJose E. Roman 
1071c3034d77SJose E. Roman       /* v += 9; */
1072c3034d77SJose E. Roman       w0 = _mm256_set1_pd((t + bs * vi[k])[8]);
1073c3034d77SJose E. Roman       a0 = _mm256_loadu_pd(&v[72]);
1074c3034d77SJose E. Roman       s0 = _mm256_fnmadd_pd(a0, w0, s0);
1075c3034d77SJose E. Roman       a1 = _mm256_loadu_pd(&v[76]);
1076c3034d77SJose E. Roman       s1 = _mm256_fnmadd_pd(a1, w0, s1);
1077c3034d77SJose E. Roman       a2 = _mm256_maskload_pd(v + 80, _mm256_set_epi64x(0LL, 0LL, 0LL, 1LL << 63));
1078c3034d77SJose E. Roman       s2 = _mm256_fnmadd_pd(a2, w0, s2);
1079c3034d77SJose E. Roman       v += bs2;
1080c3034d77SJose E. Roman     }
1081c3034d77SJose E. Roman 
1082c3034d77SJose E. Roman     _mm256_storeu_pd(&ls[0], s0);
1083c3034d77SJose E. Roman     _mm256_storeu_pd(&ls[4], s1);
1084c3034d77SJose E. Roman     _mm256_maskstore_pd(&ls[8], _mm256_set_epi64x(0LL, 0LL, 0LL, 1LL << 63), s2);
1085c3034d77SJose E. Roman 
1086c3034d77SJose E. Roman     w0 = _mm256_setzero_pd();
1087c3034d77SJose E. Roman     w1 = _mm256_setzero_pd();
1088c3034d77SJose E. Roman     w2 = _mm256_setzero_pd();
1089c3034d77SJose E. Roman 
1090c3034d77SJose E. Roman     /* first row */
1091c3034d77SJose E. Roman     v0 = _mm256_set1_pd(ls[0]);
1092c3034d77SJose E. Roman     a0 = _mm256_loadu_pd(&(aa + bs2 * adiag[i])[0]);
1093c3034d77SJose E. Roman     w0 = _mm256_fmadd_pd(a0, v0, w0);
1094c3034d77SJose E. Roman     a1 = _mm256_loadu_pd(&(aa + bs2 * adiag[i])[4]);
1095c3034d77SJose E. Roman     w1 = _mm256_fmadd_pd(a1, v0, w1);
1096c3034d77SJose E. Roman     a2 = _mm256_loadu_pd(&(aa + bs2 * adiag[i])[8]);
1097c3034d77SJose E. Roman     w2 = _mm256_fmadd_pd(a2, v0, w2);
1098c3034d77SJose E. Roman 
1099c3034d77SJose E. Roman     /* second row */
1100c3034d77SJose E. Roman     v1 = _mm256_set1_pd(ls[1]);
1101c3034d77SJose E. Roman     a3 = _mm256_loadu_pd(&(aa + bs2 * adiag[i])[9]);
1102c3034d77SJose E. Roman     w0 = _mm256_fmadd_pd(a3, v1, w0);
1103c3034d77SJose E. Roman     a4 = _mm256_loadu_pd(&(aa + bs2 * adiag[i])[13]);
1104c3034d77SJose E. Roman     w1 = _mm256_fmadd_pd(a4, v1, w1);
1105c3034d77SJose E. Roman     a5 = _mm256_loadu_pd(&(aa + bs2 * adiag[i])[17]);
1106c3034d77SJose E. Roman     w2 = _mm256_fmadd_pd(a5, v1, w2);
1107c3034d77SJose E. Roman 
1108c3034d77SJose E. Roman     /* third row */
1109c3034d77SJose E. Roman     v2 = _mm256_set1_pd(ls[2]);
1110c3034d77SJose E. Roman     a0 = _mm256_loadu_pd(&(aa + bs2 * adiag[i])[18]);
1111c3034d77SJose E. Roman     w0 = _mm256_fmadd_pd(a0, v2, w0);
1112c3034d77SJose E. Roman     a1 = _mm256_loadu_pd(&(aa + bs2 * adiag[i])[22]);
1113c3034d77SJose E. Roman     w1 = _mm256_fmadd_pd(a1, v2, w1);
1114c3034d77SJose E. Roman     a2 = _mm256_loadu_pd(&(aa + bs2 * adiag[i])[26]);
1115c3034d77SJose E. Roman     w2 = _mm256_fmadd_pd(a2, v2, w2);
1116c3034d77SJose E. Roman 
1117c3034d77SJose E. Roman     /* fourth row */
1118c3034d77SJose E. Roman     v3 = _mm256_set1_pd(ls[3]);
1119c3034d77SJose E. Roman     a3 = _mm256_loadu_pd(&(aa + bs2 * adiag[i])[27]);
1120c3034d77SJose E. Roman     w0 = _mm256_fmadd_pd(a3, v3, w0);
1121c3034d77SJose E. Roman     a4 = _mm256_loadu_pd(&(aa + bs2 * adiag[i])[31]);
1122c3034d77SJose E. Roman     w1 = _mm256_fmadd_pd(a4, v3, w1);
1123c3034d77SJose E. Roman     a5 = _mm256_loadu_pd(&(aa + bs2 * adiag[i])[35]);
1124c3034d77SJose E. Roman     w2 = _mm256_fmadd_pd(a5, v3, w2);
1125c3034d77SJose E. Roman 
1126c3034d77SJose E. Roman     /* fifth row */
1127c3034d77SJose E. Roman     v0 = _mm256_set1_pd(ls[4]);
1128c3034d77SJose E. Roman     a0 = _mm256_loadu_pd(&(aa + bs2 * adiag[i])[36]);
1129c3034d77SJose E. Roman     w0 = _mm256_fmadd_pd(a0, v0, w0);
1130c3034d77SJose E. Roman     a1 = _mm256_loadu_pd(&(aa + bs2 * adiag[i])[40]);
1131c3034d77SJose E. Roman     w1 = _mm256_fmadd_pd(a1, v0, w1);
1132c3034d77SJose E. Roman     a2 = _mm256_loadu_pd(&(aa + bs2 * adiag[i])[44]);
1133c3034d77SJose E. Roman     w2 = _mm256_fmadd_pd(a2, v0, w2);
1134c3034d77SJose E. Roman 
1135c3034d77SJose E. Roman     /* sixth row */
1136c3034d77SJose E. Roman     v1 = _mm256_set1_pd(ls[5]);
1137c3034d77SJose E. Roman     a3 = _mm256_loadu_pd(&(aa + bs2 * adiag[i])[45]);
1138c3034d77SJose E. Roman     w0 = _mm256_fmadd_pd(a3, v1, w0);
1139c3034d77SJose E. Roman     a4 = _mm256_loadu_pd(&(aa + bs2 * adiag[i])[49]);
1140c3034d77SJose E. Roman     w1 = _mm256_fmadd_pd(a4, v1, w1);
1141c3034d77SJose E. Roman     a5 = _mm256_loadu_pd(&(aa + bs2 * adiag[i])[53]);
1142c3034d77SJose E. Roman     w2 = _mm256_fmadd_pd(a5, v1, w2);
1143c3034d77SJose E. Roman 
1144c3034d77SJose E. Roman     /* seventh row */
1145c3034d77SJose E. Roman     v2 = _mm256_set1_pd(ls[6]);
1146c3034d77SJose E. Roman     a0 = _mm256_loadu_pd(&(aa + bs2 * adiag[i])[54]);
1147c3034d77SJose E. Roman     w0 = _mm256_fmadd_pd(a0, v2, w0);
1148c3034d77SJose E. Roman     a1 = _mm256_loadu_pd(&(aa + bs2 * adiag[i])[58]);
1149c3034d77SJose E. Roman     w1 = _mm256_fmadd_pd(a1, v2, w1);
1150c3034d77SJose E. Roman     a2 = _mm256_loadu_pd(&(aa + bs2 * adiag[i])[62]);
1151c3034d77SJose E. Roman     w2 = _mm256_fmadd_pd(a2, v2, w2);
1152c3034d77SJose E. Roman 
1153c3034d77SJose E. Roman     /* eighth row */
1154c3034d77SJose E. Roman     v3 = _mm256_set1_pd(ls[7]);
1155c3034d77SJose E. Roman     a3 = _mm256_loadu_pd(&(aa + bs2 * adiag[i])[63]);
1156c3034d77SJose E. Roman     w0 = _mm256_fmadd_pd(a3, v3, w0);
1157c3034d77SJose E. Roman     a4 = _mm256_loadu_pd(&(aa + bs2 * adiag[i])[67]);
1158c3034d77SJose E. Roman     w1 = _mm256_fmadd_pd(a4, v3, w1);
1159c3034d77SJose E. Roman     a5 = _mm256_loadu_pd(&(aa + bs2 * adiag[i])[71]);
1160c3034d77SJose E. Roman     w2 = _mm256_fmadd_pd(a5, v3, w2);
1161c3034d77SJose E. Roman 
1162c3034d77SJose E. Roman     /* ninth row */
1163c3034d77SJose E. Roman     v0 = _mm256_set1_pd(ls[8]);
1164c3034d77SJose E. Roman     a3 = _mm256_loadu_pd(&(aa + bs2 * adiag[i])[72]);
1165c3034d77SJose E. Roman     w0 = _mm256_fmadd_pd(a3, v0, w0);
1166c3034d77SJose E. Roman     a4 = _mm256_loadu_pd(&(aa + bs2 * adiag[i])[76]);
1167c3034d77SJose E. Roman     w1 = _mm256_fmadd_pd(a4, v0, w1);
1168c3034d77SJose E. Roman     a2 = _mm256_maskload_pd(&(aa + bs2 * adiag[i])[80], _mm256_set_epi64x(0LL, 0LL, 0LL, 1LL << 63));
1169c3034d77SJose E. Roman     w2 = _mm256_fmadd_pd(a2, v0, w2);
1170c3034d77SJose E. Roman 
1171c3034d77SJose E. Roman     _mm256_storeu_pd(&(t + i * bs)[0], w0);
1172c3034d77SJose E. Roman     _mm256_storeu_pd(&(t + i * bs)[4], w1);
1173c3034d77SJose E. Roman     _mm256_maskstore_pd(&(t + i * bs)[8], _mm256_set_epi64x(0LL, 0LL, 0LL, 1LL << 63), w2);
1174c3034d77SJose E. Roman 
1175c3034d77SJose E. Roman     PetscCall(PetscArraycpy(x + i * bs, t + i * bs, bs));
1176c3034d77SJose E. Roman   }
1177c3034d77SJose E. Roman 
1178c3034d77SJose E. Roman   PetscCall(VecRestoreArrayRead(bb, &b));
1179c3034d77SJose E. Roman   PetscCall(VecRestoreArray(xx, &x));
1180c3034d77SJose E. Roman   PetscCall(PetscLogFlops(2.0 * (a->bs2) * (a->nz) - A->rmap->bs * A->cmap->n));
1181c3034d77SJose E. Roman   PetscFunctionReturn(PETSC_SUCCESS);
1182c3034d77SJose E. Roman }
1183c3034d77SJose E. Roman #endif
1184