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 1109573ac7SBarry Smith extern PetscErrorCode MatDuplicateNoCreate_SeqBAIJ(Mat, Mat, MatDuplicateOption, PetscBool); 126bce7ff8SHong Zhang 13766f9fbaSBarry Smith /* 14766f9fbaSBarry Smith This is not much faster than MatLUFactorNumeric_SeqBAIJ_N() but the solve is faster at least sometimes 15766f9fbaSBarry Smith */ 16d71ae5a4SJacob Faibussowitsch PetscErrorCode MatLUFactorNumeric_SeqBAIJ_15_NaturalOrdering(Mat B, Mat A, const MatFactorInfo *info) 17d71ae5a4SJacob Faibussowitsch { 182b0b2ea7SShri Abhyankar Mat C = B; 192b0b2ea7SShri Abhyankar Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data; 20766f9fbaSBarry Smith PetscInt i, j, k, ipvt[15]; 21766f9fbaSBarry Smith const PetscInt n = a->mbs, *ai = a->i, *aj = a->j, *bi = b->i, *bj = b->j, *ajtmp, *bjtmp, *bdiag = b->diag, *pj; 22766f9fbaSBarry Smith PetscInt nz, nzL, row; 23766f9fbaSBarry Smith MatScalar *rtmp, *pc, *mwork, *pv, *vv, work[225]; 24766f9fbaSBarry Smith const MatScalar *v, *aa = a->a; 252b0b2ea7SShri Abhyankar PetscInt bs2 = a->bs2, bs = A->rmap->bs, flg; 260fa040f9SShri Abhyankar PetscInt sol_ver; 27a455e926SHong Zhang PetscBool allowzeropivot, zeropivotdetected; 282b0b2ea7SShri Abhyankar 292b0b2ea7SShri Abhyankar PetscFunctionBegin; 300164db54SHong Zhang allowzeropivot = PetscNot(A->erroriffailure); 319566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, ((PetscObject)A)->prefix, "-sol_ver", &sol_ver, NULL)); 320fa040f9SShri Abhyankar 332b0b2ea7SShri Abhyankar /* generate work space needed by the factorization */ 349566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(bs2 * n, &rtmp, bs2, &mwork)); 359566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(rtmp, bs2 * n)); 362b0b2ea7SShri Abhyankar 372b0b2ea7SShri Abhyankar for (i = 0; i < n; i++) { 382b0b2ea7SShri Abhyankar /* zero rtmp */ 392b0b2ea7SShri Abhyankar /* L part */ 402b0b2ea7SShri Abhyankar nz = bi[i + 1] - bi[i]; 412b0b2ea7SShri Abhyankar bjtmp = bj + bi[i]; 4248a46eb9SPierre Jolivet for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2)); 432b0b2ea7SShri Abhyankar 442b0b2ea7SShri Abhyankar /* U part */ 452b0b2ea7SShri Abhyankar nz = bdiag[i] - bdiag[i + 1]; 462b0b2ea7SShri Abhyankar bjtmp = bj + bdiag[i + 1] + 1; 4748a46eb9SPierre Jolivet for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2)); 482b0b2ea7SShri Abhyankar 492b0b2ea7SShri Abhyankar /* load in initial (unfactored row) */ 5029a97285SShri Abhyankar nz = ai[i + 1] - ai[i]; 5129a97285SShri Abhyankar ajtmp = aj + ai[i]; 5229a97285SShri Abhyankar v = aa + bs2 * ai[i]; 5348a46eb9SPierre Jolivet for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(rtmp + bs2 * ajtmp[j], v + bs2 * j, bs2)); 542b0b2ea7SShri Abhyankar 552b0b2ea7SShri Abhyankar /* elimination */ 562b0b2ea7SShri Abhyankar bjtmp = bj + bi[i]; 572b0b2ea7SShri Abhyankar nzL = bi[i + 1] - bi[i]; 582b0b2ea7SShri Abhyankar for (k = 0; k < nzL; k++) { 592b0b2ea7SShri Abhyankar row = bjtmp[k]; 602b0b2ea7SShri Abhyankar pc = rtmp + bs2 * row; 61c35f09e5SBarry Smith for (flg = 0, j = 0; j < bs2; j++) { 62c35f09e5SBarry Smith if (pc[j] != 0.0) { 63c35f09e5SBarry Smith flg = 1; 64c35f09e5SBarry Smith break; 65c35f09e5SBarry Smith } 66c35f09e5SBarry Smith } 672b0b2ea7SShri Abhyankar if (flg) { 682b0b2ea7SShri Abhyankar pv = b->a + bs2 * bdiag[row]; 6996b95a6bSBarry Smith PetscKernel_A_gets_A_times_B(bs, pc, pv, mwork); 709566063dSJacob Faibussowitsch /* PetscCall(PetscKernel_A_gets_A_times_B_15(pc,pv,mwork)); */ 71a5b23f4aSJose E. Roman pj = b->j + bdiag[row + 1] + 1; /* beginning of U(row,:) */ 722b0b2ea7SShri Abhyankar pv = b->a + bs2 * (bdiag[row + 1] + 1); 732b0b2ea7SShri Abhyankar nz = bdiag[row] - bdiag[row + 1] - 1; /* num of entries inU(row,:), excluding diag */ 742b0b2ea7SShri Abhyankar for (j = 0; j < nz; j++) { 75766f9fbaSBarry Smith vv = rtmp + bs2 * pj[j]; 7696b95a6bSBarry Smith PetscKernel_A_gets_A_minus_B_times_C(bs, vv, pc, pv); 779566063dSJacob Faibussowitsch /* PetscCall(PetscKernel_A_gets_A_minus_B_times_C_15(vv,pc,pv)); */ 782b0b2ea7SShri Abhyankar pv += bs2; 792b0b2ea7SShri Abhyankar } 809566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * bs2 * bs * (nz + 1) - bs2)); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */ 812b0b2ea7SShri Abhyankar } 822b0b2ea7SShri Abhyankar } 832b0b2ea7SShri Abhyankar 842b0b2ea7SShri Abhyankar /* finished row so stick it into b->a */ 852b0b2ea7SShri Abhyankar /* L part */ 862b0b2ea7SShri Abhyankar pv = b->a + bs2 * bi[i]; 872b0b2ea7SShri Abhyankar pj = b->j + bi[i]; 882b0b2ea7SShri Abhyankar nz = bi[i + 1] - bi[i]; 8948a46eb9SPierre Jolivet for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2)); 902b0b2ea7SShri Abhyankar 91a5b23f4aSJose E. Roman /* Mark diagonal and invert diagonal for simpler triangular solves */ 922b0b2ea7SShri Abhyankar pv = b->a + bs2 * bdiag[i]; 932b0b2ea7SShri Abhyankar pj = b->j + bdiag[i]; 949566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(pv, rtmp + bs2 * pj[0], bs2)); 959566063dSJacob Faibussowitsch PetscCall(PetscKernel_A_gets_inverse_A_15(pv, ipvt, work, info->shiftamount, allowzeropivot, &zeropivotdetected)); 967b6c816cSBarry Smith if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 972b0b2ea7SShri Abhyankar 982b0b2ea7SShri Abhyankar /* U part */ 992b0b2ea7SShri Abhyankar pv = b->a + bs2 * (bdiag[i + 1] + 1); 1002b0b2ea7SShri Abhyankar pj = b->j + bdiag[i + 1] + 1; 1012b0b2ea7SShri Abhyankar nz = bdiag[i] - bdiag[i + 1] - 1; 10248a46eb9SPierre Jolivet for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2)); 1032b0b2ea7SShri Abhyankar } 1042b0b2ea7SShri Abhyankar 1059566063dSJacob Faibussowitsch PetscCall(PetscFree2(rtmp, mwork)); 10626fbe8dcSKarl Rupp 107832cc040SShri Abhyankar C->ops->solve = MatSolve_SeqBAIJ_15_NaturalOrdering_ver1; 108766f9fbaSBarry Smith C->ops->solvetranspose = MatSolve_SeqBAIJ_N_NaturalOrdering; 1092b0b2ea7SShri Abhyankar C->assembled = PETSC_TRUE; 11026fbe8dcSKarl Rupp 1119566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(1.333333333333 * bs * bs2 * b->mbs)); /* from inverting diagonal blocks */ 1123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1132b0b2ea7SShri Abhyankar } 1142b0b2ea7SShri Abhyankar 115d71ae5a4SJacob Faibussowitsch PetscErrorCode MatLUFactorNumeric_SeqBAIJ_N(Mat B, Mat A, const MatFactorInfo *info) 116d71ae5a4SJacob Faibussowitsch { 1176bce7ff8SHong Zhang Mat C = B; 1186bce7ff8SHong Zhang Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data; 1196bce7ff8SHong Zhang IS isrow = b->row, isicol = b->icol; 1205a586d82SBarry Smith const PetscInt *r, *ic; 1216bce7ff8SHong Zhang PetscInt i, j, k, n = a->mbs, *ai = a->i, *aj = a->j, *bi = b->i, *bj = b->j; 1226bce7ff8SHong Zhang PetscInt *ajtmp, *bjtmp, nz, nzL, row, *bdiag = b->diag, *pj; 123b588c5a2SHong Zhang MatScalar *rtmp, *pc, *mwork, *v, *pv, *aa = a->a; 124914a18a2SHong Zhang PetscInt bs = A->rmap->bs, bs2 = a->bs2, *v_pivots, flg; 125914a18a2SHong Zhang MatScalar *v_work; 126ace3abfcSBarry Smith PetscBool col_identity, row_identity, both_identity; 1275f8bbccaSHong Zhang PetscBool allowzeropivot, zeropivotdetected; 1286bce7ff8SHong Zhang 1296bce7ff8SHong Zhang PetscFunctionBegin; 1309566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isrow, &r)); 1319566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isicol, &ic)); 1325f8bbccaSHong Zhang allowzeropivot = PetscNot(A->erroriffailure); 133ae3d28f0SHong Zhang 1349566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(bs2 * n, &rtmp)); 1356bce7ff8SHong Zhang 136914a18a2SHong Zhang /* generate work space needed by dense LU factorization */ 1379566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(bs, &v_work, bs2, &mwork, bs, &v_pivots)); 138914a18a2SHong Zhang 1396bce7ff8SHong Zhang for (i = 0; i < n; i++) { 1406bce7ff8SHong Zhang /* zero rtmp */ 1416bce7ff8SHong Zhang /* L part */ 1426bce7ff8SHong Zhang nz = bi[i + 1] - bi[i]; 1436bce7ff8SHong Zhang bjtmp = bj + bi[i]; 14448a46eb9SPierre Jolivet for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2)); 1456bce7ff8SHong Zhang 1466bce7ff8SHong Zhang /* U part */ 1471a83e813SShri Abhyankar nz = bdiag[i] - bdiag[i + 1]; 1481a83e813SShri Abhyankar bjtmp = bj + bdiag[i + 1] + 1; 14948a46eb9SPierre Jolivet for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2)); 1501a83e813SShri Abhyankar 1511a83e813SShri Abhyankar /* load in initial (unfactored row) */ 1521a83e813SShri Abhyankar nz = ai[r[i] + 1] - ai[r[i]]; 1531a83e813SShri Abhyankar ajtmp = aj + ai[r[i]]; 1541a83e813SShri Abhyankar v = aa + bs2 * ai[r[i]]; 15548a46eb9SPierre Jolivet for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(rtmp + bs2 * ic[ajtmp[j]], v + bs2 * j, bs2)); 1561a83e813SShri Abhyankar 1571a83e813SShri Abhyankar /* elimination */ 1581a83e813SShri Abhyankar bjtmp = bj + bi[i]; 1591a83e813SShri Abhyankar nzL = bi[i + 1] - bi[i]; 1601a83e813SShri Abhyankar for (k = 0; k < nzL; k++) { 1611a83e813SShri Abhyankar row = bjtmp[k]; 1621a83e813SShri Abhyankar pc = rtmp + bs2 * row; 163c35f09e5SBarry Smith for (flg = 0, j = 0; j < bs2; j++) { 164c35f09e5SBarry Smith if (pc[j] != 0.0) { 165c35f09e5SBarry Smith flg = 1; 166c35f09e5SBarry Smith break; 167c35f09e5SBarry Smith } 168c35f09e5SBarry Smith } 1691a83e813SShri Abhyankar if (flg) { 1701a83e813SShri Abhyankar pv = b->a + bs2 * bdiag[row]; 17196b95a6bSBarry Smith PetscKernel_A_gets_A_times_B(bs, pc, pv, mwork); /* *pc = *pc * (*pv); */ 172a5b23f4aSJose E. Roman pj = b->j + bdiag[row + 1] + 1; /* beginning of U(row,:) */ 1731a83e813SShri Abhyankar pv = b->a + bs2 * (bdiag[row + 1] + 1); 1741a83e813SShri Abhyankar nz = bdiag[row] - bdiag[row + 1] - 1; /* num of entries inU(row,:), excluding diag */ 175ad540459SPierre Jolivet for (j = 0; j < nz; j++) PetscKernel_A_gets_A_minus_B_times_C(bs, rtmp + bs2 * pj[j], pc, pv + bs2 * j); 1769566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * bs2 * bs * (nz + 1) - bs2)); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */ 1771a83e813SShri Abhyankar } 1781a83e813SShri Abhyankar } 1791a83e813SShri Abhyankar 1801a83e813SShri Abhyankar /* finished row so stick it into b->a */ 1811a83e813SShri Abhyankar /* L part */ 1821a83e813SShri Abhyankar pv = b->a + bs2 * bi[i]; 1831a83e813SShri Abhyankar pj = b->j + bi[i]; 1841a83e813SShri Abhyankar nz = bi[i + 1] - bi[i]; 18548a46eb9SPierre Jolivet for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2)); 1861a83e813SShri Abhyankar 187a5b23f4aSJose E. Roman /* Mark diagonal and invert diagonal for simpler triangular solves */ 1881a83e813SShri Abhyankar pv = b->a + bs2 * bdiag[i]; 1891a83e813SShri Abhyankar pj = b->j + bdiag[i]; 1909566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(pv, rtmp + bs2 * pj[0], bs2)); 1915f8bbccaSHong Zhang 1929566063dSJacob Faibussowitsch PetscCall(PetscKernel_A_gets_inverse_A(bs, pv, v_pivots, v_work, allowzeropivot, &zeropivotdetected)); 1937b6c816cSBarry Smith if (zeropivotdetected) B->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 1941a83e813SShri Abhyankar 1951a83e813SShri Abhyankar /* U part */ 1961a83e813SShri Abhyankar pv = b->a + bs2 * (bdiag[i + 1] + 1); 1971a83e813SShri Abhyankar pj = b->j + bdiag[i + 1] + 1; 1981a83e813SShri Abhyankar nz = bdiag[i] - bdiag[i + 1] - 1; 19948a46eb9SPierre Jolivet for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2)); 2001a83e813SShri Abhyankar } 2011a83e813SShri Abhyankar 2029566063dSJacob Faibussowitsch PetscCall(PetscFree(rtmp)); 2039566063dSJacob Faibussowitsch PetscCall(PetscFree3(v_work, mwork, v_pivots)); 2049566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isicol, &ic)); 2059566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrow, &r)); 2061a83e813SShri Abhyankar 2079566063dSJacob Faibussowitsch PetscCall(ISIdentity(isrow, &row_identity)); 2089566063dSJacob Faibussowitsch PetscCall(ISIdentity(isicol, &col_identity)); 20926fbe8dcSKarl Rupp 210ace3abfcSBarry Smith both_identity = (PetscBool)(row_identity && col_identity); 211ae3d28f0SHong Zhang if (both_identity) { 212ba7f0461SHong Zhang switch (bs) { 21396e086a2SDaniel Kokron case 9: 2145f70456aSHong 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) 21596e086a2SDaniel Kokron C->ops->solve = MatSolve_SeqBAIJ_9_NaturalOrdering; 216aee5c371SBarry Smith #else 217aee5c371SBarry Smith C->ops->solve = MatSolve_SeqBAIJ_N_NaturalOrdering; 218aee5c371SBarry Smith #endif 21996e086a2SDaniel Kokron break; 220d71ae5a4SJacob Faibussowitsch case 11: 221d71ae5a4SJacob Faibussowitsch C->ops->solve = MatSolve_SeqBAIJ_11_NaturalOrdering; 222d71ae5a4SJacob Faibussowitsch break; 223d71ae5a4SJacob Faibussowitsch case 12: 224d71ae5a4SJacob Faibussowitsch C->ops->solve = MatSolve_SeqBAIJ_12_NaturalOrdering; 225d71ae5a4SJacob Faibussowitsch break; 226d71ae5a4SJacob Faibussowitsch case 13: 227d71ae5a4SJacob Faibussowitsch C->ops->solve = MatSolve_SeqBAIJ_13_NaturalOrdering; 228d71ae5a4SJacob Faibussowitsch break; 229d71ae5a4SJacob Faibussowitsch case 14: 230d71ae5a4SJacob Faibussowitsch C->ops->solve = MatSolve_SeqBAIJ_14_NaturalOrdering; 231d71ae5a4SJacob Faibussowitsch break; 232d71ae5a4SJacob Faibussowitsch default: 233d71ae5a4SJacob Faibussowitsch C->ops->solve = MatSolve_SeqBAIJ_N_NaturalOrdering; 234d71ae5a4SJacob Faibussowitsch break; 235ba7f0461SHong Zhang } 236ae3d28f0SHong Zhang } else { 2374dd39f65SShri Abhyankar C->ops->solve = MatSolve_SeqBAIJ_N; 238ae3d28f0SHong Zhang } 2394dd39f65SShri Abhyankar C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_N; 240ae3d28f0SHong Zhang 2411a83e813SShri Abhyankar C->assembled = PETSC_TRUE; 24226fbe8dcSKarl Rupp 2439566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(1.333333333333 * bs * bs2 * b->mbs)); /* from inverting diagonal blocks */ 2443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2451a83e813SShri Abhyankar } 2461a83e813SShri Abhyankar 2476bce7ff8SHong Zhang /* 2486bce7ff8SHong Zhang ilu(0) with natural ordering under new data structure. 2494dd39f65SShri Abhyankar See MatILUFactorSymbolic_SeqAIJ_ilu0() for detailed description 2504dd39f65SShri Abhyankar because this code is almost identical to MatILUFactorSymbolic_SeqAIJ_ilu0_inplace(). 2516bce7ff8SHong Zhang */ 252c0c7eb62SShri Abhyankar 253d71ae5a4SJacob Faibussowitsch PetscErrorCode MatILUFactorSymbolic_SeqBAIJ_ilu0(Mat fact, Mat A, IS isrow, IS iscol, const MatFactorInfo *info) 254d71ae5a4SJacob Faibussowitsch { 2556bce7ff8SHong Zhang Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b; 25616a2bf60SHong Zhang PetscInt n = a->mbs, *ai = a->i, *aj, *adiag = a->diag, bs2 = a->bs2; 25735aa4fcfSShri Abhyankar PetscInt i, j, nz, *bi, *bj, *bdiag, bi_temp; 25835aa4fcfSShri Abhyankar 25935aa4fcfSShri Abhyankar PetscFunctionBegin; 2609566063dSJacob Faibussowitsch PetscCall(MatDuplicateNoCreate_SeqBAIJ(fact, A, MAT_DO_NOT_COPY_VALUES, PETSC_FALSE)); 26135aa4fcfSShri Abhyankar b = (Mat_SeqBAIJ *)(fact)->data; 26235aa4fcfSShri Abhyankar 26335aa4fcfSShri Abhyankar /* allocate matrix arrays for new data structure */ 2649566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(bs2 * ai[n] + 1, &b->a, ai[n] + 1, &b->j, n + 1, &b->i)); 26526fbe8dcSKarl Rupp 26635aa4fcfSShri Abhyankar b->singlemalloc = PETSC_TRUE; 267379be0ddSLisandro Dalcin b->free_a = PETSC_TRUE; 268379be0ddSLisandro Dalcin b->free_ij = PETSC_TRUE; 2691e40a84eSLisandro Dalcin fact->preallocated = PETSC_TRUE; 2701e40a84eSLisandro Dalcin fact->assembled = PETSC_TRUE; 271aa624791SPierre Jolivet if (!b->diag) PetscCall(PetscMalloc1(n + 1, &b->diag)); 27235aa4fcfSShri Abhyankar bdiag = b->diag; 27335aa4fcfSShri Abhyankar 27448a46eb9SPierre Jolivet if (n > 0) PetscCall(PetscArrayzero(b->a, bs2 * ai[n])); 27535aa4fcfSShri Abhyankar 27635aa4fcfSShri Abhyankar /* set bi and bj with new data structure */ 27735aa4fcfSShri Abhyankar bi = b->i; 27835aa4fcfSShri Abhyankar bj = b->j; 27935aa4fcfSShri Abhyankar 28035aa4fcfSShri Abhyankar /* L part */ 28135aa4fcfSShri Abhyankar bi[0] = 0; 28235aa4fcfSShri Abhyankar for (i = 0; i < n; i++) { 28335aa4fcfSShri Abhyankar nz = adiag[i] - ai[i]; 28435aa4fcfSShri Abhyankar bi[i + 1] = bi[i] + nz; 28535aa4fcfSShri Abhyankar aj = a->j + ai[i]; 28635aa4fcfSShri Abhyankar for (j = 0; j < nz; j++) { 2879371c9d4SSatish Balay *bj = aj[j]; 2889371c9d4SSatish Balay bj++; 28935aa4fcfSShri Abhyankar } 29035aa4fcfSShri Abhyankar } 29135aa4fcfSShri Abhyankar 29235aa4fcfSShri Abhyankar /* U part */ 29335aa4fcfSShri Abhyankar bi_temp = bi[n]; 29435aa4fcfSShri Abhyankar bdiag[n] = bi[n] - 1; 29535aa4fcfSShri Abhyankar for (i = n - 1; i >= 0; i--) { 29635aa4fcfSShri Abhyankar nz = ai[i + 1] - adiag[i] - 1; 29735aa4fcfSShri Abhyankar bi_temp = bi_temp + nz + 1; 29835aa4fcfSShri Abhyankar aj = a->j + adiag[i] + 1; 29935aa4fcfSShri Abhyankar for (j = 0; j < nz; j++) { 3009371c9d4SSatish Balay *bj = aj[j]; 3019371c9d4SSatish Balay bj++; 30235aa4fcfSShri Abhyankar } 30335aa4fcfSShri Abhyankar /* diag[i] */ 3049371c9d4SSatish Balay *bj = i; 3059371c9d4SSatish Balay bj++; 30635aa4fcfSShri Abhyankar bdiag[i] = bi_temp - 1; 30735aa4fcfSShri Abhyankar } 3083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 30935aa4fcfSShri Abhyankar } 31035aa4fcfSShri Abhyankar 311d71ae5a4SJacob Faibussowitsch PetscErrorCode MatILUFactorSymbolic_SeqBAIJ(Mat fact, Mat A, IS isrow, IS iscol, const MatFactorInfo *info) 312d71ae5a4SJacob Faibussowitsch { 31316a2bf60SHong Zhang Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b; 31416a2bf60SHong Zhang IS isicol; 31516a2bf60SHong Zhang const PetscInt *r, *ic; 3167fa3a6a0SHong Zhang PetscInt n = a->mbs, *ai = a->i, *aj = a->j, d; 31716a2bf60SHong Zhang PetscInt *bi, *cols, nnz, *cols_lvl; 31816a2bf60SHong Zhang PetscInt *bdiag, prow, fm, nzbd, reallocs = 0, dcount = 0; 31916a2bf60SHong Zhang PetscInt i, levels, diagonal_fill; 320ace3abfcSBarry Smith PetscBool col_identity, row_identity, both_identity; 32116a2bf60SHong Zhang PetscReal f; 3220298fd71SBarry Smith PetscInt nlnk, *lnk, *lnk_lvl = NULL; 32316a2bf60SHong Zhang PetscBT lnkbt; 32416a2bf60SHong Zhang PetscInt nzi, *bj, **bj_ptr, **bjlvl_ptr; 3250298fd71SBarry Smith PetscFreeSpaceList free_space = NULL, current_space = NULL; 3260298fd71SBarry Smith PetscFreeSpaceList free_space_lvl = NULL, current_space_lvl = NULL; 327ace3abfcSBarry Smith PetscBool missing; 3287fa3a6a0SHong Zhang PetscInt bs = A->rmap->bs, bs2 = a->bs2; 32916a2bf60SHong Zhang 33016a2bf60SHong Zhang PetscFunctionBegin; 33108401ef6SPierre 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); 3326ba06ab7SHong Zhang if (bs > 1) { /* check shifttype */ 333*13bcc0bdSJacob Faibussowitsch PetscCheck(info->shifttype != (PetscReal)MAT_SHIFT_NONZERO && info->shifttype != (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only MAT_SHIFT_NONE and MAT_SHIFT_INBLOCKS are supported for BAIJ matrix"); 3346ba06ab7SHong Zhang } 3356ba06ab7SHong Zhang 3369566063dSJacob Faibussowitsch PetscCall(MatMissingDiagonal(A, &missing, &d)); 33728b400f6SJacob Faibussowitsch PetscCheck(!missing, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Matrix is missing diagonal entry %" PetscInt_FMT, d); 33816a2bf60SHong Zhang 33916a2bf60SHong Zhang f = info->fill; 34016a2bf60SHong Zhang levels = (PetscInt)info->levels; 34116a2bf60SHong Zhang diagonal_fill = (PetscInt)info->diagonal_fill; 34226fbe8dcSKarl Rupp 3439566063dSJacob Faibussowitsch PetscCall(ISInvertPermutation(iscol, PETSC_DECIDE, &isicol)); 34416a2bf60SHong Zhang 3459566063dSJacob Faibussowitsch PetscCall(ISIdentity(isrow, &row_identity)); 3469566063dSJacob Faibussowitsch PetscCall(ISIdentity(iscol, &col_identity)); 34726fbe8dcSKarl Rupp 348ace3abfcSBarry Smith both_identity = (PetscBool)(row_identity && col_identity); 34916a2bf60SHong Zhang 3507fa3a6a0SHong Zhang if (!levels && both_identity) { 35116a2bf60SHong Zhang /* special case: ilu(0) with natural ordering */ 3529566063dSJacob Faibussowitsch PetscCall(MatILUFactorSymbolic_SeqBAIJ_ilu0(fact, A, isrow, iscol, info)); 3539566063dSJacob Faibussowitsch PetscCall(MatSeqBAIJSetNumericFactorization(fact, both_identity)); 35435aa4fcfSShri Abhyankar 355d5f3da31SBarry Smith fact->factortype = MAT_FACTOR_ILU; 35635aa4fcfSShri Abhyankar (fact)->info.factor_mallocs = 0; 35735aa4fcfSShri Abhyankar (fact)->info.fill_ratio_given = info->fill; 35835aa4fcfSShri Abhyankar (fact)->info.fill_ratio_needed = 1.0; 35926fbe8dcSKarl Rupp 36035aa4fcfSShri Abhyankar b = (Mat_SeqBAIJ *)(fact)->data; 36135aa4fcfSShri Abhyankar b->row = isrow; 36235aa4fcfSShri Abhyankar b->col = iscol; 36335aa4fcfSShri Abhyankar b->icol = isicol; 3649566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)isrow)); 3659566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)iscol)); 36635aa4fcfSShri Abhyankar b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE; 36726fbe8dcSKarl Rupp 3689566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((n + 1) * bs, &b->solve_work)); 3693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 37035aa4fcfSShri Abhyankar } 37135aa4fcfSShri Abhyankar 3729566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isrow, &r)); 3739566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isicol, &ic)); 37435aa4fcfSShri Abhyankar 37535aa4fcfSShri Abhyankar /* get new row pointers */ 3769566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n + 1, &bi)); 37735aa4fcfSShri Abhyankar bi[0] = 0; 37835aa4fcfSShri Abhyankar /* bdiag is location of diagonal in factor */ 3799566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n + 1, &bdiag)); 38035aa4fcfSShri Abhyankar bdiag[0] = 0; 38135aa4fcfSShri Abhyankar 3829566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(n, &bj_ptr, n, &bjlvl_ptr)); 38335aa4fcfSShri Abhyankar 38435aa4fcfSShri Abhyankar /* create a linked list for storing column indices of the active row */ 38535aa4fcfSShri Abhyankar nlnk = n + 1; 3869566063dSJacob Faibussowitsch PetscCall(PetscIncompleteLLCreate(n, n, nlnk, lnk, lnk_lvl, lnkbt)); 38735aa4fcfSShri Abhyankar 38835aa4fcfSShri Abhyankar /* initial FreeSpace size is f*(ai[n]+1) */ 3899566063dSJacob Faibussowitsch PetscCall(PetscFreeSpaceGet(PetscRealIntMultTruncate(f, ai[n] + 1), &free_space)); 39035aa4fcfSShri Abhyankar current_space = free_space; 3919566063dSJacob Faibussowitsch PetscCall(PetscFreeSpaceGet(PetscRealIntMultTruncate(f, ai[n] + 1), &free_space_lvl)); 39235aa4fcfSShri Abhyankar current_space_lvl = free_space_lvl; 39335aa4fcfSShri Abhyankar 39435aa4fcfSShri Abhyankar for (i = 0; i < n; i++) { 39535aa4fcfSShri Abhyankar nzi = 0; 39635aa4fcfSShri Abhyankar /* copy current row into linked list */ 39735aa4fcfSShri Abhyankar nnz = ai[r[i] + 1] - ai[r[i]]; 39894bad497SJacob 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); 39935aa4fcfSShri Abhyankar cols = aj + ai[r[i]]; 40035aa4fcfSShri Abhyankar lnk[i] = -1; /* marker to indicate if diagonal exists */ 4019566063dSJacob Faibussowitsch PetscCall(PetscIncompleteLLInit(nnz, cols, n, ic, &nlnk, lnk, lnk_lvl, lnkbt)); 40235aa4fcfSShri Abhyankar nzi += nlnk; 40335aa4fcfSShri Abhyankar 40435aa4fcfSShri Abhyankar /* make sure diagonal entry is included */ 40535aa4fcfSShri Abhyankar if (diagonal_fill && lnk[i] == -1) { 40635aa4fcfSShri Abhyankar fm = n; 40735aa4fcfSShri Abhyankar while (lnk[fm] < i) fm = lnk[fm]; 40835aa4fcfSShri Abhyankar lnk[i] = lnk[fm]; /* insert diagonal into linked list */ 40935aa4fcfSShri Abhyankar lnk[fm] = i; 41035aa4fcfSShri Abhyankar lnk_lvl[i] = 0; 4119371c9d4SSatish Balay nzi++; 4129371c9d4SSatish Balay dcount++; 41335aa4fcfSShri Abhyankar } 41435aa4fcfSShri Abhyankar 41535aa4fcfSShri Abhyankar /* add pivot rows into the active row */ 41635aa4fcfSShri Abhyankar nzbd = 0; 41735aa4fcfSShri Abhyankar prow = lnk[n]; 41835aa4fcfSShri Abhyankar while (prow < i) { 41935aa4fcfSShri Abhyankar nnz = bdiag[prow]; 42035aa4fcfSShri Abhyankar cols = bj_ptr[prow] + nnz + 1; 42135aa4fcfSShri Abhyankar cols_lvl = bjlvl_ptr[prow] + nnz + 1; 42235aa4fcfSShri Abhyankar nnz = bi[prow + 1] - bi[prow] - nnz - 1; 42326fbe8dcSKarl Rupp 4249566063dSJacob Faibussowitsch PetscCall(PetscILULLAddSorted(nnz, cols, levels, cols_lvl, prow, &nlnk, lnk, lnk_lvl, lnkbt, prow)); 42535aa4fcfSShri Abhyankar nzi += nlnk; 42635aa4fcfSShri Abhyankar prow = lnk[prow]; 42735aa4fcfSShri Abhyankar nzbd++; 42835aa4fcfSShri Abhyankar } 42935aa4fcfSShri Abhyankar bdiag[i] = nzbd; 43035aa4fcfSShri Abhyankar bi[i + 1] = bi[i] + nzi; 43135aa4fcfSShri Abhyankar 43235aa4fcfSShri Abhyankar /* if free space is not available, make more free space */ 43335aa4fcfSShri Abhyankar if (current_space->local_remaining < nzi) { 434f91af8c7SBarry Smith nnz = PetscIntMultTruncate(2, PetscIntMultTruncate(nzi, (n - i))); /* estimated and max additional space needed */ 4359566063dSJacob Faibussowitsch PetscCall(PetscFreeSpaceGet(nnz, ¤t_space)); 4369566063dSJacob Faibussowitsch PetscCall(PetscFreeSpaceGet(nnz, ¤t_space_lvl)); 43735aa4fcfSShri Abhyankar reallocs++; 43835aa4fcfSShri Abhyankar } 43935aa4fcfSShri Abhyankar 44035aa4fcfSShri Abhyankar /* copy data into free_space and free_space_lvl, then initialize lnk */ 4419566063dSJacob Faibussowitsch PetscCall(PetscIncompleteLLClean(n, n, nzi, lnk, lnk_lvl, current_space->array, current_space_lvl->array, lnkbt)); 44226fbe8dcSKarl Rupp 44335aa4fcfSShri Abhyankar bj_ptr[i] = current_space->array; 44435aa4fcfSShri Abhyankar bjlvl_ptr[i] = current_space_lvl->array; 44535aa4fcfSShri Abhyankar 44635aa4fcfSShri Abhyankar /* make sure the active row i has diagonal entry */ 44794bad497SJacob 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); 44835aa4fcfSShri Abhyankar 44935aa4fcfSShri Abhyankar current_space->array += nzi; 45035aa4fcfSShri Abhyankar current_space->local_used += nzi; 45135aa4fcfSShri Abhyankar current_space->local_remaining -= nzi; 45226fbe8dcSKarl Rupp 45335aa4fcfSShri Abhyankar current_space_lvl->array += nzi; 45435aa4fcfSShri Abhyankar current_space_lvl->local_used += nzi; 45535aa4fcfSShri Abhyankar current_space_lvl->local_remaining -= nzi; 45635aa4fcfSShri Abhyankar } 45735aa4fcfSShri Abhyankar 4589566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrow, &r)); 4599566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isicol, &ic)); 46035aa4fcfSShri Abhyankar 46135aa4fcfSShri Abhyankar /* copy free_space into bj and free free_space; set bi, bj, bdiag in new datastructure; */ 4629566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bi[n] + 1, &bj)); 4639566063dSJacob Faibussowitsch PetscCall(PetscFreeSpaceContiguous_LU(&free_space, bj, n, bi, bdiag)); 46435aa4fcfSShri Abhyankar 4659566063dSJacob Faibussowitsch PetscCall(PetscIncompleteLLDestroy(lnk, lnkbt)); 4669566063dSJacob Faibussowitsch PetscCall(PetscFreeSpaceDestroy(free_space_lvl)); 4679566063dSJacob Faibussowitsch PetscCall(PetscFree2(bj_ptr, bjlvl_ptr)); 46835aa4fcfSShri Abhyankar 46935aa4fcfSShri Abhyankar #if defined(PETSC_USE_INFO) 47035aa4fcfSShri Abhyankar { 471aef85c9fSShri Abhyankar PetscReal af = ((PetscReal)(bdiag[0] + 1)) / ((PetscReal)ai[n]); 4729566063dSJacob Faibussowitsch PetscCall(PetscInfo(A, "Reallocs %" PetscInt_FMT " Fill ratio:given %g needed %g\n", reallocs, (double)f, (double)af)); 4739566063dSJacob Faibussowitsch PetscCall(PetscInfo(A, "Run with -[sub_]pc_factor_fill %g or use \n", (double)af)); 4749566063dSJacob Faibussowitsch PetscCall(PetscInfo(A, "PCFactorSetFill([sub]pc,%g);\n", (double)af)); 4759566063dSJacob Faibussowitsch PetscCall(PetscInfo(A, "for best performance.\n")); 47648a46eb9SPierre Jolivet if (diagonal_fill) PetscCall(PetscInfo(A, "Detected and replaced %" PetscInt_FMT " missing diagonals\n", dcount)); 47735aa4fcfSShri Abhyankar } 47835aa4fcfSShri Abhyankar #endif 47935aa4fcfSShri Abhyankar 48035aa4fcfSShri Abhyankar /* put together the new matrix */ 4819566063dSJacob Faibussowitsch PetscCall(MatSeqBAIJSetPreallocation(fact, bs, MAT_SKIP_ALLOCATION, NULL)); 48226fbe8dcSKarl Rupp 48335aa4fcfSShri Abhyankar b = (Mat_SeqBAIJ *)(fact)->data; 48435aa4fcfSShri Abhyankar b->free_a = PETSC_TRUE; 48535aa4fcfSShri Abhyankar b->free_ij = PETSC_TRUE; 48635aa4fcfSShri Abhyankar b->singlemalloc = PETSC_FALSE; 48726fbe8dcSKarl Rupp 4889566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs2 * (bdiag[0] + 1), &b->a)); 48926fbe8dcSKarl Rupp 49035aa4fcfSShri Abhyankar b->j = bj; 49135aa4fcfSShri Abhyankar b->i = bi; 49235aa4fcfSShri Abhyankar b->diag = bdiag; 49335aa4fcfSShri Abhyankar b->free_diag = PETSC_TRUE; 494f4259b30SLisandro Dalcin b->ilen = NULL; 495f4259b30SLisandro Dalcin b->imax = NULL; 49635aa4fcfSShri Abhyankar b->row = isrow; 49735aa4fcfSShri Abhyankar b->col = iscol; 4989566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)isrow)); 4999566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)iscol)); 50035aa4fcfSShri Abhyankar b->icol = isicol; 50126fbe8dcSKarl Rupp 5029566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs * n + bs, &b->solve_work)); 50335aa4fcfSShri Abhyankar /* In b structure: Free imax, ilen, old a, old j. 50435aa4fcfSShri Abhyankar Allocate bdiag, solve_work, new a, new j */ 50535aa4fcfSShri Abhyankar b->maxnz = b->nz = bdiag[0] + 1; 50626fbe8dcSKarl Rupp 507ae3d28f0SHong Zhang fact->info.factor_mallocs = reallocs; 508ae3d28f0SHong Zhang fact->info.fill_ratio_given = f; 509ae3d28f0SHong Zhang fact->info.fill_ratio_needed = ((PetscReal)(bdiag[0] + 1)) / ((PetscReal)ai[n]); 51026fbe8dcSKarl Rupp 5119566063dSJacob Faibussowitsch PetscCall(MatSeqBAIJSetNumericFactorization(fact, both_identity)); 5123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 51335aa4fcfSShri Abhyankar } 51435aa4fcfSShri Abhyankar 515ff6a9541SJacob Faibussowitsch #if 0 516ff6a9541SJacob Faibussowitsch // unused 5174e2b4712SSatish Balay /* 5184e2b4712SSatish Balay This code is virtually identical to MatILUFactorSymbolic_SeqAIJ 5194e2b4712SSatish Balay except that the data structure of Mat_SeqAIJ is slightly different. 5204e2b4712SSatish Balay Not a good example of code reuse. 5214e2b4712SSatish Balay */ 522ff6a9541SJacob Faibussowitsch static PetscErrorCode MatILUFactorSymbolic_SeqBAIJ_inplace(Mat fact, Mat A, IS isrow, IS iscol, const MatFactorInfo *info) 523d71ae5a4SJacob Faibussowitsch { 5244e2b4712SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b; 5254e2b4712SSatish Balay IS isicol; 5265d0c19d7SBarry Smith const PetscInt *r, *ic, *ai = a->i, *aj = a->j, *xi; 5275d0c19d7SBarry Smith PetscInt prow, n = a->mbs, *ainew, *ajnew, jmax, *fill, nz, *im, *ajfill, *flev, *xitmp; 528a96a251dSBarry Smith PetscInt *dloc, idx, row, m, fm, nzf, nzi, reallocate = 0, dcount = 0; 529d0f46423SBarry Smith PetscInt incrlev, nnz, i, bs = A->rmap->bs, bs2 = a->bs2, levels, diagonal_fill, dd; 530ace3abfcSBarry Smith PetscBool col_identity, row_identity, both_identity, flg; 531329f5518SBarry Smith PetscReal f; 5324e2b4712SSatish Balay 5334e2b4712SSatish Balay PetscFunctionBegin; 5349566063dSJacob Faibussowitsch PetscCall(MatMissingDiagonal_SeqBAIJ(A, &flg, &dd)); 53528b400f6SJacob Faibussowitsch PetscCheck(!flg, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Matrix A is missing diagonal entry in row %" PetscInt_FMT, dd); 5366bce7ff8SHong Zhang 537435faa5fSBarry Smith f = info->fill; 538690b6cddSBarry Smith levels = (PetscInt)info->levels; 539690b6cddSBarry Smith diagonal_fill = (PetscInt)info->diagonal_fill; 54026fbe8dcSKarl Rupp 5419566063dSJacob Faibussowitsch PetscCall(ISInvertPermutation(iscol, PETSC_DECIDE, &isicol)); 54216a2bf60SHong Zhang 5439566063dSJacob Faibussowitsch PetscCall(ISIdentity(isrow, &row_identity)); 5449566063dSJacob Faibussowitsch PetscCall(ISIdentity(iscol, &col_identity)); 545ace3abfcSBarry Smith both_identity = (PetscBool)(row_identity && col_identity); 546309c388cSBarry Smith 54741df41f0SMatthew Knepley if (!levels && both_identity) { /* special case copy the nonzero structure */ 5489566063dSJacob Faibussowitsch PetscCall(MatDuplicateNoCreate_SeqBAIJ(fact, A, MAT_DO_NOT_COPY_VALUES, PETSC_TRUE)); 5499566063dSJacob Faibussowitsch PetscCall(MatSeqBAIJSetNumericFactorization_inplace(fact, both_identity)); 5506bce7ff8SHong Zhang 551d5f3da31SBarry Smith fact->factortype = MAT_FACTOR_ILU; 552ae3d28f0SHong Zhang b = (Mat_SeqBAIJ *)fact->data; 553bb3d539aSBarry Smith b->row = isrow; 554bb3d539aSBarry Smith b->col = iscol; 5559566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)isrow)); 5569566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)iscol)); 557bb3d539aSBarry Smith b->icol = isicol; 558bcd9e38bSBarry Smith b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE; 55926fbe8dcSKarl Rupp 5609566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((n + 1) * bs, &b->solve_work)); 5613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5626bce7ff8SHong Zhang } 5636bce7ff8SHong Zhang 5646bce7ff8SHong Zhang /* general case perform the symbolic factorization */ 5659566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isrow, &r)); 5669566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isicol, &ic)); 5674e2b4712SSatish Balay 5684e2b4712SSatish Balay /* get new row pointers */ 5699566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n + 1, &ainew)); 5704e2b4712SSatish Balay ainew[0] = 0; 5714e2b4712SSatish Balay /* don't know how many column pointers are needed so estimate */ 572690b6cddSBarry Smith jmax = (PetscInt)(f * ai[n] + 1); 5739566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jmax, &ajnew)); 5744e2b4712SSatish Balay /* ajfill is level of fill for each fill entry */ 5759566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jmax, &ajfill)); 5764e2b4712SSatish Balay /* fill is a linked list of nonzeros in active row */ 5779566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n + 1, &fill)); 5784e2b4712SSatish Balay /* im is level for each filled value */ 5799566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n + 1, &im)); 5804e2b4712SSatish Balay /* dloc is location of diagonal in factor */ 5819566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n + 1, &dloc)); 5824e2b4712SSatish Balay dloc[0] = 0; 5834e2b4712SSatish Balay for (prow = 0; prow < n; prow++) { 584435faa5fSBarry Smith /* copy prow into linked list */ 5854e2b4712SSatish Balay nzf = nz = ai[r[prow] + 1] - ai[r[prow]]; 58628b400f6SJacob 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); 5874e2b4712SSatish Balay xi = aj + ai[r[prow]]; 5884e2b4712SSatish Balay fill[n] = n; 589435faa5fSBarry Smith fill[prow] = -1; /* marker for diagonal entry */ 5904e2b4712SSatish Balay while (nz--) { 5914e2b4712SSatish Balay fm = n; 5924e2b4712SSatish Balay idx = ic[*xi++]; 5934e2b4712SSatish Balay do { 5944e2b4712SSatish Balay m = fm; 5954e2b4712SSatish Balay fm = fill[m]; 5964e2b4712SSatish Balay } while (fm < idx); 5974e2b4712SSatish Balay fill[m] = idx; 5984e2b4712SSatish Balay fill[idx] = fm; 5994e2b4712SSatish Balay im[idx] = 0; 6004e2b4712SSatish Balay } 601435faa5fSBarry Smith 602435faa5fSBarry Smith /* make sure diagonal entry is included */ 603435faa5fSBarry Smith if (diagonal_fill && fill[prow] == -1) { 604435faa5fSBarry Smith fm = n; 605435faa5fSBarry Smith while (fill[fm] < prow) fm = fill[fm]; 606435faa5fSBarry Smith fill[prow] = fill[fm]; /* insert diagonal into linked list */ 607435faa5fSBarry Smith fill[fm] = prow; 608435faa5fSBarry Smith im[prow] = 0; 609435faa5fSBarry Smith nzf++; 610335d9088SBarry Smith dcount++; 611435faa5fSBarry Smith } 612435faa5fSBarry Smith 6134e2b4712SSatish Balay nzi = 0; 6144e2b4712SSatish Balay row = fill[n]; 6154e2b4712SSatish Balay while (row < prow) { 6164e2b4712SSatish Balay incrlev = im[row] + 1; 6174e2b4712SSatish Balay nz = dloc[row]; 618435faa5fSBarry Smith xi = ajnew + ainew[row] + nz + 1; 6194e2b4712SSatish Balay flev = ajfill + ainew[row] + nz + 1; 6204e2b4712SSatish Balay nnz = ainew[row + 1] - ainew[row] - nz - 1; 6214e2b4712SSatish Balay fm = row; 6224e2b4712SSatish Balay while (nnz-- > 0) { 6234e2b4712SSatish Balay idx = *xi++; 6244e2b4712SSatish Balay if (*flev + incrlev > levels) { 6254e2b4712SSatish Balay flev++; 6264e2b4712SSatish Balay continue; 6274e2b4712SSatish Balay } 6284e2b4712SSatish Balay do { 6294e2b4712SSatish Balay m = fm; 6304e2b4712SSatish Balay fm = fill[m]; 6314e2b4712SSatish Balay } while (fm < idx); 6324e2b4712SSatish Balay if (fm != idx) { 6334e2b4712SSatish Balay im[idx] = *flev + incrlev; 6344e2b4712SSatish Balay fill[m] = idx; 6354e2b4712SSatish Balay fill[idx] = fm; 6364e2b4712SSatish Balay fm = idx; 6374e2b4712SSatish Balay nzf++; 63826fbe8dcSKarl Rupp } else if (im[idx] > *flev + incrlev) im[idx] = *flev + incrlev; 6394e2b4712SSatish Balay flev++; 6404e2b4712SSatish Balay } 6414e2b4712SSatish Balay row = fill[row]; 6424e2b4712SSatish Balay nzi++; 6434e2b4712SSatish Balay } 6444e2b4712SSatish Balay /* copy new filled row into permanent storage */ 6454e2b4712SSatish Balay ainew[prow + 1] = ainew[prow] + nzf; 6464e2b4712SSatish Balay if (ainew[prow + 1] > jmax) { 647ecf371e4SBarry Smith /* estimate how much additional space we will need */ 648ecf371e4SBarry Smith /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */ 649ecf371e4SBarry Smith /* just double the memory each time */ 650690b6cddSBarry Smith PetscInt maxadd = jmax; 651ecf371e4SBarry Smith /* maxadd = (int)(((f*ai[n]+1)*(n-prow+5))/n); */ 6524e2b4712SSatish Balay if (maxadd < nzf) maxadd = (n - prow) * (nzf + 1); 6534e2b4712SSatish Balay jmax += maxadd; 654ecf371e4SBarry Smith 655ecf371e4SBarry Smith /* allocate a longer ajnew and ajfill */ 6569566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jmax, &xitmp)); 6579566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(xitmp, ajnew, ainew[prow])); 6589566063dSJacob Faibussowitsch PetscCall(PetscFree(ajnew)); 6595d0c19d7SBarry Smith ajnew = xitmp; 6609566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jmax, &xitmp)); 6619566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(xitmp, ajfill, ainew[prow])); 6629566063dSJacob Faibussowitsch PetscCall(PetscFree(ajfill)); 6635d0c19d7SBarry Smith ajfill = xitmp; 664eb150c5cSKris Buschelman reallocate++; /* count how many reallocations are needed */ 6654e2b4712SSatish Balay } 6665d0c19d7SBarry Smith xitmp = ajnew + ainew[prow]; 6674e2b4712SSatish Balay flev = ajfill + ainew[prow]; 6684e2b4712SSatish Balay dloc[prow] = nzi; 6694e2b4712SSatish Balay fm = fill[n]; 6704e2b4712SSatish Balay while (nzf--) { 6715d0c19d7SBarry Smith *xitmp++ = fm; 6724e2b4712SSatish Balay *flev++ = im[fm]; 6734e2b4712SSatish Balay fm = fill[fm]; 6744e2b4712SSatish Balay } 675435faa5fSBarry Smith /* make sure row has diagonal entry */ 676aed4548fSBarry 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\ 6779371c9d4SSatish Balay try running with -pc_factor_nonzeros_along_diagonal or -pc_factor_diagonal_fill", 6789371c9d4SSatish Balay prow); 679435faa5fSBarry Smith } 6809566063dSJacob Faibussowitsch PetscCall(PetscFree(ajfill)); 6819566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrow, &r)); 6829566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isicol, &ic)); 6839566063dSJacob Faibussowitsch PetscCall(PetscFree(fill)); 6849566063dSJacob Faibussowitsch PetscCall(PetscFree(im)); 6854e2b4712SSatish Balay 6866cf91177SBarry Smith #if defined(PETSC_USE_INFO) 6874e2b4712SSatish Balay { 688329f5518SBarry Smith PetscReal af = ((PetscReal)ainew[n]) / ((PetscReal)ai[n]); 6899566063dSJacob Faibussowitsch PetscCall(PetscInfo(A, "Reallocs %" PetscInt_FMT " Fill ratio:given %g needed %g\n", reallocate, (double)f, (double)af)); 6909566063dSJacob Faibussowitsch PetscCall(PetscInfo(A, "Run with -pc_factor_fill %g or use \n", (double)af)); 6919566063dSJacob Faibussowitsch PetscCall(PetscInfo(A, "PCFactorSetFill(pc,%g);\n", (double)af)); 6929566063dSJacob Faibussowitsch PetscCall(PetscInfo(A, "for best performance.\n")); 69348a46eb9SPierre Jolivet if (diagonal_fill) PetscCall(PetscInfo(A, "Detected and replaced %" PetscInt_FMT " missing diagonals\n", dcount)); 6944e2b4712SSatish Balay } 69563ba0a88SBarry Smith #endif 6964e2b4712SSatish Balay 6974e2b4712SSatish Balay /* put together the new matrix */ 6989566063dSJacob Faibussowitsch PetscCall(MatSeqBAIJSetPreallocation(fact, bs, MAT_SKIP_ALLOCATION, NULL)); 699ae3d28f0SHong Zhang b = (Mat_SeqBAIJ *)fact->data; 70026fbe8dcSKarl Rupp 701e6b907acSBarry Smith b->free_a = PETSC_TRUE; 702e6b907acSBarry Smith b->free_ij = PETSC_TRUE; 7037c922b88SBarry Smith b->singlemalloc = PETSC_FALSE; 70426fbe8dcSKarl Rupp 7059566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs2 * ainew[n], &b->a)); 70626fbe8dcSKarl Rupp 7074e2b4712SSatish Balay b->j = ajnew; 7084e2b4712SSatish Balay b->i = ainew; 7094e2b4712SSatish Balay for (i = 0; i < n; i++) dloc[i] += ainew[i]; 7104e2b4712SSatish Balay b->diag = dloc; 7117f53bb6cSHong Zhang b->free_diag = PETSC_TRUE; 712f4259b30SLisandro Dalcin b->ilen = NULL; 713f4259b30SLisandro Dalcin b->imax = NULL; 7144e2b4712SSatish Balay b->row = isrow; 7154e2b4712SSatish Balay b->col = iscol; 716bcd9e38bSBarry Smith b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE; 71726fbe8dcSKarl Rupp 7189566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)isrow)); 7199566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)iscol)); 720e51c0b9cSSatish Balay b->icol = isicol; 7219566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs * n + bs, &b->solve_work)); 7224e2b4712SSatish Balay /* In b structure: Free imax, ilen, old a, old j. 7234e2b4712SSatish Balay Allocate dloc, solve_work, new a, new j */ 7244e2b4712SSatish Balay b->maxnz = b->nz = ainew[n]; 7254e2b4712SSatish Balay 726ae3d28f0SHong Zhang fact->info.factor_mallocs = reallocate; 727ae3d28f0SHong Zhang fact->info.fill_ratio_given = f; 728ae3d28f0SHong Zhang fact->info.fill_ratio_needed = ((PetscReal)ainew[n]) / ((PetscReal)ai[prow]); 7296bce7ff8SHong Zhang 7309566063dSJacob Faibussowitsch PetscCall(MatSeqBAIJSetNumericFactorization_inplace(fact, both_identity)); 7313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7328661488fSKris Buschelman } 733ff6a9541SJacob Faibussowitsch #endif 7348661488fSKris Buschelman 735d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE(Mat A) 736d71ae5a4SJacob Faibussowitsch { 73712272027SHong Zhang /* Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; */ 73812272027SHong Zhang /* int i,*AJ=a->j,nz=a->nz; */ 7395fd66863SKarl Rupp 7405a9542e3SKris Buschelman PetscFunctionBegin; 7417cf1b8d3SKris Buschelman /* Undo Column scaling */ 7427cf1b8d3SKris Buschelman /* while (nz--) { */ 7437cf1b8d3SKris Buschelman /* AJ[i] = AJ[i]/4; */ 7447cf1b8d3SKris Buschelman /* } */ 745c115a38dSKris Buschelman /* This should really invoke a push/pop logic, but we don't have that yet. */ 7460298fd71SBarry Smith A->ops->setunfactored = NULL; 7473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7487cf1b8d3SKris Buschelman } 7497cf1b8d3SKris Buschelman 750d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE_usj(Mat A) 751d71ae5a4SJacob Faibussowitsch { 7527cf1b8d3SKris Buschelman Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data; 753b24ad042SBarry Smith PetscInt *AJ = a->j, nz = a->nz; 7542aa5897fSKris Buschelman unsigned short *aj = (unsigned short *)AJ; 7555fd66863SKarl Rupp 7565a9542e3SKris Buschelman PetscFunctionBegin; 7570b9da03eSKris Buschelman /* Is this really necessary? */ 7589371c9d4SSatish Balay while (nz--) { AJ[nz] = (int)((unsigned int)aj[nz]); /* First extend, then convert to signed. */ } 7590298fd71SBarry Smith A->ops->setunfactored = NULL; 7603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7617e7071cdSKris Buschelman } 762