xref: /petsc/src/mat/utils/zerodiag.c (revision b3cc6726fbddf13cd805f591b2b0fc9b163819ab)
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 
848b35521SBarry Smith #define SWAP(a,b) {int _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 @*/
53329f5518SBarry Smith int MatReorderForNonzeroDiagonal(Mat mat,PetscReal atol,IS ris,IS cis)
54ad608de0SBarry Smith {
5505b94e36SKris Buschelman   int 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*b3cc6726SBarry Smith EXTERN int MatGetRow_SeqAIJ(Mat,int,int*,int**,PetscScalar**);
66*b3cc6726SBarry Smith EXTERN int MatRestoreRow_SeqAIJ(Mat,int,int*,int**,PetscScalar**);
67*b3cc6726SBarry Smith 
6805b94e36SKris Buschelman EXTERN_C_BEGIN
6905b94e36SKris Buschelman #undef __FUNCT__
7005b94e36SKris Buschelman #define __FUNCT__ "MatReorderForNonzeroDiagonal_SeqAIJ"
7105b94e36SKris Buschelman int MatReorderForNonzeroDiagonal_SeqAIJ(Mat mat,PetscReal atol,IS ris,IS cis)
7205b94e36SKris Buschelman {
7372af6520SSatish Balay   int         ierr,prow,k,nz,n,repl,*j,*col,*row,m,*icol,nnz,*jj,kk;
7487828ca2SBarry Smith   PetscScalar *v,*vv;
75329f5518SBarry Smith   PetscReal   repla;
76fd61f322SBarry Smith   IS          icis;
77ad608de0SBarry Smith 
783a40ed3dSBarry Smith   PetscFunctionBegin;
7948b35521SBarry Smith   ierr = ISGetIndices(ris,&row);CHKERRQ(ierr);
8048b35521SBarry Smith   ierr = ISGetIndices(cis,&col);CHKERRQ(ierr);
814c49b128SBarry Smith   ierr = ISInvertPermutation(cis,PETSC_DECIDE,&icis);CHKERRQ(ierr);
8233f51a72SBarry Smith   ierr = ISGetIndices(icis,&icol);CHKERRQ(ierr);
8348b35521SBarry Smith   ierr = MatGetSize(mat,&m,&n);CHKERRQ(ierr);
84ad608de0SBarry Smith 
85ad608de0SBarry Smith   for (prow=0; prow<n; prow++) {
86*b3cc6726SBarry Smith     ierr = MatGetRow_SeqAIJ(mat,row[prow],&nz,&j,&v);CHKERRQ(ierr);
8733f51a72SBarry Smith     for (k=0; k<nz; k++) {if (icol[j[k]] == prow) break;}
88cddf8d76SBarry Smith     if (k >= nz || PetscAbsScalar(v[k]) <= atol) {
89ad608de0SBarry Smith       /* Element too small or zero; find the best candidate */
90cddf8d76SBarry Smith       repla = (k >= nz) ? 0.0 : PetscAbsScalar(v[k]);
91d4bb536fSBarry Smith       /*
92fd61f322SBarry Smith           Look for a later column we can swap with this one
93d4bb536fSBarry Smith       */
9448b35521SBarry Smith       for (k=0; k<nz; k++) {
9533f51a72SBarry Smith 	if (icol[j[k]] > prow && PetscAbsScalar(v[k]) > repla) {
96fd61f322SBarry Smith           /* found a suitable later column */
9733f51a72SBarry Smith 	  repl  = icol[j[k]];
9833f51a72SBarry Smith           SWAP(icol[col[prow]],icol[col[repl]]);
99ad608de0SBarry Smith           SWAP(col[prow],col[repl]);
100fd61f322SBarry Smith           goto found;
101fd61f322SBarry Smith         }
102fd61f322SBarry Smith       }
103fd61f322SBarry Smith       /*
104fd61f322SBarry Smith            Did not find a suitable later column so look for an earlier column
105fd61f322SBarry Smith 	   We need to be sure that we don't introduce a zero in a previous
106fd61f322SBarry Smith 	   diagonal
107fd61f322SBarry Smith       */
108fd61f322SBarry Smith       for (k=0; k<nz; k++) {
109fd61f322SBarry Smith         if (icol[j[k]] < prow && PetscAbsScalar(v[k]) > repla) {
110fd61f322SBarry Smith           /* See if this one will work */
111fd61f322SBarry Smith           repl  = icol[j[k]];
112*b3cc6726SBarry Smith           ierr = MatGetRow_SeqAIJ(mat,row[repl],&nnz,&jj,&vv);CHKERRQ(ierr);
113fd61f322SBarry Smith           for (kk=0; kk<nnz; kk++) {
114fd61f322SBarry Smith             if (icol[jj[kk]] == prow && PetscAbsScalar(vv[kk]) > atol) {
115*b3cc6726SBarry Smith               ierr = MatRestoreRow_SeqAIJ(mat,row[repl],&nnz,&jj,&vv);CHKERRQ(ierr);
116fd61f322SBarry Smith               SWAP(icol[col[prow]],icol[col[repl]]);
117fd61f322SBarry Smith               SWAP(col[prow],col[repl]);
118fd61f322SBarry Smith               goto found;
119fd61f322SBarry Smith 	    }
120fd61f322SBarry Smith           }
121*b3cc6726SBarry Smith           ierr = MatRestoreRow_SeqAIJ(mat,row[repl],&nnz,&jj,&vv);CHKERRQ(ierr);
122fd61f322SBarry Smith         }
123fd61f322SBarry Smith       }
124fd61f322SBarry Smith       /*
125fd61f322SBarry Smith           No column  suitable; instead check all future rows
126fd61f322SBarry Smith           Note: this will be very slow
127fd61f322SBarry Smith       */
128fd61f322SBarry Smith       for (k=prow+1; k<n; k++) {
129*b3cc6726SBarry Smith         ierr = MatGetRow_SeqAIJ(mat,row[k],&nnz,&jj,&vv);CHKERRQ(ierr);
130fd61f322SBarry Smith         for (kk=0; kk<nnz; kk++) {
131fd61f322SBarry Smith           if (icol[jj[kk]] == prow && PetscAbsScalar(vv[kk]) > atol) {
132fd61f322SBarry Smith             /* found a row */
133fd61f322SBarry Smith             SWAP(row[prow],row[k]);
134fd61f322SBarry Smith             goto found;
135fd61f322SBarry Smith           }
136fd61f322SBarry Smith         }
137*b3cc6726SBarry Smith         ierr = MatRestoreRow_SeqAIJ(mat,row[k],&nnz,&jj,&vv);CHKERRQ(ierr);
138fd61f322SBarry Smith       }
139fd61f322SBarry Smith 
140fd61f322SBarry Smith       found:;
141ad608de0SBarry Smith     }
142*b3cc6726SBarry Smith     ierr = MatRestoreRow_SeqAIJ(mat,row[prow],&nz,&j,&v);CHKERRQ(ierr);
143ad608de0SBarry Smith   }
14448b35521SBarry Smith   ierr = ISRestoreIndices(ris,&row);CHKERRQ(ierr);
14548b35521SBarry Smith   ierr = ISRestoreIndices(cis,&col);CHKERRQ(ierr);
14633f51a72SBarry Smith   ierr = ISRestoreIndices(icis,&icol);CHKERRQ(ierr);
14733f51a72SBarry Smith   ierr = ISDestroy(icis);CHKERRQ(ierr);
1483a40ed3dSBarry Smith   PetscFunctionReturn(0);
149ad608de0SBarry Smith }
15005b94e36SKris Buschelman EXTERN_C_END
15148b35521SBarry Smith 
15248b35521SBarry Smith 
153