xref: /petsc/src/mat/utils/zerodiag.c (revision c877dd1b0ef103f8a5e1c6faaa97a532fd33b12d)
1 /*
2     This file contains routines to reorder a matrix so that the diagonal
3     elements are nonzero.
4  */
5 
6 #include "src/mat/matimpl.h"       /*I  "petscmat.h"  I*/
7 
8 #define SWAP(a,b) {PetscInt _t; _t = a; a = b; b = _t; }
9 
10 #undef __FUNCT__
11 #define __FUNCT__ "MatReorderForNonzeroDiagonal"
12 /*@
13     MatReorderForNonzeroDiagonal - Changes matrix ordering to remove
14     zeros from diagonal. This may help in the LU factorization to
15     prevent a zero pivot.
16 
17     Collective on Mat
18 
19     Input Parameters:
20 +   mat  - matrix to reorder
21 -   rmap,cmap - row and column permutations.  Usually obtained from
22                MatGetOrdering().
23 
24     Level: intermediate
25 
26     Notes:
27     This is not intended as a replacement for pivoting for matrices that
28     have ``bad'' structure. It is only a stop-gap measure. Should be called
29     after a call to MatGetOrdering(), this routine changes the column
30     ordering defined in cis.
31 
32     Only works for SeqAIJ matrices
33 
34     Options Database Keys (When using KSP):
35 +      -pc_ilu_nonzeros_along_diagonal
36 -      -pc_lu_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,(*f)(Mat,PetscReal,IS,IS);
56 
57   PetscFunctionBegin;
58   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatReorderForNonzeroDiagonal_C",(void (**)(void))&f);CHKERRQ(ierr);
59   if (f) {
60     ierr = (*f)(mat,abstol,ris,cis);CHKERRQ(ierr);
61   }
62   PetscFunctionReturn(0);
63 }
64 
65 EXTERN PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
66 EXTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
67 
68 EXTERN_C_BEGIN
69 #undef __FUNCT__
70 #define __FUNCT__ "MatReorderForNonzeroDiagonal_SeqAIJ"
71 PetscErrorCode MatReorderForNonzeroDiagonal_SeqAIJ(Mat mat,PetscReal abstol,IS ris,IS cis)
72 {
73   PetscErrorCode ierr;
74   PetscInt       prow,k,nz,n,repl,*j,*col,*row,m,*icol,nnz,*jj,kk;
75   PetscScalar    *v,*vv;
76   PetscReal      repla;
77   IS             icis;
78 
79   PetscFunctionBegin;
80   ierr = ISGetIndices(ris,&row);CHKERRQ(ierr);
81   ierr = ISGetIndices(cis,&col);CHKERRQ(ierr);
82   ierr = ISInvertPermutation(cis,PETSC_DECIDE,&icis);CHKERRQ(ierr);
83   ierr = ISGetIndices(icis,&icol);CHKERRQ(ierr);
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++) {if (icol[j[k]] == prow) break;}
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           ierr = MatGetRow_SeqAIJ(mat,row[repl],&nnz,&jj,&vv);CHKERRQ(ierr);
114           for (kk=0; kk<nnz; kk++) {
115             if (icol[jj[kk]] == prow && PetscAbsScalar(vv[kk]) > abstol) {
116               ierr = MatRestoreRow_SeqAIJ(mat,row[repl],&nnz,&jj,&vv);CHKERRQ(ierr);
117               SWAP(icol[col[prow]],icol[col[repl]]);
118               SWAP(col[prow],col[repl]);
119               goto found;
120 	    }
121           }
122           ierr = MatRestoreRow_SeqAIJ(mat,row[repl],&nnz,&jj,&vv);CHKERRQ(ierr);
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         ierr = MatGetRow_SeqAIJ(mat,row[k],&nnz,&jj,&vv);CHKERRQ(ierr);
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         ierr = MatRestoreRow_SeqAIJ(mat,row[k],&nnz,&jj,&vv);CHKERRQ(ierr);
139       }
140 
141       found:;
142     }
143     ierr = MatRestoreRow_SeqAIJ(mat,row[prow],&nz,&j,&v);CHKERRQ(ierr);
144   }
145   ierr = ISRestoreIndices(ris,&row);CHKERRQ(ierr);
146   ierr = ISRestoreIndices(cis,&col);CHKERRQ(ierr);
147   ierr = ISRestoreIndices(icis,&icol);CHKERRQ(ierr);
148   ierr = ISDestroy(icis);CHKERRQ(ierr);
149   PetscFunctionReturn(0);
150 }
151 EXTERN_C_END
152 
153 
154