xref: /petsc/src/mat/utils/zerodiag.c (revision c69effb2e05caa58810396c101d616d7ad9da0a7)
1 
2 /*
3     This file contains routines to reorder a matrix so that the diagonal
4     elements are nonzero.
5  */
6 
7 #include <petsc/private/matimpl.h> /*I  "petscmat.h"  I*/
8 
9 #define SWAP(a, b) \
10   { \
11     PetscInt _t; \
12     _t = a; \
13     a  = b; \
14     b  = _t; \
15   }
16 
17 /*@
18     MatReorderForNonzeroDiagonal - Changes matrix ordering to remove
19     zeros from diagonal. This may help in the `PCLU` factorization to
20     prevent a zero pivot.
21 
22     Collective
23 
24     Input Parameters:
25 +   mat  - matrix to reorder
26 -   rmap,cmap - row and column permutations.  Usually obtained from
27                `MatGetOrdering()`.
28 
29     Level: intermediate
30 
31     Options Database Key:
32 .      -pc_factor_nonzeros_along_diagonal - Reorder to remove zeros from diagonal
33 
34     Notes:
35     This is not intended as a replacement for pivoting for matrices that
36     have ``bad'' structure. It is only a stop-gap measure. Should be called
37     after a call to `MatGetOrdering()`, this routine changes the column
38     ordering defined in cis.
39 
40     Only works for `MATSEQAIJ` matrices
41 
42     Developer Notes:
43     Column pivoting is used.
44 
45     1) Choice of column is made by looking at the
46        non-zero elements in the troublesome row for columns that are not yet
47        included (moving from left to right).
48 
49     2) If (1) fails we check all the columns to the left of the current row
50        and see if one of them has could be swapped. It can be swapped if
51        its corresponding row has a non-zero in the column it is being
52        swapped with; to make sure the previous nonzero diagonal remains
53        nonzero
54 
55 .seealso: `Mat`
56 @*/
57 PetscErrorCode MatReorderForNonzeroDiagonal(Mat mat, PetscReal abstol, IS ris, IS cis)
58 {
59   PetscFunctionBegin;
60   PetscTryMethod(mat, "MatReorderForNonzeroDiagonal_C", (Mat, PetscReal, IS, IS), (mat, abstol, ris, cis));
61   PetscFunctionReturn(PETSC_SUCCESS);
62 }
63 
64 PETSC_INTERN PetscErrorCode MatGetRow_SeqAIJ(Mat, PetscInt, PetscInt *, PetscInt **, PetscScalar **);
65 PETSC_INTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat, PetscInt, PetscInt *, PetscInt **, PetscScalar **);
66 
67 #include <../src/vec/is/is/impls/general/general.h>
68 
69 PETSC_INTERN PetscErrorCode MatReorderForNonzeroDiagonal_SeqAIJ(Mat mat, PetscReal abstol, IS ris, IS cis)
70 {
71   PetscInt     prow, k, nz, n, repl, *j, *col, *row, m, *icol, nnz, *jj, kk;
72   PetscScalar *v, *vv;
73   PetscReal    repla;
74   IS           icis;
75 
76   PetscFunctionBegin;
77   /* access the indices of the IS directly, because it changes them */
78   row = ((IS_General *)ris->data)->idx;
79   col = ((IS_General *)cis->data)->idx;
80   PetscCall(ISInvertPermutation(cis, PETSC_DECIDE, &icis));
81   icol = ((IS_General *)icis->data)->idx;
82   PetscCall(MatGetSize(mat, &m, &n));
83 
84   for (prow = 0; prow < n; prow++) {
85     PetscCall(MatGetRow_SeqAIJ(mat, row[prow], &nz, &j, &v));
86     for (k = 0; k < nz; k++) {
87       if (icol[j[k]] == prow) break;
88     }
89     if (k >= nz || PetscAbsScalar(v[k]) <= abstol) {
90       /* Element too small or zero; find the best candidate */
91       repla = (k >= nz) ? 0.0 : PetscAbsScalar(v[k]);
92       /*
93           Look for a later column we can swap with this one
94       */
95       for (k = 0; k < nz; k++) {
96         if (icol[j[k]] > prow && PetscAbsScalar(v[k]) > repla) {
97           /* found a suitable later column */
98           repl = icol[j[k]];
99           SWAP(icol[col[prow]], icol[col[repl]]);
100           SWAP(col[prow], col[repl]);
101           goto found;
102         }
103       }
104       /*
105            Did not find a suitable later column so look for an earlier column
106            We need to be sure that we don't introduce a zero in a previous
107            diagonal
108       */
109       for (k = 0; k < nz; k++) {
110         if (icol[j[k]] < prow && PetscAbsScalar(v[k]) > repla) {
111           /* See if this one will work */
112           repl = icol[j[k]];
113           PetscCall(MatGetRow_SeqAIJ(mat, row[repl], &nnz, &jj, &vv));
114           for (kk = 0; kk < nnz; kk++) {
115             if (icol[jj[kk]] == prow && PetscAbsScalar(vv[kk]) > abstol) {
116               PetscCall(MatRestoreRow_SeqAIJ(mat, row[repl], &nnz, &jj, &vv));
117               SWAP(icol[col[prow]], icol[col[repl]]);
118               SWAP(col[prow], col[repl]);
119               goto found;
120             }
121           }
122           PetscCall(MatRestoreRow_SeqAIJ(mat, row[repl], &nnz, &jj, &vv));
123         }
124       }
125       /*
126           No column  suitable; instead check all future rows
127           Note: this will be very slow
128       */
129       for (k = prow + 1; k < n; k++) {
130         PetscCall(MatGetRow_SeqAIJ(mat, row[k], &nnz, &jj, &vv));
131         for (kk = 0; kk < nnz; kk++) {
132           if (icol[jj[kk]] == prow && PetscAbsScalar(vv[kk]) > abstol) {
133             /* found a row */
134             SWAP(row[prow], row[k]);
135             goto found;
136           }
137         }
138         PetscCall(MatRestoreRow_SeqAIJ(mat, row[k], &nnz, &jj, &vv));
139       }
140 
141     found:;
142     }
143     PetscCall(MatRestoreRow_SeqAIJ(mat, row[prow], &nz, &j, &v));
144   }
145   PetscCall(ISDestroy(&icis));
146   PetscFunctionReturn(PETSC_SUCCESS);
147 }
148