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 938baddfdSBarry Smith #define SWAP(a,b) {PetscInt _t; _t = a; a = b; b = _t; } 1048b35521SBarry Smith 11ad608de0SBarry Smith /*@ 1248b35521SBarry Smith MatReorderForNonzeroDiagonal - Changes matrix ordering to remove 1348b35521SBarry Smith zeros from diagonal. This may help in the LU factorization to 1448b35521SBarry Smith prevent a zero pivot. 15ad608de0SBarry Smith 16fee21e36SBarry Smith Collective on Mat 17fee21e36SBarry Smith 1898a79cdbSBarry Smith Input Parameters: 1998a79cdbSBarry Smith + mat - matrix to reorder 2098a79cdbSBarry Smith - rmap,cmap - row and column permutations. Usually obtained from 2191e9ee9fSBarry Smith MatGetOrdering(). 2298a79cdbSBarry Smith 2315091d37SBarry Smith Level: intermediate 2415091d37SBarry Smith 25ad608de0SBarry Smith Notes: 26ad608de0SBarry Smith This is not intended as a replacement for pivoting for matrices that 27d252947aSBarry Smith have ``bad'' structure. It is only a stop-gap measure. Should be called 2891e9ee9fSBarry Smith after a call to MatGetOrdering(), this routine changes the column 29d252947aSBarry Smith ordering defined in cis. 30d252947aSBarry Smith 31106f7b34SBarry Smith Only works for SeqAIJ matrices 32106f7b34SBarry Smith 3394b7f48cSBarry Smith Options Database Keys (When using KSP): 3467b8a455SSatish Balay . -pc_factor_nonzeros_along_diagonal - Reorder to remove zeros from diagonal 35ad608de0SBarry Smith 3633f51a72SBarry Smith Algorithm Notes: 3733f51a72SBarry Smith Column pivoting is used. 3833f51a72SBarry Smith 3933f51a72SBarry Smith 1) Choice of column is made by looking at the 4033f51a72SBarry Smith non-zero elements in the troublesome row for columns that are not yet 4133f51a72SBarry Smith included (moving from left to right). 4233f51a72SBarry Smith 4333f51a72SBarry Smith 2) If (1) fails we check all the columns to the left of the current row 44814823dcSBarry Smith and see if one of them has could be swapped. It can be swapped if 45814823dcSBarry Smith its corresponding row has a non-zero in the column it is being 46814823dcSBarry Smith swapped with; to make sure the previous nonzero diagonal remains 47814823dcSBarry Smith nonzero 4833f51a72SBarry Smith 49ad608de0SBarry Smith @*/ 507087cfbeSBarry Smith PetscErrorCode MatReorderForNonzeroDiagonal(Mat mat,PetscReal abstol,IS ris,IS cis) 51ad608de0SBarry Smith { 5205b94e36SKris Buschelman PetscFunctionBegin; 53*cac4c232SBarry Smith PetscTryMethod(mat,"MatReorderForNonzeroDiagonal_C",(Mat,PetscReal,IS,IS),(mat,abstol,ris,cis)); 5405b94e36SKris Buschelman PetscFunctionReturn(0); 5505b94e36SKris Buschelman } 5605b94e36SKris Buschelman 575a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 585a576424SJed Brown PETSC_INTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 59b3cc6726SBarry Smith 60d9e81513SBarry Smith #include <../src/vec/is/is/impls/general/general.h> 615d0c19d7SBarry Smith 62d8fa6a54SBarry Smith PETSC_INTERN PetscErrorCode MatReorderForNonzeroDiagonal_SeqAIJ(Mat mat,PetscReal abstol,IS ris,IS cis) 6305b94e36SKris Buschelman { 6438baddfdSBarry Smith PetscInt prow,k,nz,n,repl,*j,*col,*row,m,*icol,nnz,*jj,kk; 6587828ca2SBarry Smith PetscScalar *v,*vv; 66329f5518SBarry Smith PetscReal repla; 67fd61f322SBarry Smith IS icis; 68ad608de0SBarry Smith 693a40ed3dSBarry Smith PetscFunctionBegin; 705d0c19d7SBarry Smith /* access the indices of the IS directly, because it changes them */ 715d0c19d7SBarry Smith row = ((IS_General*)ris->data)->idx; 725d0c19d7SBarry Smith col = ((IS_General*)cis->data)->idx; 739566063dSJacob Faibussowitsch PetscCall(ISInvertPermutation(cis,PETSC_DECIDE,&icis)); 745d0c19d7SBarry Smith icol = ((IS_General*)icis->data)->idx; 759566063dSJacob Faibussowitsch PetscCall(MatGetSize(mat,&m,&n)); 76ad608de0SBarry Smith 77ad608de0SBarry Smith for (prow=0; prow<n; prow++) { 789566063dSJacob Faibussowitsch PetscCall(MatGetRow_SeqAIJ(mat,row[prow],&nz,&j,&v)); 798865f1eaSKarl Rupp for (k=0; k<nz; k++) { 808865f1eaSKarl Rupp if (icol[j[k]] == prow) break; 818865f1eaSKarl Rupp } 8270441072SBarry Smith if (k >= nz || PetscAbsScalar(v[k]) <= abstol) { 83ad608de0SBarry Smith /* Element too small or zero; find the best candidate */ 84cddf8d76SBarry Smith repla = (k >= nz) ? 0.0 : PetscAbsScalar(v[k]); 85d4bb536fSBarry Smith /* 86fd61f322SBarry Smith Look for a later column we can swap with this one 87d4bb536fSBarry Smith */ 8848b35521SBarry Smith for (k=0; k<nz; k++) { 8933f51a72SBarry Smith if (icol[j[k]] > prow && PetscAbsScalar(v[k]) > repla) { 90fd61f322SBarry Smith /* found a suitable later column */ 9133f51a72SBarry Smith repl = icol[j[k]]; 9233f51a72SBarry Smith SWAP(icol[col[prow]],icol[col[repl]]); 93ad608de0SBarry Smith SWAP(col[prow],col[repl]); 94fd61f322SBarry Smith goto found; 95fd61f322SBarry Smith } 96fd61f322SBarry Smith } 97fd61f322SBarry Smith /* 98fd61f322SBarry Smith Did not find a suitable later column so look for an earlier column 99fd61f322SBarry Smith We need to be sure that we don't introduce a zero in a previous 100fd61f322SBarry Smith diagonal 101fd61f322SBarry Smith */ 102fd61f322SBarry Smith for (k=0; k<nz; k++) { 103fd61f322SBarry Smith if (icol[j[k]] < prow && PetscAbsScalar(v[k]) > repla) { 104fd61f322SBarry Smith /* See if this one will work */ 105fd61f322SBarry Smith repl = icol[j[k]]; 1069566063dSJacob Faibussowitsch PetscCall(MatGetRow_SeqAIJ(mat,row[repl],&nnz,&jj,&vv)); 107fd61f322SBarry Smith for (kk=0; kk<nnz; kk++) { 10870441072SBarry Smith if (icol[jj[kk]] == prow && PetscAbsScalar(vv[kk]) > abstol) { 1099566063dSJacob Faibussowitsch PetscCall(MatRestoreRow_SeqAIJ(mat,row[repl],&nnz,&jj,&vv)); 110fd61f322SBarry Smith SWAP(icol[col[prow]],icol[col[repl]]); 111fd61f322SBarry Smith SWAP(col[prow],col[repl]); 112fd61f322SBarry Smith goto found; 113fd61f322SBarry Smith } 114fd61f322SBarry Smith } 1159566063dSJacob Faibussowitsch PetscCall(MatRestoreRow_SeqAIJ(mat,row[repl],&nnz,&jj,&vv)); 116fd61f322SBarry Smith } 117fd61f322SBarry Smith } 118fd61f322SBarry Smith /* 119fd61f322SBarry Smith No column suitable; instead check all future rows 120fd61f322SBarry Smith Note: this will be very slow 121fd61f322SBarry Smith */ 122fd61f322SBarry Smith for (k=prow+1; k<n; k++) { 1239566063dSJacob Faibussowitsch PetscCall(MatGetRow_SeqAIJ(mat,row[k],&nnz,&jj,&vv)); 124fd61f322SBarry Smith for (kk=0; kk<nnz; kk++) { 12570441072SBarry Smith if (icol[jj[kk]] == prow && PetscAbsScalar(vv[kk]) > abstol) { 126fd61f322SBarry Smith /* found a row */ 127fd61f322SBarry Smith SWAP(row[prow],row[k]); 128fd61f322SBarry Smith goto found; 129fd61f322SBarry Smith } 130fd61f322SBarry Smith } 1319566063dSJacob Faibussowitsch PetscCall(MatRestoreRow_SeqAIJ(mat,row[k],&nnz,&jj,&vv)); 132fd61f322SBarry Smith } 133fd61f322SBarry Smith 134fd61f322SBarry Smith found:; 135ad608de0SBarry Smith } 1369566063dSJacob Faibussowitsch PetscCall(MatRestoreRow_SeqAIJ(mat,row[prow],&nz,&j,&v)); 137ad608de0SBarry Smith } 1389566063dSJacob Faibussowitsch PetscCall(ISDestroy(&icis)); 1393a40ed3dSBarry Smith PetscFunctionReturn(0); 140ad608de0SBarry Smith } 141