xref: /petsc/src/mat/utils/zerodiag.c (revision 05b94e364f8082b74a6c3d9e85326f85fec037b9)
173f4d377SMatthew Knepley /*$Id: zerodiag.c,v 1.44 2001/08/06 21:16:10 bsmith Exp $*/
2ad608de0SBarry Smith 
3ad608de0SBarry Smith /*
4ad608de0SBarry Smith     This file contains routines to reorder a matrix so that the diagonal
548b35521SBarry Smith     elements are nonzero.
6ad608de0SBarry Smith  */
7ad608de0SBarry Smith 
8e090d566SSatish Balay #include "src/mat/matimpl.h"       /*I  "petscmat.h"  I*/
9ad608de0SBarry Smith 
1048b35521SBarry Smith #define SWAP(a,b) {int _t; _t = a; a = b; b = _t; }
1148b35521SBarry Smith 
124a2ae208SSatish Balay #undef __FUNCT__
134a2ae208SSatish Balay #define __FUNCT__ "MatReorderForNonzeroDiagonal"
14ad608de0SBarry Smith /*@
1548b35521SBarry Smith     MatReorderForNonzeroDiagonal - Changes matrix ordering to remove
1648b35521SBarry Smith     zeros from diagonal. This may help in the LU factorization to
1748b35521SBarry Smith     prevent a zero pivot.
18ad608de0SBarry Smith 
19fee21e36SBarry Smith     Collective on Mat
20fee21e36SBarry Smith 
2198a79cdbSBarry Smith     Input Parameters:
2298a79cdbSBarry Smith +   mat  - matrix to reorder
2398a79cdbSBarry Smith -   rmap,cmap - row and column permutations.  Usually obtained from
2491e9ee9fSBarry Smith                MatGetOrdering().
2598a79cdbSBarry Smith 
2615091d37SBarry Smith     Level: intermediate
2715091d37SBarry Smith 
28ad608de0SBarry Smith     Notes:
29ad608de0SBarry Smith     This is not intended as a replacement for pivoting for matrices that
30d252947aSBarry Smith     have ``bad'' structure. It is only a stop-gap measure. Should be called
3191e9ee9fSBarry Smith     after a call to MatGetOrdering(), this routine changes the column
32d252947aSBarry Smith     ordering defined in cis.
33d252947aSBarry Smith 
34106f7b34SBarry Smith     Only works for SeqAIJ matrices
35106f7b34SBarry Smith 
36b259b22eSLois Curfman McInnes     Options Database Keys (When using SLES):
3798a79cdbSBarry Smith +      -pc_ilu_nonzeros_along_diagonal
3898a79cdbSBarry Smith -      -pc_lu_nonzeros_along_diagonal
39ad608de0SBarry Smith 
4033f51a72SBarry Smith     Algorithm Notes:
4133f51a72SBarry Smith     Column pivoting is used.
4233f51a72SBarry Smith 
4333f51a72SBarry Smith     1) Choice of column is made by looking at the
4433f51a72SBarry Smith        non-zero elements in the troublesome row for columns that are not yet
4533f51a72SBarry Smith        included (moving from left to right).
4633f51a72SBarry Smith 
4733f51a72SBarry Smith     2) If (1) fails we check all the columns to the left of the current row
48814823dcSBarry Smith        and see if one of them has could be swapped. It can be swapped if
49814823dcSBarry Smith        its corresponding row has a non-zero in the column it is being
50814823dcSBarry Smith        swapped with; to make sure the previous nonzero diagonal remains
51814823dcSBarry Smith        nonzero
5233f51a72SBarry Smith 
5398a79cdbSBarry Smith 
54ad608de0SBarry Smith @*/
55329f5518SBarry Smith int MatReorderForNonzeroDiagonal(Mat mat,PetscReal atol,IS ris,IS cis)
56ad608de0SBarry Smith {
57*05b94e36SKris Buschelman   int ierr,(*f)(Mat,PetscReal,IS,IS);
58*05b94e36SKris Buschelman 
59*05b94e36SKris Buschelman   PetscFunctionBegin;
60*05b94e36SKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatReorderForNonzeroDiagonal_C",(void (**)(void))&f);CHKERRQ(ierr);
61*05b94e36SKris Buschelman   if (f) {
62*05b94e36SKris Buschelman     ierr = (*f)(mat,atol,ris,cis);CHKERRQ(ierr);
63*05b94e36SKris Buschelman   }
64*05b94e36SKris Buschelman   PetscFunctionReturn(0);
65*05b94e36SKris Buschelman }
66*05b94e36SKris Buschelman 
67*05b94e36SKris Buschelman EXTERN_C_BEGIN
68*05b94e36SKris Buschelman #undef __FUNCT__
69*05b94e36SKris Buschelman #define __FUNCT__ "MatReorderForNonzeroDiagonal_SeqAIJ"
70*05b94e36SKris Buschelman int MatReorderForNonzeroDiagonal_SeqAIJ(Mat mat,PetscReal atol,IS ris,IS cis)
71*05b94e36SKris Buschelman {
7272af6520SSatish Balay   int         ierr,prow,k,nz,n,repl,*j,*col,*row,m,*icol,nnz,*jj,kk;
7387828ca2SBarry Smith   PetscScalar *v,*vv;
74329f5518SBarry Smith   PetscReal   repla;
75fd61f322SBarry Smith   IS          icis;
76106f7b34SBarry Smith   PetscTruth  flg;
77ad608de0SBarry Smith 
783a40ed3dSBarry Smith   PetscFunctionBegin;
7990f02eecSBarry Smith   PetscValidHeaderSpecific(mat,MAT_COOKIE);
8090f02eecSBarry Smith   PetscValidHeaderSpecific(ris,IS_COOKIE);
8190f02eecSBarry Smith   PetscValidHeaderSpecific(cis,IS_COOKIE);
8290f02eecSBarry Smith 
8348b35521SBarry Smith   ierr = ISGetIndices(ris,&row);CHKERRQ(ierr);
8448b35521SBarry Smith   ierr = ISGetIndices(cis,&col);CHKERRQ(ierr);
854c49b128SBarry Smith   ierr = ISInvertPermutation(cis,PETSC_DECIDE,&icis);CHKERRQ(ierr);
8633f51a72SBarry Smith   ierr = ISGetIndices(icis,&icol);CHKERRQ(ierr);
8748b35521SBarry Smith   ierr = MatGetSize(mat,&m,&n);CHKERRQ(ierr);
88ad608de0SBarry Smith 
89ad608de0SBarry Smith   for (prow=0; prow<n; prow++) {
9090f02eecSBarry Smith     ierr = MatGetRow(mat,row[prow],&nz,&j,&v);CHKERRQ(ierr);
9133f51a72SBarry Smith     for (k=0; k<nz; k++) {if (icol[j[k]] == prow) break;}
92cddf8d76SBarry Smith     if (k >= nz || PetscAbsScalar(v[k]) <= atol) {
93ad608de0SBarry Smith       /* Element too small or zero; find the best candidate */
94cddf8d76SBarry Smith       repla = (k >= nz) ? 0.0 : PetscAbsScalar(v[k]);
95d4bb536fSBarry Smith       /*
96fd61f322SBarry Smith           Look for a later column we can swap with this one
97d4bb536fSBarry Smith       */
9848b35521SBarry Smith       for (k=0; k<nz; k++) {
9933f51a72SBarry Smith 	if (icol[j[k]] > prow && PetscAbsScalar(v[k]) > repla) {
100fd61f322SBarry Smith           /* found a suitable later column */
10133f51a72SBarry Smith 	  repl  = icol[j[k]];
10233f51a72SBarry Smith           SWAP(icol[col[prow]],icol[col[repl]]);
103ad608de0SBarry Smith           SWAP(col[prow],col[repl]);
104fd61f322SBarry Smith           goto found;
105fd61f322SBarry Smith         }
106fd61f322SBarry Smith       }
107fd61f322SBarry Smith       /*
108fd61f322SBarry Smith            Did not find a suitable later column so look for an earlier column
109fd61f322SBarry Smith 	   We need to be sure that we don't introduce a zero in a previous
110fd61f322SBarry Smith 	   diagonal
111fd61f322SBarry Smith       */
112fd61f322SBarry Smith       for (k=0; k<nz; k++) {
113fd61f322SBarry Smith         if (icol[j[k]] < prow && PetscAbsScalar(v[k]) > repla) {
114fd61f322SBarry Smith           /* See if this one will work */
115fd61f322SBarry Smith           repl  = icol[j[k]];
116fd61f322SBarry Smith           ierr = MatGetRow(mat,row[repl],&nnz,&jj,&vv);CHKERRQ(ierr);
117fd61f322SBarry Smith           for (kk=0; kk<nnz; kk++) {
118fd61f322SBarry Smith             if (icol[jj[kk]] == prow && PetscAbsScalar(vv[kk]) > atol) {
119fd61f322SBarry Smith               ierr = MatRestoreRow(mat,row[repl],&nnz,&jj,&vv);CHKERRQ(ierr);
120fd61f322SBarry Smith               SWAP(icol[col[prow]],icol[col[repl]]);
121fd61f322SBarry Smith               SWAP(col[prow],col[repl]);
122fd61f322SBarry Smith               goto found;
123fd61f322SBarry Smith 	    }
124fd61f322SBarry Smith           }
125fd61f322SBarry Smith           ierr = MatRestoreRow(mat,row[repl],&nnz,&jj,&vv);CHKERRQ(ierr);
126fd61f322SBarry Smith         }
127fd61f322SBarry Smith       }
128fd61f322SBarry Smith       /*
129fd61f322SBarry Smith           No column  suitable; instead check all future rows
130fd61f322SBarry Smith           Note: this will be very slow
131fd61f322SBarry Smith       */
132fd61f322SBarry Smith       for (k=prow+1; k<n; k++) {
133fd61f322SBarry Smith         ierr = MatGetRow(mat,row[k],&nnz,&jj,&vv);CHKERRQ(ierr);
134fd61f322SBarry Smith         for (kk=0; kk<nnz; kk++) {
135fd61f322SBarry Smith           if (icol[jj[kk]] == prow && PetscAbsScalar(vv[kk]) > atol) {
136fd61f322SBarry Smith             /* found a row */
137fd61f322SBarry Smith             SWAP(row[prow],row[k]);
138fd61f322SBarry Smith             goto found;
139fd61f322SBarry Smith           }
140fd61f322SBarry Smith         }
141fd61f322SBarry Smith         ierr = MatRestoreRow(mat,row[k],&nnz,&jj,&vv);CHKERRQ(ierr);
142fd61f322SBarry Smith       }
143fd61f322SBarry Smith 
144fd61f322SBarry Smith       found:;
145ad608de0SBarry Smith     }
14690f02eecSBarry Smith     ierr = MatRestoreRow(mat,row[prow],&nz,&j,&v);CHKERRQ(ierr);
147ad608de0SBarry Smith   }
14848b35521SBarry Smith   ierr = ISRestoreIndices(ris,&row);CHKERRQ(ierr);
14948b35521SBarry Smith   ierr = ISRestoreIndices(cis,&col);CHKERRQ(ierr);
15033f51a72SBarry Smith   ierr = ISRestoreIndices(icis,&icol);CHKERRQ(ierr);
15133f51a72SBarry Smith   ierr = ISDestroy(icis);CHKERRQ(ierr);
1523a40ed3dSBarry Smith   PetscFunctionReturn(0);
153ad608de0SBarry Smith }
154*05b94e36SKris Buschelman EXTERN_C_END
15548b35521SBarry Smith 
15648b35521SBarry Smith 
157