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