xref: /petsc/src/mat/utils/zerodiag.c (revision 8865f1ea4ea560cd84ab8db62e98b7095cdff96f)
1be1d678aSKris Buschelman 
2ad608de0SBarry Smith /*
3ad608de0SBarry Smith     This file contains routines to reorder a matrix so that the diagonal
448b35521SBarry Smith     elements are nonzero.
5ad608de0SBarry Smith  */
6ad608de0SBarry Smith 
7b45d2f2cSJed Brown #include <petsc-private/matimpl.h>       /*I  "petscmat.h"  I*/
8ad608de0SBarry Smith 
938baddfdSBarry Smith #define SWAP(a,b) {PetscInt _t; _t = a; a = b; b = _t; }
1048b35521SBarry Smith 
114a2ae208SSatish Balay #undef __FUNCT__
124a2ae208SSatish Balay #define __FUNCT__ "MatReorderForNonzeroDiagonal"
13ad608de0SBarry Smith /*@
1448b35521SBarry Smith     MatReorderForNonzeroDiagonal - Changes matrix ordering to remove
1548b35521SBarry Smith     zeros from diagonal. This may help in the LU factorization to
1648b35521SBarry Smith     prevent a zero pivot.
17ad608de0SBarry Smith 
18fee21e36SBarry Smith     Collective on Mat
19fee21e36SBarry Smith 
2098a79cdbSBarry Smith     Input Parameters:
2198a79cdbSBarry Smith +   mat  - matrix to reorder
2298a79cdbSBarry Smith -   rmap,cmap - row and column permutations.  Usually obtained from
2391e9ee9fSBarry Smith                MatGetOrdering().
2498a79cdbSBarry Smith 
2515091d37SBarry Smith     Level: intermediate
2615091d37SBarry Smith 
27ad608de0SBarry Smith     Notes:
28ad608de0SBarry Smith     This is not intended as a replacement for pivoting for matrices that
29d252947aSBarry Smith     have ``bad'' structure. It is only a stop-gap measure. Should be called
3091e9ee9fSBarry Smith     after a call to MatGetOrdering(), this routine changes the column
31d252947aSBarry Smith     ordering defined in cis.
32d252947aSBarry Smith 
33106f7b34SBarry Smith     Only works for SeqAIJ matrices
34106f7b34SBarry Smith 
3594b7f48cSBarry Smith     Options Database Keys (When using KSP):
362401956bSBarry Smith .      -pc_factor_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 @*/
537087cfbeSBarry Smith PetscErrorCode  MatReorderForNonzeroDiagonal(Mat mat,PetscReal abstol,IS ris,IS cis)
54ad608de0SBarry Smith {
554ac538c5SBarry Smith   PetscErrorCode ierr;
5605b94e36SKris Buschelman 
5705b94e36SKris Buschelman   PetscFunctionBegin;
584ac538c5SBarry Smith   ierr = PetscTryMethod(mat,"MatReorderForNonzeroDiagonal_C",(Mat,PetscReal,IS,IS),(mat,abstol,ris,cis));CHKERRQ(ierr);
5905b94e36SKris Buschelman   PetscFunctionReturn(0);
6005b94e36SKris Buschelman }
6105b94e36SKris Buschelman 
6209573ac7SBarry Smith extern PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
6309573ac7SBarry Smith extern PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
64b3cc6726SBarry Smith 
65c6db04a5SJed Brown #include <../src/vec/is/impls/general/general.h>
665d0c19d7SBarry Smith 
6705b94e36SKris Buschelman EXTERN_C_BEGIN
6805b94e36SKris Buschelman #undef __FUNCT__
6905b94e36SKris Buschelman #define __FUNCT__ "MatReorderForNonzeroDiagonal_SeqAIJ"
707087cfbeSBarry Smith PetscErrorCode  MatReorderForNonzeroDiagonal_SeqAIJ(Mat mat,PetscReal abstol,IS ris,IS cis)
7105b94e36SKris Buschelman {
72dfbe8321SBarry Smith   PetscErrorCode ierr;
7338baddfdSBarry Smith   PetscInt       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;
795d0c19d7SBarry Smith   /* access the indices of the IS directly, because it changes them */
805d0c19d7SBarry Smith   row  = ((IS_General*)ris->data)->idx;
815d0c19d7SBarry Smith   col  = ((IS_General*)cis->data)->idx;
824c49b128SBarry Smith   ierr = ISInvertPermutation(cis,PETSC_DECIDE,&icis);CHKERRQ(ierr);
835d0c19d7SBarry Smith   icol = ((IS_General*)icis->data)->idx;
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);
88*8865f1eaSKarl Rupp     for (k=0; k<nz; k++) {
89*8865f1eaSKarl Rupp       if (icol[j[k]] == prow) break;
90*8865f1eaSKarl Rupp     }
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   }
1476bf464f9SBarry Smith   ierr = ISDestroy(&icis);CHKERRQ(ierr);
1483a40ed3dSBarry Smith   PetscFunctionReturn(0);
149ad608de0SBarry Smith }
15005b94e36SKris Buschelman EXTERN_C_END
15148b35521SBarry Smith 
15248b35521SBarry Smith 
153