xref: /petsc/src/mat/utils/zerodiag.c (revision a8f51744601f91dc30d5f98153b1ecd91eebb0e5)
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 
7af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I  "petscmat.h"  I*/
8ad608de0SBarry Smith 
99371c9d4SSatish Balay #define SWAP(a, b) \
10*a8f51744SPierre Jolivet   do { \
119371c9d4SSatish Balay     PetscInt _t; \
129371c9d4SSatish Balay     _t = a; \
139371c9d4SSatish Balay     a  = b; \
149371c9d4SSatish Balay     b  = _t; \
15*a8f51744SPierre Jolivet   } while (0)
1648b35521SBarry Smith 
17ad608de0SBarry Smith /*@
1848b35521SBarry Smith   MatReorderForNonzeroDiagonal - Changes matrix ordering to remove
1911a5261eSBarry Smith   zeros from diagonal. This may help in the `PCLU` factorization to
2048b35521SBarry Smith   prevent a zero pivot.
21ad608de0SBarry Smith 
22c3339decSBarry Smith   Collective
23fee21e36SBarry Smith 
2498a79cdbSBarry Smith   Input Parameters:
2598a79cdbSBarry Smith + mat    - matrix to reorder
262ef1f0ffSBarry Smith . abstol - absolute tolerance, it attempts to move all values smaller off the diagonal
272ef1f0ffSBarry Smith . ris    - the row reordering
282ef1f0ffSBarry Smith - cis    - the column reordering; this may be changed
2998a79cdbSBarry Smith 
3015091d37SBarry Smith   Level: intermediate
3115091d37SBarry Smith 
3227430b45SBarry Smith   Options Database Key:
3327430b45SBarry Smith . -pc_factor_nonzeros_along_diagonal - Reorder to remove zeros from diagonal
3427430b45SBarry Smith 
35ad608de0SBarry Smith   Notes:
36ad608de0SBarry Smith   This is not intended as a replacement for pivoting for matrices that
372ef1f0ffSBarry Smith   have ``bad'' structure. It is only a stop-gap measure.
382ef1f0ffSBarry Smith 
392ef1f0ffSBarry Smith   Should be called
402ef1f0ffSBarry Smith   after a call to `MatGetOrdering()`.
41d252947aSBarry Smith 
4211a5261eSBarry Smith   Only works for `MATSEQAIJ` matrices
43106f7b34SBarry Smith 
4427430b45SBarry Smith   Developer Notes:
4533f51a72SBarry Smith   Column pivoting is used.
4633f51a72SBarry Smith 
4733f51a72SBarry Smith   1) Choice of column is made by looking at the
4833f51a72SBarry Smith   non-zero elements in the troublesome row for columns that are not yet
4933f51a72SBarry Smith   included (moving from left to right).
5033f51a72SBarry Smith 
5133f51a72SBarry Smith   2) If (1) fails we check all the columns to the left of the current row
52814823dcSBarry Smith   and see if one of them has could be swapped. It can be swapped if
53814823dcSBarry Smith   its corresponding row has a non-zero in the column it is being
54814823dcSBarry Smith   swapped with; to make sure the previous nonzero diagonal remains
55814823dcSBarry Smith   nonzero
5633f51a72SBarry Smith 
572ef1f0ffSBarry Smith .seealso: `Mat`, `MatGetFactor()`, `MatGetOrdering()`
58ad608de0SBarry Smith @*/
59d71ae5a4SJacob Faibussowitsch PetscErrorCode MatReorderForNonzeroDiagonal(Mat mat, PetscReal abstol, IS ris, IS cis)
60d71ae5a4SJacob Faibussowitsch {
6105b94e36SKris Buschelman   PetscFunctionBegin;
62cac4c232SBarry Smith   PetscTryMethod(mat, "MatReorderForNonzeroDiagonal_C", (Mat, PetscReal, IS, IS), (mat, abstol, ris, cis));
633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6405b94e36SKris Buschelman }
6505b94e36SKris Buschelman 
665a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetRow_SeqAIJ(Mat, PetscInt, PetscInt *, PetscInt **, PetscScalar **);
675a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat, PetscInt, PetscInt *, PetscInt **, PetscScalar **);
68b3cc6726SBarry Smith 
69d9e81513SBarry Smith #include <../src/vec/is/is/impls/general/general.h>
705d0c19d7SBarry Smith 
71d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatReorderForNonzeroDiagonal_SeqAIJ(Mat mat, PetscReal abstol, IS ris, IS cis)
72d71ae5a4SJacob Faibussowitsch {
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;
829566063dSJacob Faibussowitsch   PetscCall(ISInvertPermutation(cis, PETSC_DECIDE, &icis));
835d0c19d7SBarry Smith   icol = ((IS_General *)icis->data)->idx;
849566063dSJacob Faibussowitsch   PetscCall(MatGetSize(mat, &m, &n));
85ad608de0SBarry Smith 
86ad608de0SBarry Smith   for (prow = 0; prow < n; prow++) {
879566063dSJacob Faibussowitsch     PetscCall(MatGetRow_SeqAIJ(mat, row[prow], &nz, &j, &v));
888865f1eaSKarl Rupp     for (k = 0; k < nz; k++) {
898865f1eaSKarl Rupp       if (icol[j[k]] == prow) break;
908865f1eaSKarl 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]];
1159566063dSJacob Faibussowitsch           PetscCall(MatGetRow_SeqAIJ(mat, row[repl], &nnz, &jj, &vv));
116fd61f322SBarry Smith           for (kk = 0; kk < nnz; kk++) {
11770441072SBarry Smith             if (icol[jj[kk]] == prow && PetscAbsScalar(vv[kk]) > abstol) {
1189566063dSJacob Faibussowitsch               PetscCall(MatRestoreRow_SeqAIJ(mat, row[repl], &nnz, &jj, &vv));
119fd61f322SBarry Smith               SWAP(icol[col[prow]], icol[col[repl]]);
120fd61f322SBarry Smith               SWAP(col[prow], col[repl]);
121fd61f322SBarry Smith               goto found;
122fd61f322SBarry Smith             }
123fd61f322SBarry Smith           }
1249566063dSJacob Faibussowitsch           PetscCall(MatRestoreRow_SeqAIJ(mat, row[repl], &nnz, &jj, &vv));
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++) {
1329566063dSJacob Faibussowitsch         PetscCall(MatGetRow_SeqAIJ(mat, row[k], &nnz, &jj, &vv));
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         }
1409566063dSJacob Faibussowitsch         PetscCall(MatRestoreRow_SeqAIJ(mat, row[k], &nnz, &jj, &vv));
141fd61f322SBarry Smith       }
142fd61f322SBarry Smith 
143fd61f322SBarry Smith     found:;
144ad608de0SBarry Smith     }
1459566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow_SeqAIJ(mat, row[prow], &nz, &j, &v));
146ad608de0SBarry Smith   }
1479566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&icis));
1483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
149ad608de0SBarry Smith }
150