xref: /petsc/src/mat/utils/zerodiag.c (revision 38baddfdda9fcd85a080f141519fb97fb01bc0a4)
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 
6e090d566SSatish Balay #include "src/mat/matimpl.h"       /*I  "petscmat.h"  I*/
7ad608de0SBarry Smith 
8*38baddfdSBarry Smith #define SWAP(a,b) {PetscInt _t; _t = a; a = b; b = _t; }
948b35521SBarry Smith 
104a2ae208SSatish Balay #undef __FUNCT__
114a2ae208SSatish Balay #define __FUNCT__ "MatReorderForNonzeroDiagonal"
12ad608de0SBarry Smith /*@
1348b35521SBarry Smith     MatReorderForNonzeroDiagonal - Changes matrix ordering to remove
1448b35521SBarry Smith     zeros from diagonal. This may help in the LU factorization to
1548b35521SBarry Smith     prevent a zero pivot.
16ad608de0SBarry Smith 
17fee21e36SBarry Smith     Collective on Mat
18fee21e36SBarry Smith 
1998a79cdbSBarry Smith     Input Parameters:
2098a79cdbSBarry Smith +   mat  - matrix to reorder
2198a79cdbSBarry Smith -   rmap,cmap - row and column permutations.  Usually obtained from
2291e9ee9fSBarry Smith                MatGetOrdering().
2398a79cdbSBarry Smith 
2415091d37SBarry Smith     Level: intermediate
2515091d37SBarry Smith 
26ad608de0SBarry Smith     Notes:
27ad608de0SBarry Smith     This is not intended as a replacement for pivoting for matrices that
28d252947aSBarry Smith     have ``bad'' structure. It is only a stop-gap measure. Should be called
2991e9ee9fSBarry Smith     after a call to MatGetOrdering(), this routine changes the column
30d252947aSBarry Smith     ordering defined in cis.
31d252947aSBarry Smith 
32106f7b34SBarry Smith     Only works for SeqAIJ matrices
33106f7b34SBarry Smith 
3494b7f48cSBarry Smith     Options Database Keys (When using KSP):
3598a79cdbSBarry Smith +      -pc_ilu_nonzeros_along_diagonal
3698a79cdbSBarry Smith -      -pc_lu_nonzeros_along_diagonal
37ad608de0SBarry Smith 
3833f51a72SBarry Smith     Algorithm Notes:
3933f51a72SBarry Smith     Column pivoting is used.
4033f51a72SBarry Smith 
4133f51a72SBarry Smith     1) Choice of column is made by looking at the
4233f51a72SBarry Smith        non-zero elements in the troublesome row for columns that are not yet
4333f51a72SBarry Smith        included (moving from left to right).
4433f51a72SBarry Smith 
4533f51a72SBarry Smith     2) If (1) fails we check all the columns to the left of the current row
46814823dcSBarry Smith        and see if one of them has could be swapped. It can be swapped if
47814823dcSBarry Smith        its corresponding row has a non-zero in the column it is being
48814823dcSBarry Smith        swapped with; to make sure the previous nonzero diagonal remains
49814823dcSBarry Smith        nonzero
5033f51a72SBarry Smith 
5198a79cdbSBarry Smith 
52ad608de0SBarry Smith @*/
53dfbe8321SBarry Smith PetscErrorCode MatReorderForNonzeroDiagonal(Mat mat,PetscReal atol,IS ris,IS cis)
54ad608de0SBarry Smith {
55dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(Mat,PetscReal,IS,IS);
5605b94e36SKris Buschelman 
5705b94e36SKris Buschelman   PetscFunctionBegin;
5805b94e36SKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatReorderForNonzeroDiagonal_C",(void (**)(void))&f);CHKERRQ(ierr);
5905b94e36SKris Buschelman   if (f) {
6005b94e36SKris Buschelman     ierr = (*f)(mat,atol,ris,cis);CHKERRQ(ierr);
6105b94e36SKris Buschelman   }
6205b94e36SKris Buschelman   PetscFunctionReturn(0);
6305b94e36SKris Buschelman }
6405b94e36SKris Buschelman 
65*38baddfdSBarry Smith EXTERN PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
66*38baddfdSBarry Smith EXTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
67b3cc6726SBarry Smith 
6805b94e36SKris Buschelman EXTERN_C_BEGIN
6905b94e36SKris Buschelman #undef __FUNCT__
7005b94e36SKris Buschelman #define __FUNCT__ "MatReorderForNonzeroDiagonal_SeqAIJ"
71dfbe8321SBarry Smith PetscErrorCode MatReorderForNonzeroDiagonal_SeqAIJ(Mat mat,PetscReal atol,IS ris,IS cis)
7205b94e36SKris Buschelman {
73dfbe8321SBarry Smith   PetscErrorCode ierr;
74*38baddfdSBarry Smith   PetscInt       prow,k,nz,n,repl,*j,*col,*row,m,*icol,nnz,*jj,kk;
7587828ca2SBarry Smith   PetscScalar    *v,*vv;
76329f5518SBarry Smith   PetscReal      repla;
77fd61f322SBarry Smith   IS             icis;
78ad608de0SBarry Smith 
793a40ed3dSBarry Smith   PetscFunctionBegin;
8048b35521SBarry Smith   ierr = ISGetIndices(ris,&row);CHKERRQ(ierr);
8148b35521SBarry Smith   ierr = ISGetIndices(cis,&col);CHKERRQ(ierr);
824c49b128SBarry Smith   ierr = ISInvertPermutation(cis,PETSC_DECIDE,&icis);CHKERRQ(ierr);
8333f51a72SBarry Smith   ierr = ISGetIndices(icis,&icol);CHKERRQ(ierr);
8448b35521SBarry Smith   ierr = MatGetSize(mat,&m,&n);CHKERRQ(ierr);
85ad608de0SBarry Smith 
86ad608de0SBarry Smith   for (prow=0; prow<n; prow++) {
87b3cc6726SBarry Smith     ierr = MatGetRow_SeqAIJ(mat,row[prow],&nz,&j,&v);CHKERRQ(ierr);
8833f51a72SBarry Smith     for (k=0; k<nz; k++) {if (icol[j[k]] == prow) break;}
89cddf8d76SBarry Smith     if (k >= nz || PetscAbsScalar(v[k]) <= atol) {
90ad608de0SBarry Smith       /* Element too small or zero; find the best candidate */
91cddf8d76SBarry Smith       repla = (k >= nz) ? 0.0 : PetscAbsScalar(v[k]);
92d4bb536fSBarry Smith       /*
93fd61f322SBarry Smith           Look for a later column we can swap with this one
94d4bb536fSBarry Smith       */
9548b35521SBarry Smith       for (k=0; k<nz; k++) {
9633f51a72SBarry Smith 	if (icol[j[k]] > prow && PetscAbsScalar(v[k]) > repla) {
97fd61f322SBarry Smith           /* found a suitable later column */
9833f51a72SBarry Smith 	  repl  = icol[j[k]];
9933f51a72SBarry Smith           SWAP(icol[col[prow]],icol[col[repl]]);
100ad608de0SBarry Smith           SWAP(col[prow],col[repl]);
101fd61f322SBarry Smith           goto found;
102fd61f322SBarry Smith         }
103fd61f322SBarry Smith       }
104fd61f322SBarry Smith       /*
105fd61f322SBarry Smith            Did not find a suitable later column so look for an earlier column
106fd61f322SBarry Smith 	   We need to be sure that we don't introduce a zero in a previous
107fd61f322SBarry Smith 	   diagonal
108fd61f322SBarry Smith       */
109fd61f322SBarry Smith       for (k=0; k<nz; k++) {
110fd61f322SBarry Smith         if (icol[j[k]] < prow && PetscAbsScalar(v[k]) > repla) {
111fd61f322SBarry Smith           /* See if this one will work */
112fd61f322SBarry Smith           repl  = icol[j[k]];
113b3cc6726SBarry Smith           ierr = MatGetRow_SeqAIJ(mat,row[repl],&nnz,&jj,&vv);CHKERRQ(ierr);
114fd61f322SBarry Smith           for (kk=0; kk<nnz; kk++) {
115fd61f322SBarry Smith             if (icol[jj[kk]] == prow && PetscAbsScalar(vv[kk]) > atol) {
116b3cc6726SBarry Smith               ierr = MatRestoreRow_SeqAIJ(mat,row[repl],&nnz,&jj,&vv);CHKERRQ(ierr);
117fd61f322SBarry Smith               SWAP(icol[col[prow]],icol[col[repl]]);
118fd61f322SBarry Smith               SWAP(col[prow],col[repl]);
119fd61f322SBarry Smith               goto found;
120fd61f322SBarry Smith 	    }
121fd61f322SBarry Smith           }
122b3cc6726SBarry Smith           ierr = MatRestoreRow_SeqAIJ(mat,row[repl],&nnz,&jj,&vv);CHKERRQ(ierr);
123fd61f322SBarry Smith         }
124fd61f322SBarry Smith       }
125fd61f322SBarry Smith       /*
126fd61f322SBarry Smith           No column  suitable; instead check all future rows
127fd61f322SBarry Smith           Note: this will be very slow
128fd61f322SBarry Smith       */
129fd61f322SBarry Smith       for (k=prow+1; k<n; k++) {
130b3cc6726SBarry Smith         ierr = MatGetRow_SeqAIJ(mat,row[k],&nnz,&jj,&vv);CHKERRQ(ierr);
131fd61f322SBarry Smith         for (kk=0; kk<nnz; kk++) {
132fd61f322SBarry Smith           if (icol[jj[kk]] == prow && PetscAbsScalar(vv[kk]) > atol) {
133fd61f322SBarry Smith             /* found a row */
134fd61f322SBarry Smith             SWAP(row[prow],row[k]);
135fd61f322SBarry Smith             goto found;
136fd61f322SBarry Smith           }
137fd61f322SBarry Smith         }
138b3cc6726SBarry Smith         ierr = MatRestoreRow_SeqAIJ(mat,row[k],&nnz,&jj,&vv);CHKERRQ(ierr);
139fd61f322SBarry Smith       }
140fd61f322SBarry Smith 
141fd61f322SBarry Smith       found:;
142ad608de0SBarry Smith     }
143b3cc6726SBarry Smith     ierr = MatRestoreRow_SeqAIJ(mat,row[prow],&nz,&j,&v);CHKERRQ(ierr);
144ad608de0SBarry Smith   }
14548b35521SBarry Smith   ierr = ISRestoreIndices(ris,&row);CHKERRQ(ierr);
14648b35521SBarry Smith   ierr = ISRestoreIndices(cis,&col);CHKERRQ(ierr);
14733f51a72SBarry Smith   ierr = ISRestoreIndices(icis,&icol);CHKERRQ(ierr);
14833f51a72SBarry Smith   ierr = ISDestroy(icis);CHKERRQ(ierr);
1493a40ed3dSBarry Smith   PetscFunctionReturn(0);
150ad608de0SBarry Smith }
15105b94e36SKris Buschelman EXTERN_C_END
15248b35521SBarry Smith 
15348b35521SBarry Smith 
154