xref: /petsc/src/mat/impls/baij/seq/baijfact2.c (revision ff6a95418ff72e09afb68819ffbfc86bad111fa0)
1be1d678aSKris Buschelman 
24e2b4712SSatish Balay /*
34e2b4712SSatish Balay     Factorization code for BAIJ format.
44e2b4712SSatish Balay */
54e2b4712SSatish Balay 
6c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h>
7af0996ceSBarry Smith #include <petsc/private/kernels/blockinvert.h>
8c6db04a5SJed Brown #include <petscbt.h>
9c6db04a5SJed Brown #include <../src/mat/utils/freespace.h>
104e2b4712SSatish Balay 
114e2b4712SSatish Balay /* ----------------------------------------------------------------*/
1209573ac7SBarry Smith extern PetscErrorCode MatDuplicateNoCreate_SeqBAIJ(Mat, Mat, MatDuplicateOption, PetscBool);
136bce7ff8SHong Zhang 
14766f9fbaSBarry Smith /*
15766f9fbaSBarry Smith    This is not much faster than MatLUFactorNumeric_SeqBAIJ_N() but the solve is faster at least sometimes
16766f9fbaSBarry Smith */
17d71ae5a4SJacob Faibussowitsch PetscErrorCode MatLUFactorNumeric_SeqBAIJ_15_NaturalOrdering(Mat B, Mat A, const MatFactorInfo *info)
18d71ae5a4SJacob Faibussowitsch {
192b0b2ea7SShri Abhyankar   Mat              C = B;
202b0b2ea7SShri Abhyankar   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data;
21766f9fbaSBarry Smith   PetscInt         i, j, k, ipvt[15];
22766f9fbaSBarry Smith   const PetscInt   n = a->mbs, *ai = a->i, *aj = a->j, *bi = b->i, *bj = b->j, *ajtmp, *bjtmp, *bdiag = b->diag, *pj;
23766f9fbaSBarry Smith   PetscInt         nz, nzL, row;
24766f9fbaSBarry Smith   MatScalar       *rtmp, *pc, *mwork, *pv, *vv, work[225];
25766f9fbaSBarry Smith   const MatScalar *v, *aa = a->a;
262b0b2ea7SShri Abhyankar   PetscInt         bs2 = a->bs2, bs = A->rmap->bs, flg;
270fa040f9SShri Abhyankar   PetscInt         sol_ver;
28a455e926SHong Zhang   PetscBool        allowzeropivot, zeropivotdetected;
292b0b2ea7SShri Abhyankar 
302b0b2ea7SShri Abhyankar   PetscFunctionBegin;
310164db54SHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
329566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, ((PetscObject)A)->prefix, "-sol_ver", &sol_ver, NULL));
330fa040f9SShri Abhyankar 
342b0b2ea7SShri Abhyankar   /* generate work space needed by the factorization */
359566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(bs2 * n, &rtmp, bs2, &mwork));
369566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(rtmp, bs2 * n));
372b0b2ea7SShri Abhyankar 
382b0b2ea7SShri Abhyankar   for (i = 0; i < n; i++) {
392b0b2ea7SShri Abhyankar     /* zero rtmp */
402b0b2ea7SShri Abhyankar     /* L part */
412b0b2ea7SShri Abhyankar     nz    = bi[i + 1] - bi[i];
422b0b2ea7SShri Abhyankar     bjtmp = bj + bi[i];
4348a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2));
442b0b2ea7SShri Abhyankar 
452b0b2ea7SShri Abhyankar     /* U part */
462b0b2ea7SShri Abhyankar     nz    = bdiag[i] - bdiag[i + 1];
472b0b2ea7SShri Abhyankar     bjtmp = bj + bdiag[i + 1] + 1;
4848a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2));
492b0b2ea7SShri Abhyankar 
502b0b2ea7SShri Abhyankar     /* load in initial (unfactored row) */
5129a97285SShri Abhyankar     nz    = ai[i + 1] - ai[i];
5229a97285SShri Abhyankar     ajtmp = aj + ai[i];
5329a97285SShri Abhyankar     v     = aa + bs2 * ai[i];
5448a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(rtmp + bs2 * ajtmp[j], v + bs2 * j, bs2));
552b0b2ea7SShri Abhyankar 
562b0b2ea7SShri Abhyankar     /* elimination */
572b0b2ea7SShri Abhyankar     bjtmp = bj + bi[i];
582b0b2ea7SShri Abhyankar     nzL   = bi[i + 1] - bi[i];
592b0b2ea7SShri Abhyankar     for (k = 0; k < nzL; k++) {
602b0b2ea7SShri Abhyankar       row = bjtmp[k];
612b0b2ea7SShri Abhyankar       pc  = rtmp + bs2 * row;
62c35f09e5SBarry Smith       for (flg = 0, j = 0; j < bs2; j++) {
63c35f09e5SBarry Smith         if (pc[j] != 0.0) {
64c35f09e5SBarry Smith           flg = 1;
65c35f09e5SBarry Smith           break;
66c35f09e5SBarry Smith         }
67c35f09e5SBarry Smith       }
682b0b2ea7SShri Abhyankar       if (flg) {
692b0b2ea7SShri Abhyankar         pv = b->a + bs2 * bdiag[row];
7096b95a6bSBarry Smith         PetscKernel_A_gets_A_times_B(bs, pc, pv, mwork);
719566063dSJacob Faibussowitsch         /* PetscCall(PetscKernel_A_gets_A_times_B_15(pc,pv,mwork)); */
72a5b23f4aSJose E. Roman         pj = b->j + bdiag[row + 1] + 1; /* beginning of U(row,:) */
732b0b2ea7SShri Abhyankar         pv = b->a + bs2 * (bdiag[row + 1] + 1);
742b0b2ea7SShri Abhyankar         nz = bdiag[row] - bdiag[row + 1] - 1; /* num of entries inU(row,:), excluding diag */
752b0b2ea7SShri Abhyankar         for (j = 0; j < nz; j++) {
76766f9fbaSBarry Smith           vv = rtmp + bs2 * pj[j];
7796b95a6bSBarry Smith           PetscKernel_A_gets_A_minus_B_times_C(bs, vv, pc, pv);
789566063dSJacob Faibussowitsch           /* PetscCall(PetscKernel_A_gets_A_minus_B_times_C_15(vv,pc,pv)); */
792b0b2ea7SShri Abhyankar           pv += bs2;
802b0b2ea7SShri Abhyankar         }
819566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(2.0 * bs2 * bs * (nz + 1) - bs2)); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
822b0b2ea7SShri Abhyankar       }
832b0b2ea7SShri Abhyankar     }
842b0b2ea7SShri Abhyankar 
852b0b2ea7SShri Abhyankar     /* finished row so stick it into b->a */
862b0b2ea7SShri Abhyankar     /* L part */
872b0b2ea7SShri Abhyankar     pv = b->a + bs2 * bi[i];
882b0b2ea7SShri Abhyankar     pj = b->j + bi[i];
892b0b2ea7SShri Abhyankar     nz = bi[i + 1] - bi[i];
9048a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2));
912b0b2ea7SShri Abhyankar 
92a5b23f4aSJose E. Roman     /* Mark diagonal and invert diagonal for simpler triangular solves */
932b0b2ea7SShri Abhyankar     pv = b->a + bs2 * bdiag[i];
942b0b2ea7SShri Abhyankar     pj = b->j + bdiag[i];
959566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(pv, rtmp + bs2 * pj[0], bs2));
969566063dSJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A_15(pv, ipvt, work, info->shiftamount, allowzeropivot, &zeropivotdetected));
977b6c816cSBarry Smith     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
982b0b2ea7SShri Abhyankar 
992b0b2ea7SShri Abhyankar     /* U part */
1002b0b2ea7SShri Abhyankar     pv = b->a + bs2 * (bdiag[i + 1] + 1);
1012b0b2ea7SShri Abhyankar     pj = b->j + bdiag[i + 1] + 1;
1022b0b2ea7SShri Abhyankar     nz = bdiag[i] - bdiag[i + 1] - 1;
10348a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2));
1042b0b2ea7SShri Abhyankar   }
1052b0b2ea7SShri Abhyankar 
1069566063dSJacob Faibussowitsch   PetscCall(PetscFree2(rtmp, mwork));
10726fbe8dcSKarl Rupp 
108832cc040SShri Abhyankar   C->ops->solve          = MatSolve_SeqBAIJ_15_NaturalOrdering_ver1;
109766f9fbaSBarry Smith   C->ops->solvetranspose = MatSolve_SeqBAIJ_N_NaturalOrdering;
1102b0b2ea7SShri Abhyankar   C->assembled           = PETSC_TRUE;
11126fbe8dcSKarl Rupp 
1129566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(1.333333333333 * bs * bs2 * b->mbs)); /* from inverting diagonal blocks */
1133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1142b0b2ea7SShri Abhyankar }
1152b0b2ea7SShri Abhyankar 
116d71ae5a4SJacob Faibussowitsch PetscErrorCode MatLUFactorNumeric_SeqBAIJ_N(Mat B, Mat A, const MatFactorInfo *info)
117d71ae5a4SJacob Faibussowitsch {
1186bce7ff8SHong Zhang   Mat             C = B;
1196bce7ff8SHong Zhang   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data;
1206bce7ff8SHong Zhang   IS              isrow = b->row, isicol = b->icol;
1215a586d82SBarry Smith   const PetscInt *r, *ic;
1226bce7ff8SHong Zhang   PetscInt        i, j, k, n = a->mbs, *ai = a->i, *aj = a->j, *bi = b->i, *bj = b->j;
1236bce7ff8SHong Zhang   PetscInt       *ajtmp, *bjtmp, nz, nzL, row, *bdiag = b->diag, *pj;
124b588c5a2SHong Zhang   MatScalar      *rtmp, *pc, *mwork, *v, *pv, *aa     = a->a;
125914a18a2SHong Zhang   PetscInt        bs = A->rmap->bs, bs2 = a->bs2, *v_pivots, flg;
126914a18a2SHong Zhang   MatScalar      *v_work;
127ace3abfcSBarry Smith   PetscBool       col_identity, row_identity, both_identity;
1285f8bbccaSHong Zhang   PetscBool       allowzeropivot, zeropivotdetected;
1296bce7ff8SHong Zhang 
1306bce7ff8SHong Zhang   PetscFunctionBegin;
1319566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isrow, &r));
1329566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isicol, &ic));
1335f8bbccaSHong Zhang   allowzeropivot = PetscNot(A->erroriffailure);
134ae3d28f0SHong Zhang 
1359566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(bs2 * n, &rtmp));
1366bce7ff8SHong Zhang 
137914a18a2SHong Zhang   /* generate work space needed by dense LU factorization */
1389566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(bs, &v_work, bs2, &mwork, bs, &v_pivots));
139914a18a2SHong Zhang 
1406bce7ff8SHong Zhang   for (i = 0; i < n; i++) {
1416bce7ff8SHong Zhang     /* zero rtmp */
1426bce7ff8SHong Zhang     /* L part */
1436bce7ff8SHong Zhang     nz    = bi[i + 1] - bi[i];
1446bce7ff8SHong Zhang     bjtmp = bj + bi[i];
14548a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2));
1466bce7ff8SHong Zhang 
1476bce7ff8SHong Zhang     /* U part */
1481a83e813SShri Abhyankar     nz    = bdiag[i] - bdiag[i + 1];
1491a83e813SShri Abhyankar     bjtmp = bj + bdiag[i + 1] + 1;
15048a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2));
1511a83e813SShri Abhyankar 
1521a83e813SShri Abhyankar     /* load in initial (unfactored row) */
1531a83e813SShri Abhyankar     nz    = ai[r[i] + 1] - ai[r[i]];
1541a83e813SShri Abhyankar     ajtmp = aj + ai[r[i]];
1551a83e813SShri Abhyankar     v     = aa + bs2 * ai[r[i]];
15648a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(rtmp + bs2 * ic[ajtmp[j]], v + bs2 * j, bs2));
1571a83e813SShri Abhyankar 
1581a83e813SShri Abhyankar     /* elimination */
1591a83e813SShri Abhyankar     bjtmp = bj + bi[i];
1601a83e813SShri Abhyankar     nzL   = bi[i + 1] - bi[i];
1611a83e813SShri Abhyankar     for (k = 0; k < nzL; k++) {
1621a83e813SShri Abhyankar       row = bjtmp[k];
1631a83e813SShri Abhyankar       pc  = rtmp + bs2 * row;
164c35f09e5SBarry Smith       for (flg = 0, j = 0; j < bs2; j++) {
165c35f09e5SBarry Smith         if (pc[j] != 0.0) {
166c35f09e5SBarry Smith           flg = 1;
167c35f09e5SBarry Smith           break;
168c35f09e5SBarry Smith         }
169c35f09e5SBarry Smith       }
1701a83e813SShri Abhyankar       if (flg) {
1711a83e813SShri Abhyankar         pv = b->a + bs2 * bdiag[row];
17296b95a6bSBarry Smith         PetscKernel_A_gets_A_times_B(bs, pc, pv, mwork); /* *pc = *pc * (*pv); */
173a5b23f4aSJose E. Roman         pj = b->j + bdiag[row + 1] + 1;                  /* beginning of U(row,:) */
1741a83e813SShri Abhyankar         pv = b->a + bs2 * (bdiag[row + 1] + 1);
1751a83e813SShri Abhyankar         nz = bdiag[row] - bdiag[row + 1] - 1; /* num of entries inU(row,:), excluding diag */
176ad540459SPierre Jolivet         for (j = 0; j < nz; j++) PetscKernel_A_gets_A_minus_B_times_C(bs, rtmp + bs2 * pj[j], pc, pv + bs2 * j);
1779566063dSJacob Faibussowitsch         PetscCall(PetscLogFlops(2.0 * bs2 * bs * (nz + 1) - bs2)); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
1781a83e813SShri Abhyankar       }
1791a83e813SShri Abhyankar     }
1801a83e813SShri Abhyankar 
1811a83e813SShri Abhyankar     /* finished row so stick it into b->a */
1821a83e813SShri Abhyankar     /* L part */
1831a83e813SShri Abhyankar     pv = b->a + bs2 * bi[i];
1841a83e813SShri Abhyankar     pj = b->j + bi[i];
1851a83e813SShri Abhyankar     nz = bi[i + 1] - bi[i];
18648a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2));
1871a83e813SShri Abhyankar 
188a5b23f4aSJose E. Roman     /* Mark diagonal and invert diagonal for simpler triangular solves */
1891a83e813SShri Abhyankar     pv = b->a + bs2 * bdiag[i];
1901a83e813SShri Abhyankar     pj = b->j + bdiag[i];
1919566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(pv, rtmp + bs2 * pj[0], bs2));
1925f8bbccaSHong Zhang 
1939566063dSJacob Faibussowitsch     PetscCall(PetscKernel_A_gets_inverse_A(bs, pv, v_pivots, v_work, allowzeropivot, &zeropivotdetected));
1947b6c816cSBarry Smith     if (zeropivotdetected) B->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1951a83e813SShri Abhyankar 
1961a83e813SShri Abhyankar     /* U part */
1971a83e813SShri Abhyankar     pv = b->a + bs2 * (bdiag[i + 1] + 1);
1981a83e813SShri Abhyankar     pj = b->j + bdiag[i + 1] + 1;
1991a83e813SShri Abhyankar     nz = bdiag[i] - bdiag[i + 1] - 1;
20048a46eb9SPierre Jolivet     for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2));
2011a83e813SShri Abhyankar   }
2021a83e813SShri Abhyankar 
2039566063dSJacob Faibussowitsch   PetscCall(PetscFree(rtmp));
2049566063dSJacob Faibussowitsch   PetscCall(PetscFree3(v_work, mwork, v_pivots));
2059566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isicol, &ic));
2069566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow, &r));
2071a83e813SShri Abhyankar 
2089566063dSJacob Faibussowitsch   PetscCall(ISIdentity(isrow, &row_identity));
2099566063dSJacob Faibussowitsch   PetscCall(ISIdentity(isicol, &col_identity));
21026fbe8dcSKarl Rupp 
211ace3abfcSBarry Smith   both_identity = (PetscBool)(row_identity && col_identity);
212ae3d28f0SHong Zhang   if (both_identity) {
213ba7f0461SHong Zhang     switch (bs) {
21496e086a2SDaniel Kokron     case 9:
2155f70456aSHong Zhang #if defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX2__) && defined(__FMA__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_64BIT_INDICES)
21696e086a2SDaniel Kokron       C->ops->solve = MatSolve_SeqBAIJ_9_NaturalOrdering;
217aee5c371SBarry Smith #else
218aee5c371SBarry Smith       C->ops->solve = MatSolve_SeqBAIJ_N_NaturalOrdering;
219aee5c371SBarry Smith #endif
22096e086a2SDaniel Kokron       break;
221d71ae5a4SJacob Faibussowitsch     case 11:
222d71ae5a4SJacob Faibussowitsch       C->ops->solve = MatSolve_SeqBAIJ_11_NaturalOrdering;
223d71ae5a4SJacob Faibussowitsch       break;
224d71ae5a4SJacob Faibussowitsch     case 12:
225d71ae5a4SJacob Faibussowitsch       C->ops->solve = MatSolve_SeqBAIJ_12_NaturalOrdering;
226d71ae5a4SJacob Faibussowitsch       break;
227d71ae5a4SJacob Faibussowitsch     case 13:
228d71ae5a4SJacob Faibussowitsch       C->ops->solve = MatSolve_SeqBAIJ_13_NaturalOrdering;
229d71ae5a4SJacob Faibussowitsch       break;
230d71ae5a4SJacob Faibussowitsch     case 14:
231d71ae5a4SJacob Faibussowitsch       C->ops->solve = MatSolve_SeqBAIJ_14_NaturalOrdering;
232d71ae5a4SJacob Faibussowitsch       break;
233d71ae5a4SJacob Faibussowitsch     default:
234d71ae5a4SJacob Faibussowitsch       C->ops->solve = MatSolve_SeqBAIJ_N_NaturalOrdering;
235d71ae5a4SJacob Faibussowitsch       break;
236ba7f0461SHong Zhang     }
237ae3d28f0SHong Zhang   } else {
2384dd39f65SShri Abhyankar     C->ops->solve = MatSolve_SeqBAIJ_N;
239ae3d28f0SHong Zhang   }
2404dd39f65SShri Abhyankar   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_N;
241ae3d28f0SHong Zhang 
2421a83e813SShri Abhyankar   C->assembled = PETSC_TRUE;
24326fbe8dcSKarl Rupp 
2449566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(1.333333333333 * bs * bs2 * b->mbs)); /* from inverting diagonal blocks */
2453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2461a83e813SShri Abhyankar }
2471a83e813SShri Abhyankar 
2486bce7ff8SHong Zhang /*
2496bce7ff8SHong Zhang    ilu(0) with natural ordering under new data structure.
2504dd39f65SShri Abhyankar    See MatILUFactorSymbolic_SeqAIJ_ilu0() for detailed description
2514dd39f65SShri Abhyankar    because this code is almost identical to MatILUFactorSymbolic_SeqAIJ_ilu0_inplace().
2526bce7ff8SHong Zhang */
253c0c7eb62SShri Abhyankar 
254d71ae5a4SJacob Faibussowitsch PetscErrorCode MatILUFactorSymbolic_SeqBAIJ_ilu0(Mat fact, Mat A, IS isrow, IS iscol, const MatFactorInfo *info)
255d71ae5a4SJacob Faibussowitsch {
2566bce7ff8SHong Zhang   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b;
25716a2bf60SHong Zhang   PetscInt     n = a->mbs, *ai = a->i, *aj, *adiag = a->diag, bs2 = a->bs2;
25835aa4fcfSShri Abhyankar   PetscInt     i, j, nz, *bi, *bj, *bdiag, bi_temp;
25935aa4fcfSShri Abhyankar 
26035aa4fcfSShri Abhyankar   PetscFunctionBegin;
2619566063dSJacob Faibussowitsch   PetscCall(MatDuplicateNoCreate_SeqBAIJ(fact, A, MAT_DO_NOT_COPY_VALUES, PETSC_FALSE));
26235aa4fcfSShri Abhyankar   b = (Mat_SeqBAIJ *)(fact)->data;
26335aa4fcfSShri Abhyankar 
26435aa4fcfSShri Abhyankar   /* allocate matrix arrays for new data structure */
2659566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(bs2 * ai[n] + 1, &b->a, ai[n] + 1, &b->j, n + 1, &b->i));
26626fbe8dcSKarl Rupp 
26735aa4fcfSShri Abhyankar   b->singlemalloc    = PETSC_TRUE;
268379be0ddSLisandro Dalcin   b->free_a          = PETSC_TRUE;
269379be0ddSLisandro Dalcin   b->free_ij         = PETSC_TRUE;
2701e40a84eSLisandro Dalcin   fact->preallocated = PETSC_TRUE;
2711e40a84eSLisandro Dalcin   fact->assembled    = PETSC_TRUE;
2724dfa11a4SJacob Faibussowitsch   if (!b->diag) { PetscCall(PetscMalloc1(n + 1, &b->diag)); }
27335aa4fcfSShri Abhyankar   bdiag = b->diag;
27435aa4fcfSShri Abhyankar 
27548a46eb9SPierre Jolivet   if (n > 0) PetscCall(PetscArrayzero(b->a, bs2 * ai[n]));
27635aa4fcfSShri Abhyankar 
27735aa4fcfSShri Abhyankar   /* set bi and bj with new data structure */
27835aa4fcfSShri Abhyankar   bi = b->i;
27935aa4fcfSShri Abhyankar   bj = b->j;
28035aa4fcfSShri Abhyankar 
28135aa4fcfSShri Abhyankar   /* L part */
28235aa4fcfSShri Abhyankar   bi[0] = 0;
28335aa4fcfSShri Abhyankar   for (i = 0; i < n; i++) {
28435aa4fcfSShri Abhyankar     nz        = adiag[i] - ai[i];
28535aa4fcfSShri Abhyankar     bi[i + 1] = bi[i] + nz;
28635aa4fcfSShri Abhyankar     aj        = a->j + ai[i];
28735aa4fcfSShri Abhyankar     for (j = 0; j < nz; j++) {
2889371c9d4SSatish Balay       *bj = aj[j];
2899371c9d4SSatish Balay       bj++;
29035aa4fcfSShri Abhyankar     }
29135aa4fcfSShri Abhyankar   }
29235aa4fcfSShri Abhyankar 
29335aa4fcfSShri Abhyankar   /* U part */
29435aa4fcfSShri Abhyankar   bi_temp  = bi[n];
29535aa4fcfSShri Abhyankar   bdiag[n] = bi[n] - 1;
29635aa4fcfSShri Abhyankar   for (i = n - 1; i >= 0; i--) {
29735aa4fcfSShri Abhyankar     nz      = ai[i + 1] - adiag[i] - 1;
29835aa4fcfSShri Abhyankar     bi_temp = bi_temp + nz + 1;
29935aa4fcfSShri Abhyankar     aj      = a->j + adiag[i] + 1;
30035aa4fcfSShri Abhyankar     for (j = 0; j < nz; j++) {
3019371c9d4SSatish Balay       *bj = aj[j];
3029371c9d4SSatish Balay       bj++;
30335aa4fcfSShri Abhyankar     }
30435aa4fcfSShri Abhyankar     /* diag[i] */
3059371c9d4SSatish Balay     *bj = i;
3069371c9d4SSatish Balay     bj++;
30735aa4fcfSShri Abhyankar     bdiag[i] = bi_temp - 1;
30835aa4fcfSShri Abhyankar   }
3093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
31035aa4fcfSShri Abhyankar }
31135aa4fcfSShri Abhyankar 
312d71ae5a4SJacob Faibussowitsch PetscErrorCode MatILUFactorSymbolic_SeqBAIJ(Mat fact, Mat A, IS isrow, IS iscol, const MatFactorInfo *info)
313d71ae5a4SJacob Faibussowitsch {
31416a2bf60SHong Zhang   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ *)A->data, *b;
31516a2bf60SHong Zhang   IS                 isicol;
31616a2bf60SHong Zhang   const PetscInt    *r, *ic;
3177fa3a6a0SHong Zhang   PetscInt           n = a->mbs, *ai = a->i, *aj = a->j, d;
31816a2bf60SHong Zhang   PetscInt          *bi, *cols, nnz, *cols_lvl;
31916a2bf60SHong Zhang   PetscInt          *bdiag, prow, fm, nzbd, reallocs = 0, dcount = 0;
32016a2bf60SHong Zhang   PetscInt           i, levels, diagonal_fill;
321ace3abfcSBarry Smith   PetscBool          col_identity, row_identity, both_identity;
32216a2bf60SHong Zhang   PetscReal          f;
3230298fd71SBarry Smith   PetscInt           nlnk, *lnk, *lnk_lvl = NULL;
32416a2bf60SHong Zhang   PetscBT            lnkbt;
32516a2bf60SHong Zhang   PetscInt           nzi, *bj, **bj_ptr, **bjlvl_ptr;
3260298fd71SBarry Smith   PetscFreeSpaceList free_space = NULL, current_space = NULL;
3270298fd71SBarry Smith   PetscFreeSpaceList free_space_lvl = NULL, current_space_lvl = NULL;
328ace3abfcSBarry Smith   PetscBool          missing;
3297fa3a6a0SHong Zhang   PetscInt           bs = A->rmap->bs, bs2 = a->bs2;
33016a2bf60SHong Zhang 
33116a2bf60SHong Zhang   PetscFunctionBegin;
33208401ef6SPierre Jolivet   PetscCheck(A->rmap->n == A->cmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Must be square matrix, rows %" PetscInt_FMT " columns %" PetscInt_FMT, A->rmap->n, A->cmap->n);
3336ba06ab7SHong Zhang   if (bs > 1) { /* check shifttype */
334aed4548fSBarry Smith     PetscCheck(info->shifttype != MAT_SHIFT_NONZERO && info->shifttype != MAT_SHIFT_POSITIVE_DEFINITE, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only MAT_SHIFT_NONE and MAT_SHIFT_INBLOCKS are supported for BAIJ matrix");
3356ba06ab7SHong Zhang   }
3366ba06ab7SHong Zhang 
3379566063dSJacob Faibussowitsch   PetscCall(MatMissingDiagonal(A, &missing, &d));
33828b400f6SJacob Faibussowitsch   PetscCheck(!missing, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Matrix is missing diagonal entry %" PetscInt_FMT, d);
33916a2bf60SHong Zhang 
34016a2bf60SHong Zhang   f             = info->fill;
34116a2bf60SHong Zhang   levels        = (PetscInt)info->levels;
34216a2bf60SHong Zhang   diagonal_fill = (PetscInt)info->diagonal_fill;
34326fbe8dcSKarl Rupp 
3449566063dSJacob Faibussowitsch   PetscCall(ISInvertPermutation(iscol, PETSC_DECIDE, &isicol));
34516a2bf60SHong Zhang 
3469566063dSJacob Faibussowitsch   PetscCall(ISIdentity(isrow, &row_identity));
3479566063dSJacob Faibussowitsch   PetscCall(ISIdentity(iscol, &col_identity));
34826fbe8dcSKarl Rupp 
349ace3abfcSBarry Smith   both_identity = (PetscBool)(row_identity && col_identity);
35016a2bf60SHong Zhang 
3517fa3a6a0SHong Zhang   if (!levels && both_identity) {
35216a2bf60SHong Zhang     /* special case: ilu(0) with natural ordering */
3539566063dSJacob Faibussowitsch     PetscCall(MatILUFactorSymbolic_SeqBAIJ_ilu0(fact, A, isrow, iscol, info));
3549566063dSJacob Faibussowitsch     PetscCall(MatSeqBAIJSetNumericFactorization(fact, both_identity));
35535aa4fcfSShri Abhyankar 
356d5f3da31SBarry Smith     fact->factortype               = MAT_FACTOR_ILU;
35735aa4fcfSShri Abhyankar     (fact)->info.factor_mallocs    = 0;
35835aa4fcfSShri Abhyankar     (fact)->info.fill_ratio_given  = info->fill;
35935aa4fcfSShri Abhyankar     (fact)->info.fill_ratio_needed = 1.0;
36026fbe8dcSKarl Rupp 
36135aa4fcfSShri Abhyankar     b       = (Mat_SeqBAIJ *)(fact)->data;
36235aa4fcfSShri Abhyankar     b->row  = isrow;
36335aa4fcfSShri Abhyankar     b->col  = iscol;
36435aa4fcfSShri Abhyankar     b->icol = isicol;
3659566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)isrow));
3669566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)iscol));
36735aa4fcfSShri Abhyankar     b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE;
36826fbe8dcSKarl Rupp 
3699566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1((n + 1) * bs, &b->solve_work));
3703ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
37135aa4fcfSShri Abhyankar   }
37235aa4fcfSShri Abhyankar 
3739566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isrow, &r));
3749566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isicol, &ic));
37535aa4fcfSShri Abhyankar 
37635aa4fcfSShri Abhyankar   /* get new row pointers */
3779566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n + 1, &bi));
37835aa4fcfSShri Abhyankar   bi[0] = 0;
37935aa4fcfSShri Abhyankar   /* bdiag is location of diagonal in factor */
3809566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n + 1, &bdiag));
38135aa4fcfSShri Abhyankar   bdiag[0] = 0;
38235aa4fcfSShri Abhyankar 
3839566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(n, &bj_ptr, n, &bjlvl_ptr));
38435aa4fcfSShri Abhyankar 
38535aa4fcfSShri Abhyankar   /* create a linked list for storing column indices of the active row */
38635aa4fcfSShri Abhyankar   nlnk = n + 1;
3879566063dSJacob Faibussowitsch   PetscCall(PetscIncompleteLLCreate(n, n, nlnk, lnk, lnk_lvl, lnkbt));
38835aa4fcfSShri Abhyankar 
38935aa4fcfSShri Abhyankar   /* initial FreeSpace size is f*(ai[n]+1) */
3909566063dSJacob Faibussowitsch   PetscCall(PetscFreeSpaceGet(PetscRealIntMultTruncate(f, ai[n] + 1), &free_space));
39135aa4fcfSShri Abhyankar   current_space = free_space;
3929566063dSJacob Faibussowitsch   PetscCall(PetscFreeSpaceGet(PetscRealIntMultTruncate(f, ai[n] + 1), &free_space_lvl));
39335aa4fcfSShri Abhyankar   current_space_lvl = free_space_lvl;
39435aa4fcfSShri Abhyankar 
39535aa4fcfSShri Abhyankar   for (i = 0; i < n; i++) {
39635aa4fcfSShri Abhyankar     nzi = 0;
39735aa4fcfSShri Abhyankar     /* copy current row into linked list */
39835aa4fcfSShri Abhyankar     nnz = ai[r[i] + 1] - ai[r[i]];
39994bad497SJacob Faibussowitsch     PetscCheck(nnz, PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Empty row in matrix: row in original ordering %" PetscInt_FMT " in permuted ordering %" PetscInt_FMT, r[i], i);
40035aa4fcfSShri Abhyankar     cols   = aj + ai[r[i]];
40135aa4fcfSShri Abhyankar     lnk[i] = -1; /* marker to indicate if diagonal exists */
4029566063dSJacob Faibussowitsch     PetscCall(PetscIncompleteLLInit(nnz, cols, n, ic, &nlnk, lnk, lnk_lvl, lnkbt));
40335aa4fcfSShri Abhyankar     nzi += nlnk;
40435aa4fcfSShri Abhyankar 
40535aa4fcfSShri Abhyankar     /* make sure diagonal entry is included */
40635aa4fcfSShri Abhyankar     if (diagonal_fill && lnk[i] == -1) {
40735aa4fcfSShri Abhyankar       fm = n;
40835aa4fcfSShri Abhyankar       while (lnk[fm] < i) fm = lnk[fm];
40935aa4fcfSShri Abhyankar       lnk[i]     = lnk[fm]; /* insert diagonal into linked list */
41035aa4fcfSShri Abhyankar       lnk[fm]    = i;
41135aa4fcfSShri Abhyankar       lnk_lvl[i] = 0;
4129371c9d4SSatish Balay       nzi++;
4139371c9d4SSatish Balay       dcount++;
41435aa4fcfSShri Abhyankar     }
41535aa4fcfSShri Abhyankar 
41635aa4fcfSShri Abhyankar     /* add pivot rows into the active row */
41735aa4fcfSShri Abhyankar     nzbd = 0;
41835aa4fcfSShri Abhyankar     prow = lnk[n];
41935aa4fcfSShri Abhyankar     while (prow < i) {
42035aa4fcfSShri Abhyankar       nnz      = bdiag[prow];
42135aa4fcfSShri Abhyankar       cols     = bj_ptr[prow] + nnz + 1;
42235aa4fcfSShri Abhyankar       cols_lvl = bjlvl_ptr[prow] + nnz + 1;
42335aa4fcfSShri Abhyankar       nnz      = bi[prow + 1] - bi[prow] - nnz - 1;
42426fbe8dcSKarl Rupp 
4259566063dSJacob Faibussowitsch       PetscCall(PetscILULLAddSorted(nnz, cols, levels, cols_lvl, prow, &nlnk, lnk, lnk_lvl, lnkbt, prow));
42635aa4fcfSShri Abhyankar       nzi += nlnk;
42735aa4fcfSShri Abhyankar       prow = lnk[prow];
42835aa4fcfSShri Abhyankar       nzbd++;
42935aa4fcfSShri Abhyankar     }
43035aa4fcfSShri Abhyankar     bdiag[i]  = nzbd;
43135aa4fcfSShri Abhyankar     bi[i + 1] = bi[i] + nzi;
43235aa4fcfSShri Abhyankar 
43335aa4fcfSShri Abhyankar     /* if free space is not available, make more free space */
43435aa4fcfSShri Abhyankar     if (current_space->local_remaining < nzi) {
435f91af8c7SBarry Smith       nnz = PetscIntMultTruncate(2, PetscIntMultTruncate(nzi, (n - i))); /* estimated and max additional space needed */
4369566063dSJacob Faibussowitsch       PetscCall(PetscFreeSpaceGet(nnz, &current_space));
4379566063dSJacob Faibussowitsch       PetscCall(PetscFreeSpaceGet(nnz, &current_space_lvl));
43835aa4fcfSShri Abhyankar       reallocs++;
43935aa4fcfSShri Abhyankar     }
44035aa4fcfSShri Abhyankar 
44135aa4fcfSShri Abhyankar     /* copy data into free_space and free_space_lvl, then initialize lnk */
4429566063dSJacob Faibussowitsch     PetscCall(PetscIncompleteLLClean(n, n, nzi, lnk, lnk_lvl, current_space->array, current_space_lvl->array, lnkbt));
44326fbe8dcSKarl Rupp 
44435aa4fcfSShri Abhyankar     bj_ptr[i]    = current_space->array;
44535aa4fcfSShri Abhyankar     bjlvl_ptr[i] = current_space_lvl->array;
44635aa4fcfSShri Abhyankar 
44735aa4fcfSShri Abhyankar     /* make sure the active row i has diagonal entry */
44894bad497SJacob Faibussowitsch     PetscCheck(*(bj_ptr[i] + bdiag[i]) == i, PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Row %" PetscInt_FMT " has missing diagonal in factored matrix\ntry running with -pc_factor_nonzeros_along_diagonal or -pc_factor_diagonal_fill", i);
44935aa4fcfSShri Abhyankar 
45035aa4fcfSShri Abhyankar     current_space->array += nzi;
45135aa4fcfSShri Abhyankar     current_space->local_used += nzi;
45235aa4fcfSShri Abhyankar     current_space->local_remaining -= nzi;
45326fbe8dcSKarl Rupp 
45435aa4fcfSShri Abhyankar     current_space_lvl->array += nzi;
45535aa4fcfSShri Abhyankar     current_space_lvl->local_used += nzi;
45635aa4fcfSShri Abhyankar     current_space_lvl->local_remaining -= nzi;
45735aa4fcfSShri Abhyankar   }
45835aa4fcfSShri Abhyankar 
4599566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow, &r));
4609566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isicol, &ic));
46135aa4fcfSShri Abhyankar 
46235aa4fcfSShri Abhyankar   /* copy free_space into bj and free free_space; set bi, bj, bdiag in new datastructure; */
4639566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bi[n] + 1, &bj));
4649566063dSJacob Faibussowitsch   PetscCall(PetscFreeSpaceContiguous_LU(&free_space, bj, n, bi, bdiag));
46535aa4fcfSShri Abhyankar 
4669566063dSJacob Faibussowitsch   PetscCall(PetscIncompleteLLDestroy(lnk, lnkbt));
4679566063dSJacob Faibussowitsch   PetscCall(PetscFreeSpaceDestroy(free_space_lvl));
4689566063dSJacob Faibussowitsch   PetscCall(PetscFree2(bj_ptr, bjlvl_ptr));
46935aa4fcfSShri Abhyankar 
47035aa4fcfSShri Abhyankar #if defined(PETSC_USE_INFO)
47135aa4fcfSShri Abhyankar   {
472aef85c9fSShri Abhyankar     PetscReal af = ((PetscReal)(bdiag[0] + 1)) / ((PetscReal)ai[n]);
4739566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "Reallocs %" PetscInt_FMT " Fill ratio:given %g needed %g\n", reallocs, (double)f, (double)af));
4749566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "Run with -[sub_]pc_factor_fill %g or use \n", (double)af));
4759566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "PCFactorSetFill([sub]pc,%g);\n", (double)af));
4769566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "for best performance.\n"));
47748a46eb9SPierre Jolivet     if (diagonal_fill) PetscCall(PetscInfo(A, "Detected and replaced %" PetscInt_FMT " missing diagonals\n", dcount));
47835aa4fcfSShri Abhyankar   }
47935aa4fcfSShri Abhyankar #endif
48035aa4fcfSShri Abhyankar 
48135aa4fcfSShri Abhyankar   /* put together the new matrix */
4829566063dSJacob Faibussowitsch   PetscCall(MatSeqBAIJSetPreallocation(fact, bs, MAT_SKIP_ALLOCATION, NULL));
48326fbe8dcSKarl Rupp 
48435aa4fcfSShri Abhyankar   b               = (Mat_SeqBAIJ *)(fact)->data;
48535aa4fcfSShri Abhyankar   b->free_a       = PETSC_TRUE;
48635aa4fcfSShri Abhyankar   b->free_ij      = PETSC_TRUE;
48735aa4fcfSShri Abhyankar   b->singlemalloc = PETSC_FALSE;
48826fbe8dcSKarl Rupp 
4899566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs2 * (bdiag[0] + 1), &b->a));
49026fbe8dcSKarl Rupp 
49135aa4fcfSShri Abhyankar   b->j         = bj;
49235aa4fcfSShri Abhyankar   b->i         = bi;
49335aa4fcfSShri Abhyankar   b->diag      = bdiag;
49435aa4fcfSShri Abhyankar   b->free_diag = PETSC_TRUE;
495f4259b30SLisandro Dalcin   b->ilen      = NULL;
496f4259b30SLisandro Dalcin   b->imax      = NULL;
49735aa4fcfSShri Abhyankar   b->row       = isrow;
49835aa4fcfSShri Abhyankar   b->col       = iscol;
4999566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)isrow));
5009566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)iscol));
50135aa4fcfSShri Abhyankar   b->icol = isicol;
50226fbe8dcSKarl Rupp 
5039566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs * n + bs, &b->solve_work));
50435aa4fcfSShri Abhyankar   /* In b structure:  Free imax, ilen, old a, old j.
50535aa4fcfSShri Abhyankar      Allocate bdiag, solve_work, new a, new j */
50635aa4fcfSShri Abhyankar   b->maxnz = b->nz = bdiag[0] + 1;
50726fbe8dcSKarl Rupp 
508ae3d28f0SHong Zhang   fact->info.factor_mallocs    = reallocs;
509ae3d28f0SHong Zhang   fact->info.fill_ratio_given  = f;
510ae3d28f0SHong Zhang   fact->info.fill_ratio_needed = ((PetscReal)(bdiag[0] + 1)) / ((PetscReal)ai[n]);
51126fbe8dcSKarl Rupp 
5129566063dSJacob Faibussowitsch   PetscCall(MatSeqBAIJSetNumericFactorization(fact, both_identity));
5133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
51435aa4fcfSShri Abhyankar }
51535aa4fcfSShri Abhyankar 
516*ff6a9541SJacob Faibussowitsch #if 0
517*ff6a9541SJacob Faibussowitsch // unused
5184e2b4712SSatish Balay /*
5194e2b4712SSatish Balay      This code is virtually identical to MatILUFactorSymbolic_SeqAIJ
5204e2b4712SSatish Balay    except that the data structure of Mat_SeqAIJ is slightly different.
5214e2b4712SSatish Balay    Not a good example of code reuse.
5224e2b4712SSatish Balay */
523*ff6a9541SJacob Faibussowitsch static PetscErrorCode MatILUFactorSymbolic_SeqBAIJ_inplace(Mat fact, Mat A, IS isrow, IS iscol, const MatFactorInfo *info)
524d71ae5a4SJacob Faibussowitsch {
5254e2b4712SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *)A->data, *b;
5264e2b4712SSatish Balay   IS              isicol;
5275d0c19d7SBarry Smith   const PetscInt *r, *ic, *ai = a->i, *aj = a->j, *xi;
5285d0c19d7SBarry Smith   PetscInt        prow, n = a->mbs, *ainew, *ajnew, jmax, *fill, nz, *im, *ajfill, *flev, *xitmp;
529a96a251dSBarry Smith   PetscInt       *dloc, idx, row, m, fm, nzf, nzi, reallocate = 0, dcount = 0;
530d0f46423SBarry Smith   PetscInt        incrlev, nnz, i, bs = A->rmap->bs, bs2 = a->bs2, levels, diagonal_fill, dd;
531ace3abfcSBarry Smith   PetscBool       col_identity, row_identity, both_identity, flg;
532329f5518SBarry Smith   PetscReal       f;
5334e2b4712SSatish Balay 
5344e2b4712SSatish Balay   PetscFunctionBegin;
5359566063dSJacob Faibussowitsch   PetscCall(MatMissingDiagonal_SeqBAIJ(A, &flg, &dd));
53628b400f6SJacob Faibussowitsch   PetscCheck(!flg, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Matrix A is missing diagonal entry in row %" PetscInt_FMT, dd);
5376bce7ff8SHong Zhang 
538435faa5fSBarry Smith   f             = info->fill;
539690b6cddSBarry Smith   levels        = (PetscInt)info->levels;
540690b6cddSBarry Smith   diagonal_fill = (PetscInt)info->diagonal_fill;
54126fbe8dcSKarl Rupp 
5429566063dSJacob Faibussowitsch   PetscCall(ISInvertPermutation(iscol, PETSC_DECIDE, &isicol));
54316a2bf60SHong Zhang 
5449566063dSJacob Faibussowitsch   PetscCall(ISIdentity(isrow, &row_identity));
5459566063dSJacob Faibussowitsch   PetscCall(ISIdentity(iscol, &col_identity));
546ace3abfcSBarry Smith   both_identity = (PetscBool)(row_identity && col_identity);
547309c388cSBarry Smith 
54841df41f0SMatthew Knepley   if (!levels && both_identity) { /* special case copy the nonzero structure */
5499566063dSJacob Faibussowitsch     PetscCall(MatDuplicateNoCreate_SeqBAIJ(fact, A, MAT_DO_NOT_COPY_VALUES, PETSC_TRUE));
5509566063dSJacob Faibussowitsch     PetscCall(MatSeqBAIJSetNumericFactorization_inplace(fact, both_identity));
5516bce7ff8SHong Zhang 
552d5f3da31SBarry Smith     fact->factortype = MAT_FACTOR_ILU;
553ae3d28f0SHong Zhang     b                = (Mat_SeqBAIJ *)fact->data;
554bb3d539aSBarry Smith     b->row           = isrow;
555bb3d539aSBarry Smith     b->col           = iscol;
5569566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)isrow));
5579566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)iscol));
558bb3d539aSBarry Smith     b->icol          = isicol;
559bcd9e38bSBarry Smith     b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE;
56026fbe8dcSKarl Rupp 
5619566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1((n + 1) * bs, &b->solve_work));
5623ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
5636bce7ff8SHong Zhang   }
5646bce7ff8SHong Zhang 
5656bce7ff8SHong Zhang   /* general case perform the symbolic factorization */
5669566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isrow, &r));
5679566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isicol, &ic));
5684e2b4712SSatish Balay 
5694e2b4712SSatish Balay   /* get new row pointers */
5709566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n + 1, &ainew));
5714e2b4712SSatish Balay   ainew[0] = 0;
5724e2b4712SSatish Balay   /* don't know how many column pointers are needed so estimate */
573690b6cddSBarry Smith   jmax = (PetscInt)(f * ai[n] + 1);
5749566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(jmax, &ajnew));
5754e2b4712SSatish Balay   /* ajfill is level of fill for each fill entry */
5769566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(jmax, &ajfill));
5774e2b4712SSatish Balay   /* fill is a linked list of nonzeros in active row */
5789566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n + 1, &fill));
5794e2b4712SSatish Balay   /* im is level for each filled value */
5809566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n + 1, &im));
5814e2b4712SSatish Balay   /* dloc is location of diagonal in factor */
5829566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n + 1, &dloc));
5834e2b4712SSatish Balay   dloc[0] = 0;
5844e2b4712SSatish Balay   for (prow = 0; prow < n; prow++) {
585435faa5fSBarry Smith     /* copy prow into linked list */
5864e2b4712SSatish Balay     nzf = nz = ai[r[prow] + 1] - ai[r[prow]];
58728b400f6SJacob Faibussowitsch     PetscCheck(nz, PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Empty row in matrix: row in original ordering %" PetscInt_FMT " in permuted ordering %" PetscInt_FMT, r[prow], prow);
5884e2b4712SSatish Balay     xi         = aj + ai[r[prow]];
5894e2b4712SSatish Balay     fill[n]    = n;
590435faa5fSBarry Smith     fill[prow] = -1; /* marker for diagonal entry */
5914e2b4712SSatish Balay     while (nz--) {
5924e2b4712SSatish Balay       fm  = n;
5934e2b4712SSatish Balay       idx = ic[*xi++];
5944e2b4712SSatish Balay       do {
5954e2b4712SSatish Balay         m  = fm;
5964e2b4712SSatish Balay         fm = fill[m];
5974e2b4712SSatish Balay       } while (fm < idx);
5984e2b4712SSatish Balay       fill[m]   = idx;
5994e2b4712SSatish Balay       fill[idx] = fm;
6004e2b4712SSatish Balay       im[idx]   = 0;
6014e2b4712SSatish Balay     }
602435faa5fSBarry Smith 
603435faa5fSBarry Smith     /* make sure diagonal entry is included */
604435faa5fSBarry Smith     if (diagonal_fill && fill[prow] == -1) {
605435faa5fSBarry Smith       fm = n;
606435faa5fSBarry Smith       while (fill[fm] < prow) fm = fill[fm];
607435faa5fSBarry Smith       fill[prow] = fill[fm]; /* insert diagonal into linked list */
608435faa5fSBarry Smith       fill[fm]   = prow;
609435faa5fSBarry Smith       im[prow]   = 0;
610435faa5fSBarry Smith       nzf++;
611335d9088SBarry Smith       dcount++;
612435faa5fSBarry Smith     }
613435faa5fSBarry Smith 
6144e2b4712SSatish Balay     nzi = 0;
6154e2b4712SSatish Balay     row = fill[n];
6164e2b4712SSatish Balay     while (row < prow) {
6174e2b4712SSatish Balay       incrlev = im[row] + 1;
6184e2b4712SSatish Balay       nz      = dloc[row];
619435faa5fSBarry Smith       xi      = ajnew + ainew[row] + nz + 1;
6204e2b4712SSatish Balay       flev    = ajfill + ainew[row] + nz + 1;
6214e2b4712SSatish Balay       nnz     = ainew[row + 1] - ainew[row] - nz - 1;
6224e2b4712SSatish Balay       fm      = row;
6234e2b4712SSatish Balay       while (nnz-- > 0) {
6244e2b4712SSatish Balay         idx = *xi++;
6254e2b4712SSatish Balay         if (*flev + incrlev > levels) {
6264e2b4712SSatish Balay           flev++;
6274e2b4712SSatish Balay           continue;
6284e2b4712SSatish Balay         }
6294e2b4712SSatish Balay         do {
6304e2b4712SSatish Balay           m  = fm;
6314e2b4712SSatish Balay           fm = fill[m];
6324e2b4712SSatish Balay         } while (fm < idx);
6334e2b4712SSatish Balay         if (fm != idx) {
6344e2b4712SSatish Balay           im[idx]   = *flev + incrlev;
6354e2b4712SSatish Balay           fill[m]   = idx;
6364e2b4712SSatish Balay           fill[idx] = fm;
6374e2b4712SSatish Balay           fm        = idx;
6384e2b4712SSatish Balay           nzf++;
63926fbe8dcSKarl Rupp         } else if (im[idx] > *flev + incrlev) im[idx] = *flev + incrlev;
6404e2b4712SSatish Balay         flev++;
6414e2b4712SSatish Balay       }
6424e2b4712SSatish Balay       row = fill[row];
6434e2b4712SSatish Balay       nzi++;
6444e2b4712SSatish Balay     }
6454e2b4712SSatish Balay     /* copy new filled row into permanent storage */
6464e2b4712SSatish Balay     ainew[prow + 1] = ainew[prow] + nzf;
6474e2b4712SSatish Balay     if (ainew[prow + 1] > jmax) {
648ecf371e4SBarry Smith       /* estimate how much additional space we will need */
649ecf371e4SBarry Smith       /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */
650ecf371e4SBarry Smith       /* just double the memory each time */
651690b6cddSBarry Smith       PetscInt maxadd = jmax;
652ecf371e4SBarry Smith       /* maxadd = (int)(((f*ai[n]+1)*(n-prow+5))/n); */
6534e2b4712SSatish Balay       if (maxadd < nzf) maxadd = (n - prow) * (nzf + 1);
6544e2b4712SSatish Balay       jmax += maxadd;
655ecf371e4SBarry Smith 
656ecf371e4SBarry Smith       /* allocate a longer ajnew and ajfill */
6579566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(jmax, &xitmp));
6589566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(xitmp, ajnew, ainew[prow]));
6599566063dSJacob Faibussowitsch       PetscCall(PetscFree(ajnew));
6605d0c19d7SBarry Smith       ajnew = xitmp;
6619566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(jmax, &xitmp));
6629566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(xitmp, ajfill, ainew[prow]));
6639566063dSJacob Faibussowitsch       PetscCall(PetscFree(ajfill));
6645d0c19d7SBarry Smith       ajfill = xitmp;
665eb150c5cSKris Buschelman       reallocate++; /* count how many reallocations are needed */
6664e2b4712SSatish Balay     }
6675d0c19d7SBarry Smith     xitmp      = ajnew + ainew[prow];
6684e2b4712SSatish Balay     flev       = ajfill + ainew[prow];
6694e2b4712SSatish Balay     dloc[prow] = nzi;
6704e2b4712SSatish Balay     fm         = fill[n];
6714e2b4712SSatish Balay     while (nzf--) {
6725d0c19d7SBarry Smith       *xitmp++ = fm;
6734e2b4712SSatish Balay       *flev++  = im[fm];
6744e2b4712SSatish Balay       fm       = fill[fm];
6754e2b4712SSatish Balay     }
676435faa5fSBarry Smith     /* make sure row has diagonal entry */
677aed4548fSBarry Smith     PetscCheck(ajnew[ainew[prow] + dloc[prow]] == prow, PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Row %" PetscInt_FMT " has missing diagonal in factored matrix\n\
6789371c9d4SSatish Balay                                                         try running with -pc_factor_nonzeros_along_diagonal or -pc_factor_diagonal_fill",
6799371c9d4SSatish Balay                prow);
680435faa5fSBarry Smith   }
6819566063dSJacob Faibussowitsch   PetscCall(PetscFree(ajfill));
6829566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrow, &r));
6839566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isicol, &ic));
6849566063dSJacob Faibussowitsch   PetscCall(PetscFree(fill));
6859566063dSJacob Faibussowitsch   PetscCall(PetscFree(im));
6864e2b4712SSatish Balay 
6876cf91177SBarry Smith   #if defined(PETSC_USE_INFO)
6884e2b4712SSatish Balay   {
689329f5518SBarry Smith     PetscReal af = ((PetscReal)ainew[n]) / ((PetscReal)ai[n]);
6909566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "Reallocs %" PetscInt_FMT " Fill ratio:given %g needed %g\n", reallocate, (double)f, (double)af));
6919566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "Run with -pc_factor_fill %g or use \n", (double)af));
6929566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "PCFactorSetFill(pc,%g);\n", (double)af));
6939566063dSJacob Faibussowitsch     PetscCall(PetscInfo(A, "for best performance.\n"));
69448a46eb9SPierre Jolivet     if (diagonal_fill) PetscCall(PetscInfo(A, "Detected and replaced %" PetscInt_FMT " missing diagonals\n", dcount));
6954e2b4712SSatish Balay   }
69663ba0a88SBarry Smith   #endif
6974e2b4712SSatish Balay 
6984e2b4712SSatish Balay   /* put together the new matrix */
6999566063dSJacob Faibussowitsch   PetscCall(MatSeqBAIJSetPreallocation(fact, bs, MAT_SKIP_ALLOCATION, NULL));
700ae3d28f0SHong Zhang   b = (Mat_SeqBAIJ *)fact->data;
70126fbe8dcSKarl Rupp 
702e6b907acSBarry Smith   b->free_a       = PETSC_TRUE;
703e6b907acSBarry Smith   b->free_ij      = PETSC_TRUE;
7047c922b88SBarry Smith   b->singlemalloc = PETSC_FALSE;
70526fbe8dcSKarl Rupp 
7069566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs2 * ainew[n], &b->a));
70726fbe8dcSKarl Rupp 
7084e2b4712SSatish Balay   b->j = ajnew;
7094e2b4712SSatish Balay   b->i = ainew;
7104e2b4712SSatish Balay   for (i = 0; i < n; i++) dloc[i] += ainew[i];
7114e2b4712SSatish Balay   b->diag          = dloc;
7127f53bb6cSHong Zhang   b->free_diag     = PETSC_TRUE;
713f4259b30SLisandro Dalcin   b->ilen          = NULL;
714f4259b30SLisandro Dalcin   b->imax          = NULL;
7154e2b4712SSatish Balay   b->row           = isrow;
7164e2b4712SSatish Balay   b->col           = iscol;
717bcd9e38bSBarry Smith   b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE;
71826fbe8dcSKarl Rupp 
7199566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)isrow));
7209566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)iscol));
721e51c0b9cSSatish Balay   b->icol = isicol;
7229566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bs * n + bs, &b->solve_work));
7234e2b4712SSatish Balay   /* In b structure:  Free imax, ilen, old a, old j.
7244e2b4712SSatish Balay      Allocate dloc, solve_work, new a, new j */
7254e2b4712SSatish Balay   b->maxnz = b->nz = ainew[n];
7264e2b4712SSatish Balay 
727ae3d28f0SHong Zhang   fact->info.factor_mallocs    = reallocate;
728ae3d28f0SHong Zhang   fact->info.fill_ratio_given  = f;
729ae3d28f0SHong Zhang   fact->info.fill_ratio_needed = ((PetscReal)ainew[n]) / ((PetscReal)ai[prow]);
7306bce7ff8SHong Zhang 
7319566063dSJacob Faibussowitsch   PetscCall(MatSeqBAIJSetNumericFactorization_inplace(fact, both_identity));
7323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7338661488fSKris Buschelman }
734*ff6a9541SJacob Faibussowitsch #endif
7358661488fSKris Buschelman 
736d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE(Mat A)
737d71ae5a4SJacob Faibussowitsch {
73812272027SHong Zhang   /* Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; */
73912272027SHong Zhang   /* int i,*AJ=a->j,nz=a->nz; */
7405fd66863SKarl Rupp 
7415a9542e3SKris Buschelman   PetscFunctionBegin;
7427cf1b8d3SKris Buschelman   /* Undo Column scaling */
7437cf1b8d3SKris Buschelman   /*    while (nz--) { */
7447cf1b8d3SKris Buschelman   /*      AJ[i] = AJ[i]/4; */
7457cf1b8d3SKris Buschelman   /*    } */
746c115a38dSKris Buschelman   /* This should really invoke a push/pop logic, but we don't have that yet. */
7470298fd71SBarry Smith   A->ops->setunfactored = NULL;
7483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7497cf1b8d3SKris Buschelman }
7507cf1b8d3SKris Buschelman 
751d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE_usj(Mat A)
752d71ae5a4SJacob Faibussowitsch {
7537cf1b8d3SKris Buschelman   Mat_SeqBAIJ    *a  = (Mat_SeqBAIJ *)A->data;
754b24ad042SBarry Smith   PetscInt       *AJ = a->j, nz = a->nz;
7552aa5897fSKris Buschelman   unsigned short *aj = (unsigned short *)AJ;
7565fd66863SKarl Rupp 
7575a9542e3SKris Buschelman   PetscFunctionBegin;
7580b9da03eSKris Buschelman   /* Is this really necessary? */
7599371c9d4SSatish Balay   while (nz--) { AJ[nz] = (int)((unsigned int)aj[nz]); /* First extend, then convert to signed. */ }
7600298fd71SBarry Smith   A->ops->setunfactored = NULL;
7613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7627e7071cdSKris Buschelman }
763