xref: /petsc/src/mat/utils/zerodiag.c (revision e8e8640d1cb9a3a2f50c0c0d7b26e5c4d521e2e4)
1ad608de0SBarry Smith /*
2ad608de0SBarry Smith     This file contains routines to reorder a matrix so that the diagonal
348b35521SBarry Smith     elements are nonzero.
4ad608de0SBarry Smith  */
5ad608de0SBarry Smith 
6af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I  "petscmat.h"  I*/
7ad608de0SBarry Smith 
89371c9d4SSatish Balay #define SWAP(a, b) \
9*a8f51744SPierre Jolivet   do { \
109371c9d4SSatish Balay     PetscInt _t; \
119371c9d4SSatish Balay     _t = a; \
129371c9d4SSatish Balay     a  = b; \
139371c9d4SSatish Balay     b  = _t; \
14*a8f51744SPierre Jolivet   } while (0)
1548b35521SBarry Smith 
16ad608de0SBarry Smith /*@
1748b35521SBarry Smith   MatReorderForNonzeroDiagonal - Changes matrix ordering to remove
1811a5261eSBarry Smith   zeros from diagonal. This may help in the `PCLU` factorization to
1948b35521SBarry Smith   prevent a zero pivot.
20ad608de0SBarry Smith 
21c3339decSBarry Smith   Collective
22fee21e36SBarry Smith 
2398a79cdbSBarry Smith   Input Parameters:
2498a79cdbSBarry Smith + mat    - matrix to reorder
252ef1f0ffSBarry Smith . abstol - absolute tolerance, it attempts to move all values smaller off the diagonal
262ef1f0ffSBarry Smith . ris    - the row reordering
272ef1f0ffSBarry Smith - cis    - the column reordering; this may be changed
2898a79cdbSBarry Smith 
2915091d37SBarry Smith   Level: intermediate
3015091d37SBarry Smith 
3127430b45SBarry Smith   Options Database Key:
3227430b45SBarry Smith . -pc_factor_nonzeros_along_diagonal - Reorder to remove zeros from diagonal
3327430b45SBarry Smith 
34ad608de0SBarry Smith   Notes:
35ad608de0SBarry Smith   This is not intended as a replacement for pivoting for matrices that
362ef1f0ffSBarry Smith   have ``bad'' structure. It is only a stop-gap measure.
372ef1f0ffSBarry Smith 
382ef1f0ffSBarry Smith   Should be called
392ef1f0ffSBarry Smith   after a call to `MatGetOrdering()`.
40d252947aSBarry Smith 
4111a5261eSBarry Smith   Only works for `MATSEQAIJ` matrices
42106f7b34SBarry Smith 
4327430b45SBarry Smith   Developer Notes:
4433f51a72SBarry Smith   Column pivoting is used.
4533f51a72SBarry Smith 
4633f51a72SBarry Smith   1) Choice of column is made by looking at the
4733f51a72SBarry Smith   non-zero elements in the troublesome row for columns that are not yet
4833f51a72SBarry Smith   included (moving from left to right).
4933f51a72SBarry Smith 
5033f51a72SBarry Smith   2) If (1) fails we check all the columns to the left of the current row
51814823dcSBarry Smith   and see if one of them has could be swapped. It can be swapped if
52814823dcSBarry Smith   its corresponding row has a non-zero in the column it is being
53814823dcSBarry Smith   swapped with; to make sure the previous nonzero diagonal remains
54814823dcSBarry Smith   nonzero
5533f51a72SBarry Smith 
562ef1f0ffSBarry Smith .seealso: `Mat`, `MatGetFactor()`, `MatGetOrdering()`
57ad608de0SBarry Smith @*/
MatReorderForNonzeroDiagonal(Mat mat,PetscReal abstol,IS ris,IS cis)58d71ae5a4SJacob Faibussowitsch PetscErrorCode MatReorderForNonzeroDiagonal(Mat mat, PetscReal abstol, IS ris, IS cis)
59d71ae5a4SJacob Faibussowitsch {
6005b94e36SKris Buschelman   PetscFunctionBegin;
61cac4c232SBarry Smith   PetscTryMethod(mat, "MatReorderForNonzeroDiagonal_C", (Mat, PetscReal, IS, IS), (mat, abstol, ris, cis));
623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6305b94e36SKris Buschelman }
6405b94e36SKris Buschelman 
655a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetRow_SeqAIJ(Mat, PetscInt, PetscInt *, PetscInt **, PetscScalar **);
665a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat, PetscInt, PetscInt *, PetscInt **, PetscScalar **);
67b3cc6726SBarry Smith 
68d9e81513SBarry Smith #include <../src/vec/is/is/impls/general/general.h>
695d0c19d7SBarry Smith 
MatReorderForNonzeroDiagonal_SeqAIJ(Mat mat,PetscReal abstol,IS ris,IS cis)70d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatReorderForNonzeroDiagonal_SeqAIJ(Mat mat, PetscReal abstol, IS ris, IS cis)
71d71ae5a4SJacob Faibussowitsch {
7238baddfdSBarry Smith   PetscInt     prow, k, nz, n, repl, *j, *col, *row, m, *icol, nnz, *jj, kk;
7387828ca2SBarry Smith   PetscScalar *v, *vv;
74329f5518SBarry Smith   PetscReal    repla;
75fd61f322SBarry Smith   IS           icis;
76ad608de0SBarry Smith 
773a40ed3dSBarry Smith   PetscFunctionBegin;
785d0c19d7SBarry Smith   /* access the indices of the IS directly, because it changes them */
795d0c19d7SBarry Smith   row = ((IS_General *)ris->data)->idx;
805d0c19d7SBarry Smith   col = ((IS_General *)cis->data)->idx;
819566063dSJacob Faibussowitsch   PetscCall(ISInvertPermutation(cis, PETSC_DECIDE, &icis));
825d0c19d7SBarry Smith   icol = ((IS_General *)icis->data)->idx;
839566063dSJacob Faibussowitsch   PetscCall(MatGetSize(mat, &m, &n));
84ad608de0SBarry Smith 
85ad608de0SBarry Smith   for (prow = 0; prow < n; prow++) {
869566063dSJacob Faibussowitsch     PetscCall(MatGetRow_SeqAIJ(mat, row[prow], &nz, &j, &v));
878865f1eaSKarl Rupp     for (k = 0; k < nz; k++) {
888865f1eaSKarl Rupp       if (icol[j[k]] == prow) break;
898865f1eaSKarl Rupp     }
9070441072SBarry Smith     if (k >= nz || PetscAbsScalar(v[k]) <= abstol) {
91ad608de0SBarry Smith       /* Element too small or zero; find the best candidate */
92cddf8d76SBarry Smith       repla = (k >= nz) ? 0.0 : PetscAbsScalar(v[k]);
93d4bb536fSBarry Smith       /*
94fd61f322SBarry Smith           Look for a later column we can swap with this one
95d4bb536fSBarry Smith       */
9648b35521SBarry Smith       for (k = 0; k < nz; k++) {
9733f51a72SBarry Smith         if (icol[j[k]] > prow && PetscAbsScalar(v[k]) > repla) {
98fd61f322SBarry Smith           /* found a suitable later column */
9933f51a72SBarry Smith           repl = icol[j[k]];
10033f51a72SBarry Smith           SWAP(icol[col[prow]], icol[col[repl]]);
101ad608de0SBarry Smith           SWAP(col[prow], col[repl]);
102fd61f322SBarry Smith           goto found;
103fd61f322SBarry Smith         }
104fd61f322SBarry Smith       }
105fd61f322SBarry Smith       /*
106fd61f322SBarry Smith            Did not find a suitable later column so look for an earlier column
107fd61f322SBarry Smith            We need to be sure that we don't introduce a zero in a previous
108fd61f322SBarry Smith            diagonal
109fd61f322SBarry Smith       */
110fd61f322SBarry Smith       for (k = 0; k < nz; k++) {
111fd61f322SBarry Smith         if (icol[j[k]] < prow && PetscAbsScalar(v[k]) > repla) {
112fd61f322SBarry Smith           /* See if this one will work */
113fd61f322SBarry Smith           repl = icol[j[k]];
1149566063dSJacob Faibussowitsch           PetscCall(MatGetRow_SeqAIJ(mat, row[repl], &nnz, &jj, &vv));
115fd61f322SBarry Smith           for (kk = 0; kk < nnz; kk++) {
11670441072SBarry Smith             if (icol[jj[kk]] == prow && PetscAbsScalar(vv[kk]) > abstol) {
1179566063dSJacob Faibussowitsch               PetscCall(MatRestoreRow_SeqAIJ(mat, row[repl], &nnz, &jj, &vv));
118fd61f322SBarry Smith               SWAP(icol[col[prow]], icol[col[repl]]);
119fd61f322SBarry Smith               SWAP(col[prow], col[repl]);
120fd61f322SBarry Smith               goto found;
121fd61f322SBarry Smith             }
122fd61f322SBarry Smith           }
1239566063dSJacob Faibussowitsch           PetscCall(MatRestoreRow_SeqAIJ(mat, row[repl], &nnz, &jj, &vv));
124fd61f322SBarry Smith         }
125fd61f322SBarry Smith       }
126fd61f322SBarry Smith       /*
127fd61f322SBarry Smith           No column  suitable; instead check all future rows
128fd61f322SBarry Smith           Note: this will be very slow
129fd61f322SBarry Smith       */
130fd61f322SBarry Smith       for (k = prow + 1; k < n; k++) {
1319566063dSJacob Faibussowitsch         PetscCall(MatGetRow_SeqAIJ(mat, row[k], &nnz, &jj, &vv));
132fd61f322SBarry Smith         for (kk = 0; kk < nnz; kk++) {
13370441072SBarry Smith           if (icol[jj[kk]] == prow && PetscAbsScalar(vv[kk]) > abstol) {
134fd61f322SBarry Smith             /* found a row */
135fd61f322SBarry Smith             SWAP(row[prow], row[k]);
136fd61f322SBarry Smith             goto found;
137fd61f322SBarry Smith           }
138fd61f322SBarry Smith         }
1399566063dSJacob Faibussowitsch         PetscCall(MatRestoreRow_SeqAIJ(mat, row[k], &nnz, &jj, &vv));
140fd61f322SBarry Smith       }
141fd61f322SBarry Smith 
142fd61f322SBarry Smith     found:;
143ad608de0SBarry Smith     }
1449566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow_SeqAIJ(mat, row[prow], &nz, &j, &v));
145ad608de0SBarry Smith   }
1469566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&icis));
1473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
148ad608de0SBarry Smith }
149