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
MatICCFactorSymbolic_SeqAIJ_Bas(Mat fact,Mat A,IS perm,const MatFactorInfo * info)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;
9*421480d9SBarry Smith PetscBool perm_identity, diagDense;
10b2f1ab58SBarry Smith PetscInt reallocs = 0, i, *ai = a->i, *aj = a->j, am = A->rmap->n, *ui;
11*421480d9SBarry Smith const PetscInt *rip, *riip, *adiag;
12b2f1ab58SBarry Smith PetscInt j;
13b2f1ab58SBarry Smith PetscInt ncols, *cols, *uj;
14b2f1ab58SBarry Smith PetscReal fill = info->fill, levels = info->levels;
15b2f1ab58SBarry Smith IS iperm;
16b2f1ab58SBarry Smith spbas_matrix Pattern_0, Pattern_P;
17b2f1ab58SBarry Smith
18b2f1ab58SBarry Smith PetscFunctionBegin;
195f80ce2aSJacob 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);
20*421480d9SBarry Smith PetscCall(MatGetDiagonalMarkers_SeqAIJ(A, &adiag, &diagDense));
21*421480d9SBarry Smith PetscCheck(diagDense, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Matrix is missing diagonal entries");
229566063dSJacob Faibussowitsch PetscCall(ISIdentity(perm, &perm_identity));
239566063dSJacob Faibussowitsch PetscCall(ISInvertPermutation(perm, PETSC_DECIDE, &iperm));
24b2f1ab58SBarry Smith
25b2f1ab58SBarry Smith /* ICC(0) without matrix ordering: simply copies fill pattern */
26b2f1ab58SBarry Smith if (!levels && perm_identity) {
279566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(am + 1, &ui));
28b2f1ab58SBarry Smith ui[0] = 0;
29b2f1ab58SBarry Smith
30*421480d9SBarry Smith for (i = 0; i < am; i++) ui[i + 1] = ui[i] + ai[i + 1] - adiag[i];
319566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ui[am] + 1, &uj));
32b2f1ab58SBarry Smith cols = uj;
33b2f1ab58SBarry Smith for (i = 0; i < am; i++) {
34*421480d9SBarry Smith aj = a->j + adiag[i];
35b2f1ab58SBarry Smith ncols = ui[i + 1] - ui[i];
36b2f1ab58SBarry Smith for (j = 0; j < ncols; j++) *cols++ = *aj++;
37b2f1ab58SBarry Smith }
38b2f1ab58SBarry Smith } else { /* case: levels>0 || (levels=0 && !perm_identity) */
399566063dSJacob Faibussowitsch PetscCall(ISGetIndices(iperm, &riip));
409566063dSJacob Faibussowitsch PetscCall(ISGetIndices(perm, &rip));
41b2f1ab58SBarry Smith
4271d55d18SBarry Smith /* Create spbas_matrix for pattern */
439566063dSJacob Faibussowitsch PetscCall(spbas_pattern_only(am, am, ai, aj, &Pattern_0));
44b2f1ab58SBarry Smith
4571d55d18SBarry Smith /* Apply the permutation */
469566063dSJacob Faibussowitsch PetscCall(spbas_apply_reordering(&Pattern_0, rip, riip));
47b2f1ab58SBarry Smith
4871d55d18SBarry Smith /* Raise the power */
499566063dSJacob Faibussowitsch PetscCall(spbas_power(Pattern_0, (int)levels + 1, &Pattern_P));
509566063dSJacob Faibussowitsch PetscCall(spbas_delete(Pattern_0));
51b2f1ab58SBarry Smith
5271d55d18SBarry Smith /* Keep only upper triangle of pattern */
539566063dSJacob Faibussowitsch PetscCall(spbas_keep_upper(&Pattern_P));
54b2f1ab58SBarry Smith
5571d55d18SBarry Smith /* Convert to Sparse Row Storage */
569566063dSJacob Faibussowitsch PetscCall(spbas_matrix_to_crs(Pattern_P, NULL, &ui, &uj));
579566063dSJacob Faibussowitsch PetscCall(spbas_delete(Pattern_P));
58b2f1ab58SBarry Smith } /* end of case: levels>0 || (levels=0 && !perm_identity) */
59b2f1ab58SBarry Smith
60b2f1ab58SBarry Smith /* put together the new matrix in MATSEQSBAIJ format */
61b2f1ab58SBarry Smith
6257508eceSPierre Jolivet b = (Mat_SeqSBAIJ *)fact->data;
6384648c2dSPierre Jolivet PetscCall(PetscMalloc1(ui[am], &b->a));
642205254eSKarl Rupp
65b2f1ab58SBarry Smith b->j = uj;
66b2f1ab58SBarry Smith b->i = ui;
67f4259b30SLisandro Dalcin b->diag = NULL;
68f4259b30SLisandro Dalcin b->ilen = NULL;
69f4259b30SLisandro Dalcin b->imax = NULL;
70b2f1ab58SBarry Smith b->row = perm;
71b2f1ab58SBarry Smith b->col = perm;
722205254eSKarl Rupp
739566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)perm));
749566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)perm));
752205254eSKarl Rupp
76b2f1ab58SBarry Smith b->icol = iperm;
77b2f1ab58SBarry Smith b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
7884648c2dSPierre Jolivet PetscCall(PetscMalloc1(am, &b->solve_work));
79b2f1ab58SBarry Smith b->maxnz = b->nz = ui[am];
80b2f1ab58SBarry Smith b->free_a = PETSC_TRUE;
81b2f1ab58SBarry Smith b->free_ij = PETSC_TRUE;
82b2f1ab58SBarry Smith
8357508eceSPierre Jolivet fact->info.factor_mallocs = reallocs;
8457508eceSPierre Jolivet fact->info.fill_ratio_given = fill;
85b2f1ab58SBarry Smith if (ai[am] != 0) {
8657508eceSPierre Jolivet fact->info.fill_ratio_needed = (PetscReal)ui[am] / (PetscReal)ai[am];
87b2f1ab58SBarry Smith } else {
8857508eceSPierre Jolivet fact->info.fill_ratio_needed = 0.0;
89b2f1ab58SBarry Smith }
9057508eceSPierre Jolivet /* fact->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ_inplace; */
913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
92b2f1ab58SBarry Smith }
93b2f1ab58SBarry Smith
MatCholeskyFactorNumeric_SeqAIJ_Bas(Mat B,Mat A,const MatFactorInfo * info)9466976f2fSJacob Faibussowitsch static PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_Bas(Mat B, Mat A, const MatFactorInfo *info)
95d71ae5a4SJacob Faibussowitsch {
96b2f1ab58SBarry Smith Mat C = B;
97b2f1ab58SBarry Smith Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ *)C->data;
98b2f1ab58SBarry Smith IS ip = b->row, iip = b->icol;
99b2f1ab58SBarry Smith const PetscInt *rip, *riip;
100b2f1ab58SBarry Smith PetscInt mbs = A->rmap->n, *bi = b->i, *bj = b->j;
101b2f1ab58SBarry Smith MatScalar *ba = b->a;
102c88783f4SHong Zhang PetscReal shiftnz = info->shiftamount;
103d36aa34eSBarry Smith PetscReal droptol = -1;
104ace3abfcSBarry Smith PetscBool perm_identity;
105b2f1ab58SBarry Smith spbas_matrix Pattern, matrix_L, matrix_LT;
1069767453cSBarry Smith PetscReal mem_reduction;
107b2f1ab58SBarry Smith
108b2f1ab58SBarry Smith PetscFunctionBegin;
10971d55d18SBarry Smith /* Reduce memory requirements: erase values of B-matrix */
1109566063dSJacob Faibussowitsch PetscCall(PetscFree(ba));
11171d55d18SBarry Smith /* Compress (maximum) sparseness pattern of B-matrix */
1129566063dSJacob Faibussowitsch PetscCall(spbas_compress_pattern(bi, bj, mbs, mbs, SPBAS_DIAGONAL_OFFSETS, &Pattern, &mem_reduction));
1139566063dSJacob Faibussowitsch PetscCall(PetscFree(bi));
1149566063dSJacob Faibussowitsch PetscCall(PetscFree(bj));
115b2f1ab58SBarry Smith
1169566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, " compression rate for spbas_compress_pattern %g \n", (double)mem_reduction));
117b2f1ab58SBarry Smith
1184e1ff37aSBarry Smith /* Make Cholesky decompositions with larger Manteuffel shifts until no more negative diagonals are found. */
1199566063dSJacob Faibussowitsch PetscCall(ISGetIndices(ip, &rip));
1209566063dSJacob Faibussowitsch PetscCall(ISGetIndices(iip, &riip));
121b2f1ab58SBarry Smith
1225f80ce2aSJacob Faibussowitsch if (info->usedt) droptol = info->dt;
1235f80ce2aSJacob Faibussowitsch
1243ba16761SJacob Faibussowitsch for (int ierr = NEGATIVE_DIAGONAL; ierr == NEGATIVE_DIAGONAL;) {
125cc442ecdSBarry Smith PetscBool success;
1265f80ce2aSJacob Faibussowitsch
1273ba16761SJacob Faibussowitsch ierr = (int)spbas_incomplete_cholesky(A, rip, riip, Pattern, droptol, shiftnz, &matrix_LT, &success);
128cc442ecdSBarry Smith if (!success) {
129b2f1ab58SBarry Smith shiftnz *= 1.5;
1309767453cSBarry Smith if (shiftnz < 1e-5) shiftnz = 1e-5;
1319566063dSJacob Faibussowitsch PetscCall(PetscInfo(NULL, "spbas_incomplete_cholesky found a negative diagonal. Trying again with Manteuffel shift=%g\n", (double)shiftnz));
132b2f1ab58SBarry Smith }
133b2f1ab58SBarry Smith }
1349566063dSJacob Faibussowitsch PetscCall(spbas_delete(Pattern));
135b2f1ab58SBarry Smith
136835f2295SStefano Zampini PetscCall(PetscInfo(NULL, " memory_usage for spbas_incomplete_cholesky %g bytes per row\n", (double)(spbas_memory_requirement(matrix_LT) / (PetscReal)mbs)));
137b2f1ab58SBarry Smith
1389566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(ip, &rip));
1399566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(iip, &riip));
140b2f1ab58SBarry Smith
14171d55d18SBarry Smith /* Convert spbas_matrix to compressed row storage */
1429566063dSJacob Faibussowitsch PetscCall(spbas_transpose(matrix_LT, &matrix_L));
1439566063dSJacob Faibussowitsch PetscCall(spbas_delete(matrix_LT));
1449566063dSJacob Faibussowitsch PetscCall(spbas_matrix_to_crs(matrix_L, &ba, &bi, &bj));
1459371c9d4SSatish Balay b->i = bi;
1469371c9d4SSatish Balay b->j = bj;
1479371c9d4SSatish Balay b->a = ba;
1489566063dSJacob Faibussowitsch PetscCall(spbas_delete(matrix_L));
149b2f1ab58SBarry Smith
15071d55d18SBarry Smith /* Set the appropriate solution functions */
1519566063dSJacob Faibussowitsch PetscCall(ISIdentity(ip, &perm_identity));
152b2f1ab58SBarry Smith if (perm_identity) {
15357508eceSPierre Jolivet B->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
15457508eceSPierre Jolivet B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
15557508eceSPierre Jolivet B->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
15657508eceSPierre Jolivet B->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
157b2f1ab58SBarry Smith } else {
15857508eceSPierre Jolivet B->ops->solve = MatSolve_SeqSBAIJ_1_inplace;
15957508eceSPierre Jolivet B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_inplace;
16057508eceSPierre Jolivet B->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_1_inplace;
16157508eceSPierre Jolivet B->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_1_inplace;
162b2f1ab58SBarry Smith }
163b2f1ab58SBarry Smith
164b2f1ab58SBarry Smith C->assembled = PETSC_TRUE;
165b2f1ab58SBarry Smith C->preallocated = PETSC_TRUE;
1662205254eSKarl Rupp
1679566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(C->rmap->n));
1683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
169b2f1ab58SBarry Smith }
170b2f1ab58SBarry Smith
MatFactorGetSolverType_seqaij_bas(Mat A,MatSolverType * type)17166976f2fSJacob Faibussowitsch static PetscErrorCode MatFactorGetSolverType_seqaij_bas(Mat A, MatSolverType *type)
172d71ae5a4SJacob Faibussowitsch {
1733dfa2556SBarry Smith PetscFunctionBegin;
1743dfa2556SBarry Smith *type = MATSOLVERBAS;
1753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1763dfa2556SBarry Smith }
1773dfa2556SBarry Smith
MatGetFactor_seqaij_bas(Mat A,MatFactorType ftype,Mat * B)178d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_bas(Mat A, MatFactorType ftype, Mat *B)
179d71ae5a4SJacob Faibussowitsch {
180b2f1ab58SBarry Smith PetscInt n = A->rmap->n;
181b2f1ab58SBarry Smith
182b2f1ab58SBarry Smith PetscFunctionBegin;
1839566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), B));
1849566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*B, n, n, n, n));
185966bd95aSPierre Jolivet PetscCheck(ftype == MAT_FACTOR_ICC, PETSC_COMM_SELF, PETSC_ERR_SUP, "Factor type not supported");
1869566063dSJacob Faibussowitsch PetscCall(MatSetType(*B, MATSEQSBAIJ));
1879566063dSJacob Faibussowitsch PetscCall(MatSeqSBAIJSetPreallocation(*B, 1, MAT_SKIP_ALLOCATION, NULL));
1882205254eSKarl Rupp
189b2f1ab58SBarry Smith (*B)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqAIJ_Bas;
190b2f1ab58SBarry Smith (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ_Bas;
1919566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)*B, "MatFactorGetSolverType_C", MatFactorGetSolverType_seqaij_bas));
1929566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATORDERINGND, (char **)&(*B)->preferredordering[MAT_FACTOR_LU]));
1939566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATORDERINGND, (char **)&(*B)->preferredordering[MAT_FACTOR_CHOLESKY]));
194d5f3da31SBarry Smith (*B)->factortype = ftype;
19500c67f3bSHong Zhang
1969566063dSJacob Faibussowitsch PetscCall(PetscFree((*B)->solvertype));
1979566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATSOLVERBAS, &(*B)->solvertype));
198f73b0415SBarry Smith (*B)->canuseordering = PETSC_TRUE;
1999566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATORDERINGNATURAL, (char **)&(*B)->preferredordering[MAT_FACTOR_ICC]));
2003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
201b2f1ab58SBarry Smith }
202