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 */ 179371c9d4SSatish Balay PetscErrorCode MatLUFactorNumeric_SeqBAIJ_15_NaturalOrdering(Mat B, Mat A, const MatFactorInfo *info) { 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]; 42*48a46eb9SPierre 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; 47*48a46eb9SPierre 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]; 53*48a46eb9SPierre 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]; 89*48a46eb9SPierre 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; 102*48a46eb9SPierre 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 */ 1122b0b2ea7SShri Abhyankar PetscFunctionReturn(0); 1132b0b2ea7SShri Abhyankar } 1142b0b2ea7SShri Abhyankar 1159371c9d4SSatish Balay PetscErrorCode MatLUFactorNumeric_SeqBAIJ_N(Mat B, Mat A, const MatFactorInfo *info) { 1166bce7ff8SHong Zhang Mat C = B; 1176bce7ff8SHong Zhang Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data; 1186bce7ff8SHong Zhang IS isrow = b->row, isicol = b->icol; 1195a586d82SBarry Smith const PetscInt *r, *ic; 1206bce7ff8SHong Zhang PetscInt i, j, k, n = a->mbs, *ai = a->i, *aj = a->j, *bi = b->i, *bj = b->j; 1216bce7ff8SHong Zhang PetscInt *ajtmp, *bjtmp, nz, nzL, row, *bdiag = b->diag, *pj; 122b588c5a2SHong Zhang MatScalar *rtmp, *pc, *mwork, *v, *pv, *aa = a->a; 123914a18a2SHong Zhang PetscInt bs = A->rmap->bs, bs2 = a->bs2, *v_pivots, flg; 124914a18a2SHong Zhang MatScalar *v_work; 125ace3abfcSBarry Smith PetscBool col_identity, row_identity, both_identity; 1265f8bbccaSHong Zhang PetscBool allowzeropivot, zeropivotdetected; 1276bce7ff8SHong Zhang 1286bce7ff8SHong Zhang PetscFunctionBegin; 1299566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isrow, &r)); 1309566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isicol, &ic)); 1315f8bbccaSHong Zhang allowzeropivot = PetscNot(A->erroriffailure); 132ae3d28f0SHong Zhang 1339566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(bs2 * n, &rtmp)); 1346bce7ff8SHong Zhang 135914a18a2SHong Zhang /* generate work space needed by dense LU factorization */ 1369566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(bs, &v_work, bs2, &mwork, bs, &v_pivots)); 137914a18a2SHong Zhang 1386bce7ff8SHong Zhang for (i = 0; i < n; i++) { 1396bce7ff8SHong Zhang /* zero rtmp */ 1406bce7ff8SHong Zhang /* L part */ 1416bce7ff8SHong Zhang nz = bi[i + 1] - bi[i]; 1426bce7ff8SHong Zhang bjtmp = bj + bi[i]; 143*48a46eb9SPierre Jolivet for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2)); 1446bce7ff8SHong Zhang 1456bce7ff8SHong Zhang /* U part */ 1461a83e813SShri Abhyankar nz = bdiag[i] - bdiag[i + 1]; 1471a83e813SShri Abhyankar bjtmp = bj + bdiag[i + 1] + 1; 148*48a46eb9SPierre Jolivet for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2)); 1491a83e813SShri Abhyankar 1501a83e813SShri Abhyankar /* load in initial (unfactored row) */ 1511a83e813SShri Abhyankar nz = ai[r[i] + 1] - ai[r[i]]; 1521a83e813SShri Abhyankar ajtmp = aj + ai[r[i]]; 1531a83e813SShri Abhyankar v = aa + bs2 * ai[r[i]]; 154*48a46eb9SPierre Jolivet for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(rtmp + bs2 * ic[ajtmp[j]], v + bs2 * j, bs2)); 1551a83e813SShri Abhyankar 1561a83e813SShri Abhyankar /* elimination */ 1571a83e813SShri Abhyankar bjtmp = bj + bi[i]; 1581a83e813SShri Abhyankar nzL = bi[i + 1] - bi[i]; 1591a83e813SShri Abhyankar for (k = 0; k < nzL; k++) { 1601a83e813SShri Abhyankar row = bjtmp[k]; 1611a83e813SShri Abhyankar pc = rtmp + bs2 * row; 162c35f09e5SBarry Smith for (flg = 0, j = 0; j < bs2; j++) { 163c35f09e5SBarry Smith if (pc[j] != 0.0) { 164c35f09e5SBarry Smith flg = 1; 165c35f09e5SBarry Smith break; 166c35f09e5SBarry Smith } 167c35f09e5SBarry Smith } 1681a83e813SShri Abhyankar if (flg) { 1691a83e813SShri Abhyankar pv = b->a + bs2 * bdiag[row]; 17096b95a6bSBarry Smith PetscKernel_A_gets_A_times_B(bs, pc, pv, mwork); /* *pc = *pc * (*pv); */ 171a5b23f4aSJose E. Roman pj = b->j + bdiag[row + 1] + 1; /* beginning of U(row,:) */ 1721a83e813SShri Abhyankar pv = b->a + bs2 * (bdiag[row + 1] + 1); 1731a83e813SShri Abhyankar nz = bdiag[row] - bdiag[row + 1] - 1; /* num of entries inU(row,:), excluding diag */ 1749371c9d4SSatish Balay for (j = 0; j < nz; j++) { PetscKernel_A_gets_A_minus_B_times_C(bs, rtmp + bs2 * pj[j], pc, pv + bs2 * j); } 1759566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * bs2 * bs * (nz + 1) - bs2)); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */ 1761a83e813SShri Abhyankar } 1771a83e813SShri Abhyankar } 1781a83e813SShri Abhyankar 1791a83e813SShri Abhyankar /* finished row so stick it into b->a */ 1801a83e813SShri Abhyankar /* L part */ 1811a83e813SShri Abhyankar pv = b->a + bs2 * bi[i]; 1821a83e813SShri Abhyankar pj = b->j + bi[i]; 1831a83e813SShri Abhyankar nz = bi[i + 1] - bi[i]; 184*48a46eb9SPierre Jolivet for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2)); 1851a83e813SShri Abhyankar 186a5b23f4aSJose E. Roman /* Mark diagonal and invert diagonal for simpler triangular solves */ 1871a83e813SShri Abhyankar pv = b->a + bs2 * bdiag[i]; 1881a83e813SShri Abhyankar pj = b->j + bdiag[i]; 1899566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(pv, rtmp + bs2 * pj[0], bs2)); 1905f8bbccaSHong Zhang 1919566063dSJacob Faibussowitsch PetscCall(PetscKernel_A_gets_inverse_A(bs, pv, v_pivots, v_work, allowzeropivot, &zeropivotdetected)); 1927b6c816cSBarry Smith if (zeropivotdetected) B->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 1931a83e813SShri Abhyankar 1941a83e813SShri Abhyankar /* U part */ 1951a83e813SShri Abhyankar pv = b->a + bs2 * (bdiag[i + 1] + 1); 1961a83e813SShri Abhyankar pj = b->j + bdiag[i + 1] + 1; 1971a83e813SShri Abhyankar nz = bdiag[i] - bdiag[i + 1] - 1; 198*48a46eb9SPierre Jolivet for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2)); 1991a83e813SShri Abhyankar } 2001a83e813SShri Abhyankar 2019566063dSJacob Faibussowitsch PetscCall(PetscFree(rtmp)); 2029566063dSJacob Faibussowitsch PetscCall(PetscFree3(v_work, mwork, v_pivots)); 2039566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isicol, &ic)); 2049566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrow, &r)); 2051a83e813SShri Abhyankar 2069566063dSJacob Faibussowitsch PetscCall(ISIdentity(isrow, &row_identity)); 2079566063dSJacob Faibussowitsch PetscCall(ISIdentity(isicol, &col_identity)); 20826fbe8dcSKarl Rupp 209ace3abfcSBarry Smith both_identity = (PetscBool)(row_identity && col_identity); 210ae3d28f0SHong Zhang if (both_identity) { 211ba7f0461SHong Zhang switch (bs) { 21296e086a2SDaniel Kokron case 9: 2135f70456aSHong 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) 21496e086a2SDaniel Kokron C->ops->solve = MatSolve_SeqBAIJ_9_NaturalOrdering; 215aee5c371SBarry Smith #else 216aee5c371SBarry Smith C->ops->solve = MatSolve_SeqBAIJ_N_NaturalOrdering; 217aee5c371SBarry Smith #endif 21896e086a2SDaniel Kokron break; 2199371c9d4SSatish Balay case 11: C->ops->solve = MatSolve_SeqBAIJ_11_NaturalOrdering; break; 2209371c9d4SSatish Balay case 12: C->ops->solve = MatSolve_SeqBAIJ_12_NaturalOrdering; break; 2219371c9d4SSatish Balay case 13: C->ops->solve = MatSolve_SeqBAIJ_13_NaturalOrdering; break; 2229371c9d4SSatish Balay case 14: C->ops->solve = MatSolve_SeqBAIJ_14_NaturalOrdering; break; 2239371c9d4SSatish Balay default: C->ops->solve = MatSolve_SeqBAIJ_N_NaturalOrdering; break; 224ba7f0461SHong Zhang } 225ae3d28f0SHong Zhang } else { 2264dd39f65SShri Abhyankar C->ops->solve = MatSolve_SeqBAIJ_N; 227ae3d28f0SHong Zhang } 2284dd39f65SShri Abhyankar C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_N; 229ae3d28f0SHong Zhang 2301a83e813SShri Abhyankar C->assembled = PETSC_TRUE; 23126fbe8dcSKarl Rupp 2329566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(1.333333333333 * bs * bs2 * b->mbs)); /* from inverting diagonal blocks */ 2331a83e813SShri Abhyankar PetscFunctionReturn(0); 2341a83e813SShri Abhyankar } 2351a83e813SShri Abhyankar 2366bce7ff8SHong Zhang /* 2376bce7ff8SHong Zhang ilu(0) with natural ordering under new data structure. 2384dd39f65SShri Abhyankar See MatILUFactorSymbolic_SeqAIJ_ilu0() for detailed description 2394dd39f65SShri Abhyankar because this code is almost identical to MatILUFactorSymbolic_SeqAIJ_ilu0_inplace(). 2406bce7ff8SHong Zhang */ 241c0c7eb62SShri Abhyankar 2429371c9d4SSatish Balay PetscErrorCode MatILUFactorSymbolic_SeqBAIJ_ilu0(Mat fact, Mat A, IS isrow, IS iscol, const MatFactorInfo *info) { 2436bce7ff8SHong Zhang Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b; 24416a2bf60SHong Zhang PetscInt n = a->mbs, *ai = a->i, *aj, *adiag = a->diag, bs2 = a->bs2; 24535aa4fcfSShri Abhyankar PetscInt i, j, nz, *bi, *bj, *bdiag, bi_temp; 24635aa4fcfSShri Abhyankar 24735aa4fcfSShri Abhyankar PetscFunctionBegin; 2489566063dSJacob Faibussowitsch PetscCall(MatDuplicateNoCreate_SeqBAIJ(fact, A, MAT_DO_NOT_COPY_VALUES, PETSC_FALSE)); 24935aa4fcfSShri Abhyankar b = (Mat_SeqBAIJ *)(fact)->data; 25035aa4fcfSShri Abhyankar 25135aa4fcfSShri Abhyankar /* allocate matrix arrays for new data structure */ 2529566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(bs2 * ai[n] + 1, &b->a, ai[n] + 1, &b->j, n + 1, &b->i)); 2539566063dSJacob Faibussowitsch PetscCall(PetscLogObjectMemory((PetscObject)fact, ai[n] * (bs2 * sizeof(PetscScalar) + sizeof(PetscInt)) + (n + 1) * sizeof(PetscInt))); 25426fbe8dcSKarl Rupp 25535aa4fcfSShri Abhyankar b->singlemalloc = PETSC_TRUE; 256379be0ddSLisandro Dalcin b->free_a = PETSC_TRUE; 257379be0ddSLisandro Dalcin b->free_ij = PETSC_TRUE; 2581e40a84eSLisandro Dalcin fact->preallocated = PETSC_TRUE; 2591e40a84eSLisandro Dalcin fact->assembled = PETSC_TRUE; 26035aa4fcfSShri Abhyankar if (!b->diag) { 2619566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n + 1, &b->diag)); 2629566063dSJacob Faibussowitsch PetscCall(PetscLogObjectMemory((PetscObject)fact, (n + 1) * sizeof(PetscInt))); 26335aa4fcfSShri Abhyankar } 26435aa4fcfSShri Abhyankar bdiag = b->diag; 26535aa4fcfSShri Abhyankar 266*48a46eb9SPierre Jolivet if (n > 0) PetscCall(PetscArrayzero(b->a, bs2 * ai[n])); 26735aa4fcfSShri Abhyankar 26835aa4fcfSShri Abhyankar /* set bi and bj with new data structure */ 26935aa4fcfSShri Abhyankar bi = b->i; 27035aa4fcfSShri Abhyankar bj = b->j; 27135aa4fcfSShri Abhyankar 27235aa4fcfSShri Abhyankar /* L part */ 27335aa4fcfSShri Abhyankar bi[0] = 0; 27435aa4fcfSShri Abhyankar for (i = 0; i < n; i++) { 27535aa4fcfSShri Abhyankar nz = adiag[i] - ai[i]; 27635aa4fcfSShri Abhyankar bi[i + 1] = bi[i] + nz; 27735aa4fcfSShri Abhyankar aj = a->j + ai[i]; 27835aa4fcfSShri Abhyankar for (j = 0; j < nz; j++) { 2799371c9d4SSatish Balay *bj = aj[j]; 2809371c9d4SSatish Balay bj++; 28135aa4fcfSShri Abhyankar } 28235aa4fcfSShri Abhyankar } 28335aa4fcfSShri Abhyankar 28435aa4fcfSShri Abhyankar /* U part */ 28535aa4fcfSShri Abhyankar bi_temp = bi[n]; 28635aa4fcfSShri Abhyankar bdiag[n] = bi[n] - 1; 28735aa4fcfSShri Abhyankar for (i = n - 1; i >= 0; i--) { 28835aa4fcfSShri Abhyankar nz = ai[i + 1] - adiag[i] - 1; 28935aa4fcfSShri Abhyankar bi_temp = bi_temp + nz + 1; 29035aa4fcfSShri Abhyankar aj = a->j + adiag[i] + 1; 29135aa4fcfSShri Abhyankar for (j = 0; j < nz; j++) { 2929371c9d4SSatish Balay *bj = aj[j]; 2939371c9d4SSatish Balay bj++; 29435aa4fcfSShri Abhyankar } 29535aa4fcfSShri Abhyankar /* diag[i] */ 2969371c9d4SSatish Balay *bj = i; 2979371c9d4SSatish Balay bj++; 29835aa4fcfSShri Abhyankar bdiag[i] = bi_temp - 1; 29935aa4fcfSShri Abhyankar } 30035aa4fcfSShri Abhyankar PetscFunctionReturn(0); 30135aa4fcfSShri Abhyankar } 30235aa4fcfSShri Abhyankar 3039371c9d4SSatish Balay PetscErrorCode MatILUFactorSymbolic_SeqBAIJ(Mat fact, Mat A, IS isrow, IS iscol, const MatFactorInfo *info) { 30416a2bf60SHong Zhang Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b; 30516a2bf60SHong Zhang IS isicol; 30616a2bf60SHong Zhang const PetscInt *r, *ic; 3077fa3a6a0SHong Zhang PetscInt n = a->mbs, *ai = a->i, *aj = a->j, d; 30816a2bf60SHong Zhang PetscInt *bi, *cols, nnz, *cols_lvl; 30916a2bf60SHong Zhang PetscInt *bdiag, prow, fm, nzbd, reallocs = 0, dcount = 0; 31016a2bf60SHong Zhang PetscInt i, levels, diagonal_fill; 311ace3abfcSBarry Smith PetscBool col_identity, row_identity, both_identity; 31216a2bf60SHong Zhang PetscReal f; 3130298fd71SBarry Smith PetscInt nlnk, *lnk, *lnk_lvl = NULL; 31416a2bf60SHong Zhang PetscBT lnkbt; 31516a2bf60SHong Zhang PetscInt nzi, *bj, **bj_ptr, **bjlvl_ptr; 3160298fd71SBarry Smith PetscFreeSpaceList free_space = NULL, current_space = NULL; 3170298fd71SBarry Smith PetscFreeSpaceList free_space_lvl = NULL, current_space_lvl = NULL; 318ace3abfcSBarry Smith PetscBool missing; 3197fa3a6a0SHong Zhang PetscInt bs = A->rmap->bs, bs2 = a->bs2; 32016a2bf60SHong Zhang 32116a2bf60SHong Zhang PetscFunctionBegin; 32208401ef6SPierre 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); 3236ba06ab7SHong Zhang if (bs > 1) { /* check shifttype */ 324aed4548fSBarry 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"); 3256ba06ab7SHong Zhang } 3266ba06ab7SHong Zhang 3279566063dSJacob Faibussowitsch PetscCall(MatMissingDiagonal(A, &missing, &d)); 32828b400f6SJacob Faibussowitsch PetscCheck(!missing, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Matrix is missing diagonal entry %" PetscInt_FMT, d); 32916a2bf60SHong Zhang 33016a2bf60SHong Zhang f = info->fill; 33116a2bf60SHong Zhang levels = (PetscInt)info->levels; 33216a2bf60SHong Zhang diagonal_fill = (PetscInt)info->diagonal_fill; 33326fbe8dcSKarl Rupp 3349566063dSJacob Faibussowitsch PetscCall(ISInvertPermutation(iscol, PETSC_DECIDE, &isicol)); 33516a2bf60SHong Zhang 3369566063dSJacob Faibussowitsch PetscCall(ISIdentity(isrow, &row_identity)); 3379566063dSJacob Faibussowitsch PetscCall(ISIdentity(iscol, &col_identity)); 33826fbe8dcSKarl Rupp 339ace3abfcSBarry Smith both_identity = (PetscBool)(row_identity && col_identity); 34016a2bf60SHong Zhang 3417fa3a6a0SHong Zhang if (!levels && both_identity) { 34216a2bf60SHong Zhang /* special case: ilu(0) with natural ordering */ 3439566063dSJacob Faibussowitsch PetscCall(MatILUFactorSymbolic_SeqBAIJ_ilu0(fact, A, isrow, iscol, info)); 3449566063dSJacob Faibussowitsch PetscCall(MatSeqBAIJSetNumericFactorization(fact, both_identity)); 34535aa4fcfSShri Abhyankar 346d5f3da31SBarry Smith fact->factortype = MAT_FACTOR_ILU; 34735aa4fcfSShri Abhyankar (fact)->info.factor_mallocs = 0; 34835aa4fcfSShri Abhyankar (fact)->info.fill_ratio_given = info->fill; 34935aa4fcfSShri Abhyankar (fact)->info.fill_ratio_needed = 1.0; 35026fbe8dcSKarl Rupp 35135aa4fcfSShri Abhyankar b = (Mat_SeqBAIJ *)(fact)->data; 35235aa4fcfSShri Abhyankar b->row = isrow; 35335aa4fcfSShri Abhyankar b->col = iscol; 35435aa4fcfSShri Abhyankar b->icol = isicol; 3559566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)isrow)); 3569566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)iscol)); 35735aa4fcfSShri Abhyankar b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE; 35826fbe8dcSKarl Rupp 3599566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((n + 1) * bs, &b->solve_work)); 36035aa4fcfSShri Abhyankar PetscFunctionReturn(0); 36135aa4fcfSShri Abhyankar } 36235aa4fcfSShri Abhyankar 3639566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isrow, &r)); 3649566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isicol, &ic)); 36535aa4fcfSShri Abhyankar 36635aa4fcfSShri Abhyankar /* get new row pointers */ 3679566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n + 1, &bi)); 36835aa4fcfSShri Abhyankar bi[0] = 0; 36935aa4fcfSShri Abhyankar /* bdiag is location of diagonal in factor */ 3709566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n + 1, &bdiag)); 37135aa4fcfSShri Abhyankar bdiag[0] = 0; 37235aa4fcfSShri Abhyankar 3739566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(n, &bj_ptr, n, &bjlvl_ptr)); 37435aa4fcfSShri Abhyankar 37535aa4fcfSShri Abhyankar /* create a linked list for storing column indices of the active row */ 37635aa4fcfSShri Abhyankar nlnk = n + 1; 3779566063dSJacob Faibussowitsch PetscCall(PetscIncompleteLLCreate(n, n, nlnk, lnk, lnk_lvl, lnkbt)); 37835aa4fcfSShri Abhyankar 37935aa4fcfSShri Abhyankar /* initial FreeSpace size is f*(ai[n]+1) */ 3809566063dSJacob Faibussowitsch PetscCall(PetscFreeSpaceGet(PetscRealIntMultTruncate(f, ai[n] + 1), &free_space)); 38135aa4fcfSShri Abhyankar current_space = free_space; 3829566063dSJacob Faibussowitsch PetscCall(PetscFreeSpaceGet(PetscRealIntMultTruncate(f, ai[n] + 1), &free_space_lvl)); 38335aa4fcfSShri Abhyankar current_space_lvl = free_space_lvl; 38435aa4fcfSShri Abhyankar 38535aa4fcfSShri Abhyankar for (i = 0; i < n; i++) { 38635aa4fcfSShri Abhyankar nzi = 0; 38735aa4fcfSShri Abhyankar /* copy current row into linked list */ 38835aa4fcfSShri Abhyankar nnz = ai[r[i] + 1] - ai[r[i]]; 38994bad497SJacob 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); 39035aa4fcfSShri Abhyankar cols = aj + ai[r[i]]; 39135aa4fcfSShri Abhyankar lnk[i] = -1; /* marker to indicate if diagonal exists */ 3929566063dSJacob Faibussowitsch PetscCall(PetscIncompleteLLInit(nnz, cols, n, ic, &nlnk, lnk, lnk_lvl, lnkbt)); 39335aa4fcfSShri Abhyankar nzi += nlnk; 39435aa4fcfSShri Abhyankar 39535aa4fcfSShri Abhyankar /* make sure diagonal entry is included */ 39635aa4fcfSShri Abhyankar if (diagonal_fill && lnk[i] == -1) { 39735aa4fcfSShri Abhyankar fm = n; 39835aa4fcfSShri Abhyankar while (lnk[fm] < i) fm = lnk[fm]; 39935aa4fcfSShri Abhyankar lnk[i] = lnk[fm]; /* insert diagonal into linked list */ 40035aa4fcfSShri Abhyankar lnk[fm] = i; 40135aa4fcfSShri Abhyankar lnk_lvl[i] = 0; 4029371c9d4SSatish Balay nzi++; 4039371c9d4SSatish Balay dcount++; 40435aa4fcfSShri Abhyankar } 40535aa4fcfSShri Abhyankar 40635aa4fcfSShri Abhyankar /* add pivot rows into the active row */ 40735aa4fcfSShri Abhyankar nzbd = 0; 40835aa4fcfSShri Abhyankar prow = lnk[n]; 40935aa4fcfSShri Abhyankar while (prow < i) { 41035aa4fcfSShri Abhyankar nnz = bdiag[prow]; 41135aa4fcfSShri Abhyankar cols = bj_ptr[prow] + nnz + 1; 41235aa4fcfSShri Abhyankar cols_lvl = bjlvl_ptr[prow] + nnz + 1; 41335aa4fcfSShri Abhyankar nnz = bi[prow + 1] - bi[prow] - nnz - 1; 41426fbe8dcSKarl Rupp 4159566063dSJacob Faibussowitsch PetscCall(PetscILULLAddSorted(nnz, cols, levels, cols_lvl, prow, &nlnk, lnk, lnk_lvl, lnkbt, prow)); 41635aa4fcfSShri Abhyankar nzi += nlnk; 41735aa4fcfSShri Abhyankar prow = lnk[prow]; 41835aa4fcfSShri Abhyankar nzbd++; 41935aa4fcfSShri Abhyankar } 42035aa4fcfSShri Abhyankar bdiag[i] = nzbd; 42135aa4fcfSShri Abhyankar bi[i + 1] = bi[i] + nzi; 42235aa4fcfSShri Abhyankar 42335aa4fcfSShri Abhyankar /* if free space is not available, make more free space */ 42435aa4fcfSShri Abhyankar if (current_space->local_remaining < nzi) { 425f91af8c7SBarry Smith nnz = PetscIntMultTruncate(2, PetscIntMultTruncate(nzi, (n - i))); /* estimated and max additional space needed */ 4269566063dSJacob Faibussowitsch PetscCall(PetscFreeSpaceGet(nnz, ¤t_space)); 4279566063dSJacob Faibussowitsch PetscCall(PetscFreeSpaceGet(nnz, ¤t_space_lvl)); 42835aa4fcfSShri Abhyankar reallocs++; 42935aa4fcfSShri Abhyankar } 43035aa4fcfSShri Abhyankar 43135aa4fcfSShri Abhyankar /* copy data into free_space and free_space_lvl, then initialize lnk */ 4329566063dSJacob Faibussowitsch PetscCall(PetscIncompleteLLClean(n, n, nzi, lnk, lnk_lvl, current_space->array, current_space_lvl->array, lnkbt)); 43326fbe8dcSKarl Rupp 43435aa4fcfSShri Abhyankar bj_ptr[i] = current_space->array; 43535aa4fcfSShri Abhyankar bjlvl_ptr[i] = current_space_lvl->array; 43635aa4fcfSShri Abhyankar 43735aa4fcfSShri Abhyankar /* make sure the active row i has diagonal entry */ 43894bad497SJacob 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); 43935aa4fcfSShri Abhyankar 44035aa4fcfSShri Abhyankar current_space->array += nzi; 44135aa4fcfSShri Abhyankar current_space->local_used += nzi; 44235aa4fcfSShri Abhyankar current_space->local_remaining -= nzi; 44326fbe8dcSKarl Rupp 44435aa4fcfSShri Abhyankar current_space_lvl->array += nzi; 44535aa4fcfSShri Abhyankar current_space_lvl->local_used += nzi; 44635aa4fcfSShri Abhyankar current_space_lvl->local_remaining -= nzi; 44735aa4fcfSShri Abhyankar } 44835aa4fcfSShri Abhyankar 4499566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrow, &r)); 4509566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isicol, &ic)); 45135aa4fcfSShri Abhyankar 45235aa4fcfSShri Abhyankar /* copy free_space into bj and free free_space; set bi, bj, bdiag in new datastructure; */ 4539566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bi[n] + 1, &bj)); 4549566063dSJacob Faibussowitsch PetscCall(PetscFreeSpaceContiguous_LU(&free_space, bj, n, bi, bdiag)); 45535aa4fcfSShri Abhyankar 4569566063dSJacob Faibussowitsch PetscCall(PetscIncompleteLLDestroy(lnk, lnkbt)); 4579566063dSJacob Faibussowitsch PetscCall(PetscFreeSpaceDestroy(free_space_lvl)); 4589566063dSJacob Faibussowitsch PetscCall(PetscFree2(bj_ptr, bjlvl_ptr)); 45935aa4fcfSShri Abhyankar 46035aa4fcfSShri Abhyankar #if defined(PETSC_USE_INFO) 46135aa4fcfSShri Abhyankar { 462aef85c9fSShri Abhyankar PetscReal af = ((PetscReal)(bdiag[0] + 1)) / ((PetscReal)ai[n]); 4639566063dSJacob Faibussowitsch PetscCall(PetscInfo(A, "Reallocs %" PetscInt_FMT " Fill ratio:given %g needed %g\n", reallocs, (double)f, (double)af)); 4649566063dSJacob Faibussowitsch PetscCall(PetscInfo(A, "Run with -[sub_]pc_factor_fill %g or use \n", (double)af)); 4659566063dSJacob Faibussowitsch PetscCall(PetscInfo(A, "PCFactorSetFill([sub]pc,%g);\n", (double)af)); 4669566063dSJacob Faibussowitsch PetscCall(PetscInfo(A, "for best performance.\n")); 467*48a46eb9SPierre Jolivet if (diagonal_fill) PetscCall(PetscInfo(A, "Detected and replaced %" PetscInt_FMT " missing diagonals\n", dcount)); 46835aa4fcfSShri Abhyankar } 46935aa4fcfSShri Abhyankar #endif 47035aa4fcfSShri Abhyankar 47135aa4fcfSShri Abhyankar /* put together the new matrix */ 4729566063dSJacob Faibussowitsch PetscCall(MatSeqBAIJSetPreallocation(fact, bs, MAT_SKIP_ALLOCATION, NULL)); 4739566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)fact, (PetscObject)isicol)); 47426fbe8dcSKarl Rupp 47535aa4fcfSShri Abhyankar b = (Mat_SeqBAIJ *)(fact)->data; 47635aa4fcfSShri Abhyankar b->free_a = PETSC_TRUE; 47735aa4fcfSShri Abhyankar b->free_ij = PETSC_TRUE; 47835aa4fcfSShri Abhyankar b->singlemalloc = PETSC_FALSE; 47926fbe8dcSKarl Rupp 4809566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs2 * (bdiag[0] + 1), &b->a)); 48126fbe8dcSKarl Rupp 48235aa4fcfSShri Abhyankar b->j = bj; 48335aa4fcfSShri Abhyankar b->i = bi; 48435aa4fcfSShri Abhyankar b->diag = bdiag; 48535aa4fcfSShri Abhyankar b->free_diag = PETSC_TRUE; 486f4259b30SLisandro Dalcin b->ilen = NULL; 487f4259b30SLisandro Dalcin b->imax = NULL; 48835aa4fcfSShri Abhyankar b->row = isrow; 48935aa4fcfSShri Abhyankar b->col = iscol; 4909566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)isrow)); 4919566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)iscol)); 49235aa4fcfSShri Abhyankar b->icol = isicol; 49326fbe8dcSKarl Rupp 4949566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs * n + bs, &b->solve_work)); 49535aa4fcfSShri Abhyankar /* In b structure: Free imax, ilen, old a, old j. 49635aa4fcfSShri Abhyankar Allocate bdiag, solve_work, new a, new j */ 4979566063dSJacob Faibussowitsch PetscCall(PetscLogObjectMemory((PetscObject)fact, (bdiag[0] + 1) * (sizeof(PetscInt) + bs2 * sizeof(PetscScalar)))); 49835aa4fcfSShri Abhyankar b->maxnz = b->nz = bdiag[0] + 1; 49926fbe8dcSKarl Rupp 500ae3d28f0SHong Zhang fact->info.factor_mallocs = reallocs; 501ae3d28f0SHong Zhang fact->info.fill_ratio_given = f; 502ae3d28f0SHong Zhang fact->info.fill_ratio_needed = ((PetscReal)(bdiag[0] + 1)) / ((PetscReal)ai[n]); 50326fbe8dcSKarl Rupp 5049566063dSJacob Faibussowitsch PetscCall(MatSeqBAIJSetNumericFactorization(fact, both_identity)); 50535aa4fcfSShri Abhyankar PetscFunctionReturn(0); 50635aa4fcfSShri Abhyankar } 50735aa4fcfSShri Abhyankar 5084e2b4712SSatish Balay /* 5094e2b4712SSatish Balay This code is virtually identical to MatILUFactorSymbolic_SeqAIJ 5104e2b4712SSatish Balay except that the data structure of Mat_SeqAIJ is slightly different. 5114e2b4712SSatish Balay Not a good example of code reuse. 5124e2b4712SSatish Balay */ 5139371c9d4SSatish Balay PetscErrorCode MatILUFactorSymbolic_SeqBAIJ_inplace(Mat fact, Mat A, IS isrow, IS iscol, const MatFactorInfo *info) { 5144e2b4712SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b; 5154e2b4712SSatish Balay IS isicol; 5165d0c19d7SBarry Smith const PetscInt *r, *ic, *ai = a->i, *aj = a->j, *xi; 5175d0c19d7SBarry Smith PetscInt prow, n = a->mbs, *ainew, *ajnew, jmax, *fill, nz, *im, *ajfill, *flev, *xitmp; 518a96a251dSBarry Smith PetscInt *dloc, idx, row, m, fm, nzf, nzi, reallocate = 0, dcount = 0; 519d0f46423SBarry Smith PetscInt incrlev, nnz, i, bs = A->rmap->bs, bs2 = a->bs2, levels, diagonal_fill, dd; 520ace3abfcSBarry Smith PetscBool col_identity, row_identity, both_identity, flg; 521329f5518SBarry Smith PetscReal f; 5224e2b4712SSatish Balay 5234e2b4712SSatish Balay PetscFunctionBegin; 5249566063dSJacob Faibussowitsch PetscCall(MatMissingDiagonal_SeqBAIJ(A, &flg, &dd)); 52528b400f6SJacob Faibussowitsch PetscCheck(!flg, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Matrix A is missing diagonal entry in row %" PetscInt_FMT, dd); 5266bce7ff8SHong Zhang 527435faa5fSBarry Smith f = info->fill; 528690b6cddSBarry Smith levels = (PetscInt)info->levels; 529690b6cddSBarry Smith diagonal_fill = (PetscInt)info->diagonal_fill; 53026fbe8dcSKarl Rupp 5319566063dSJacob Faibussowitsch PetscCall(ISInvertPermutation(iscol, PETSC_DECIDE, &isicol)); 53216a2bf60SHong Zhang 5339566063dSJacob Faibussowitsch PetscCall(ISIdentity(isrow, &row_identity)); 5349566063dSJacob Faibussowitsch PetscCall(ISIdentity(iscol, &col_identity)); 535ace3abfcSBarry Smith both_identity = (PetscBool)(row_identity && col_identity); 536309c388cSBarry Smith 53741df41f0SMatthew Knepley if (!levels && both_identity) { /* special case copy the nonzero structure */ 5389566063dSJacob Faibussowitsch PetscCall(MatDuplicateNoCreate_SeqBAIJ(fact, A, MAT_DO_NOT_COPY_VALUES, PETSC_TRUE)); 5399566063dSJacob Faibussowitsch PetscCall(MatSeqBAIJSetNumericFactorization_inplace(fact, both_identity)); 5406bce7ff8SHong Zhang 541d5f3da31SBarry Smith fact->factortype = MAT_FACTOR_ILU; 542ae3d28f0SHong Zhang b = (Mat_SeqBAIJ *)fact->data; 543bb3d539aSBarry Smith b->row = isrow; 544bb3d539aSBarry Smith b->col = iscol; 5459566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)isrow)); 5469566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)iscol)); 547bb3d539aSBarry Smith b->icol = isicol; 548bcd9e38bSBarry Smith b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE; 54926fbe8dcSKarl Rupp 5509566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((n + 1) * bs, &b->solve_work)); 5516bce7ff8SHong Zhang PetscFunctionReturn(0); 5526bce7ff8SHong Zhang } 5536bce7ff8SHong Zhang 5546bce7ff8SHong Zhang /* general case perform the symbolic factorization */ 5559566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isrow, &r)); 5569566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isicol, &ic)); 5574e2b4712SSatish Balay 5584e2b4712SSatish Balay /* get new row pointers */ 5599566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n + 1, &ainew)); 5604e2b4712SSatish Balay ainew[0] = 0; 5614e2b4712SSatish Balay /* don't know how many column pointers are needed so estimate */ 562690b6cddSBarry Smith jmax = (PetscInt)(f * ai[n] + 1); 5639566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jmax, &ajnew)); 5644e2b4712SSatish Balay /* ajfill is level of fill for each fill entry */ 5659566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jmax, &ajfill)); 5664e2b4712SSatish Balay /* fill is a linked list of nonzeros in active row */ 5679566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n + 1, &fill)); 5684e2b4712SSatish Balay /* im is level for each filled value */ 5699566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n + 1, &im)); 5704e2b4712SSatish Balay /* dloc is location of diagonal in factor */ 5719566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n + 1, &dloc)); 5724e2b4712SSatish Balay dloc[0] = 0; 5734e2b4712SSatish Balay for (prow = 0; prow < n; prow++) { 574435faa5fSBarry Smith /* copy prow into linked list */ 5754e2b4712SSatish Balay nzf = nz = ai[r[prow] + 1] - ai[r[prow]]; 57628b400f6SJacob 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); 5774e2b4712SSatish Balay xi = aj + ai[r[prow]]; 5784e2b4712SSatish Balay fill[n] = n; 579435faa5fSBarry Smith fill[prow] = -1; /* marker for diagonal entry */ 5804e2b4712SSatish Balay while (nz--) { 5814e2b4712SSatish Balay fm = n; 5824e2b4712SSatish Balay idx = ic[*xi++]; 5834e2b4712SSatish Balay do { 5844e2b4712SSatish Balay m = fm; 5854e2b4712SSatish Balay fm = fill[m]; 5864e2b4712SSatish Balay } while (fm < idx); 5874e2b4712SSatish Balay fill[m] = idx; 5884e2b4712SSatish Balay fill[idx] = fm; 5894e2b4712SSatish Balay im[idx] = 0; 5904e2b4712SSatish Balay } 591435faa5fSBarry Smith 592435faa5fSBarry Smith /* make sure diagonal entry is included */ 593435faa5fSBarry Smith if (diagonal_fill && fill[prow] == -1) { 594435faa5fSBarry Smith fm = n; 595435faa5fSBarry Smith while (fill[fm] < prow) fm = fill[fm]; 596435faa5fSBarry Smith fill[prow] = fill[fm]; /* insert diagonal into linked list */ 597435faa5fSBarry Smith fill[fm] = prow; 598435faa5fSBarry Smith im[prow] = 0; 599435faa5fSBarry Smith nzf++; 600335d9088SBarry Smith dcount++; 601435faa5fSBarry Smith } 602435faa5fSBarry Smith 6034e2b4712SSatish Balay nzi = 0; 6044e2b4712SSatish Balay row = fill[n]; 6054e2b4712SSatish Balay while (row < prow) { 6064e2b4712SSatish Balay incrlev = im[row] + 1; 6074e2b4712SSatish Balay nz = dloc[row]; 608435faa5fSBarry Smith xi = ajnew + ainew[row] + nz + 1; 6094e2b4712SSatish Balay flev = ajfill + ainew[row] + nz + 1; 6104e2b4712SSatish Balay nnz = ainew[row + 1] - ainew[row] - nz - 1; 6114e2b4712SSatish Balay fm = row; 6124e2b4712SSatish Balay while (nnz-- > 0) { 6134e2b4712SSatish Balay idx = *xi++; 6144e2b4712SSatish Balay if (*flev + incrlev > levels) { 6154e2b4712SSatish Balay flev++; 6164e2b4712SSatish Balay continue; 6174e2b4712SSatish Balay } 6184e2b4712SSatish Balay do { 6194e2b4712SSatish Balay m = fm; 6204e2b4712SSatish Balay fm = fill[m]; 6214e2b4712SSatish Balay } while (fm < idx); 6224e2b4712SSatish Balay if (fm != idx) { 6234e2b4712SSatish Balay im[idx] = *flev + incrlev; 6244e2b4712SSatish Balay fill[m] = idx; 6254e2b4712SSatish Balay fill[idx] = fm; 6264e2b4712SSatish Balay fm = idx; 6274e2b4712SSatish Balay nzf++; 62826fbe8dcSKarl Rupp } else if (im[idx] > *flev + incrlev) im[idx] = *flev + incrlev; 6294e2b4712SSatish Balay flev++; 6304e2b4712SSatish Balay } 6314e2b4712SSatish Balay row = fill[row]; 6324e2b4712SSatish Balay nzi++; 6334e2b4712SSatish Balay } 6344e2b4712SSatish Balay /* copy new filled row into permanent storage */ 6354e2b4712SSatish Balay ainew[prow + 1] = ainew[prow] + nzf; 6364e2b4712SSatish Balay if (ainew[prow + 1] > jmax) { 637ecf371e4SBarry Smith /* estimate how much additional space we will need */ 638ecf371e4SBarry Smith /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */ 639ecf371e4SBarry Smith /* just double the memory each time */ 640690b6cddSBarry Smith PetscInt maxadd = jmax; 641ecf371e4SBarry Smith /* maxadd = (int)(((f*ai[n]+1)*(n-prow+5))/n); */ 6424e2b4712SSatish Balay if (maxadd < nzf) maxadd = (n - prow) * (nzf + 1); 6434e2b4712SSatish Balay jmax += maxadd; 644ecf371e4SBarry Smith 645ecf371e4SBarry Smith /* allocate a longer ajnew and ajfill */ 6469566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jmax, &xitmp)); 6479566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(xitmp, ajnew, ainew[prow])); 6489566063dSJacob Faibussowitsch PetscCall(PetscFree(ajnew)); 6495d0c19d7SBarry Smith ajnew = xitmp; 6509566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(jmax, &xitmp)); 6519566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(xitmp, ajfill, ainew[prow])); 6529566063dSJacob Faibussowitsch PetscCall(PetscFree(ajfill)); 6535d0c19d7SBarry Smith ajfill = xitmp; 654eb150c5cSKris Buschelman reallocate++; /* count how many reallocations are needed */ 6554e2b4712SSatish Balay } 6565d0c19d7SBarry Smith xitmp = ajnew + ainew[prow]; 6574e2b4712SSatish Balay flev = ajfill + ainew[prow]; 6584e2b4712SSatish Balay dloc[prow] = nzi; 6594e2b4712SSatish Balay fm = fill[n]; 6604e2b4712SSatish Balay while (nzf--) { 6615d0c19d7SBarry Smith *xitmp++ = fm; 6624e2b4712SSatish Balay *flev++ = im[fm]; 6634e2b4712SSatish Balay fm = fill[fm]; 6644e2b4712SSatish Balay } 665435faa5fSBarry Smith /* make sure row has diagonal entry */ 666aed4548fSBarry 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\ 6679371c9d4SSatish Balay try running with -pc_factor_nonzeros_along_diagonal or -pc_factor_diagonal_fill", 6689371c9d4SSatish Balay prow); 669435faa5fSBarry Smith } 6709566063dSJacob Faibussowitsch PetscCall(PetscFree(ajfill)); 6719566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrow, &r)); 6729566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isicol, &ic)); 6739566063dSJacob Faibussowitsch PetscCall(PetscFree(fill)); 6749566063dSJacob Faibussowitsch PetscCall(PetscFree(im)); 6754e2b4712SSatish Balay 6766cf91177SBarry Smith #if defined(PETSC_USE_INFO) 6774e2b4712SSatish Balay { 678329f5518SBarry Smith PetscReal af = ((PetscReal)ainew[n]) / ((PetscReal)ai[n]); 6799566063dSJacob Faibussowitsch PetscCall(PetscInfo(A, "Reallocs %" PetscInt_FMT " Fill ratio:given %g needed %g\n", reallocate, (double)f, (double)af)); 6809566063dSJacob Faibussowitsch PetscCall(PetscInfo(A, "Run with -pc_factor_fill %g or use \n", (double)af)); 6819566063dSJacob Faibussowitsch PetscCall(PetscInfo(A, "PCFactorSetFill(pc,%g);\n", (double)af)); 6829566063dSJacob Faibussowitsch PetscCall(PetscInfo(A, "for best performance.\n")); 683*48a46eb9SPierre Jolivet if (diagonal_fill) PetscCall(PetscInfo(A, "Detected and replaced %" PetscInt_FMT " missing diagonals\n", dcount)); 6844e2b4712SSatish Balay } 68563ba0a88SBarry Smith #endif 6864e2b4712SSatish Balay 6874e2b4712SSatish Balay /* put together the new matrix */ 6889566063dSJacob Faibussowitsch PetscCall(MatSeqBAIJSetPreallocation(fact, bs, MAT_SKIP_ALLOCATION, NULL)); 6899566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)fact, (PetscObject)isicol)); 690ae3d28f0SHong Zhang b = (Mat_SeqBAIJ *)fact->data; 69126fbe8dcSKarl Rupp 692e6b907acSBarry Smith b->free_a = PETSC_TRUE; 693e6b907acSBarry Smith b->free_ij = PETSC_TRUE; 6947c922b88SBarry Smith b->singlemalloc = PETSC_FALSE; 69526fbe8dcSKarl Rupp 6969566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs2 * ainew[n], &b->a)); 69726fbe8dcSKarl Rupp 6984e2b4712SSatish Balay b->j = ajnew; 6994e2b4712SSatish Balay b->i = ainew; 7004e2b4712SSatish Balay for (i = 0; i < n; i++) dloc[i] += ainew[i]; 7014e2b4712SSatish Balay b->diag = dloc; 7027f53bb6cSHong Zhang b->free_diag = PETSC_TRUE; 703f4259b30SLisandro Dalcin b->ilen = NULL; 704f4259b30SLisandro Dalcin b->imax = NULL; 7054e2b4712SSatish Balay b->row = isrow; 7064e2b4712SSatish Balay b->col = iscol; 707bcd9e38bSBarry Smith b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE; 70826fbe8dcSKarl Rupp 7099566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)isrow)); 7109566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)iscol)); 711e51c0b9cSSatish Balay b->icol = isicol; 7129566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs * n + bs, &b->solve_work)); 7134e2b4712SSatish Balay /* In b structure: Free imax, ilen, old a, old j. 7144e2b4712SSatish Balay Allocate dloc, solve_work, new a, new j */ 7159566063dSJacob Faibussowitsch PetscCall(PetscLogObjectMemory((PetscObject)fact, (ainew[n] - n) * (sizeof(PetscInt)) + bs2 * ainew[n] * sizeof(PetscScalar))); 7164e2b4712SSatish Balay b->maxnz = b->nz = ainew[n]; 7174e2b4712SSatish Balay 718ae3d28f0SHong Zhang fact->info.factor_mallocs = reallocate; 719ae3d28f0SHong Zhang fact->info.fill_ratio_given = f; 720ae3d28f0SHong Zhang fact->info.fill_ratio_needed = ((PetscReal)ainew[n]) / ((PetscReal)ai[prow]); 7216bce7ff8SHong Zhang 7229566063dSJacob Faibussowitsch PetscCall(MatSeqBAIJSetNumericFactorization_inplace(fact, both_identity)); 7238661488fSKris Buschelman PetscFunctionReturn(0); 7248661488fSKris Buschelman } 7258661488fSKris Buschelman 7269371c9d4SSatish Balay PetscErrorCode MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE(Mat A) { 72712272027SHong Zhang /* Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; */ 72812272027SHong Zhang /* int i,*AJ=a->j,nz=a->nz; */ 7295fd66863SKarl Rupp 7305a9542e3SKris Buschelman PetscFunctionBegin; 7317cf1b8d3SKris Buschelman /* Undo Column scaling */ 7327cf1b8d3SKris Buschelman /* while (nz--) { */ 7337cf1b8d3SKris Buschelman /* AJ[i] = AJ[i]/4; */ 7347cf1b8d3SKris Buschelman /* } */ 735c115a38dSKris Buschelman /* This should really invoke a push/pop logic, but we don't have that yet. */ 7360298fd71SBarry Smith A->ops->setunfactored = NULL; 7377cf1b8d3SKris Buschelman PetscFunctionReturn(0); 7387cf1b8d3SKris Buschelman } 7397cf1b8d3SKris Buschelman 7409371c9d4SSatish Balay PetscErrorCode MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE_usj(Mat A) { 7417cf1b8d3SKris Buschelman Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data; 742b24ad042SBarry Smith PetscInt *AJ = a->j, nz = a->nz; 7432aa5897fSKris Buschelman unsigned short *aj = (unsigned short *)AJ; 7445fd66863SKarl Rupp 7455a9542e3SKris Buschelman PetscFunctionBegin; 7460b9da03eSKris Buschelman /* Is this really necessary? */ 7479371c9d4SSatish Balay while (nz--) { AJ[nz] = (int)((unsigned int)aj[nz]); /* First extend, then convert to signed. */ } 7480298fd71SBarry Smith A->ops->setunfactored = NULL; 7497e7071cdSKris Buschelman PetscFunctionReturn(0); 7507e7071cdSKris Buschelman } 751