xref: /petsc/src/mat/utils/zerodiag.c (revision dba47a550923b04c7c4ebbb735eb62a1b3e4e9ae)
1 #define PETSCMAT_DLL
2 
3 /*
4     This file contains routines to reorder a matrix so that the diagonal
5     elements are nonzero.
6  */
7 
8 #include "src/mat/matimpl.h"       /*I  "petscmat.h"  I*/
9 
10 #define SWAP(a,b) {PetscInt _t; _t = a; a = b; b = _t; }
11 
12 #undef __FUNCT__
13 #define __FUNCT__ "MatReorderForNonzeroDiagonal"
14 /*@
15     MatReorderForNonzeroDiagonal - Changes matrix ordering to remove
16     zeros from diagonal. This may help in the LU factorization to
17     prevent a zero pivot.
18 
19     Collective on Mat
20 
21     Input Parameters:
22 +   mat  - matrix to reorder
23 -   rmap,cmap - row and column permutations.  Usually obtained from
24                MatGetOrdering().
25 
26     Level: intermediate
27 
28     Notes:
29     This is not intended as a replacement for pivoting for matrices that
30     have ``bad'' structure. It is only a stop-gap measure. Should be called
31     after a call to MatGetOrdering(), this routine changes the column
32     ordering defined in cis.
33 
34     Only works for SeqAIJ matrices
35 
36     Options Database Keys (When using KSP):
37 +      -pc_ilu_nonzeros_along_diagonal
38 -      -pc_lu_nonzeros_along_diagonal
39 
40     Algorithm Notes:
41     Column pivoting is used.
42 
43     1) Choice of column is made by looking at the
44        non-zero elements in the troublesome row for columns that are not yet
45        included (moving from left to right).
46 
47     2) If (1) fails we check all the columns to the left of the current row
48        and see if one of them has could be swapped. It can be swapped if
49        its corresponding row has a non-zero in the column it is being
50        swapped with; to make sure the previous nonzero diagonal remains
51        nonzero
52 
53 
54 @*/
55 PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal(Mat mat,PetscReal abstol,IS ris,IS cis)
56 {
57   PetscErrorCode ierr,(*f)(Mat,PetscReal,IS,IS);
58 
59   PetscFunctionBegin;
60   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatReorderForNonzeroDiagonal_C",(void (**)(void))&f);CHKERRQ(ierr);
61   if (f) {
62     ierr = (*f)(mat,abstol,ris,cis);CHKERRQ(ierr);
63   }
64   PetscFunctionReturn(0);
65 }
66 
67 EXTERN PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
68 EXTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
69 
70 EXTERN_C_BEGIN
71 #undef __FUNCT__
72 #define __FUNCT__ "MatReorderForNonzeroDiagonal_SeqAIJ"
73 PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal_SeqAIJ(Mat mat,PetscReal abstol,IS ris,IS cis)
74 {
75   PetscErrorCode ierr;
76   PetscInt       prow,k,nz,n,repl,*j,*col,*row,m,*icol,nnz,*jj,kk;
77   PetscScalar    *v,*vv;
78   PetscReal      repla;
79   IS             icis;
80 
81   PetscFunctionBegin;
82   ierr = ISGetIndices(ris,&row);CHKERRQ(ierr);
83   ierr = ISGetIndices(cis,&col);CHKERRQ(ierr);
84   ierr = ISInvertPermutation(cis,PETSC_DECIDE,&icis);CHKERRQ(ierr);
85   ierr = ISGetIndices(icis,&icol);CHKERRQ(ierr);
86   ierr = MatGetSize(mat,&m,&n);CHKERRQ(ierr);
87 
88   for (prow=0; prow<n; prow++) {
89     ierr = MatGetRow_SeqAIJ(mat,row[prow],&nz,&j,&v);CHKERRQ(ierr);
90     for (k=0; k<nz; k++) {if (icol[j[k]] == prow) break;}
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 = ISRestoreIndices(ris,&row);CHKERRQ(ierr);
148   ierr = ISRestoreIndices(cis,&col);CHKERRQ(ierr);
149   ierr = ISRestoreIndices(icis,&icol);CHKERRQ(ierr);
150   ierr = ISDestroy(icis);CHKERRQ(ierr);
151   PetscFunctionReturn(0);
152 }
153 EXTERN_C_END
154 
155 
156