xref: /petsc/src/mat/utils/zerodiag.c (revision be1d678a52e6eff2808b2fa31ae986cdbf03c9fe)
1*be1d678aSKris Buschelman #define PETSCMAT_DLL
2*be1d678aSKris Buschelman 
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 
1038baddfdSBarry Smith #define SWAP(a,b) {PetscInt _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 
3694b7f48cSBarry Smith     Options Database Keys (When using KSP):
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 @*/
55*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal(Mat mat,PetscReal abstol,IS ris,IS cis)
56ad608de0SBarry Smith {
57dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(Mat,PetscReal,IS,IS);
5805b94e36SKris Buschelman 
5905b94e36SKris Buschelman   PetscFunctionBegin;
6005b94e36SKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatReorderForNonzeroDiagonal_C",(void (**)(void))&f);CHKERRQ(ierr);
6105b94e36SKris Buschelman   if (f) {
6270441072SBarry Smith     ierr = (*f)(mat,abstol,ris,cis);CHKERRQ(ierr);
6305b94e36SKris Buschelman   }
6405b94e36SKris Buschelman   PetscFunctionReturn(0);
6505b94e36SKris Buschelman }
6605b94e36SKris Buschelman 
6738baddfdSBarry Smith EXTERN PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
6838baddfdSBarry Smith EXTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
69b3cc6726SBarry Smith 
7005b94e36SKris Buschelman EXTERN_C_BEGIN
7105b94e36SKris Buschelman #undef __FUNCT__
7205b94e36SKris Buschelman #define __FUNCT__ "MatReorderForNonzeroDiagonal_SeqAIJ"
73*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal_SeqAIJ(Mat mat,PetscReal abstol,IS ris,IS cis)
7405b94e36SKris Buschelman {
75dfbe8321SBarry Smith   PetscErrorCode ierr;
7638baddfdSBarry Smith   PetscInt       prow,k,nz,n,repl,*j,*col,*row,m,*icol,nnz,*jj,kk;
7787828ca2SBarry Smith   PetscScalar    *v,*vv;
78329f5518SBarry Smith   PetscReal      repla;
79fd61f322SBarry Smith   IS             icis;
80ad608de0SBarry Smith 
813a40ed3dSBarry Smith   PetscFunctionBegin;
8248b35521SBarry Smith   ierr = ISGetIndices(ris,&row);CHKERRQ(ierr);
8348b35521SBarry Smith   ierr = ISGetIndices(cis,&col);CHKERRQ(ierr);
844c49b128SBarry Smith   ierr = ISInvertPermutation(cis,PETSC_DECIDE,&icis);CHKERRQ(ierr);
8533f51a72SBarry Smith   ierr = ISGetIndices(icis,&icol);CHKERRQ(ierr);
8648b35521SBarry Smith   ierr = MatGetSize(mat,&m,&n);CHKERRQ(ierr);
87ad608de0SBarry Smith 
88ad608de0SBarry Smith   for (prow=0; prow<n; prow++) {
89b3cc6726SBarry Smith     ierr = MatGetRow_SeqAIJ(mat,row[prow],&nz,&j,&v);CHKERRQ(ierr);
9033f51a72SBarry Smith     for (k=0; k<nz; k++) {if (icol[j[k]] == prow) break;}
9170441072SBarry Smith     if (k >= nz || PetscAbsScalar(v[k]) <= abstol) {
92ad608de0SBarry Smith       /* Element too small or zero; find the best candidate */
93cddf8d76SBarry Smith       repla = (k >= nz) ? 0.0 : PetscAbsScalar(v[k]);
94d4bb536fSBarry Smith       /*
95fd61f322SBarry Smith           Look for a later column we can swap with this one
96d4bb536fSBarry Smith       */
9748b35521SBarry Smith       for (k=0; k<nz; k++) {
9833f51a72SBarry Smith 	if (icol[j[k]] > prow && PetscAbsScalar(v[k]) > repla) {
99fd61f322SBarry Smith           /* found a suitable later column */
10033f51a72SBarry Smith 	  repl  = icol[j[k]];
10133f51a72SBarry Smith           SWAP(icol[col[prow]],icol[col[repl]]);
102ad608de0SBarry Smith           SWAP(col[prow],col[repl]);
103fd61f322SBarry Smith           goto found;
104fd61f322SBarry Smith         }
105fd61f322SBarry Smith       }
106fd61f322SBarry Smith       /*
107fd61f322SBarry Smith            Did not find a suitable later column so look for an earlier column
108fd61f322SBarry Smith 	   We need to be sure that we don't introduce a zero in a previous
109fd61f322SBarry Smith 	   diagonal
110fd61f322SBarry Smith       */
111fd61f322SBarry Smith       for (k=0; k<nz; k++) {
112fd61f322SBarry Smith         if (icol[j[k]] < prow && PetscAbsScalar(v[k]) > repla) {
113fd61f322SBarry Smith           /* See if this one will work */
114fd61f322SBarry Smith           repl  = icol[j[k]];
115b3cc6726SBarry Smith           ierr = MatGetRow_SeqAIJ(mat,row[repl],&nnz,&jj,&vv);CHKERRQ(ierr);
116fd61f322SBarry Smith           for (kk=0; kk<nnz; kk++) {
11770441072SBarry Smith             if (icol[jj[kk]] == prow && PetscAbsScalar(vv[kk]) > abstol) {
118b3cc6726SBarry Smith               ierr = MatRestoreRow_SeqAIJ(mat,row[repl],&nnz,&jj,&vv);CHKERRQ(ierr);
119fd61f322SBarry Smith               SWAP(icol[col[prow]],icol[col[repl]]);
120fd61f322SBarry Smith               SWAP(col[prow],col[repl]);
121fd61f322SBarry Smith               goto found;
122fd61f322SBarry Smith 	    }
123fd61f322SBarry Smith           }
124b3cc6726SBarry Smith           ierr = MatRestoreRow_SeqAIJ(mat,row[repl],&nnz,&jj,&vv);CHKERRQ(ierr);
125fd61f322SBarry Smith         }
126fd61f322SBarry Smith       }
127fd61f322SBarry Smith       /*
128fd61f322SBarry Smith           No column  suitable; instead check all future rows
129fd61f322SBarry Smith           Note: this will be very slow
130fd61f322SBarry Smith       */
131fd61f322SBarry Smith       for (k=prow+1; k<n; k++) {
132b3cc6726SBarry Smith         ierr = MatGetRow_SeqAIJ(mat,row[k],&nnz,&jj,&vv);CHKERRQ(ierr);
133fd61f322SBarry Smith         for (kk=0; kk<nnz; kk++) {
13470441072SBarry Smith           if (icol[jj[kk]] == prow && PetscAbsScalar(vv[kk]) > abstol) {
135fd61f322SBarry Smith             /* found a row */
136fd61f322SBarry Smith             SWAP(row[prow],row[k]);
137fd61f322SBarry Smith             goto found;
138fd61f322SBarry Smith           }
139fd61f322SBarry Smith         }
140b3cc6726SBarry Smith         ierr = MatRestoreRow_SeqAIJ(mat,row[k],&nnz,&jj,&vv);CHKERRQ(ierr);
141fd61f322SBarry Smith       }
142fd61f322SBarry Smith 
143fd61f322SBarry Smith       found:;
144ad608de0SBarry Smith     }
145b3cc6726SBarry Smith     ierr = MatRestoreRow_SeqAIJ(mat,row[prow],&nz,&j,&v);CHKERRQ(ierr);
146ad608de0SBarry Smith   }
14748b35521SBarry Smith   ierr = ISRestoreIndices(ris,&row);CHKERRQ(ierr);
14848b35521SBarry Smith   ierr = ISRestoreIndices(cis,&col);CHKERRQ(ierr);
14933f51a72SBarry Smith   ierr = ISRestoreIndices(icis,&icol);CHKERRQ(ierr);
15033f51a72SBarry Smith   ierr = ISDestroy(icis);CHKERRQ(ierr);
1513a40ed3dSBarry Smith   PetscFunctionReturn(0);
152ad608de0SBarry Smith }
15305b94e36SKris Buschelman EXTERN_C_END
15448b35521SBarry Smith 
15548b35521SBarry Smith 
156