1be1d678aSKris Buschelman #define PETSCMAT_DLL 2be1d678aSKris 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 8b9147fbbSdalcinl #include "include/private/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): 372401956bSBarry Smith . -pc_factor_nonzeros_along_diagonal 38ad608de0SBarry Smith 3933f51a72SBarry Smith Algorithm Notes: 4033f51a72SBarry Smith Column pivoting is used. 4133f51a72SBarry Smith 4233f51a72SBarry Smith 1) Choice of column is made by looking at the 4333f51a72SBarry Smith non-zero elements in the troublesome row for columns that are not yet 4433f51a72SBarry Smith included (moving from left to right). 4533f51a72SBarry Smith 4633f51a72SBarry Smith 2) If (1) fails we check all the columns to the left of the current row 47814823dcSBarry Smith and see if one of them has could be swapped. It can be swapped if 48814823dcSBarry Smith its corresponding row has a non-zero in the column it is being 49814823dcSBarry Smith swapped with; to make sure the previous nonzero diagonal remains 50814823dcSBarry Smith nonzero 5133f51a72SBarry Smith 5298a79cdbSBarry Smith 53ad608de0SBarry Smith @*/ 54be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal(Mat mat,PetscReal abstol,IS ris,IS cis) 55ad608de0SBarry Smith { 56dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscReal,IS,IS); 5705b94e36SKris Buschelman 5805b94e36SKris Buschelman PetscFunctionBegin; 5905b94e36SKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)mat,"MatReorderForNonzeroDiagonal_C",(void (**)(void))&f);CHKERRQ(ierr); 6005b94e36SKris Buschelman if (f) { 6170441072SBarry Smith ierr = (*f)(mat,abstol,ris,cis);CHKERRQ(ierr); 6205b94e36SKris Buschelman } 6305b94e36SKris Buschelman PetscFunctionReturn(0); 6405b94e36SKris Buschelman } 6505b94e36SKris Buschelman 6638baddfdSBarry Smith EXTERN PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 6738baddfdSBarry Smith EXTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**); 68b3cc6726SBarry Smith 69*5d0c19d7SBarry Smith #include "src/vec/is/impls/general/general.h" 70*5d0c19d7SBarry Smith 7105b94e36SKris Buschelman EXTERN_C_BEGIN 7205b94e36SKris Buschelman #undef __FUNCT__ 7305b94e36SKris Buschelman #define __FUNCT__ "MatReorderForNonzeroDiagonal_SeqAIJ" 74be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal_SeqAIJ(Mat mat,PetscReal abstol,IS ris,IS cis) 7505b94e36SKris Buschelman { 76dfbe8321SBarry Smith PetscErrorCode ierr; 7738baddfdSBarry Smith PetscInt prow,k,nz,n,repl,*j,*col,*row,m,*icol,nnz,*jj,kk; 7887828ca2SBarry Smith PetscScalar *v,*vv; 79329f5518SBarry Smith PetscReal repla; 80fd61f322SBarry Smith IS icis; 81ad608de0SBarry Smith 823a40ed3dSBarry Smith PetscFunctionBegin; 83*5d0c19d7SBarry Smith /* access the indices of the IS directly, because it changes them */ 84*5d0c19d7SBarry Smith row = ((IS_General*)ris->data)->idx; 85*5d0c19d7SBarry Smith col = ((IS_General*)cis->data)->idx; 864c49b128SBarry Smith ierr = ISInvertPermutation(cis,PETSC_DECIDE,&icis);CHKERRQ(ierr); 87*5d0c19d7SBarry Smith icol = ((IS_General*)icis->data)->idx; 8848b35521SBarry Smith ierr = MatGetSize(mat,&m,&n);CHKERRQ(ierr); 89ad608de0SBarry Smith 90ad608de0SBarry Smith for (prow=0; prow<n; prow++) { 91b3cc6726SBarry Smith ierr = MatGetRow_SeqAIJ(mat,row[prow],&nz,&j,&v);CHKERRQ(ierr); 9233f51a72SBarry Smith for (k=0; k<nz; k++) {if (icol[j[k]] == prow) break;} 9370441072SBarry Smith if (k >= nz || PetscAbsScalar(v[k]) <= abstol) { 94ad608de0SBarry Smith /* Element too small or zero; find the best candidate */ 95cddf8d76SBarry Smith repla = (k >= nz) ? 0.0 : PetscAbsScalar(v[k]); 96d4bb536fSBarry Smith /* 97fd61f322SBarry Smith Look for a later column we can swap with this one 98d4bb536fSBarry Smith */ 9948b35521SBarry Smith for (k=0; k<nz; k++) { 10033f51a72SBarry Smith if (icol[j[k]] > prow && PetscAbsScalar(v[k]) > repla) { 101fd61f322SBarry Smith /* found a suitable later column */ 10233f51a72SBarry Smith repl = icol[j[k]]; 10333f51a72SBarry Smith SWAP(icol[col[prow]],icol[col[repl]]); 104ad608de0SBarry Smith SWAP(col[prow],col[repl]); 105fd61f322SBarry Smith goto found; 106fd61f322SBarry Smith } 107fd61f322SBarry Smith } 108fd61f322SBarry Smith /* 109fd61f322SBarry Smith Did not find a suitable later column so look for an earlier column 110fd61f322SBarry Smith We need to be sure that we don't introduce a zero in a previous 111fd61f322SBarry Smith diagonal 112fd61f322SBarry Smith */ 113fd61f322SBarry Smith for (k=0; k<nz; k++) { 114fd61f322SBarry Smith if (icol[j[k]] < prow && PetscAbsScalar(v[k]) > repla) { 115fd61f322SBarry Smith /* See if this one will work */ 116fd61f322SBarry Smith repl = icol[j[k]]; 117b3cc6726SBarry Smith ierr = MatGetRow_SeqAIJ(mat,row[repl],&nnz,&jj,&vv);CHKERRQ(ierr); 118fd61f322SBarry Smith for (kk=0; kk<nnz; kk++) { 11970441072SBarry Smith if (icol[jj[kk]] == prow && PetscAbsScalar(vv[kk]) > abstol) { 120b3cc6726SBarry Smith ierr = MatRestoreRow_SeqAIJ(mat,row[repl],&nnz,&jj,&vv);CHKERRQ(ierr); 121fd61f322SBarry Smith SWAP(icol[col[prow]],icol[col[repl]]); 122fd61f322SBarry Smith SWAP(col[prow],col[repl]); 123fd61f322SBarry Smith goto found; 124fd61f322SBarry Smith } 125fd61f322SBarry Smith } 126b3cc6726SBarry Smith ierr = MatRestoreRow_SeqAIJ(mat,row[repl],&nnz,&jj,&vv);CHKERRQ(ierr); 127fd61f322SBarry Smith } 128fd61f322SBarry Smith } 129fd61f322SBarry Smith /* 130fd61f322SBarry Smith No column suitable; instead check all future rows 131fd61f322SBarry Smith Note: this will be very slow 132fd61f322SBarry Smith */ 133fd61f322SBarry Smith for (k=prow+1; k<n; k++) { 134b3cc6726SBarry Smith ierr = MatGetRow_SeqAIJ(mat,row[k],&nnz,&jj,&vv);CHKERRQ(ierr); 135fd61f322SBarry Smith for (kk=0; kk<nnz; kk++) { 13670441072SBarry Smith if (icol[jj[kk]] == prow && PetscAbsScalar(vv[kk]) > abstol) { 137fd61f322SBarry Smith /* found a row */ 138fd61f322SBarry Smith SWAP(row[prow],row[k]); 139fd61f322SBarry Smith goto found; 140fd61f322SBarry Smith } 141fd61f322SBarry Smith } 142b3cc6726SBarry Smith ierr = MatRestoreRow_SeqAIJ(mat,row[k],&nnz,&jj,&vv);CHKERRQ(ierr); 143fd61f322SBarry Smith } 144fd61f322SBarry Smith 145fd61f322SBarry Smith found:; 146ad608de0SBarry Smith } 147b3cc6726SBarry Smith ierr = MatRestoreRow_SeqAIJ(mat,row[prow],&nz,&j,&v);CHKERRQ(ierr); 148ad608de0SBarry Smith } 14933f51a72SBarry Smith ierr = ISDestroy(icis);CHKERRQ(ierr); 1503a40ed3dSBarry Smith PetscFunctionReturn(0); 151ad608de0SBarry Smith } 15205b94e36SKris Buschelman EXTERN_C_END 15348b35521SBarry Smith 15448b35521SBarry Smith 155