xref: /petsc/src/mat/impls/baij/seq/baijfact11.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 4 by 4
983287d42SBarry Smith */
MatILUFactorNumeric_SeqBAIJ_4_inplace(Mat C,Mat A,const MatFactorInfo * info)10*421480d9SBarry Smith PetscErrorCode MatILUFactorNumeric_SeqBAIJ_4_inplace(Mat C, Mat A, const MatFactorInfo *info)
11d71ae5a4SJacob Faibussowitsch {
1283287d42SBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data;
1383287d42SBarry Smith   IS              isrow = b->row, isicol = b->icol;
145d0c19d7SBarry Smith   const PetscInt *r, *ic;
155d0c19d7SBarry Smith   PetscInt        i, j, n = a->mbs, *bi = b->i, *bj = b->j;
16690b6cddSBarry Smith   PetscInt       *ajtmpold, *ajtmp, nz, row;
17*421480d9SBarry Smith   const PetscInt *diag_offset;
18*421480d9SBarry Smith   PetscInt        idx, *ai = a->i, *aj = a->j, *pj;
1983287d42SBarry Smith   MatScalar      *pv, *v, *rtmp, *pc, *w, *x;
2083287d42SBarry Smith   MatScalar       p1, p2, p3, p4, m1, m2, m3, m4, m5, m6, m7, m8, m9, x1, x2, x3, x4;
2183287d42SBarry Smith   MatScalar       p5, p6, p7, p8, p9, x5, x6, x7, x8, x9, x10, x11, x12, x13, x14, x15, x16;
2283287d42SBarry Smith   MatScalar       p10, p11, p12, p13, p14, p15, p16, m10, m11, m12;
2383287d42SBarry Smith   MatScalar       m13, m14, m15, m16;
2483287d42SBarry Smith   MatScalar      *ba = b->a, *aa = a->a;
25a455e926SHong Zhang   PetscBool       pivotinblocks = b->pivotinblocks;
26182b8fbaSHong Zhang   PetscReal       shift         = info->shiftamount;
270164db54SHong Zhang   PetscBool       allowzeropivot, zeropivotdetected = PETSC_FALSE;
2883287d42SBarry Smith 
2983287d42SBarry Smith   PetscFunctionBegin;
30*421480d9SBarry Smith   /* Since A is C and C is labeled as a factored matrix we need to lie to MatGetDiagonalMarkers_SeqBAIJ() to get it to compute the diagonals */
31*421480d9SBarry Smith   A->factortype = MAT_FACTOR_NONE;
32*421480d9SBarry Smith   PetscCall(MatGetDiagonalMarkers_SeqBAIJ(A, &diag_offset, NULL));
33*421480d9SBarry Smith   A->factortype = MAT_FACTOR_ILU;
349566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isrow, &r));
359566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isicol, &ic));
369566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(16 * (n + 1), &rtmp));
370164db54SHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
3883287d42SBarry Smith 
3983287d42SBarry Smith   for (i = 0; i < n; i++) {
4083287d42SBarry Smith     nz    = bi[i + 1] - bi[i];
4183287d42SBarry Smith     ajtmp = bj + bi[i];
4283287d42SBarry Smith     for (j = 0; j < nz; j++) {
4383287d42SBarry Smith       x    = rtmp + 16 * ajtmp[j];
4483287d42SBarry Smith       x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = x[9] = 0.0;
4583287d42SBarry Smith       x[10] = x[11] = x[12] = x[13] = x[14] = x[15] = 0.0;
4683287d42SBarry Smith     }
4783287d42SBarry Smith     /* load in initial (unfactored row) */
4883287d42SBarry Smith     idx      = r[i];
4983287d42SBarry Smith     nz       = ai[idx + 1] - ai[idx];
5083287d42SBarry Smith     ajtmpold = aj + ai[idx];
5183287d42SBarry Smith     v        = aa + 16 * ai[idx];
5283287d42SBarry Smith     for (j = 0; j < nz; j++) {
5383287d42SBarry Smith       x     = rtmp + 16 * ic[ajtmpold[j]];
549371c9d4SSatish Balay       x[0]  = v[0];
559371c9d4SSatish Balay       x[1]  = v[1];
569371c9d4SSatish Balay       x[2]  = v[2];
579371c9d4SSatish Balay       x[3]  = v[3];
589371c9d4SSatish Balay       x[4]  = v[4];
599371c9d4SSatish Balay       x[5]  = v[5];
609371c9d4SSatish Balay       x[6]  = v[6];
619371c9d4SSatish Balay       x[7]  = v[7];
629371c9d4SSatish Balay       x[8]  = v[8];
639371c9d4SSatish Balay       x[9]  = v[9];
649371c9d4SSatish Balay       x[10] = v[10];
659371c9d4SSatish Balay       x[11] = v[11];
669371c9d4SSatish Balay       x[12] = v[12];
679371c9d4SSatish Balay       x[13] = v[13];
689371c9d4SSatish Balay       x[14] = v[14];
699371c9d4SSatish Balay       x[15] = v[15];
7083287d42SBarry Smith       v += 16;
7183287d42SBarry Smith     }
7283287d42SBarry Smith     row = *ajtmp++;
7383287d42SBarry Smith     while (row < i) {
7483287d42SBarry Smith       pc  = rtmp + 16 * row;
759371c9d4SSatish Balay       p1  = pc[0];
769371c9d4SSatish Balay       p2  = pc[1];
779371c9d4SSatish Balay       p3  = pc[2];
789371c9d4SSatish Balay       p4  = pc[3];
799371c9d4SSatish Balay       p5  = pc[4];
809371c9d4SSatish Balay       p6  = pc[5];
819371c9d4SSatish Balay       p7  = pc[6];
829371c9d4SSatish Balay       p8  = pc[7];
839371c9d4SSatish Balay       p9  = pc[8];
849371c9d4SSatish Balay       p10 = pc[9];
859371c9d4SSatish Balay       p11 = pc[10];
869371c9d4SSatish Balay       p12 = pc[11];
879371c9d4SSatish Balay       p13 = pc[12];
889371c9d4SSatish Balay       p14 = pc[13];
899371c9d4SSatish Balay       p15 = pc[14];
909371c9d4SSatish Balay       p16 = pc[15];
919371c9d4SSatish 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) {
9283287d42SBarry Smith         pv    = ba + 16 * diag_offset[row];
9383287d42SBarry Smith         pj    = bj + diag_offset[row] + 1;
949371c9d4SSatish Balay         x1    = pv[0];
959371c9d4SSatish Balay         x2    = pv[1];
969371c9d4SSatish Balay         x3    = pv[2];
979371c9d4SSatish Balay         x4    = pv[3];
989371c9d4SSatish Balay         x5    = pv[4];
999371c9d4SSatish Balay         x6    = pv[5];
1009371c9d4SSatish Balay         x7    = pv[6];
1019371c9d4SSatish Balay         x8    = pv[7];
1029371c9d4SSatish Balay         x9    = pv[8];
1039371c9d4SSatish Balay         x10   = pv[9];
1049371c9d4SSatish Balay         x11   = pv[10];
1059371c9d4SSatish Balay         x12   = pv[11];
1069371c9d4SSatish Balay         x13   = pv[12];
1079371c9d4SSatish Balay         x14   = pv[13];
1089371c9d4SSatish Balay         x15   = pv[14];
1099371c9d4SSatish Balay         x16   = pv[15];
11083287d42SBarry Smith         pc[0] = m1 = p1 * x1 + p5 * x2 + p9 * x3 + p13 * x4;
11183287d42SBarry Smith         pc[1] = m2 = p2 * x1 + p6 * x2 + p10 * x3 + p14 * x4;
11283287d42SBarry Smith         pc[2] = m3 = p3 * x1 + p7 * x2 + p11 * x3 + p15 * x4;
11383287d42SBarry Smith         pc[3] = m4 = p4 * x1 + p8 * x2 + p12 * x3 + p16 * x4;
11483287d42SBarry Smith 
11583287d42SBarry Smith         pc[4] = m5 = p1 * x5 + p5 * x6 + p9 * x7 + p13 * x8;
11683287d42SBarry Smith         pc[5] = m6 = p2 * x5 + p6 * x6 + p10 * x7 + p14 * x8;
11783287d42SBarry Smith         pc[6] = m7 = p3 * x5 + p7 * x6 + p11 * x7 + p15 * x8;
11883287d42SBarry Smith         pc[7] = m8 = p4 * x5 + p8 * x6 + p12 * x7 + p16 * x8;
11983287d42SBarry Smith 
12083287d42SBarry Smith         pc[8] = m9 = p1 * x9 + p5 * x10 + p9 * x11 + p13 * x12;
12183287d42SBarry Smith         pc[9] = m10 = p2 * x9 + p6 * x10 + p10 * x11 + p14 * x12;
12283287d42SBarry Smith         pc[10] = m11 = p3 * x9 + p7 * x10 + p11 * x11 + p15 * x12;
12383287d42SBarry Smith         pc[11] = m12 = p4 * x9 + p8 * x10 + p12 * x11 + p16 * x12;
12483287d42SBarry Smith 
12583287d42SBarry Smith         pc[12] = m13 = p1 * x13 + p5 * x14 + p9 * x15 + p13 * x16;
12683287d42SBarry Smith         pc[13] = m14 = p2 * x13 + p6 * x14 + p10 * x15 + p14 * x16;
12783287d42SBarry Smith         pc[14] = m15 = p3 * x13 + p7 * x14 + p11 * x15 + p15 * x16;
12883287d42SBarry Smith         pc[15] = m16 = p4 * x13 + p8 * x14 + p12 * x15 + p16 * x16;
12983287d42SBarry Smith 
13083287d42SBarry Smith         nz = bi[row + 1] - diag_offset[row] - 1;
13183287d42SBarry Smith         pv += 16;
13283287d42SBarry Smith         for (j = 0; j < nz; j++) {
1339371c9d4SSatish Balay           x1  = pv[0];
1349371c9d4SSatish Balay           x2  = pv[1];
1359371c9d4SSatish Balay           x3  = pv[2];
1369371c9d4SSatish Balay           x4  = pv[3];
1379371c9d4SSatish Balay           x5  = pv[4];
1389371c9d4SSatish Balay           x6  = pv[5];
1399371c9d4SSatish Balay           x7  = pv[6];
1409371c9d4SSatish Balay           x8  = pv[7];
1419371c9d4SSatish Balay           x9  = pv[8];
1429371c9d4SSatish Balay           x10 = pv[9];
1439371c9d4SSatish Balay           x11 = pv[10];
1449371c9d4SSatish Balay           x12 = pv[11];
1459371c9d4SSatish Balay           x13 = pv[12];
1469371c9d4SSatish Balay           x14 = pv[13];
1479371c9d4SSatish Balay           x15 = pv[14];
1489371c9d4SSatish Balay           x16 = pv[15];
14983287d42SBarry Smith           x   = rtmp + 16 * pj[j];
15083287d42SBarry Smith           x[0] -= m1 * x1 + m5 * x2 + m9 * x3 + m13 * x4;
15183287d42SBarry Smith           x[1] -= m2 * x1 + m6 * x2 + m10 * x3 + m14 * x4;
15283287d42SBarry Smith           x[2] -= m3 * x1 + m7 * x2 + m11 * x3 + m15 * x4;
15383287d42SBarry Smith           x[3] -= m4 * x1 + m8 * x2 + m12 * x3 + m16 * x4;
15483287d42SBarry Smith 
15583287d42SBarry Smith           x[4] -= m1 * x5 + m5 * x6 + m9 * x7 + m13 * x8;
15683287d42SBarry Smith           x[5] -= m2 * x5 + m6 * x6 + m10 * x7 + m14 * x8;
15783287d42SBarry Smith           x[6] -= m3 * x5 + m7 * x6 + m11 * x7 + m15 * x8;
15883287d42SBarry Smith           x[7] -= m4 * x5 + m8 * x6 + m12 * x7 + m16 * x8;
15983287d42SBarry Smith 
16083287d42SBarry Smith           x[8] -= m1 * x9 + m5 * x10 + m9 * x11 + m13 * x12;
16183287d42SBarry Smith           x[9] -= m2 * x9 + m6 * x10 + m10 * x11 + m14 * x12;
16283287d42SBarry Smith           x[10] -= m3 * x9 + m7 * x10 + m11 * x11 + m15 * x12;
16383287d42SBarry Smith           x[11] -= m4 * x9 + m8 * x10 + m12 * x11 + m16 * x12;
16483287d42SBarry Smith 
16583287d42SBarry Smith           x[12] -= m1 * x13 + m5 * x14 + m9 * x15 + m13 * x16;
16683287d42SBarry Smith           x[13] -= m2 * x13 + m6 * x14 + m10 * x15 + m14 * x16;
16783287d42SBarry Smith           x[14] -= m3 * x13 + m7 * x14 + m11 * x15 + m15 * x16;
16883287d42SBarry Smith           x[15] -= m4 * x13 + m8 * x14 + m12 * x15 + m16 * x16;
16983287d42SBarry Smith 
17083287d42SBarry Smith           pv += 16;
17183287d42SBarry Smith         }
1729566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(128.0 * nz + 112.0));
17383287d42SBarry Smith       }
17483287d42SBarry Smith       row = *ajtmp++;
17583287d42SBarry Smith     }
17683287d42SBarry Smith     /* finished row so stick it into b->a */
17783287d42SBarry Smith     pv = ba + 16 * bi[i];
17883287d42SBarry Smith     pj = bj + bi[i];
17983287d42SBarry Smith     nz = bi[i + 1] - bi[i];
18083287d42SBarry Smith     for (j = 0; j < nz; j++) {
18183287d42SBarry Smith       x      = rtmp + 16 * pj[j];
1829371c9d4SSatish Balay       pv[0]  = x[0];
1839371c9d4SSatish Balay       pv[1]  = x[1];
1849371c9d4SSatish Balay       pv[2]  = x[2];
1859371c9d4SSatish Balay       pv[3]  = x[3];
1869371c9d4SSatish Balay       pv[4]  = x[4];
1879371c9d4SSatish Balay       pv[5]  = x[5];
1889371c9d4SSatish Balay       pv[6]  = x[6];
1899371c9d4SSatish Balay       pv[7]  = x[7];
1909371c9d4SSatish Balay       pv[8]  = x[8];
1919371c9d4SSatish Balay       pv[9]  = x[9];
1929371c9d4SSatish Balay       pv[10] = x[10];
1939371c9d4SSatish Balay       pv[11] = x[11];
1949371c9d4SSatish Balay       pv[12] = x[12];
1959371c9d4SSatish Balay       pv[13] = x[13];
1969371c9d4SSatish Balay       pv[14] = x[14];
1979371c9d4SSatish Balay       pv[15] = x[15];
19883287d42SBarry Smith       pv += 16;
19983287d42SBarry Smith     }
20083287d42SBarry Smith     /* invert diagonal block */
20183287d42SBarry Smith     w = ba + 16 * diag_offset[i];
202bcd9e38bSBarry Smith     if (pivotinblocks) {
2039566063dSJacob Faibussowitsch       PetscCall(PetscKernel_A_gets_inverse_A_4(w, shift, allowzeropivot, &zeropivotdetected));
2047b6c816cSBarry Smith       if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
205bcd9e38bSBarry Smith     } else {
2069566063dSJacob Faibussowitsch       PetscCall(PetscKernel_A_gets_inverse_A_4_nopivot(w));
207bcd9e38bSBarry Smith     }
20883287d42SBarry Smith   }
20983287d42SBarry Smith 
2109566063dSJacob Faibussowitsch   PetscCall(PetscFree(rtmp));
2119566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isicol, &ic));
2129566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow, &r));
21326fbe8dcSKarl Rupp 
21406e38f1dSHong Zhang   C->ops->solve          = MatSolve_SeqBAIJ_4_inplace;
21506e38f1dSHong Zhang   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_4_inplace;
21683287d42SBarry Smith   C->assembled           = PETSC_TRUE;
21726fbe8dcSKarl Rupp 
2189566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(1.333333333333 * 4 * 4 * 4 * b->mbs)); /* from inverting diagonal blocks */
2193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
22083287d42SBarry Smith }
221e6580ceeSShri Abhyankar 
2224dd39f65SShri Abhyankar /* MatLUFactorNumeric_SeqBAIJ_4 -
2234dd39f65SShri Abhyankar      copied from MatLUFactorNumeric_SeqBAIJ_N_inplace() and manually re-implemented
22496b95a6bSBarry Smith        PetscKernel_A_gets_A_times_B()
22596b95a6bSBarry Smith        PetscKernel_A_gets_A_minus_B_times_C()
22696b95a6bSBarry Smith        PetscKernel_A_gets_inverse_A()
227209027a4SShri Abhyankar */
228c0c7eb62SShri Abhyankar 
MatLUFactorNumeric_SeqBAIJ_4(Mat B,Mat A,const MatFactorInfo * info)229d71ae5a4SJacob Faibussowitsch PetscErrorCode MatLUFactorNumeric_SeqBAIJ_4(Mat B, Mat A, const MatFactorInfo *info)
230d71ae5a4SJacob Faibussowitsch {
231209027a4SShri Abhyankar   Mat             C = B;
232209027a4SShri Abhyankar   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data;
233209027a4SShri Abhyankar   IS              isrow = b->row, isicol = b->icol;
2345a586d82SBarry Smith   const PetscInt *r, *ic;
235bbd65245SShri Abhyankar   PetscInt        i, j, k, nz, nzL, row;
236bbd65245SShri Abhyankar   const PetscInt  n = a->mbs, *ai = a->i, *aj = a->j, *bi = b->i, *bj = b->j;
237bbd65245SShri Abhyankar   const PetscInt *ajtmp, *bjtmp, *bdiag = b->diag, *pj, bs2 = a->bs2;
238209027a4SShri Abhyankar   MatScalar      *rtmp, *pc, *mwork, *v, *pv, *aa = a->a;
239bbd65245SShri Abhyankar   PetscInt        flg;
2406ba06ab7SHong Zhang   PetscReal       shift;
241a455e926SHong Zhang   PetscBool       allowzeropivot, zeropivotdetected;
242209027a4SShri Abhyankar 
243209027a4SShri Abhyankar   PetscFunctionBegin;
2440164db54SHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
2459566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isrow, &r));
2469566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isicol, &ic));
247209027a4SShri Abhyankar 
24813bcc0bdSJacob Faibussowitsch   if (info->shifttype == (PetscReal)MAT_SHIFT_NONE) {
2496ba06ab7SHong Zhang     shift = 0;
2506ba06ab7SHong Zhang   } else { /* info->shifttype == MAT_SHIFT_INBLOCKS */
2516ba06ab7SHong Zhang     shift = info->shiftamount;
2526ba06ab7SHong Zhang   }
2536ba06ab7SHong Zhang 
254209027a4SShri Abhyankar   /* generate work space needed by the factorization */
2559566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(bs2 * n, &rtmp, bs2, &mwork));
2569566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(rtmp, bs2 * n));
257209027a4SShri Abhyankar 
258209027a4SShri Abhyankar   for (i = 0; i < n; i++) {
259209027a4SShri Abhyankar     /* zero rtmp */
260209027a4SShri Abhyankar     /* L part */
261209027a4SShri Abhyankar     nz    = bi[i + 1] - bi[i];
262209027a4SShri Abhyankar     bjtmp = bj + bi[i];
26348a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2));
264209027a4SShri Abhyankar 
265209027a4SShri Abhyankar     /* U part */
26678bb4007SShri Abhyankar     nz    = bdiag[i] - bdiag[i + 1];
26778bb4007SShri Abhyankar     bjtmp = bj + bdiag[i + 1] + 1;
26848a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2));
26978bb4007SShri Abhyankar 
27078bb4007SShri Abhyankar     /* load in initial (unfactored row) */
27178bb4007SShri Abhyankar     nz    = ai[r[i] + 1] - ai[r[i]];
27278bb4007SShri Abhyankar     ajtmp = aj + ai[r[i]];
27378bb4007SShri Abhyankar     v     = aa + bs2 * ai[r[i]];
27448a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(rtmp + bs2 * ic[ajtmp[j]], v + bs2 * j, bs2));
27578bb4007SShri Abhyankar 
27678bb4007SShri Abhyankar     /* elimination */
27778bb4007SShri Abhyankar     bjtmp = bj + bi[i];
27878bb4007SShri Abhyankar     nzL   = bi[i + 1] - bi[i];
27978bb4007SShri Abhyankar     for (k = 0; k < nzL; k++) {
28078bb4007SShri Abhyankar       row = bjtmp[k];
28178bb4007SShri Abhyankar       pc  = rtmp + bs2 * row;
282c35f09e5SBarry Smith       for (flg = 0, j = 0; j < bs2; j++) {
283c35f09e5SBarry Smith         if (pc[j] != 0.0) {
284c35f09e5SBarry Smith           flg = 1;
285c35f09e5SBarry Smith           break;
286c35f09e5SBarry Smith         }
287c35f09e5SBarry Smith       }
28878bb4007SShri Abhyankar       if (flg) {
28978bb4007SShri Abhyankar         pv = b->a + bs2 * bdiag[row];
29096b95a6bSBarry Smith         /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
2919566063dSJacob Faibussowitsch         PetscCall(PetscKernel_A_gets_A_times_B_4(pc, pv, mwork));
29278bb4007SShri Abhyankar 
293a5b23f4aSJose E. Roman         pj = b->j + bdiag[row + 1] + 1; /* beginning of U(row,:) */
29478bb4007SShri Abhyankar         pv = b->a + bs2 * (bdiag[row + 1] + 1);
29578bb4007SShri Abhyankar         nz = bdiag[row] - bdiag[row + 1] - 1; /* num of entries inU(row,:), excluding diag */
29678bb4007SShri Abhyankar         for (j = 0; j < nz; j++) {
29796b95a6bSBarry Smith           /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
29878bb4007SShri Abhyankar           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
29978bb4007SShri Abhyankar           v = rtmp + bs2 * pj[j];
3009566063dSJacob Faibussowitsch           PetscCall(PetscKernel_A_gets_A_minus_B_times_C_4(v, pc, pv));
30178bb4007SShri Abhyankar           pv += bs2;
30278bb4007SShri Abhyankar         }
3039566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(128.0 * nz + 112)); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
30478bb4007SShri Abhyankar       }
30578bb4007SShri Abhyankar     }
30678bb4007SShri Abhyankar 
30778bb4007SShri Abhyankar     /* finished row so stick it into b->a */
30878bb4007SShri Abhyankar     /* L part */
30978bb4007SShri Abhyankar     pv = b->a + bs2 * bi[i];
31078bb4007SShri Abhyankar     pj = b->j + bi[i];
31178bb4007SShri Abhyankar     nz = bi[i + 1] - bi[i];
31248a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2));
31378bb4007SShri Abhyankar 
314a5b23f4aSJose E. Roman     /* Mark diagonal and invert diagonal for simpler triangular solves */
31578bb4007SShri Abhyankar     pv = b->a + bs2 * bdiag[i];
31678bb4007SShri Abhyankar     pj = b->j + bdiag[i];
3179566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(pv, rtmp + bs2 * pj[0], bs2));
3189566063dSJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_4(pv, shift, allowzeropivot, &zeropivotdetected));
3197b6c816cSBarry Smith     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
32078bb4007SShri Abhyankar 
32178bb4007SShri Abhyankar     /* U part */
32278bb4007SShri Abhyankar     pv = b->a + bs2 * (bdiag[i + 1] + 1);
32378bb4007SShri Abhyankar     pj = b->j + bdiag[i + 1] + 1;
32478bb4007SShri Abhyankar     nz = bdiag[i] - bdiag[i + 1] - 1;
32548a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2));
32678bb4007SShri Abhyankar   }
32778bb4007SShri Abhyankar 
3289566063dSJacob Faibussowitsch   PetscCall(PetscFree2(rtmp, mwork));
3299566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isicol, &ic));
3309566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow, &r));
33126fbe8dcSKarl Rupp 
3324dd39f65SShri Abhyankar   C->ops->solve          = MatSolve_SeqBAIJ_4;
3334dd39f65SShri Abhyankar   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_4;
33478bb4007SShri Abhyankar   C->assembled           = PETSC_TRUE;
33526fbe8dcSKarl Rupp 
3369566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(1.333333333333 * 4 * 4 * 4 * n)); /* from inverting diagonal blocks */
3373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
33878bb4007SShri Abhyankar }
33978bb4007SShri Abhyankar 
MatILUFactorNumeric_SeqBAIJ_4_NaturalOrdering_inplace(Mat C,Mat A,const MatFactorInfo * info)340*421480d9SBarry Smith PetscErrorCode MatILUFactorNumeric_SeqBAIJ_4_NaturalOrdering_inplace(Mat C, Mat A, const MatFactorInfo *info)
341d71ae5a4SJacob Faibussowitsch {
342e6580ceeSShri Abhyankar   /*
343e6580ceeSShri Abhyankar     Default Version for when blocks are 4 by 4 Using natural ordering
344e6580ceeSShri Abhyankar */
345e6580ceeSShri Abhyankar   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data;
346e6580ceeSShri Abhyankar   PetscInt        i, j, n = a->mbs, *bi = b->i, *bj = b->j;
347e6580ceeSShri Abhyankar   PetscInt       *ajtmpold, *ajtmp, nz, row;
348*421480d9SBarry Smith   const PetscInt *diag_offset;
349*421480d9SBarry Smith   PetscInt       *ai = a->i, *aj = a->j, *pj;
350e6580ceeSShri Abhyankar   MatScalar      *pv, *v, *rtmp, *pc, *w, *x;
351e6580ceeSShri Abhyankar   MatScalar       p1, p2, p3, p4, m1, m2, m3, m4, m5, m6, m7, m8, m9, x1, x2, x3, x4;
352e6580ceeSShri Abhyankar   MatScalar       p5, p6, p7, p8, p9, x5, x6, x7, x8, x9, x10, x11, x12, x13, x14, x15, x16;
353e6580ceeSShri Abhyankar   MatScalar       p10, p11, p12, p13, p14, p15, p16, m10, m11, m12;
354e6580ceeSShri Abhyankar   MatScalar       m13, m14, m15, m16;
355e6580ceeSShri Abhyankar   MatScalar      *ba = b->a, *aa = a->a;
356a455e926SHong Zhang   PetscBool       pivotinblocks = b->pivotinblocks;
357182b8fbaSHong Zhang   PetscReal       shift         = info->shiftamount;
3580164db54SHong Zhang   PetscBool       allowzeropivot, zeropivotdetected = PETSC_FALSE;
359e6580ceeSShri Abhyankar 
360e6580ceeSShri Abhyankar   PetscFunctionBegin;
361*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 */
362*421480d9SBarry Smith   A->factortype = MAT_FACTOR_NONE;
363*421480d9SBarry Smith   PetscCall(MatGetDiagonalMarkers_SeqBAIJ(A, &diag_offset, NULL));
364*421480d9SBarry Smith   A->factortype  = MAT_FACTOR_ILU;
3650164db54SHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
3669566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(16 * (n + 1), &rtmp));
367e6580ceeSShri Abhyankar 
368e6580ceeSShri Abhyankar   for (i = 0; i < n; i++) {
369e6580ceeSShri Abhyankar     nz    = bi[i + 1] - bi[i];
370e6580ceeSShri Abhyankar     ajtmp = bj + bi[i];
371e6580ceeSShri Abhyankar     for (j = 0; j < nz; j++) {
372e6580ceeSShri Abhyankar       x    = rtmp + 16 * ajtmp[j];
373e6580ceeSShri Abhyankar       x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = x[9] = 0.0;
374e6580ceeSShri Abhyankar       x[10] = x[11] = x[12] = x[13] = x[14] = x[15] = 0.0;
375e6580ceeSShri Abhyankar     }
376e6580ceeSShri Abhyankar     /* load in initial (unfactored row) */
377e6580ceeSShri Abhyankar     nz       = ai[i + 1] - ai[i];
378e6580ceeSShri Abhyankar     ajtmpold = aj + ai[i];
379e6580ceeSShri Abhyankar     v        = aa + 16 * ai[i];
380e6580ceeSShri Abhyankar     for (j = 0; j < nz; j++) {
381e6580ceeSShri Abhyankar       x     = rtmp + 16 * ajtmpold[j];
3829371c9d4SSatish Balay       x[0]  = v[0];
3839371c9d4SSatish Balay       x[1]  = v[1];
3849371c9d4SSatish Balay       x[2]  = v[2];
3859371c9d4SSatish Balay       x[3]  = v[3];
3869371c9d4SSatish Balay       x[4]  = v[4];
3879371c9d4SSatish Balay       x[5]  = v[5];
3889371c9d4SSatish Balay       x[6]  = v[6];
3899371c9d4SSatish Balay       x[7]  = v[7];
3909371c9d4SSatish Balay       x[8]  = v[8];
3919371c9d4SSatish Balay       x[9]  = v[9];
3929371c9d4SSatish Balay       x[10] = v[10];
3939371c9d4SSatish Balay       x[11] = v[11];
3949371c9d4SSatish Balay       x[12] = v[12];
3959371c9d4SSatish Balay       x[13] = v[13];
3969371c9d4SSatish Balay       x[14] = v[14];
3979371c9d4SSatish Balay       x[15] = v[15];
398e6580ceeSShri Abhyankar       v += 16;
399e6580ceeSShri Abhyankar     }
400e6580ceeSShri Abhyankar     row = *ajtmp++;
401e6580ceeSShri Abhyankar     while (row < i) {
402e6580ceeSShri Abhyankar       pc  = rtmp + 16 * row;
4039371c9d4SSatish Balay       p1  = pc[0];
4049371c9d4SSatish Balay       p2  = pc[1];
4059371c9d4SSatish Balay       p3  = pc[2];
4069371c9d4SSatish Balay       p4  = pc[3];
4079371c9d4SSatish Balay       p5  = pc[4];
4089371c9d4SSatish Balay       p6  = pc[5];
4099371c9d4SSatish Balay       p7  = pc[6];
4109371c9d4SSatish Balay       p8  = pc[7];
4119371c9d4SSatish Balay       p9  = pc[8];
4129371c9d4SSatish Balay       p10 = pc[9];
4139371c9d4SSatish Balay       p11 = pc[10];
4149371c9d4SSatish Balay       p12 = pc[11];
4159371c9d4SSatish Balay       p13 = pc[12];
4169371c9d4SSatish Balay       p14 = pc[13];
4179371c9d4SSatish Balay       p15 = pc[14];
4189371c9d4SSatish Balay       p16 = pc[15];
4199371c9d4SSatish 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) {
420e6580ceeSShri Abhyankar         pv    = ba + 16 * diag_offset[row];
421e6580ceeSShri Abhyankar         pj    = bj + diag_offset[row] + 1;
4229371c9d4SSatish Balay         x1    = pv[0];
4239371c9d4SSatish Balay         x2    = pv[1];
4249371c9d4SSatish Balay         x3    = pv[2];
4259371c9d4SSatish Balay         x4    = pv[3];
4269371c9d4SSatish Balay         x5    = pv[4];
4279371c9d4SSatish Balay         x6    = pv[5];
4289371c9d4SSatish Balay         x7    = pv[6];
4299371c9d4SSatish Balay         x8    = pv[7];
4309371c9d4SSatish Balay         x9    = pv[8];
4319371c9d4SSatish Balay         x10   = pv[9];
4329371c9d4SSatish Balay         x11   = pv[10];
4339371c9d4SSatish Balay         x12   = pv[11];
4349371c9d4SSatish Balay         x13   = pv[12];
4359371c9d4SSatish Balay         x14   = pv[13];
4369371c9d4SSatish Balay         x15   = pv[14];
4379371c9d4SSatish Balay         x16   = pv[15];
438e6580ceeSShri Abhyankar         pc[0] = m1 = p1 * x1 + p5 * x2 + p9 * x3 + p13 * x4;
439e6580ceeSShri Abhyankar         pc[1] = m2 = p2 * x1 + p6 * x2 + p10 * x3 + p14 * x4;
440e6580ceeSShri Abhyankar         pc[2] = m3 = p3 * x1 + p7 * x2 + p11 * x3 + p15 * x4;
441e6580ceeSShri Abhyankar         pc[3] = m4 = p4 * x1 + p8 * x2 + p12 * x3 + p16 * x4;
442e6580ceeSShri Abhyankar 
443e6580ceeSShri Abhyankar         pc[4] = m5 = p1 * x5 + p5 * x6 + p9 * x7 + p13 * x8;
444e6580ceeSShri Abhyankar         pc[5] = m6 = p2 * x5 + p6 * x6 + p10 * x7 + p14 * x8;
445e6580ceeSShri Abhyankar         pc[6] = m7 = p3 * x5 + p7 * x6 + p11 * x7 + p15 * x8;
446e6580ceeSShri Abhyankar         pc[7] = m8 = p4 * x5 + p8 * x6 + p12 * x7 + p16 * x8;
447e6580ceeSShri Abhyankar 
448e6580ceeSShri Abhyankar         pc[8] = m9 = p1 * x9 + p5 * x10 + p9 * x11 + p13 * x12;
449e6580ceeSShri Abhyankar         pc[9] = m10 = p2 * x9 + p6 * x10 + p10 * x11 + p14 * x12;
450e6580ceeSShri Abhyankar         pc[10] = m11 = p3 * x9 + p7 * x10 + p11 * x11 + p15 * x12;
451e6580ceeSShri Abhyankar         pc[11] = m12 = p4 * x9 + p8 * x10 + p12 * x11 + p16 * x12;
452e6580ceeSShri Abhyankar 
453e6580ceeSShri Abhyankar         pc[12] = m13 = p1 * x13 + p5 * x14 + p9 * x15 + p13 * x16;
454e6580ceeSShri Abhyankar         pc[13] = m14 = p2 * x13 + p6 * x14 + p10 * x15 + p14 * x16;
455e6580ceeSShri Abhyankar         pc[14] = m15 = p3 * x13 + p7 * x14 + p11 * x15 + p15 * x16;
456e6580ceeSShri Abhyankar         pc[15] = m16 = p4 * x13 + p8 * x14 + p12 * x15 + p16 * x16;
457e6580ceeSShri Abhyankar         nz           = bi[row + 1] - diag_offset[row] - 1;
458e6580ceeSShri Abhyankar         pv += 16;
459e6580ceeSShri Abhyankar         for (j = 0; j < nz; j++) {
4609371c9d4SSatish Balay           x1  = pv[0];
4619371c9d4SSatish Balay           x2  = pv[1];
4629371c9d4SSatish Balay           x3  = pv[2];
4639371c9d4SSatish Balay           x4  = pv[3];
4649371c9d4SSatish Balay           x5  = pv[4];
4659371c9d4SSatish Balay           x6  = pv[5];
4669371c9d4SSatish Balay           x7  = pv[6];
4679371c9d4SSatish Balay           x8  = pv[7];
4689371c9d4SSatish Balay           x9  = pv[8];
4699371c9d4SSatish Balay           x10 = pv[9];
4709371c9d4SSatish Balay           x11 = pv[10];
4719371c9d4SSatish Balay           x12 = pv[11];
4729371c9d4SSatish Balay           x13 = pv[12];
4739371c9d4SSatish Balay           x14 = pv[13];
4749371c9d4SSatish Balay           x15 = pv[14];
4759371c9d4SSatish Balay           x16 = pv[15];
476e6580ceeSShri Abhyankar           x   = rtmp + 16 * pj[j];
477e6580ceeSShri Abhyankar           x[0] -= m1 * x1 + m5 * x2 + m9 * x3 + m13 * x4;
478e6580ceeSShri Abhyankar           x[1] -= m2 * x1 + m6 * x2 + m10 * x3 + m14 * x4;
479e6580ceeSShri Abhyankar           x[2] -= m3 * x1 + m7 * x2 + m11 * x3 + m15 * x4;
480e6580ceeSShri Abhyankar           x[3] -= m4 * x1 + m8 * x2 + m12 * x3 + m16 * x4;
481e6580ceeSShri Abhyankar 
482e6580ceeSShri Abhyankar           x[4] -= m1 * x5 + m5 * x6 + m9 * x7 + m13 * x8;
483e6580ceeSShri Abhyankar           x[5] -= m2 * x5 + m6 * x6 + m10 * x7 + m14 * x8;
484e6580ceeSShri Abhyankar           x[6] -= m3 * x5 + m7 * x6 + m11 * x7 + m15 * x8;
485e6580ceeSShri Abhyankar           x[7] -= m4 * x5 + m8 * x6 + m12 * x7 + m16 * x8;
486e6580ceeSShri Abhyankar 
487e6580ceeSShri Abhyankar           x[8] -= m1 * x9 + m5 * x10 + m9 * x11 + m13 * x12;
488e6580ceeSShri Abhyankar           x[9] -= m2 * x9 + m6 * x10 + m10 * x11 + m14 * x12;
489e6580ceeSShri Abhyankar           x[10] -= m3 * x9 + m7 * x10 + m11 * x11 + m15 * x12;
490e6580ceeSShri Abhyankar           x[11] -= m4 * x9 + m8 * x10 + m12 * x11 + m16 * x12;
491e6580ceeSShri Abhyankar 
492e6580ceeSShri Abhyankar           x[12] -= m1 * x13 + m5 * x14 + m9 * x15 + m13 * x16;
493e6580ceeSShri Abhyankar           x[13] -= m2 * x13 + m6 * x14 + m10 * x15 + m14 * x16;
494e6580ceeSShri Abhyankar           x[14] -= m3 * x13 + m7 * x14 + m11 * x15 + m15 * x16;
495e6580ceeSShri Abhyankar           x[15] -= m4 * x13 + m8 * x14 + m12 * x15 + m16 * x16;
496e6580ceeSShri Abhyankar 
497e6580ceeSShri Abhyankar           pv += 16;
498e6580ceeSShri Abhyankar         }
4999566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(128.0 * nz + 112.0));
500e6580ceeSShri Abhyankar       }
501e6580ceeSShri Abhyankar       row = *ajtmp++;
502e6580ceeSShri Abhyankar     }
503e6580ceeSShri Abhyankar     /* finished row so stick it into b->a */
504e6580ceeSShri Abhyankar     pv = ba + 16 * bi[i];
505e6580ceeSShri Abhyankar     pj = bj + bi[i];
506e6580ceeSShri Abhyankar     nz = bi[i + 1] - bi[i];
507e6580ceeSShri Abhyankar     for (j = 0; j < nz; j++) {
508e6580ceeSShri Abhyankar       x      = rtmp + 16 * pj[j];
5099371c9d4SSatish Balay       pv[0]  = x[0];
5109371c9d4SSatish Balay       pv[1]  = x[1];
5119371c9d4SSatish Balay       pv[2]  = x[2];
5129371c9d4SSatish Balay       pv[3]  = x[3];
5139371c9d4SSatish Balay       pv[4]  = x[4];
5149371c9d4SSatish Balay       pv[5]  = x[5];
5159371c9d4SSatish Balay       pv[6]  = x[6];
5169371c9d4SSatish Balay       pv[7]  = x[7];
5179371c9d4SSatish Balay       pv[8]  = x[8];
5189371c9d4SSatish Balay       pv[9]  = x[9];
5199371c9d4SSatish Balay       pv[10] = x[10];
5209371c9d4SSatish Balay       pv[11] = x[11];
5219371c9d4SSatish Balay       pv[12] = x[12];
5229371c9d4SSatish Balay       pv[13] = x[13];
5239371c9d4SSatish Balay       pv[14] = x[14];
5249371c9d4SSatish Balay       pv[15] = x[15];
525e6580ceeSShri Abhyankar       pv += 16;
526e6580ceeSShri Abhyankar     }
527e6580ceeSShri Abhyankar     /* invert diagonal block */
528e6580ceeSShri Abhyankar     w = ba + 16 * diag_offset[i];
529e6580ceeSShri Abhyankar     if (pivotinblocks) {
5309566063dSJacob Faibussowitsch       PetscCall(PetscKernel_A_gets_inverse_A_4(w, shift, allowzeropivot, &zeropivotdetected));
5317b6c816cSBarry Smith       if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
532e6580ceeSShri Abhyankar     } else {
5339566063dSJacob Faibussowitsch       PetscCall(PetscKernel_A_gets_inverse_A_4_nopivot(w));
534e6580ceeSShri Abhyankar     }
535e6580ceeSShri Abhyankar   }
536e6580ceeSShri Abhyankar 
5379566063dSJacob Faibussowitsch   PetscCall(PetscFree(rtmp));
53826fbe8dcSKarl Rupp 
53906e38f1dSHong Zhang   C->ops->solve          = MatSolve_SeqBAIJ_4_NaturalOrdering_inplace;
54006e38f1dSHong Zhang   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_4_NaturalOrdering_inplace;
541e6580ceeSShri Abhyankar   C->assembled           = PETSC_TRUE;
54226fbe8dcSKarl Rupp 
5439566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(1.333333333333 * 4 * 4 * 4 * b->mbs)); /* from inverting diagonal blocks */
5443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
545e6580ceeSShri Abhyankar }
546c0c7eb62SShri Abhyankar 
547209027a4SShri Abhyankar /*
5484dd39f65SShri Abhyankar   MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering -
5494dd39f65SShri Abhyankar     copied from MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering_inplace()
550209027a4SShri Abhyankar */
MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering(Mat B,Mat A,const MatFactorInfo * info)551d71ae5a4SJacob Faibussowitsch PetscErrorCode MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering(Mat B, Mat A, const MatFactorInfo *info)
552d71ae5a4SJacob Faibussowitsch {
553209027a4SShri Abhyankar   Mat             C = B;
554209027a4SShri Abhyankar   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data;
555bbd65245SShri Abhyankar   PetscInt        i, j, k, nz, nzL, row;
556bbd65245SShri Abhyankar   const PetscInt  n = a->mbs, *ai = a->i, *aj = a->j, *bi = b->i, *bj = b->j;
557bbd65245SShri Abhyankar   const PetscInt *ajtmp, *bjtmp, *bdiag = b->diag, *pj, bs2 = a->bs2;
558209027a4SShri Abhyankar   MatScalar      *rtmp, *pc, *mwork, *v, *pv, *aa = a->a;
559bbd65245SShri Abhyankar   PetscInt        flg;
5606ba06ab7SHong Zhang   PetscReal       shift;
561a455e926SHong Zhang   PetscBool       allowzeropivot, zeropivotdetected;
562e6580ceeSShri Abhyankar 
563209027a4SShri Abhyankar   PetscFunctionBegin;
5640164db54SHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
5650164db54SHong Zhang 
566209027a4SShri Abhyankar   /* generate work space needed by the factorization */
5679566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(bs2 * n, &rtmp, bs2, &mwork));
5689566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(rtmp, bs2 * n));
569209027a4SShri Abhyankar 
57013bcc0bdSJacob Faibussowitsch   if (info->shifttype == (PetscReal)MAT_SHIFT_NONE) {
5716ba06ab7SHong Zhang     shift = 0;
5726ba06ab7SHong Zhang   } else { /* info->shifttype == MAT_SHIFT_INBLOCKS */
5736ba06ab7SHong Zhang     shift = info->shiftamount;
5746ba06ab7SHong Zhang   }
5756ba06ab7SHong Zhang 
576209027a4SShri Abhyankar   for (i = 0; i < n; i++) {
577209027a4SShri Abhyankar     /* zero rtmp */
578209027a4SShri Abhyankar     /* L part */
579209027a4SShri Abhyankar     nz    = bi[i + 1] - bi[i];
580209027a4SShri Abhyankar     bjtmp = bj + bi[i];
58148a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2));
582209027a4SShri Abhyankar 
583209027a4SShri Abhyankar     /* U part */
584b2b2dd24SShri Abhyankar     nz    = bdiag[i] - bdiag[i + 1];
585b2b2dd24SShri Abhyankar     bjtmp = bj + bdiag[i + 1] + 1;
58648a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2));
587b2b2dd24SShri Abhyankar 
588b2b2dd24SShri Abhyankar     /* load in initial (unfactored row) */
589b2b2dd24SShri Abhyankar     nz    = ai[i + 1] - ai[i];
590b2b2dd24SShri Abhyankar     ajtmp = aj + ai[i];
591b2b2dd24SShri Abhyankar     v     = aa + bs2 * ai[i];
59248a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(rtmp + bs2 * ajtmp[j], v + bs2 * j, bs2));
593b2b2dd24SShri Abhyankar 
594b2b2dd24SShri Abhyankar     /* elimination */
595b2b2dd24SShri Abhyankar     bjtmp = bj + bi[i];
596b2b2dd24SShri Abhyankar     nzL   = bi[i + 1] - bi[i];
597b2b2dd24SShri Abhyankar     for (k = 0; k < nzL; k++) {
598b2b2dd24SShri Abhyankar       row = bjtmp[k];
599b2b2dd24SShri Abhyankar       pc  = rtmp + bs2 * row;
600c35f09e5SBarry Smith       for (flg = 0, j = 0; j < bs2; j++) {
601c35f09e5SBarry Smith         if (pc[j] != 0.0) {
602c35f09e5SBarry Smith           flg = 1;
603c35f09e5SBarry Smith           break;
604c35f09e5SBarry Smith         }
605c35f09e5SBarry Smith       }
606b2b2dd24SShri Abhyankar       if (flg) {
607b2b2dd24SShri Abhyankar         pv = b->a + bs2 * bdiag[row];
60896b95a6bSBarry Smith         /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
6099566063dSJacob Faibussowitsch         PetscCall(PetscKernel_A_gets_A_times_B_4(pc, pv, mwork));
610b2b2dd24SShri Abhyankar 
611a5b23f4aSJose E. Roman         pj = b->j + bdiag[row + 1] + 1; /* beginning of U(row,:) */
612b2b2dd24SShri Abhyankar         pv = b->a + bs2 * (bdiag[row + 1] + 1);
613b2b2dd24SShri Abhyankar         nz = bdiag[row] - bdiag[row + 1] - 1; /* num of entries inU(row,:), excluding diag */
614b2b2dd24SShri Abhyankar         for (j = 0; j < nz; j++) {
61596b95a6bSBarry Smith           /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
616b2b2dd24SShri Abhyankar           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
617b2b2dd24SShri Abhyankar           v = rtmp + bs2 * pj[j];
6189566063dSJacob Faibussowitsch           PetscCall(PetscKernel_A_gets_A_minus_B_times_C_4(v, pc, pv));
619b2b2dd24SShri Abhyankar           pv += bs2;
620b2b2dd24SShri Abhyankar         }
6219566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(128.0 * nz + 112)); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
622b2b2dd24SShri Abhyankar       }
623b2b2dd24SShri Abhyankar     }
624b2b2dd24SShri Abhyankar 
625b2b2dd24SShri Abhyankar     /* finished row so stick it into b->a */
626b2b2dd24SShri Abhyankar     /* L part */
627b2b2dd24SShri Abhyankar     pv = b->a + bs2 * bi[i];
628b2b2dd24SShri Abhyankar     pj = b->j + bi[i];
629b2b2dd24SShri Abhyankar     nz = bi[i + 1] - bi[i];
63048a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2));
631b2b2dd24SShri Abhyankar 
632a5b23f4aSJose E. Roman     /* Mark diagonal and invert diagonal for simpler triangular solves */
633b2b2dd24SShri Abhyankar     pv = b->a + bs2 * bdiag[i];
634b2b2dd24SShri Abhyankar     pj = b->j + bdiag[i];
6359566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(pv, rtmp + bs2 * pj[0], bs2));
6369566063dSJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_4(pv, shift, allowzeropivot, &zeropivotdetected));
6377b6c816cSBarry Smith     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
638b2b2dd24SShri Abhyankar 
639b2b2dd24SShri Abhyankar     /* U part */
640b2b2dd24SShri Abhyankar     pv = b->a + bs2 * (bdiag[i + 1] + 1);
641b2b2dd24SShri Abhyankar     pj = b->j + bdiag[i + 1] + 1;
642b2b2dd24SShri Abhyankar     nz = bdiag[i] - bdiag[i + 1] - 1;
64348a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2));
644b2b2dd24SShri Abhyankar   }
6459566063dSJacob Faibussowitsch   PetscCall(PetscFree2(rtmp, mwork));
64626fbe8dcSKarl Rupp 
6474dd39f65SShri Abhyankar   C->ops->solve          = MatSolve_SeqBAIJ_4_NaturalOrdering;
6484dd39f65SShri Abhyankar   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_4_NaturalOrdering;
649b2b2dd24SShri Abhyankar   C->assembled           = PETSC_TRUE;
65026fbe8dcSKarl Rupp 
6519566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(1.333333333333 * 4 * 4 * 4 * n)); /* from inverting diagonal blocks */
6523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
653b2b2dd24SShri Abhyankar }
654