1b2f1ab58SBarry Smith 2c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h> 3c6db04a5SJed Brown #include <../src/mat/impls/sbaij/seq/sbaij.h> 4c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/bas/spbas.h> 5b2f1ab58SBarry Smith 6d71ae5a4SJacob Faibussowitsch PetscErrorCode MatICCFactorSymbolic_SeqAIJ_Bas(Mat fact, Mat A, IS perm, const MatFactorInfo *info) 7d71ae5a4SJacob Faibussowitsch { 8b2f1ab58SBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 9b2f1ab58SBarry Smith Mat_SeqSBAIJ *b; 10ace3abfcSBarry Smith PetscBool perm_identity, missing; 11b2f1ab58SBarry Smith PetscInt reallocs = 0, i, *ai = a->i, *aj = a->j, am = A->rmap->n, *ui; 12b2f1ab58SBarry Smith const PetscInt *rip, *riip; 13b2f1ab58SBarry Smith PetscInt j; 14b2f1ab58SBarry Smith PetscInt d; 15b2f1ab58SBarry Smith PetscInt ncols, *cols, *uj; 16b2f1ab58SBarry Smith PetscReal fill = info->fill, levels = info->levels; 17b2f1ab58SBarry Smith IS iperm; 18b2f1ab58SBarry Smith spbas_matrix Pattern_0, Pattern_P; 19b2f1ab58SBarry Smith 20b2f1ab58SBarry Smith PetscFunctionBegin; 215f80ce2aSJacob 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); 229566063dSJacob Faibussowitsch PetscCall(MatMissingDiagonal(A, &missing, &d)); 235f80ce2aSJacob Faibussowitsch PetscCheck(!missing, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Matrix is missing diagonal entry %" PetscInt_FMT, d); 249566063dSJacob Faibussowitsch PetscCall(ISIdentity(perm, &perm_identity)); 259566063dSJacob Faibussowitsch PetscCall(ISInvertPermutation(perm, PETSC_DECIDE, &iperm)); 26b2f1ab58SBarry Smith 27b2f1ab58SBarry Smith /* ICC(0) without matrix ordering: simply copies fill pattern */ 28b2f1ab58SBarry Smith if (!levels && perm_identity) { 299566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(am + 1, &ui)); 30b2f1ab58SBarry Smith ui[0] = 0; 31b2f1ab58SBarry Smith 32ad540459SPierre Jolivet for (i = 0; i < am; i++) ui[i + 1] = ui[i] + ai[i + 1] - a->diag[i]; 339566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ui[am] + 1, &uj)); 34b2f1ab58SBarry Smith cols = uj; 35b2f1ab58SBarry Smith for (i = 0; i < am; i++) { 36b2f1ab58SBarry Smith aj = a->j + a->diag[i]; 37b2f1ab58SBarry Smith ncols = ui[i + 1] - ui[i]; 38b2f1ab58SBarry Smith for (j = 0; j < ncols; j++) *cols++ = *aj++; 39b2f1ab58SBarry Smith } 40b2f1ab58SBarry Smith } else { /* case: levels>0 || (levels=0 && !perm_identity) */ 419566063dSJacob Faibussowitsch PetscCall(ISGetIndices(iperm, &riip)); 429566063dSJacob Faibussowitsch PetscCall(ISGetIndices(perm, &rip)); 43b2f1ab58SBarry Smith 4471d55d18SBarry Smith /* Create spbas_matrix for pattern */ 459566063dSJacob Faibussowitsch PetscCall(spbas_pattern_only(am, am, ai, aj, &Pattern_0)); 46b2f1ab58SBarry Smith 4771d55d18SBarry Smith /* Apply the permutation */ 489566063dSJacob Faibussowitsch PetscCall(spbas_apply_reordering(&Pattern_0, rip, riip)); 49b2f1ab58SBarry Smith 5071d55d18SBarry Smith /* Raise the power */ 519566063dSJacob Faibussowitsch PetscCall(spbas_power(Pattern_0, (int)levels + 1, &Pattern_P)); 529566063dSJacob Faibussowitsch PetscCall(spbas_delete(Pattern_0)); 53b2f1ab58SBarry Smith 5471d55d18SBarry Smith /* Keep only upper triangle of pattern */ 559566063dSJacob Faibussowitsch PetscCall(spbas_keep_upper(&Pattern_P)); 56b2f1ab58SBarry Smith 5771d55d18SBarry Smith /* Convert to Sparse Row Storage */ 589566063dSJacob Faibussowitsch PetscCall(spbas_matrix_to_crs(Pattern_P, NULL, &ui, &uj)); 599566063dSJacob Faibussowitsch PetscCall(spbas_delete(Pattern_P)); 60b2f1ab58SBarry Smith } /* end of case: levels>0 || (levels=0 && !perm_identity) */ 61b2f1ab58SBarry Smith 62b2f1ab58SBarry Smith /* put together the new matrix in MATSEQSBAIJ format */ 63b2f1ab58SBarry Smith 64b2f1ab58SBarry Smith b = (Mat_SeqSBAIJ *)(fact)->data; 65b2f1ab58SBarry Smith b->singlemalloc = PETSC_FALSE; 662205254eSKarl Rupp 679566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ui[am] + 1, &b->a)); 682205254eSKarl Rupp 69b2f1ab58SBarry Smith b->j = uj; 70b2f1ab58SBarry Smith b->i = ui; 71f4259b30SLisandro Dalcin b->diag = NULL; 72f4259b30SLisandro Dalcin b->ilen = NULL; 73f4259b30SLisandro Dalcin b->imax = NULL; 74b2f1ab58SBarry Smith b->row = perm; 75b2f1ab58SBarry Smith b->col = perm; 762205254eSKarl Rupp 779566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)perm)); 789566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)perm)); 792205254eSKarl Rupp 80b2f1ab58SBarry Smith b->icol = iperm; 81b2f1ab58SBarry Smith b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 829566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(am + 1, &b->solve_work)); 83b2f1ab58SBarry Smith b->maxnz = b->nz = ui[am]; 84b2f1ab58SBarry Smith b->free_a = PETSC_TRUE; 85b2f1ab58SBarry Smith b->free_ij = PETSC_TRUE; 86b2f1ab58SBarry Smith 87b2f1ab58SBarry Smith (fact)->info.factor_mallocs = reallocs; 88b2f1ab58SBarry Smith (fact)->info.fill_ratio_given = fill; 89b2f1ab58SBarry Smith if (ai[am] != 0) { 90b2f1ab58SBarry Smith (fact)->info.fill_ratio_needed = ((PetscReal)ui[am]) / ((PetscReal)ai[am]); 91b2f1ab58SBarry Smith } else { 92b2f1ab58SBarry Smith (fact)->info.fill_ratio_needed = 0.0; 93b2f1ab58SBarry Smith } 94b2f1ab58SBarry Smith /* (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ_inplace; */ 95*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 96b2f1ab58SBarry Smith } 97b2f1ab58SBarry Smith 98d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_Bas(Mat B, Mat A, const MatFactorInfo *info) 99d71ae5a4SJacob Faibussowitsch { 100b2f1ab58SBarry Smith Mat C = B; 101b2f1ab58SBarry Smith Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ *)C->data; 102b2f1ab58SBarry Smith IS ip = b->row, iip = b->icol; 103b2f1ab58SBarry Smith const PetscInt *rip, *riip; 104b2f1ab58SBarry Smith PetscInt mbs = A->rmap->n, *bi = b->i, *bj = b->j; 105b2f1ab58SBarry Smith MatScalar *ba = b->a; 106c88783f4SHong Zhang PetscReal shiftnz = info->shiftamount; 107d36aa34eSBarry Smith PetscReal droptol = -1; 108ace3abfcSBarry Smith PetscBool perm_identity; 109b2f1ab58SBarry Smith spbas_matrix Pattern, matrix_L, matrix_LT; 1109767453cSBarry Smith PetscReal mem_reduction; 111b2f1ab58SBarry Smith 112b2f1ab58SBarry Smith PetscFunctionBegin; 11371d55d18SBarry Smith /* Reduce memory requirements: erase values of B-matrix */ 1149566063dSJacob Faibussowitsch PetscCall(PetscFree(ba)); 11571d55d18SBarry Smith /* Compress (maximum) sparseness pattern of B-matrix */ 1169566063dSJacob Faibussowitsch PetscCall(spbas_compress_pattern(bi, bj, mbs, mbs, SPBAS_DIAGONAL_OFFSETS, &Pattern, &mem_reduction)); 1179566063dSJacob Faibussowitsch PetscCall(PetscFree(bi)); 1189566063dSJacob Faibussowitsch PetscCall(PetscFree(bj)); 119b2f1ab58SBarry Smith 1209566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, " compression rate for spbas_compress_pattern %g \n", (double)mem_reduction)); 121b2f1ab58SBarry Smith 1224e1ff37aSBarry Smith /* Make Cholesky decompositions with larger Manteuffel shifts until no more negative diagonals are found. */ 1239566063dSJacob Faibussowitsch PetscCall(ISGetIndices(ip, &rip)); 1249566063dSJacob Faibussowitsch PetscCall(ISGetIndices(iip, &riip)); 125b2f1ab58SBarry Smith 1265f80ce2aSJacob Faibussowitsch if (info->usedt) droptol = info->dt; 1275f80ce2aSJacob Faibussowitsch 128*3ba16761SJacob Faibussowitsch for (int ierr = NEGATIVE_DIAGONAL; ierr == NEGATIVE_DIAGONAL;) { 129cc442ecdSBarry Smith PetscBool success; 1305f80ce2aSJacob Faibussowitsch 131*3ba16761SJacob Faibussowitsch ierr = (int)spbas_incomplete_cholesky(A, rip, riip, Pattern, droptol, shiftnz, &matrix_LT, &success); 132cc442ecdSBarry Smith if (!success) { 133b2f1ab58SBarry Smith shiftnz *= 1.5; 1349767453cSBarry Smith if (shiftnz < 1e-5) shiftnz = 1e-5; 1359566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "spbas_incomplete_cholesky found a negative diagonal. Trying again with Manteuffel shift=%g\n", (double)shiftnz)); 136b2f1ab58SBarry Smith } 137b2f1ab58SBarry Smith } 1389566063dSJacob Faibussowitsch PetscCall(spbas_delete(Pattern)); 139b2f1ab58SBarry Smith 1409566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, " memory_usage for spbas_incomplete_cholesky %g bytes per row\n", (double)(PetscReal)(spbas_memory_requirement(matrix_LT) / (PetscReal)mbs))); 141b2f1ab58SBarry Smith 1429566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(ip, &rip)); 1439566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(iip, &riip)); 144b2f1ab58SBarry Smith 14571d55d18SBarry Smith /* Convert spbas_matrix to compressed row storage */ 1469566063dSJacob Faibussowitsch PetscCall(spbas_transpose(matrix_LT, &matrix_L)); 1479566063dSJacob Faibussowitsch PetscCall(spbas_delete(matrix_LT)); 1489566063dSJacob Faibussowitsch PetscCall(spbas_matrix_to_crs(matrix_L, &ba, &bi, &bj)); 1499371c9d4SSatish Balay b->i = bi; 1509371c9d4SSatish Balay b->j = bj; 1519371c9d4SSatish Balay b->a = ba; 1529566063dSJacob Faibussowitsch PetscCall(spbas_delete(matrix_L)); 153b2f1ab58SBarry Smith 15471d55d18SBarry Smith /* Set the appropriate solution functions */ 1559566063dSJacob Faibussowitsch PetscCall(ISIdentity(ip, &perm_identity)); 156b2f1ab58SBarry Smith if (perm_identity) { 157b2f1ab58SBarry Smith (B)->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace; 158b2f1ab58SBarry Smith (B)->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace; 159b2f1ab58SBarry Smith (B)->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace; 160b2f1ab58SBarry Smith (B)->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace; 161b2f1ab58SBarry Smith } else { 162b2f1ab58SBarry Smith (B)->ops->solve = MatSolve_SeqSBAIJ_1_inplace; 163b2f1ab58SBarry Smith (B)->ops->solvetranspose = MatSolve_SeqSBAIJ_1_inplace; 164b2f1ab58SBarry Smith (B)->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_1_inplace; 165b2f1ab58SBarry Smith (B)->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_1_inplace; 166b2f1ab58SBarry Smith } 167b2f1ab58SBarry Smith 168b2f1ab58SBarry Smith C->assembled = PETSC_TRUE; 169b2f1ab58SBarry Smith C->preallocated = PETSC_TRUE; 1702205254eSKarl Rupp 1719566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(C->rmap->n)); 172*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 173b2f1ab58SBarry Smith } 174b2f1ab58SBarry Smith 175d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFactorGetSolverType_seqaij_bas(Mat A, MatSolverType *type) 176d71ae5a4SJacob Faibussowitsch { 1773dfa2556SBarry Smith PetscFunctionBegin; 1783dfa2556SBarry Smith *type = MATSOLVERBAS; 179*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1803dfa2556SBarry Smith } 1813dfa2556SBarry Smith 182d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_bas(Mat A, MatFactorType ftype, Mat *B) 183d71ae5a4SJacob Faibussowitsch { 184b2f1ab58SBarry Smith PetscInt n = A->rmap->n; 185b2f1ab58SBarry Smith 186b2f1ab58SBarry Smith PetscFunctionBegin; 1879566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), B)); 1889566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*B, n, n, n, n)); 189b2f1ab58SBarry Smith if (ftype == MAT_FACTOR_ICC) { 1909566063dSJacob Faibussowitsch PetscCall(MatSetType(*B, MATSEQSBAIJ)); 1919566063dSJacob Faibussowitsch PetscCall(MatSeqSBAIJSetPreallocation(*B, 1, MAT_SKIP_ALLOCATION, NULL)); 1922205254eSKarl Rupp 193b2f1ab58SBarry Smith (*B)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqAIJ_Bas; 194b2f1ab58SBarry Smith (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ_Bas; 1959566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)*B, "MatFactorGetSolverType_C", MatFactorGetSolverType_seqaij_bas)); 1969566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATORDERINGND, (char **)&(*B)->preferredordering[MAT_FACTOR_LU])); 1979566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATORDERINGND, (char **)&(*B)->preferredordering[MAT_FACTOR_CHOLESKY])); 198e32f2f54SBarry Smith } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Factor type not supported"); 199d5f3da31SBarry Smith (*B)->factortype = ftype; 20000c67f3bSHong Zhang 2019566063dSJacob Faibussowitsch PetscCall(PetscFree((*B)->solvertype)); 2029566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATSOLVERBAS, &(*B)->solvertype)); 203f73b0415SBarry Smith (*B)->canuseordering = PETSC_TRUE; 2049566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATORDERINGNATURAL, (char **)&(*B)->preferredordering[MAT_FACTOR_ICC])); 205*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 206b2f1ab58SBarry Smith } 207