1be1d678aSKris Buschelman 2ad608de0SBarry Smith /* 3ad608de0SBarry Smith This file contains routines to reorder a matrix so that the diagonal 448b35521SBarry Smith elements are nonzero. 5ad608de0SBarry Smith */ 6ad608de0SBarry Smith 7af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 8ad608de0SBarry Smith 99371c9d4SSatish Balay #define SWAP(a, b) \ 109371c9d4SSatish Balay { \ 119371c9d4SSatish Balay PetscInt _t; \ 129371c9d4SSatish Balay _t = a; \ 139371c9d4SSatish Balay a = b; \ 149371c9d4SSatish Balay b = _t; \ 159371c9d4SSatish Balay } 1648b35521SBarry Smith 17ad608de0SBarry Smith /*@ 1848b35521SBarry Smith MatReorderForNonzeroDiagonal - Changes matrix ordering to remove 1911a5261eSBarry Smith zeros from diagonal. This may help in the `PCLU` factorization to 2048b35521SBarry Smith prevent a zero pivot. 21ad608de0SBarry Smith 22*c3339decSBarry Smith Collective 23fee21e36SBarry Smith 2498a79cdbSBarry Smith Input Parameters: 2598a79cdbSBarry Smith + mat - matrix to reorder 2698a79cdbSBarry Smith - rmap,cmap - row and column permutations. Usually obtained from 2711a5261eSBarry Smith `MatGetOrdering()`. 2898a79cdbSBarry Smith 2915091d37SBarry Smith Level: intermediate 3015091d37SBarry Smith 31ad608de0SBarry Smith Notes: 32ad608de0SBarry Smith This is not intended as a replacement for pivoting for matrices that 33d252947aSBarry Smith have ``bad'' structure. It is only a stop-gap measure. Should be called 3411a5261eSBarry Smith after a call to `MatGetOrdering()`, this routine changes the column 35d252947aSBarry Smith ordering defined in cis. 36d252947aSBarry Smith 3711a5261eSBarry Smith Only works for `MATSEQAIJ` matrices 38106f7b34SBarry Smith 3911a5261eSBarry Smith Options Database Keys (When using `KSP`): 4067b8a455SSatish Balay . -pc_factor_nonzeros_along_diagonal - Reorder to remove zeros from diagonal 41ad608de0SBarry Smith 4233f51a72SBarry Smith Algorithm Notes: 4333f51a72SBarry Smith Column pivoting is used. 4433f51a72SBarry Smith 4533f51a72SBarry Smith 1) Choice of column is made by looking at the 4633f51a72SBarry Smith non-zero elements in the troublesome row for columns that are not yet 4733f51a72SBarry Smith included (moving from left to right). 4833f51a72SBarry Smith 4933f51a72SBarry Smith 2) If (1) fails we check all the columns to the left of the current row 50814823dcSBarry Smith and see if one of them has could be swapped. It can be swapped if 51814823dcSBarry Smith its corresponding row has a non-zero in the column it is being 52814823dcSBarry Smith swapped with; to make sure the previous nonzero diagonal remains 53814823dcSBarry Smith nonzero 5433f51a72SBarry Smith 55ad608de0SBarry Smith @*/ 56d71ae5a4SJacob Faibussowitsch PetscErrorCode MatReorderForNonzeroDiagonal(Mat mat, PetscReal abstol, IS ris, IS cis) 57d71ae5a4SJacob Faibussowitsch { 5805b94e36SKris Buschelman PetscFunctionBegin; 59cac4c232SBarry Smith PetscTryMethod(mat, "MatReorderForNonzeroDiagonal_C", (Mat, PetscReal, IS, IS), (mat, abstol, ris, cis)); 6005b94e36SKris Buschelman PetscFunctionReturn(0); 6105b94e36SKris Buschelman } 6205b94e36SKris Buschelman 635a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetRow_SeqAIJ(Mat, PetscInt, PetscInt *, PetscInt **, PetscScalar **); 645a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat, PetscInt, PetscInt *, PetscInt **, PetscScalar **); 65b3cc6726SBarry Smith 66d9e81513SBarry Smith #include <../src/vec/is/is/impls/general/general.h> 675d0c19d7SBarry Smith 68d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatReorderForNonzeroDiagonal_SeqAIJ(Mat mat, PetscReal abstol, IS ris, IS cis) 69d71ae5a4SJacob Faibussowitsch { 7038baddfdSBarry Smith PetscInt prow, k, nz, n, repl, *j, *col, *row, m, *icol, nnz, *jj, kk; 7187828ca2SBarry Smith PetscScalar *v, *vv; 72329f5518SBarry Smith PetscReal repla; 73fd61f322SBarry Smith IS icis; 74ad608de0SBarry Smith 753a40ed3dSBarry Smith PetscFunctionBegin; 765d0c19d7SBarry Smith /* access the indices of the IS directly, because it changes them */ 775d0c19d7SBarry Smith row = ((IS_General *)ris->data)->idx; 785d0c19d7SBarry Smith col = ((IS_General *)cis->data)->idx; 799566063dSJacob Faibussowitsch PetscCall(ISInvertPermutation(cis, PETSC_DECIDE, &icis)); 805d0c19d7SBarry Smith icol = ((IS_General *)icis->data)->idx; 819566063dSJacob Faibussowitsch PetscCall(MatGetSize(mat, &m, &n)); 82ad608de0SBarry Smith 83ad608de0SBarry Smith for (prow = 0; prow < n; prow++) { 849566063dSJacob Faibussowitsch PetscCall(MatGetRow_SeqAIJ(mat, row[prow], &nz, &j, &v)); 858865f1eaSKarl Rupp for (k = 0; k < nz; k++) { 868865f1eaSKarl Rupp if (icol[j[k]] == prow) break; 878865f1eaSKarl Rupp } 8870441072SBarry Smith if (k >= nz || PetscAbsScalar(v[k]) <= abstol) { 89ad608de0SBarry Smith /* Element too small or zero; find the best candidate */ 90cddf8d76SBarry Smith repla = (k >= nz) ? 0.0 : PetscAbsScalar(v[k]); 91d4bb536fSBarry Smith /* 92fd61f322SBarry Smith Look for a later column we can swap with this one 93d4bb536fSBarry Smith */ 9448b35521SBarry Smith for (k = 0; k < nz; k++) { 9533f51a72SBarry Smith if (icol[j[k]] > prow && PetscAbsScalar(v[k]) > repla) { 96fd61f322SBarry Smith /* found a suitable later column */ 9733f51a72SBarry Smith repl = icol[j[k]]; 9833f51a72SBarry Smith SWAP(icol[col[prow]], icol[col[repl]]); 99ad608de0SBarry Smith SWAP(col[prow], col[repl]); 100fd61f322SBarry Smith goto found; 101fd61f322SBarry Smith } 102fd61f322SBarry Smith } 103fd61f322SBarry Smith /* 104fd61f322SBarry Smith Did not find a suitable later column so look for an earlier column 105fd61f322SBarry Smith We need to be sure that we don't introduce a zero in a previous 106fd61f322SBarry Smith diagonal 107fd61f322SBarry Smith */ 108fd61f322SBarry Smith for (k = 0; k < nz; k++) { 109fd61f322SBarry Smith if (icol[j[k]] < prow && PetscAbsScalar(v[k]) > repla) { 110fd61f322SBarry Smith /* See if this one will work */ 111fd61f322SBarry Smith repl = icol[j[k]]; 1129566063dSJacob Faibussowitsch PetscCall(MatGetRow_SeqAIJ(mat, row[repl], &nnz, &jj, &vv)); 113fd61f322SBarry Smith for (kk = 0; kk < nnz; kk++) { 11470441072SBarry Smith if (icol[jj[kk]] == prow && PetscAbsScalar(vv[kk]) > abstol) { 1159566063dSJacob Faibussowitsch PetscCall(MatRestoreRow_SeqAIJ(mat, row[repl], &nnz, &jj, &vv)); 116fd61f322SBarry Smith SWAP(icol[col[prow]], icol[col[repl]]); 117fd61f322SBarry Smith SWAP(col[prow], col[repl]); 118fd61f322SBarry Smith goto found; 119fd61f322SBarry Smith } 120fd61f322SBarry Smith } 1219566063dSJacob Faibussowitsch PetscCall(MatRestoreRow_SeqAIJ(mat, row[repl], &nnz, &jj, &vv)); 122fd61f322SBarry Smith } 123fd61f322SBarry Smith } 124fd61f322SBarry Smith /* 125fd61f322SBarry Smith No column suitable; instead check all future rows 126fd61f322SBarry Smith Note: this will be very slow 127fd61f322SBarry Smith */ 128fd61f322SBarry Smith for (k = prow + 1; k < n; k++) { 1299566063dSJacob Faibussowitsch PetscCall(MatGetRow_SeqAIJ(mat, row[k], &nnz, &jj, &vv)); 130fd61f322SBarry Smith for (kk = 0; kk < nnz; kk++) { 13170441072SBarry Smith if (icol[jj[kk]] == prow && PetscAbsScalar(vv[kk]) > abstol) { 132fd61f322SBarry Smith /* found a row */ 133fd61f322SBarry Smith SWAP(row[prow], row[k]); 134fd61f322SBarry Smith goto found; 135fd61f322SBarry Smith } 136fd61f322SBarry Smith } 1379566063dSJacob Faibussowitsch PetscCall(MatRestoreRow_SeqAIJ(mat, row[k], &nnz, &jj, &vv)); 138fd61f322SBarry Smith } 139fd61f322SBarry Smith 140fd61f322SBarry Smith found:; 141ad608de0SBarry Smith } 1429566063dSJacob Faibussowitsch PetscCall(MatRestoreRow_SeqAIJ(mat, row[prow], &nz, &j, &v)); 143ad608de0SBarry Smith } 1449566063dSJacob Faibussowitsch PetscCall(ISDestroy(&icis)); 1453a40ed3dSBarry Smith PetscFunctionReturn(0); 146ad608de0SBarry Smith } 147