1c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h> 2c6db04a5SJed Brown #include <../src/mat/impls/sbaij/seq/sbaij.h> 3c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/bas/spbas.h> 4b2f1ab58SBarry Smith 566976f2fSJacob Faibussowitsch static PetscErrorCode MatICCFactorSymbolic_SeqAIJ_Bas(Mat fact, Mat A, IS perm, const MatFactorInfo *info) 6d71ae5a4SJacob Faibussowitsch { 7b2f1ab58SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 8b2f1ab58SBarry Smith Mat_SeqSBAIJ *b; 9ace3abfcSBarry Smith PetscBool perm_identity, missing; 10b2f1ab58SBarry Smith PetscInt reallocs = 0, i, *ai = a->i, *aj = a->j, am = A->rmap->n, *ui; 11b2f1ab58SBarry Smith const PetscInt *rip, *riip; 12b2f1ab58SBarry Smith PetscInt j; 13b2f1ab58SBarry Smith PetscInt d; 14b2f1ab58SBarry Smith PetscInt ncols, *cols, *uj; 15b2f1ab58SBarry Smith PetscReal fill = info->fill, levels = info->levels; 16b2f1ab58SBarry Smith IS iperm; 17b2f1ab58SBarry Smith spbas_matrix Pattern_0, Pattern_P; 18b2f1ab58SBarry Smith 19b2f1ab58SBarry Smith PetscFunctionBegin; 205f80ce2aSJacob Faibussowitsch 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); 219566063dSJacob Faibussowitsch PetscCall(MatMissingDiagonal(A, &missing, &d)); 225f80ce2aSJacob Faibussowitsch PetscCheck(!missing, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Matrix is missing diagonal entry %" PetscInt_FMT, d); 239566063dSJacob Faibussowitsch PetscCall(ISIdentity(perm, &perm_identity)); 249566063dSJacob Faibussowitsch PetscCall(ISInvertPermutation(perm, PETSC_DECIDE, &iperm)); 25b2f1ab58SBarry Smith 26b2f1ab58SBarry Smith /* ICC(0) without matrix ordering: simply copies fill pattern */ 27b2f1ab58SBarry Smith if (!levels && perm_identity) { 289566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(am + 1, &ui)); 29b2f1ab58SBarry Smith ui[0] = 0; 30b2f1ab58SBarry Smith 31ad540459SPierre Jolivet for (i = 0; i < am; i++) ui[i + 1] = ui[i] + ai[i + 1] - a->diag[i]; 329566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ui[am] + 1, &uj)); 33b2f1ab58SBarry Smith cols = uj; 34b2f1ab58SBarry Smith for (i = 0; i < am; i++) { 35b2f1ab58SBarry Smith aj = a->j + a->diag[i]; 36b2f1ab58SBarry Smith ncols = ui[i + 1] - ui[i]; 37b2f1ab58SBarry Smith for (j = 0; j < ncols; j++) *cols++ = *aj++; 38b2f1ab58SBarry Smith } 39b2f1ab58SBarry Smith } else { /* case: levels>0 || (levels=0 && !perm_identity) */ 409566063dSJacob Faibussowitsch PetscCall(ISGetIndices(iperm, &riip)); 419566063dSJacob Faibussowitsch PetscCall(ISGetIndices(perm, &rip)); 42b2f1ab58SBarry Smith 4371d55d18SBarry Smith /* Create spbas_matrix for pattern */ 449566063dSJacob Faibussowitsch PetscCall(spbas_pattern_only(am, am, ai, aj, &Pattern_0)); 45b2f1ab58SBarry Smith 4671d55d18SBarry Smith /* Apply the permutation */ 479566063dSJacob Faibussowitsch PetscCall(spbas_apply_reordering(&Pattern_0, rip, riip)); 48b2f1ab58SBarry Smith 4971d55d18SBarry Smith /* Raise the power */ 509566063dSJacob Faibussowitsch PetscCall(spbas_power(Pattern_0, (int)levels + 1, &Pattern_P)); 519566063dSJacob Faibussowitsch PetscCall(spbas_delete(Pattern_0)); 52b2f1ab58SBarry Smith 5371d55d18SBarry Smith /* Keep only upper triangle of pattern */ 549566063dSJacob Faibussowitsch PetscCall(spbas_keep_upper(&Pattern_P)); 55b2f1ab58SBarry Smith 5671d55d18SBarry Smith /* Convert to Sparse Row Storage */ 579566063dSJacob Faibussowitsch PetscCall(spbas_matrix_to_crs(Pattern_P, NULL, &ui, &uj)); 589566063dSJacob Faibussowitsch PetscCall(spbas_delete(Pattern_P)); 59b2f1ab58SBarry Smith } /* end of case: levels>0 || (levels=0 && !perm_identity) */ 60b2f1ab58SBarry Smith 61b2f1ab58SBarry Smith /* put together the new matrix in MATSEQSBAIJ format */ 62b2f1ab58SBarry Smith 6357508eceSPierre Jolivet b = (Mat_SeqSBAIJ *)fact->data; 64*84648c2dSPierre Jolivet PetscCall(PetscMalloc1(ui[am], &b->a)); 652205254eSKarl Rupp 66b2f1ab58SBarry Smith b->j = uj; 67b2f1ab58SBarry Smith b->i = ui; 68f4259b30SLisandro Dalcin b->diag = NULL; 69f4259b30SLisandro Dalcin b->ilen = NULL; 70f4259b30SLisandro Dalcin b->imax = NULL; 71b2f1ab58SBarry Smith b->row = perm; 72b2f1ab58SBarry Smith b->col = perm; 732205254eSKarl Rupp 749566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)perm)); 759566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)perm)); 762205254eSKarl Rupp 77b2f1ab58SBarry Smith b->icol = iperm; 78b2f1ab58SBarry Smith b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 79*84648c2dSPierre Jolivet PetscCall(PetscMalloc1(am, &b->solve_work)); 80b2f1ab58SBarry Smith b->maxnz = b->nz = ui[am]; 81b2f1ab58SBarry Smith b->free_a = PETSC_TRUE; 82b2f1ab58SBarry Smith b->free_ij = PETSC_TRUE; 83b2f1ab58SBarry Smith 8457508eceSPierre Jolivet fact->info.factor_mallocs = reallocs; 8557508eceSPierre Jolivet fact->info.fill_ratio_given = fill; 86b2f1ab58SBarry Smith if (ai[am] != 0) { 8757508eceSPierre Jolivet fact->info.fill_ratio_needed = (PetscReal)ui[am] / (PetscReal)ai[am]; 88b2f1ab58SBarry Smith } else { 8957508eceSPierre Jolivet fact->info.fill_ratio_needed = 0.0; 90b2f1ab58SBarry Smith } 9157508eceSPierre Jolivet /* fact->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ_inplace; */ 923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 93b2f1ab58SBarry Smith } 94b2f1ab58SBarry Smith 9566976f2fSJacob Faibussowitsch static PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_Bas(Mat B, Mat A, const MatFactorInfo *info) 96d71ae5a4SJacob Faibussowitsch { 97b2f1ab58SBarry Smith Mat C = B; 98b2f1ab58SBarry Smith Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ *)C->data; 99b2f1ab58SBarry Smith IS ip = b->row, iip = b->icol; 100b2f1ab58SBarry Smith const PetscInt *rip, *riip; 101b2f1ab58SBarry Smith PetscInt mbs = A->rmap->n, *bi = b->i, *bj = b->j; 102b2f1ab58SBarry Smith MatScalar *ba = b->a; 103c88783f4SHong Zhang PetscReal shiftnz = info->shiftamount; 104d36aa34eSBarry Smith PetscReal droptol = -1; 105ace3abfcSBarry Smith PetscBool perm_identity; 106b2f1ab58SBarry Smith spbas_matrix Pattern, matrix_L, matrix_LT; 1079767453cSBarry Smith PetscReal mem_reduction; 108b2f1ab58SBarry Smith 109b2f1ab58SBarry Smith PetscFunctionBegin; 11071d55d18SBarry Smith /* Reduce memory requirements: erase values of B-matrix */ 1119566063dSJacob Faibussowitsch PetscCall(PetscFree(ba)); 11271d55d18SBarry Smith /* Compress (maximum) sparseness pattern of B-matrix */ 1139566063dSJacob Faibussowitsch PetscCall(spbas_compress_pattern(bi, bj, mbs, mbs, SPBAS_DIAGONAL_OFFSETS, &Pattern, &mem_reduction)); 1149566063dSJacob Faibussowitsch PetscCall(PetscFree(bi)); 1159566063dSJacob Faibussowitsch PetscCall(PetscFree(bj)); 116b2f1ab58SBarry Smith 1179566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, " compression rate for spbas_compress_pattern %g \n", (double)mem_reduction)); 118b2f1ab58SBarry Smith 1194e1ff37aSBarry Smith /* Make Cholesky decompositions with larger Manteuffel shifts until no more negative diagonals are found. */ 1209566063dSJacob Faibussowitsch PetscCall(ISGetIndices(ip, &rip)); 1219566063dSJacob Faibussowitsch PetscCall(ISGetIndices(iip, &riip)); 122b2f1ab58SBarry Smith 1235f80ce2aSJacob Faibussowitsch if (info->usedt) droptol = info->dt; 1245f80ce2aSJacob Faibussowitsch 1253ba16761SJacob Faibussowitsch for (int ierr = NEGATIVE_DIAGONAL; ierr == NEGATIVE_DIAGONAL;) { 126cc442ecdSBarry Smith PetscBool success; 1275f80ce2aSJacob Faibussowitsch 1283ba16761SJacob Faibussowitsch ierr = (int)spbas_incomplete_cholesky(A, rip, riip, Pattern, droptol, shiftnz, &matrix_LT, &success); 129cc442ecdSBarry Smith if (!success) { 130b2f1ab58SBarry Smith shiftnz *= 1.5; 1319767453cSBarry Smith if (shiftnz < 1e-5) shiftnz = 1e-5; 1329566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "spbas_incomplete_cholesky found a negative diagonal. Trying again with Manteuffel shift=%g\n", (double)shiftnz)); 133b2f1ab58SBarry Smith } 134b2f1ab58SBarry Smith } 1359566063dSJacob Faibussowitsch PetscCall(spbas_delete(Pattern)); 136b2f1ab58SBarry Smith 137835f2295SStefano Zampini PetscCall(PetscInfo(NULL, " memory_usage for spbas_incomplete_cholesky %g bytes per row\n", (double)(spbas_memory_requirement(matrix_LT) / (PetscReal)mbs))); 138b2f1ab58SBarry Smith 1399566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(ip, &rip)); 1409566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(iip, &riip)); 141b2f1ab58SBarry Smith 14271d55d18SBarry Smith /* Convert spbas_matrix to compressed row storage */ 1439566063dSJacob Faibussowitsch PetscCall(spbas_transpose(matrix_LT, &matrix_L)); 1449566063dSJacob Faibussowitsch PetscCall(spbas_delete(matrix_LT)); 1459566063dSJacob Faibussowitsch PetscCall(spbas_matrix_to_crs(matrix_L, &ba, &bi, &bj)); 1469371c9d4SSatish Balay b->i = bi; 1479371c9d4SSatish Balay b->j = bj; 1489371c9d4SSatish Balay b->a = ba; 1499566063dSJacob Faibussowitsch PetscCall(spbas_delete(matrix_L)); 150b2f1ab58SBarry Smith 15171d55d18SBarry Smith /* Set the appropriate solution functions */ 1529566063dSJacob Faibussowitsch PetscCall(ISIdentity(ip, &perm_identity)); 153b2f1ab58SBarry Smith if (perm_identity) { 15457508eceSPierre Jolivet B->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace; 15557508eceSPierre Jolivet B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace; 15657508eceSPierre Jolivet B->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace; 15757508eceSPierre Jolivet B->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace; 158b2f1ab58SBarry Smith } else { 15957508eceSPierre Jolivet B->ops->solve = MatSolve_SeqSBAIJ_1_inplace; 16057508eceSPierre Jolivet B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_inplace; 16157508eceSPierre Jolivet B->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_1_inplace; 16257508eceSPierre Jolivet B->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_1_inplace; 163b2f1ab58SBarry Smith } 164b2f1ab58SBarry Smith 165b2f1ab58SBarry Smith C->assembled = PETSC_TRUE; 166b2f1ab58SBarry Smith C->preallocated = PETSC_TRUE; 1672205254eSKarl Rupp 1689566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(C->rmap->n)); 1693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 170b2f1ab58SBarry Smith } 171b2f1ab58SBarry Smith 17266976f2fSJacob Faibussowitsch static PetscErrorCode MatFactorGetSolverType_seqaij_bas(Mat A, MatSolverType *type) 173d71ae5a4SJacob Faibussowitsch { 1743dfa2556SBarry Smith PetscFunctionBegin; 1753dfa2556SBarry Smith *type = MATSOLVERBAS; 1763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1773dfa2556SBarry Smith } 1783dfa2556SBarry Smith 179d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_bas(Mat A, MatFactorType ftype, Mat *B) 180d71ae5a4SJacob Faibussowitsch { 181b2f1ab58SBarry Smith PetscInt n = A->rmap->n; 182b2f1ab58SBarry Smith 183b2f1ab58SBarry Smith PetscFunctionBegin; 1849566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), B)); 1859566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*B, n, n, n, n)); 186966bd95aSPierre Jolivet PetscCheck(ftype == MAT_FACTOR_ICC, PETSC_COMM_SELF, PETSC_ERR_SUP, "Factor type not supported"); 1879566063dSJacob Faibussowitsch PetscCall(MatSetType(*B, MATSEQSBAIJ)); 1889566063dSJacob Faibussowitsch PetscCall(MatSeqSBAIJSetPreallocation(*B, 1, MAT_SKIP_ALLOCATION, NULL)); 1892205254eSKarl Rupp 190b2f1ab58SBarry Smith (*B)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqAIJ_Bas; 191b2f1ab58SBarry Smith (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ_Bas; 1929566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)*B, "MatFactorGetSolverType_C", MatFactorGetSolverType_seqaij_bas)); 1939566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATORDERINGND, (char **)&(*B)->preferredordering[MAT_FACTOR_LU])); 1949566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATORDERINGND, (char **)&(*B)->preferredordering[MAT_FACTOR_CHOLESKY])); 195d5f3da31SBarry Smith (*B)->factortype = ftype; 19600c67f3bSHong Zhang 1979566063dSJacob Faibussowitsch PetscCall(PetscFree((*B)->solvertype)); 1989566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATSOLVERBAS, &(*B)->solvertype)); 199f73b0415SBarry Smith (*B)->canuseordering = PETSC_TRUE; 2009566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATORDERINGNATURAL, (char **)&(*B)->preferredordering[MAT_FACTOR_ICC])); 2013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 202b2f1ab58SBarry Smith } 203