11c2a3de1SBarry Smith 2397b6df1SKris Buschelman /* 3c2b5dc30SHong Zhang Provides an interface to the MUMPS sparse solver 4397b6df1SKris Buschelman */ 551d5961aSHong Zhang 6c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h> /*I "petscmat.h" I*/ 7c6db04a5SJed Brown #include <../src/mat/impls/sbaij/mpi/mpisbaij.h> 8397b6df1SKris Buschelman 9397b6df1SKris Buschelman EXTERN_C_BEGIN 10397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 112907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 122907cef9SHong Zhang #include <cmumps_c.h> 132907cef9SHong Zhang #else 14c6db04a5SJed Brown #include <zmumps_c.h> 152907cef9SHong Zhang #endif 162907cef9SHong Zhang #else 172907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 182907cef9SHong Zhang #include <smumps_c.h> 19397b6df1SKris Buschelman #else 20c6db04a5SJed Brown #include <dmumps_c.h> 21397b6df1SKris Buschelman #endif 222907cef9SHong Zhang #endif 23397b6df1SKris Buschelman EXTERN_C_END 24397b6df1SKris Buschelman #define JOB_INIT -1 253d472b54SHong Zhang #define JOB_FACTSYMBOLIC 1 263d472b54SHong Zhang #define JOB_FACTNUMERIC 2 273d472b54SHong Zhang #define JOB_SOLVE 3 28397b6df1SKris Buschelman #define JOB_END -2 293d472b54SHong Zhang 302907cef9SHong Zhang /* calls to MUMPS */ 312907cef9SHong Zhang #if defined(PETSC_USE_COMPLEX) 322907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 332907cef9SHong Zhang #define PetscMUMPS_c cmumps_c 342907cef9SHong Zhang #else 352907cef9SHong Zhang #define PetscMUMPS_c zmumps_c 362907cef9SHong Zhang #endif 372907cef9SHong Zhang #else 382907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 392907cef9SHong Zhang #define PetscMUMPS_c smumps_c 402907cef9SHong Zhang #else 412907cef9SHong Zhang #define PetscMUMPS_c dmumps_c 422907cef9SHong Zhang #endif 432907cef9SHong Zhang #endif 442907cef9SHong Zhang 453d472b54SHong Zhang 46397b6df1SKris Buschelman /* macros s.t. indices match MUMPS documentation */ 47397b6df1SKris Buschelman #define ICNTL(I) icntl[(I)-1] 48397b6df1SKris Buschelman #define CNTL(I) cntl[(I)-1] 49397b6df1SKris Buschelman #define INFOG(I) infog[(I)-1] 50a7aca84bSHong Zhang #define INFO(I) info[(I)-1] 51397b6df1SKris Buschelman #define RINFOG(I) rinfog[(I)-1] 52adc1d99fSHong Zhang #define RINFO(I) rinfo[(I)-1] 53397b6df1SKris Buschelman 54397b6df1SKris Buschelman typedef struct { 55397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 562907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 572907cef9SHong Zhang CMUMPS_STRUC_C id; 582907cef9SHong Zhang #else 59397b6df1SKris Buschelman ZMUMPS_STRUC_C id; 602907cef9SHong Zhang #endif 612907cef9SHong Zhang #else 622907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 632907cef9SHong Zhang SMUMPS_STRUC_C id; 64397b6df1SKris Buschelman #else 65397b6df1SKris Buschelman DMUMPS_STRUC_C id; 66397b6df1SKris Buschelman #endif 672907cef9SHong Zhang #endif 682907cef9SHong Zhang 69397b6df1SKris Buschelman MatStructure matstruc; 70c1490034SHong Zhang PetscMPIInt myid,size; 71a5e57a09SHong Zhang PetscInt *irn,*jcn,nz,sym; 72397b6df1SKris Buschelman PetscScalar *val; 73397b6df1SKris Buschelman MPI_Comm comm_mumps; 74329ec9b3SHong Zhang VecScatter scat_rhs, scat_sol; 7564e6c443SBarry Smith PetscBool isAIJ,CleanUpMUMPS; 76329ec9b3SHong Zhang Vec b_seq,x_seq; 77a5e57a09SHong Zhang PetscInt ICNTL9_pre; /* check if ICNTL(9) is changed from previous MatSolve */ 782205254eSKarl Rupp 79bf0cc555SLisandro Dalcin PetscErrorCode (*Destroy)(Mat); 80bccb9932SShri Abhyankar PetscErrorCode (*ConvertToTriples)(Mat, int, MatReuse, int*, int**, int**, PetscScalar**); 81f0c56d0fSKris Buschelman } Mat_MUMPS; 82f0c56d0fSKris Buschelman 8309573ac7SBarry Smith extern PetscErrorCode MatDuplicate_MUMPS(Mat,MatDuplicateOption,Mat*); 84b24902e0SBarry Smith 85397b6df1SKris Buschelman /* 86d341cd04SHong Zhang MatConvertToTriples_A_B - convert Petsc matrix to triples: row[nz], col[nz], val[nz] 87d341cd04SHong Zhang 88397b6df1SKris Buschelman input: 8967877ebaSShri Abhyankar A - matrix in aij,baij or sbaij (bs=1) format 90397b6df1SKris Buschelman shift - 0: C style output triple; 1: Fortran style output triple. 91bccb9932SShri Abhyankar reuse - MAT_INITIAL_MATRIX: spaces are allocated and values are set for the triple 92bccb9932SShri Abhyankar MAT_REUSE_MATRIX: only the values in v array are updated 93397b6df1SKris Buschelman output: 94397b6df1SKris Buschelman nnz - dim of r, c, and v (number of local nonzero entries of A) 95397b6df1SKris Buschelman r, c, v - row and col index, matrix values (matrix triples) 96eb9baa12SBarry Smith 97eb9baa12SBarry Smith The returned values r, c, and sometimes v are obtained in a single PetscMalloc(). Then in MatDestroy_MUMPS() it is 98eb9baa12SBarry Smith freed with PetscFree((mumps->irn); This is not ideal code, the fact that v is ONLY sometimes part of mumps->irn means 99eb9baa12SBarry Smith that the PetscMalloc() cannot easily be replaced with a PetscMalloc3(). 100eb9baa12SBarry Smith 101397b6df1SKris Buschelman */ 10216ebf90aSShri Abhyankar 10316ebf90aSShri Abhyankar #undef __FUNCT__ 10416ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqaij_seqaij" 105bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 106b24902e0SBarry Smith { 107185f6596SHong Zhang const PetscInt *ai,*aj,*ajj,M=A->rmap->n; 10867877ebaSShri Abhyankar PetscInt nz,rnz,i,j; 109dfbe8321SBarry Smith PetscErrorCode ierr; 110c1490034SHong Zhang PetscInt *row,*col; 11116ebf90aSShri Abhyankar Mat_SeqAIJ *aa=(Mat_SeqAIJ*)A->data; 112397b6df1SKris Buschelman 113397b6df1SKris Buschelman PetscFunctionBegin; 11416ebf90aSShri Abhyankar *v=aa->a; 115bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 1162205254eSKarl Rupp nz = aa->nz; 1172205254eSKarl Rupp ai = aa->i; 1182205254eSKarl Rupp aj = aa->j; 11916ebf90aSShri Abhyankar *nnz = nz; 120785e854fSJed Brown ierr = PetscMalloc1(2*nz, &row);CHKERRQ(ierr); 121185f6596SHong Zhang col = row + nz; 122185f6596SHong Zhang 12316ebf90aSShri Abhyankar nz = 0; 12416ebf90aSShri Abhyankar for (i=0; i<M; i++) { 12516ebf90aSShri Abhyankar rnz = ai[i+1] - ai[i]; 12667877ebaSShri Abhyankar ajj = aj + ai[i]; 12767877ebaSShri Abhyankar for (j=0; j<rnz; j++) { 12867877ebaSShri Abhyankar row[nz] = i+shift; col[nz++] = ajj[j] + shift; 12916ebf90aSShri Abhyankar } 13016ebf90aSShri Abhyankar } 13116ebf90aSShri Abhyankar *r = row; *c = col; 13216ebf90aSShri Abhyankar } 13316ebf90aSShri Abhyankar PetscFunctionReturn(0); 13416ebf90aSShri Abhyankar } 135397b6df1SKris Buschelman 13616ebf90aSShri Abhyankar #undef __FUNCT__ 13767877ebaSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqbaij_seqaij" 138bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqbaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 13967877ebaSShri Abhyankar { 14067877ebaSShri Abhyankar Mat_SeqBAIJ *aa=(Mat_SeqBAIJ*)A->data; 14133d57670SJed Brown const PetscInt *ai,*aj,*ajj,bs2 = aa->bs2; 14233d57670SJed Brown PetscInt bs,M,nz,idx=0,rnz,i,j,k,m; 14367877ebaSShri Abhyankar PetscErrorCode ierr; 14467877ebaSShri Abhyankar PetscInt *row,*col; 14567877ebaSShri Abhyankar 14667877ebaSShri Abhyankar PetscFunctionBegin; 14733d57670SJed Brown ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 14833d57670SJed Brown M = A->rmap->N/bs; 149cf3759fdSShri Abhyankar *v = aa->a; 150bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 151cf3759fdSShri Abhyankar ai = aa->i; aj = aa->j; 15267877ebaSShri Abhyankar nz = bs2*aa->nz; 15367877ebaSShri Abhyankar *nnz = nz; 154785e854fSJed Brown ierr = PetscMalloc1(2*nz, &row);CHKERRQ(ierr); 155185f6596SHong Zhang col = row + nz; 156185f6596SHong Zhang 15767877ebaSShri Abhyankar for (i=0; i<M; i++) { 15867877ebaSShri Abhyankar ajj = aj + ai[i]; 15967877ebaSShri Abhyankar rnz = ai[i+1] - ai[i]; 16067877ebaSShri Abhyankar for (k=0; k<rnz; k++) { 16167877ebaSShri Abhyankar for (j=0; j<bs; j++) { 16267877ebaSShri Abhyankar for (m=0; m<bs; m++) { 16367877ebaSShri Abhyankar row[idx] = i*bs + m + shift; 164cf3759fdSShri Abhyankar col[idx++] = bs*(ajj[k]) + j + shift; 16567877ebaSShri Abhyankar } 16667877ebaSShri Abhyankar } 16767877ebaSShri Abhyankar } 16867877ebaSShri Abhyankar } 169cf3759fdSShri Abhyankar *r = row; *c = col; 17067877ebaSShri Abhyankar } 17167877ebaSShri Abhyankar PetscFunctionReturn(0); 17267877ebaSShri Abhyankar } 17367877ebaSShri Abhyankar 17467877ebaSShri Abhyankar #undef __FUNCT__ 17516ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqsbaij_seqsbaij" 176bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqsbaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 17716ebf90aSShri Abhyankar { 17867877ebaSShri Abhyankar const PetscInt *ai, *aj,*ajj,M=A->rmap->n; 17967877ebaSShri Abhyankar PetscInt nz,rnz,i,j; 18016ebf90aSShri Abhyankar PetscErrorCode ierr; 18116ebf90aSShri Abhyankar PetscInt *row,*col; 18216ebf90aSShri Abhyankar Mat_SeqSBAIJ *aa=(Mat_SeqSBAIJ*)A->data; 18316ebf90aSShri Abhyankar 18416ebf90aSShri Abhyankar PetscFunctionBegin; 185882afa5aSHong Zhang *v = aa->a; 186bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 1872205254eSKarl Rupp nz = aa->nz; 1882205254eSKarl Rupp ai = aa->i; 1892205254eSKarl Rupp aj = aa->j; 1902205254eSKarl Rupp *v = aa->a; 19116ebf90aSShri Abhyankar *nnz = nz; 192785e854fSJed Brown ierr = PetscMalloc1(2*nz, &row);CHKERRQ(ierr); 193185f6596SHong Zhang col = row + nz; 194185f6596SHong Zhang 19516ebf90aSShri Abhyankar nz = 0; 19616ebf90aSShri Abhyankar for (i=0; i<M; i++) { 19716ebf90aSShri Abhyankar rnz = ai[i+1] - ai[i]; 19867877ebaSShri Abhyankar ajj = aj + ai[i]; 19967877ebaSShri Abhyankar for (j=0; j<rnz; j++) { 20067877ebaSShri Abhyankar row[nz] = i+shift; col[nz++] = ajj[j] + shift; 20116ebf90aSShri Abhyankar } 20216ebf90aSShri Abhyankar } 20316ebf90aSShri Abhyankar *r = row; *c = col; 20416ebf90aSShri Abhyankar } 20516ebf90aSShri Abhyankar PetscFunctionReturn(0); 20616ebf90aSShri Abhyankar } 20716ebf90aSShri Abhyankar 20816ebf90aSShri Abhyankar #undef __FUNCT__ 20916ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqaij_seqsbaij" 210bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 21116ebf90aSShri Abhyankar { 21267877ebaSShri Abhyankar const PetscInt *ai,*aj,*ajj,*adiag,M=A->rmap->n; 21367877ebaSShri Abhyankar PetscInt nz,rnz,i,j; 21467877ebaSShri Abhyankar const PetscScalar *av,*v1; 21516ebf90aSShri Abhyankar PetscScalar *val; 21616ebf90aSShri Abhyankar PetscErrorCode ierr; 21716ebf90aSShri Abhyankar PetscInt *row,*col; 218829b1710SHong Zhang Mat_SeqAIJ *aa=(Mat_SeqAIJ*)A->data; 21916ebf90aSShri Abhyankar 22016ebf90aSShri Abhyankar PetscFunctionBegin; 22116ebf90aSShri Abhyankar ai =aa->i; aj=aa->j;av=aa->a; 22216ebf90aSShri Abhyankar adiag=aa->diag; 223bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 224829b1710SHong Zhang /* count nz in the uppper triangular part of A */ 225829b1710SHong Zhang nz = 0; 226829b1710SHong Zhang for (i=0; i<M; i++) nz += ai[i+1] - adiag[i]; 22716ebf90aSShri Abhyankar *nnz = nz; 228829b1710SHong Zhang 229185f6596SHong Zhang ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr); 230185f6596SHong Zhang col = row + nz; 231185f6596SHong Zhang val = (PetscScalar*)(col + nz); 232185f6596SHong Zhang 23316ebf90aSShri Abhyankar nz = 0; 23416ebf90aSShri Abhyankar for (i=0; i<M; i++) { 23516ebf90aSShri Abhyankar rnz = ai[i+1] - adiag[i]; 23667877ebaSShri Abhyankar ajj = aj + adiag[i]; 237cf3759fdSShri Abhyankar v1 = av + adiag[i]; 23867877ebaSShri Abhyankar for (j=0; j<rnz; j++) { 23967877ebaSShri Abhyankar row[nz] = i+shift; col[nz] = ajj[j] + shift; val[nz++] = v1[j]; 24016ebf90aSShri Abhyankar } 24116ebf90aSShri Abhyankar } 24216ebf90aSShri Abhyankar *r = row; *c = col; *v = val; 243397b6df1SKris Buschelman } else { 24416ebf90aSShri Abhyankar nz = 0; val = *v; 24516ebf90aSShri Abhyankar for (i=0; i <M; i++) { 24616ebf90aSShri Abhyankar rnz = ai[i+1] - adiag[i]; 24767877ebaSShri Abhyankar ajj = aj + adiag[i]; 24867877ebaSShri Abhyankar v1 = av + adiag[i]; 24967877ebaSShri Abhyankar for (j=0; j<rnz; j++) { 25067877ebaSShri Abhyankar val[nz++] = v1[j]; 25116ebf90aSShri Abhyankar } 25216ebf90aSShri Abhyankar } 25316ebf90aSShri Abhyankar } 25416ebf90aSShri Abhyankar PetscFunctionReturn(0); 25516ebf90aSShri Abhyankar } 25616ebf90aSShri Abhyankar 25716ebf90aSShri Abhyankar #undef __FUNCT__ 25816ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpisbaij_mpisbaij" 259bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpisbaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 26016ebf90aSShri Abhyankar { 26116ebf90aSShri Abhyankar const PetscInt *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj; 26216ebf90aSShri Abhyankar PetscErrorCode ierr; 26316ebf90aSShri Abhyankar PetscInt rstart,nz,i,j,jj,irow,countA,countB; 26416ebf90aSShri Abhyankar PetscInt *row,*col; 26516ebf90aSShri Abhyankar const PetscScalar *av, *bv,*v1,*v2; 26616ebf90aSShri Abhyankar PetscScalar *val; 267397b6df1SKris Buschelman Mat_MPISBAIJ *mat = (Mat_MPISBAIJ*)A->data; 268397b6df1SKris Buschelman Mat_SeqSBAIJ *aa = (Mat_SeqSBAIJ*)(mat->A)->data; 269397b6df1SKris Buschelman Mat_SeqBAIJ *bb = (Mat_SeqBAIJ*)(mat->B)->data; 27016ebf90aSShri Abhyankar 27116ebf90aSShri Abhyankar PetscFunctionBegin; 272d0f46423SBarry Smith ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart; 273397b6df1SKris Buschelman av=aa->a; bv=bb->a; 274397b6df1SKris Buschelman 2752205254eSKarl Rupp garray = mat->garray; 2762205254eSKarl Rupp 277bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 27816ebf90aSShri Abhyankar nz = aa->nz + bb->nz; 27916ebf90aSShri Abhyankar *nnz = nz; 280185f6596SHong Zhang ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr); 281185f6596SHong Zhang col = row + nz; 282185f6596SHong Zhang val = (PetscScalar*)(col + nz); 283185f6596SHong Zhang 284397b6df1SKris Buschelman *r = row; *c = col; *v = val; 285397b6df1SKris Buschelman } else { 286397b6df1SKris Buschelman row = *r; col = *c; val = *v; 287397b6df1SKris Buschelman } 288397b6df1SKris Buschelman 289028e57e8SHong Zhang jj = 0; irow = rstart; 290397b6df1SKris Buschelman for (i=0; i<m; i++) { 291397b6df1SKris Buschelman ajj = aj + ai[i]; /* ptr to the beginning of this row */ 292397b6df1SKris Buschelman countA = ai[i+1] - ai[i]; 293397b6df1SKris Buschelman countB = bi[i+1] - bi[i]; 294397b6df1SKris Buschelman bjj = bj + bi[i]; 29516ebf90aSShri Abhyankar v1 = av + ai[i]; 29616ebf90aSShri Abhyankar v2 = bv + bi[i]; 297397b6df1SKris Buschelman 298397b6df1SKris Buschelman /* A-part */ 299397b6df1SKris Buschelman for (j=0; j<countA; j++) { 300bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 301397b6df1SKris Buschelman row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift; 302397b6df1SKris Buschelman } 30316ebf90aSShri Abhyankar val[jj++] = v1[j]; 304397b6df1SKris Buschelman } 30516ebf90aSShri Abhyankar 30616ebf90aSShri Abhyankar /* B-part */ 30716ebf90aSShri Abhyankar for (j=0; j < countB; j++) { 308bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 309397b6df1SKris Buschelman row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift; 310397b6df1SKris Buschelman } 31116ebf90aSShri Abhyankar val[jj++] = v2[j]; 31216ebf90aSShri Abhyankar } 31316ebf90aSShri Abhyankar irow++; 31416ebf90aSShri Abhyankar } 31516ebf90aSShri Abhyankar PetscFunctionReturn(0); 31616ebf90aSShri Abhyankar } 31716ebf90aSShri Abhyankar 31816ebf90aSShri Abhyankar #undef __FUNCT__ 31916ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpiaij_mpiaij" 320bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpiaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 32116ebf90aSShri Abhyankar { 32216ebf90aSShri Abhyankar const PetscInt *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj; 32316ebf90aSShri Abhyankar PetscErrorCode ierr; 32416ebf90aSShri Abhyankar PetscInt rstart,nz,i,j,jj,irow,countA,countB; 32516ebf90aSShri Abhyankar PetscInt *row,*col; 32616ebf90aSShri Abhyankar const PetscScalar *av, *bv,*v1,*v2; 32716ebf90aSShri Abhyankar PetscScalar *val; 32816ebf90aSShri Abhyankar Mat_MPIAIJ *mat = (Mat_MPIAIJ*)A->data; 32916ebf90aSShri Abhyankar Mat_SeqAIJ *aa = (Mat_SeqAIJ*)(mat->A)->data; 33016ebf90aSShri Abhyankar Mat_SeqAIJ *bb = (Mat_SeqAIJ*)(mat->B)->data; 33116ebf90aSShri Abhyankar 33216ebf90aSShri Abhyankar PetscFunctionBegin; 33316ebf90aSShri Abhyankar ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart; 33416ebf90aSShri Abhyankar av=aa->a; bv=bb->a; 33516ebf90aSShri Abhyankar 3362205254eSKarl Rupp garray = mat->garray; 3372205254eSKarl Rupp 338bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 33916ebf90aSShri Abhyankar nz = aa->nz + bb->nz; 34016ebf90aSShri Abhyankar *nnz = nz; 341185f6596SHong Zhang ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr); 342185f6596SHong Zhang col = row + nz; 343185f6596SHong Zhang val = (PetscScalar*)(col + nz); 344185f6596SHong Zhang 34516ebf90aSShri Abhyankar *r = row; *c = col; *v = val; 34616ebf90aSShri Abhyankar } else { 34716ebf90aSShri Abhyankar row = *r; col = *c; val = *v; 34816ebf90aSShri Abhyankar } 34916ebf90aSShri Abhyankar 35016ebf90aSShri Abhyankar jj = 0; irow = rstart; 35116ebf90aSShri Abhyankar for (i=0; i<m; i++) { 35216ebf90aSShri Abhyankar ajj = aj + ai[i]; /* ptr to the beginning of this row */ 35316ebf90aSShri Abhyankar countA = ai[i+1] - ai[i]; 35416ebf90aSShri Abhyankar countB = bi[i+1] - bi[i]; 35516ebf90aSShri Abhyankar bjj = bj + bi[i]; 35616ebf90aSShri Abhyankar v1 = av + ai[i]; 35716ebf90aSShri Abhyankar v2 = bv + bi[i]; 35816ebf90aSShri Abhyankar 35916ebf90aSShri Abhyankar /* A-part */ 36016ebf90aSShri Abhyankar for (j=0; j<countA; j++) { 361bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 36216ebf90aSShri Abhyankar row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift; 36316ebf90aSShri Abhyankar } 36416ebf90aSShri Abhyankar val[jj++] = v1[j]; 36516ebf90aSShri Abhyankar } 36616ebf90aSShri Abhyankar 36716ebf90aSShri Abhyankar /* B-part */ 36816ebf90aSShri Abhyankar for (j=0; j < countB; j++) { 369bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 37016ebf90aSShri Abhyankar row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift; 37116ebf90aSShri Abhyankar } 37216ebf90aSShri Abhyankar val[jj++] = v2[j]; 37316ebf90aSShri Abhyankar } 37416ebf90aSShri Abhyankar irow++; 37516ebf90aSShri Abhyankar } 37616ebf90aSShri Abhyankar PetscFunctionReturn(0); 37716ebf90aSShri Abhyankar } 37816ebf90aSShri Abhyankar 37916ebf90aSShri Abhyankar #undef __FUNCT__ 38067877ebaSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpibaij_mpiaij" 381bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpibaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 38267877ebaSShri Abhyankar { 38367877ebaSShri Abhyankar Mat_MPIBAIJ *mat = (Mat_MPIBAIJ*)A->data; 38467877ebaSShri Abhyankar Mat_SeqBAIJ *aa = (Mat_SeqBAIJ*)(mat->A)->data; 38567877ebaSShri Abhyankar Mat_SeqBAIJ *bb = (Mat_SeqBAIJ*)(mat->B)->data; 38667877ebaSShri Abhyankar const PetscInt *ai = aa->i, *bi = bb->i, *aj = aa->j, *bj = bb->j,*ajj, *bjj; 387d985c460SShri Abhyankar const PetscInt *garray = mat->garray,mbs=mat->mbs,rstart=A->rmap->rstart; 38833d57670SJed Brown const PetscInt bs2=mat->bs2; 38967877ebaSShri Abhyankar PetscErrorCode ierr; 39033d57670SJed Brown PetscInt bs,nz,i,j,k,n,jj,irow,countA,countB,idx; 39167877ebaSShri Abhyankar PetscInt *row,*col; 39267877ebaSShri Abhyankar const PetscScalar *av=aa->a, *bv=bb->a,*v1,*v2; 39367877ebaSShri Abhyankar PetscScalar *val; 39467877ebaSShri Abhyankar 39567877ebaSShri Abhyankar PetscFunctionBegin; 39633d57670SJed Brown ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 397bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 39867877ebaSShri Abhyankar nz = bs2*(aa->nz + bb->nz); 39967877ebaSShri Abhyankar *nnz = nz; 400185f6596SHong Zhang ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr); 401185f6596SHong Zhang col = row + nz; 402185f6596SHong Zhang val = (PetscScalar*)(col + nz); 403185f6596SHong Zhang 40467877ebaSShri Abhyankar *r = row; *c = col; *v = val; 40567877ebaSShri Abhyankar } else { 40667877ebaSShri Abhyankar row = *r; col = *c; val = *v; 40767877ebaSShri Abhyankar } 40867877ebaSShri Abhyankar 409d985c460SShri Abhyankar jj = 0; irow = rstart; 41067877ebaSShri Abhyankar for (i=0; i<mbs; i++) { 41167877ebaSShri Abhyankar countA = ai[i+1] - ai[i]; 41267877ebaSShri Abhyankar countB = bi[i+1] - bi[i]; 41367877ebaSShri Abhyankar ajj = aj + ai[i]; 41467877ebaSShri Abhyankar bjj = bj + bi[i]; 41567877ebaSShri Abhyankar v1 = av + bs2*ai[i]; 41667877ebaSShri Abhyankar v2 = bv + bs2*bi[i]; 41767877ebaSShri Abhyankar 41867877ebaSShri Abhyankar idx = 0; 41967877ebaSShri Abhyankar /* A-part */ 42067877ebaSShri Abhyankar for (k=0; k<countA; k++) { 42167877ebaSShri Abhyankar for (j=0; j<bs; j++) { 42267877ebaSShri Abhyankar for (n=0; n<bs; n++) { 423bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 424d985c460SShri Abhyankar row[jj] = irow + n + shift; 425d985c460SShri Abhyankar col[jj] = rstart + bs*ajj[k] + j + shift; 42667877ebaSShri Abhyankar } 42767877ebaSShri Abhyankar val[jj++] = v1[idx++]; 42867877ebaSShri Abhyankar } 42967877ebaSShri Abhyankar } 43067877ebaSShri Abhyankar } 43167877ebaSShri Abhyankar 43267877ebaSShri Abhyankar idx = 0; 43367877ebaSShri Abhyankar /* B-part */ 43467877ebaSShri Abhyankar for (k=0; k<countB; k++) { 43567877ebaSShri Abhyankar for (j=0; j<bs; j++) { 43667877ebaSShri Abhyankar for (n=0; n<bs; n++) { 437bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 438d985c460SShri Abhyankar row[jj] = irow + n + shift; 439d985c460SShri Abhyankar col[jj] = bs*garray[bjj[k]] + j + shift; 44067877ebaSShri Abhyankar } 441d985c460SShri Abhyankar val[jj++] = v2[idx++]; 44267877ebaSShri Abhyankar } 44367877ebaSShri Abhyankar } 44467877ebaSShri Abhyankar } 445d985c460SShri Abhyankar irow += bs; 44667877ebaSShri Abhyankar } 44767877ebaSShri Abhyankar PetscFunctionReturn(0); 44867877ebaSShri Abhyankar } 44967877ebaSShri Abhyankar 45067877ebaSShri Abhyankar #undef __FUNCT__ 45116ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpiaij_mpisbaij" 452bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpiaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 45316ebf90aSShri Abhyankar { 45416ebf90aSShri Abhyankar const PetscInt *ai, *aj,*adiag, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj; 45516ebf90aSShri Abhyankar PetscErrorCode ierr; 456e0bace9bSHong Zhang PetscInt rstart,nz,nza,nzb,i,j,jj,irow,countA,countB; 45716ebf90aSShri Abhyankar PetscInt *row,*col; 45816ebf90aSShri Abhyankar const PetscScalar *av, *bv,*v1,*v2; 45916ebf90aSShri Abhyankar PetscScalar *val; 46016ebf90aSShri Abhyankar Mat_MPIAIJ *mat = (Mat_MPIAIJ*)A->data; 46116ebf90aSShri Abhyankar Mat_SeqAIJ *aa =(Mat_SeqAIJ*)(mat->A)->data; 46216ebf90aSShri Abhyankar Mat_SeqAIJ *bb =(Mat_SeqAIJ*)(mat->B)->data; 46316ebf90aSShri Abhyankar 46416ebf90aSShri Abhyankar PetscFunctionBegin; 46516ebf90aSShri Abhyankar ai=aa->i; aj=aa->j; adiag=aa->diag; 46616ebf90aSShri Abhyankar bi=bb->i; bj=bb->j; garray = mat->garray; 46716ebf90aSShri Abhyankar av=aa->a; bv=bb->a; 4682205254eSKarl Rupp 46916ebf90aSShri Abhyankar rstart = A->rmap->rstart; 47016ebf90aSShri Abhyankar 471bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 472e0bace9bSHong Zhang nza = 0; /* num of upper triangular entries in mat->A, including diagonals */ 473e0bace9bSHong Zhang nzb = 0; /* num of upper triangular entries in mat->B */ 47416ebf90aSShri Abhyankar for (i=0; i<m; i++) { 475e0bace9bSHong Zhang nza += (ai[i+1] - adiag[i]); 47616ebf90aSShri Abhyankar countB = bi[i+1] - bi[i]; 47716ebf90aSShri Abhyankar bjj = bj + bi[i]; 478e0bace9bSHong Zhang for (j=0; j<countB; j++) { 479e0bace9bSHong Zhang if (garray[bjj[j]] > rstart) nzb++; 480e0bace9bSHong Zhang } 481e0bace9bSHong Zhang } 48216ebf90aSShri Abhyankar 483e0bace9bSHong Zhang nz = nza + nzb; /* total nz of upper triangular part of mat */ 48416ebf90aSShri Abhyankar *nnz = nz; 485185f6596SHong Zhang ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr); 486185f6596SHong Zhang col = row + nz; 487185f6596SHong Zhang val = (PetscScalar*)(col + nz); 488185f6596SHong Zhang 48916ebf90aSShri Abhyankar *r = row; *c = col; *v = val; 49016ebf90aSShri Abhyankar } else { 49116ebf90aSShri Abhyankar row = *r; col = *c; val = *v; 49216ebf90aSShri Abhyankar } 49316ebf90aSShri Abhyankar 49416ebf90aSShri Abhyankar jj = 0; irow = rstart; 49516ebf90aSShri Abhyankar for (i=0; i<m; i++) { 49616ebf90aSShri Abhyankar ajj = aj + adiag[i]; /* ptr to the beginning of the diagonal of this row */ 49716ebf90aSShri Abhyankar v1 = av + adiag[i]; 49816ebf90aSShri Abhyankar countA = ai[i+1] - adiag[i]; 49916ebf90aSShri Abhyankar countB = bi[i+1] - bi[i]; 50016ebf90aSShri Abhyankar bjj = bj + bi[i]; 50116ebf90aSShri Abhyankar v2 = bv + bi[i]; 50216ebf90aSShri Abhyankar 50316ebf90aSShri Abhyankar /* A-part */ 50416ebf90aSShri Abhyankar for (j=0; j<countA; j++) { 505bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 50616ebf90aSShri Abhyankar row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift; 50716ebf90aSShri Abhyankar } 50816ebf90aSShri Abhyankar val[jj++] = v1[j]; 50916ebf90aSShri Abhyankar } 51016ebf90aSShri Abhyankar 51116ebf90aSShri Abhyankar /* B-part */ 51216ebf90aSShri Abhyankar for (j=0; j < countB; j++) { 51316ebf90aSShri Abhyankar if (garray[bjj[j]] > rstart) { 514bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 51516ebf90aSShri Abhyankar row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift; 51616ebf90aSShri Abhyankar } 51716ebf90aSShri Abhyankar val[jj++] = v2[j]; 51816ebf90aSShri Abhyankar } 519397b6df1SKris Buschelman } 520397b6df1SKris Buschelman irow++; 521397b6df1SKris Buschelman } 522397b6df1SKris Buschelman PetscFunctionReturn(0); 523397b6df1SKris Buschelman } 524397b6df1SKris Buschelman 525397b6df1SKris Buschelman #undef __FUNCT__ 52620be8e61SHong Zhang #define __FUNCT__ "MatGetDiagonal_MUMPS" 52720be8e61SHong Zhang PetscErrorCode MatGetDiagonal_MUMPS(Mat A,Vec v) 52820be8e61SHong Zhang { 52920be8e61SHong Zhang PetscFunctionBegin; 53020be8e61SHong Zhang SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Mat type: MUMPS factor"); 53120be8e61SHong Zhang PetscFunctionReturn(0); 53220be8e61SHong Zhang } 53320be8e61SHong Zhang 53420be8e61SHong Zhang #undef __FUNCT__ 5353924e44cSKris Buschelman #define __FUNCT__ "MatDestroy_MUMPS" 536dfbe8321SBarry Smith PetscErrorCode MatDestroy_MUMPS(Mat A) 537dfbe8321SBarry Smith { 538a5e57a09SHong Zhang Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 539dfbe8321SBarry Smith PetscErrorCode ierr; 540b24902e0SBarry Smith 541397b6df1SKris Buschelman PetscFunctionBegin; 542a5e57a09SHong Zhang if (mumps->CleanUpMUMPS) { 543397b6df1SKris Buschelman /* Terminate instance, deallocate memories */ 544a5e57a09SHong Zhang ierr = PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);CHKERRQ(ierr); 545a5e57a09SHong Zhang ierr = VecScatterDestroy(&mumps->scat_rhs);CHKERRQ(ierr); 546a5e57a09SHong Zhang ierr = VecDestroy(&mumps->b_seq);CHKERRQ(ierr); 547a5e57a09SHong Zhang ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr); 548a5e57a09SHong Zhang ierr = VecDestroy(&mumps->x_seq);CHKERRQ(ierr); 549a5e57a09SHong Zhang ierr = PetscFree(mumps->id.perm_in);CHKERRQ(ierr); 550a5e57a09SHong Zhang ierr = PetscFree(mumps->irn);CHKERRQ(ierr); 5512205254eSKarl Rupp 552a5e57a09SHong Zhang mumps->id.job = JOB_END; 553a5e57a09SHong Zhang PetscMUMPS_c(&mumps->id); 554a5e57a09SHong Zhang ierr = MPI_Comm_free(&(mumps->comm_mumps));CHKERRQ(ierr); 555397b6df1SKris Buschelman } 556a5e57a09SHong Zhang if (mumps->Destroy) { 557a5e57a09SHong Zhang ierr = (mumps->Destroy)(A);CHKERRQ(ierr); 558bf0cc555SLisandro Dalcin } 559bf0cc555SLisandro Dalcin ierr = PetscFree(A->spptr);CHKERRQ(ierr); 560bf0cc555SLisandro Dalcin 56197969023SHong Zhang /* clear composed functions */ 562bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverPackage_C",NULL);CHKERRQ(ierr); 563bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetIcntl_C",NULL);CHKERRQ(ierr); 564bc6112feSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetIcntl_C",NULL);CHKERRQ(ierr); 565bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetCntl_C",NULL);CHKERRQ(ierr); 566bc6112feSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetCntl_C",NULL);CHKERRQ(ierr); 567bc6112feSHong Zhang 568ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfo_C",NULL);CHKERRQ(ierr); 569ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfog_C",NULL);CHKERRQ(ierr); 570ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfo_C",NULL);CHKERRQ(ierr); 571ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfog_C",NULL);CHKERRQ(ierr); 572397b6df1SKris Buschelman PetscFunctionReturn(0); 573397b6df1SKris Buschelman } 574397b6df1SKris Buschelman 575397b6df1SKris Buschelman #undef __FUNCT__ 576f6c57405SHong Zhang #define __FUNCT__ "MatSolve_MUMPS" 577b24902e0SBarry Smith PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x) 578b24902e0SBarry Smith { 579a5e57a09SHong Zhang Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 580d54de34fSKris Buschelman PetscScalar *array; 58167877ebaSShri Abhyankar Vec b_seq; 582329ec9b3SHong Zhang IS is_iden,is_petsc; 583dfbe8321SBarry Smith PetscErrorCode ierr; 584329ec9b3SHong Zhang PetscInt i; 585883f2eb9SBarry Smith static PetscBool cite1 = PETSC_FALSE,cite2 = PETSC_FALSE; 586397b6df1SKris Buschelman 587397b6df1SKris Buschelman PetscFunctionBegin; 588883f2eb9SBarry Smith ierr = PetscCitationsRegister("@article{MUMPS01,\n author = {P.~R. Amestoy and I.~S. Duff and J.-Y. L'Excellent and J. Koster},\n title = {A fully asynchronous multifrontal solver using distributed dynamic scheduling},\n journal = {SIAM Journal on Matrix Analysis and Applications},\n volume = {23},\n number = {1},\n pages = {15--41},\n year = {2001}\n}\n",&cite1);CHKERRQ(ierr); 589883f2eb9SBarry Smith ierr = PetscCitationsRegister("@article{MUMPS02,\n author = {P.~R. Amestoy and A. Guermouche and J.-Y. L'Excellent and S. Pralet},\n title = {Hybrid scheduling for the parallel solution of linear systems},\n journal = {Parallel Computing},\n volume = {32},\n number = {2},\n pages = {136--156},\n year = {2006}\n}\n",&cite2);CHKERRQ(ierr); 590a5e57a09SHong Zhang mumps->id.nrhs = 1; 591a5e57a09SHong Zhang b_seq = mumps->b_seq; 592a5e57a09SHong Zhang if (mumps->size > 1) { 593329ec9b3SHong Zhang /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */ 594a5e57a09SHong Zhang ierr = VecScatterBegin(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 595a5e57a09SHong Zhang ierr = VecScatterEnd(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 596a5e57a09SHong Zhang if (!mumps->myid) {ierr = VecGetArray(b_seq,&array);CHKERRQ(ierr);} 597397b6df1SKris Buschelman } else { /* size == 1 */ 598397b6df1SKris Buschelman ierr = VecCopy(b,x);CHKERRQ(ierr); 599397b6df1SKris Buschelman ierr = VecGetArray(x,&array);CHKERRQ(ierr); 600397b6df1SKris Buschelman } 601a5e57a09SHong Zhang if (!mumps->myid) { /* define rhs on the host */ 602a5e57a09SHong Zhang mumps->id.nrhs = 1; 603397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 6042907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 605a5e57a09SHong Zhang mumps->id.rhs = (mumps_complex*)array; 6062907cef9SHong Zhang #else 607a5e57a09SHong Zhang mumps->id.rhs = (mumps_double_complex*)array; 6082907cef9SHong Zhang #endif 609397b6df1SKris Buschelman #else 610a5e57a09SHong Zhang mumps->id.rhs = array; 611397b6df1SKris Buschelman #endif 612397b6df1SKris Buschelman } 613397b6df1SKris Buschelman 614397b6df1SKris Buschelman /* solve phase */ 615329ec9b3SHong Zhang /*-------------*/ 616a5e57a09SHong Zhang mumps->id.job = JOB_SOLVE; 617a5e57a09SHong Zhang PetscMUMPS_c(&mumps->id); 618a5e57a09SHong Zhang if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",mumps->id.INFOG(1)); 619397b6df1SKris Buschelman 620a5e57a09SHong Zhang if (mumps->size > 1) { /* convert mumps distributed solution to petsc mpi x */ 621a5e57a09SHong Zhang if (mumps->scat_sol && mumps->ICNTL9_pre != mumps->id.ICNTL(9)) { 622a5e57a09SHong Zhang /* when id.ICNTL(9) changes, the contents of lsol_loc may change (not its size, lsol_loc), recreates scat_sol */ 623a5e57a09SHong Zhang ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr); 624397b6df1SKris Buschelman } 625a5e57a09SHong Zhang if (!mumps->scat_sol) { /* create scatter scat_sol */ 626a5e57a09SHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,mumps->id.lsol_loc,0,1,&is_iden);CHKERRQ(ierr); /* from */ 627a5e57a09SHong Zhang for (i=0; i<mumps->id.lsol_loc; i++) { 628a5e57a09SHong Zhang mumps->id.isol_loc[i] -= 1; /* change Fortran style to C style */ 629a5e57a09SHong Zhang } 630a5e57a09SHong Zhang ierr = ISCreateGeneral(PETSC_COMM_SELF,mumps->id.lsol_loc,mumps->id.isol_loc,PETSC_COPY_VALUES,&is_petsc);CHKERRQ(ierr); /* to */ 631a5e57a09SHong Zhang ierr = VecScatterCreate(mumps->x_seq,is_iden,x,is_petsc,&mumps->scat_sol);CHKERRQ(ierr); 6326bf464f9SBarry Smith ierr = ISDestroy(&is_iden);CHKERRQ(ierr); 6336bf464f9SBarry Smith ierr = ISDestroy(&is_petsc);CHKERRQ(ierr); 6342205254eSKarl Rupp 635a5e57a09SHong Zhang mumps->ICNTL9_pre = mumps->id.ICNTL(9); /* save current value of id.ICNTL(9) */ 636397b6df1SKris Buschelman } 637a5e57a09SHong Zhang 638a5e57a09SHong Zhang ierr = VecScatterBegin(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 639a5e57a09SHong Zhang ierr = VecScatterEnd(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 640329ec9b3SHong Zhang } 641397b6df1SKris Buschelman PetscFunctionReturn(0); 642397b6df1SKris Buschelman } 643397b6df1SKris Buschelman 64451d5961aSHong Zhang #undef __FUNCT__ 64551d5961aSHong Zhang #define __FUNCT__ "MatSolveTranspose_MUMPS" 64651d5961aSHong Zhang PetscErrorCode MatSolveTranspose_MUMPS(Mat A,Vec b,Vec x) 64751d5961aSHong Zhang { 648a5e57a09SHong Zhang Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 64951d5961aSHong Zhang PetscErrorCode ierr; 65051d5961aSHong Zhang 65151d5961aSHong Zhang PetscFunctionBegin; 652a5e57a09SHong Zhang mumps->id.ICNTL(9) = 0; 6532205254eSKarl Rupp 6540ad0caddSJed Brown ierr = MatSolve_MUMPS(A,b,x);CHKERRQ(ierr); 6552205254eSKarl Rupp 656a5e57a09SHong Zhang mumps->id.ICNTL(9) = 1; 65751d5961aSHong Zhang PetscFunctionReturn(0); 65851d5961aSHong Zhang } 65951d5961aSHong Zhang 660e0b74bf9SHong Zhang #undef __FUNCT__ 661e0b74bf9SHong Zhang #define __FUNCT__ "MatMatSolve_MUMPS" 662e0b74bf9SHong Zhang PetscErrorCode MatMatSolve_MUMPS(Mat A,Mat B,Mat X) 663e0b74bf9SHong Zhang { 664bda8bf91SBarry Smith PetscErrorCode ierr; 665bda8bf91SBarry Smith PetscBool flg; 666*4e34a73bSHong Zhang Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 667*4e34a73bSHong Zhang PetscInt nrhs; 668bda8bf91SBarry Smith 669e0b74bf9SHong Zhang PetscFunctionBegin; 6700298fd71SBarry Smith ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 671ce94432eSBarry Smith if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix"); 6720298fd71SBarry Smith ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 673ce94432eSBarry Smith if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix"); 674*4e34a73bSHong Zhang 675*4e34a73bSHong Zhang ierr = MatGetSize(B,NULL,&nrhs);CHKERRQ(ierr); 676*4e34a73bSHong Zhang mumps->id.ICNTL(27) = nrhs; 677*4e34a73bSHong Zhang printf("mumps_icntl_27 = %d\n",mumps->id.ICNTL(27)); 678*4e34a73bSHong Zhang 6792205254eSKarl Rupp SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatSolve_MUMPS() is not implemented yet"); 680e0b74bf9SHong Zhang PetscFunctionReturn(0); 681e0b74bf9SHong Zhang } 682e0b74bf9SHong Zhang 683ace3df97SHong Zhang #if !defined(PETSC_USE_COMPLEX) 684a58c3f20SHong Zhang /* 685a58c3f20SHong Zhang input: 686a58c3f20SHong Zhang F: numeric factor 687a58c3f20SHong Zhang output: 688a58c3f20SHong Zhang nneg: total number of negative pivots 689a58c3f20SHong Zhang nzero: 0 690a58c3f20SHong Zhang npos: (global dimension of F) - nneg 691a58c3f20SHong Zhang */ 692a58c3f20SHong Zhang 693a58c3f20SHong Zhang #undef __FUNCT__ 694a58c3f20SHong Zhang #define __FUNCT__ "MatGetInertia_SBAIJMUMPS" 695dfbe8321SBarry Smith PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos) 696a58c3f20SHong Zhang { 697a5e57a09SHong Zhang Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 698dfbe8321SBarry Smith PetscErrorCode ierr; 699c1490034SHong Zhang PetscMPIInt size; 700a58c3f20SHong Zhang 701a58c3f20SHong Zhang PetscFunctionBegin; 702ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)F),&size);CHKERRQ(ierr); 703bcb30aebSHong Zhang /* MUMPS 4.3.1 calls ScaLAPACK when ICNTL(13)=0 (default), which does not offer the possibility to compute the inertia of a dense matrix. Set ICNTL(13)=1 to skip ScaLAPACK */ 704a5e57a09SHong Zhang if (size > 1 && mumps->id.ICNTL(13) != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"ICNTL(13)=%d. -mat_mumps_icntl_13 must be set as 1 for correct global matrix inertia\n",mumps->id.INFOG(13)); 705ed85ac9fSHong Zhang 706710ac8efSHong Zhang if (nneg) *nneg = mumps->id.INFOG(12); 707ed85ac9fSHong Zhang if (nzero || npos) { 708ed85ac9fSHong Zhang if (mumps->id.ICNTL(24) != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"-mat_mumps_icntl_24 must be set as 1 for null pivot row detection"); 709710ac8efSHong Zhang if (nzero) *nzero = mumps->id.INFOG(28); 710710ac8efSHong Zhang if (npos) *npos = F->rmap->N - (mumps->id.INFOG(12) + mumps->id.INFOG(28)); 711a58c3f20SHong Zhang } 712a58c3f20SHong Zhang PetscFunctionReturn(0); 713a58c3f20SHong Zhang } 714ace3df97SHong Zhang #endif /* !defined(PETSC_USE_COMPLEX) */ 715a58c3f20SHong Zhang 716397b6df1SKris Buschelman #undef __FUNCT__ 717f6c57405SHong Zhang #define __FUNCT__ "MatFactorNumeric_MUMPS" 7180481f469SBarry Smith PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info) 719af281ebdSHong Zhang { 720a5e57a09SHong Zhang Mat_MUMPS *mumps =(Mat_MUMPS*)(F)->spptr; 7216849ba73SBarry Smith PetscErrorCode ierr; 722e09efc27SHong Zhang Mat F_diag; 723ace3abfcSBarry Smith PetscBool isMPIAIJ; 724397b6df1SKris Buschelman 725397b6df1SKris Buschelman PetscFunctionBegin; 726a5e57a09SHong Zhang ierr = (*mumps->ConvertToTriples)(A, 1, MAT_REUSE_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr); 727397b6df1SKris Buschelman 728397b6df1SKris Buschelman /* numerical factorization phase */ 729329ec9b3SHong Zhang /*-------------------------------*/ 730a5e57a09SHong Zhang mumps->id.job = JOB_FACTNUMERIC; 731*4e34a73bSHong Zhang if (!mumps->id.ICNTL(18)) { /* A is centralized */ 732a5e57a09SHong Zhang if (!mumps->myid) { 733397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 7342907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 735a5e57a09SHong Zhang mumps->id.a = (mumps_complex*)mumps->val; 7362907cef9SHong Zhang #else 737a5e57a09SHong Zhang mumps->id.a = (mumps_double_complex*)mumps->val; 7382907cef9SHong Zhang #endif 739397b6df1SKris Buschelman #else 740a5e57a09SHong Zhang mumps->id.a = mumps->val; 741397b6df1SKris Buschelman #endif 742397b6df1SKris Buschelman } 743397b6df1SKris Buschelman } else { 744397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 7452907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 746a5e57a09SHong Zhang mumps->id.a_loc = (mumps_complex*)mumps->val; 7472907cef9SHong Zhang #else 748a5e57a09SHong Zhang mumps->id.a_loc = (mumps_double_complex*)mumps->val; 7492907cef9SHong Zhang #endif 750397b6df1SKris Buschelman #else 751a5e57a09SHong Zhang mumps->id.a_loc = mumps->val; 752397b6df1SKris Buschelman #endif 753397b6df1SKris Buschelman } 754a5e57a09SHong Zhang PetscMUMPS_c(&mumps->id); 755a5e57a09SHong Zhang if (mumps->id.INFOG(1) < 0) { 756151787a6SHong Zhang if (mumps->id.INFO(1) == -13) { 757151787a6SHong Zhang if (mumps->id.INFO(2) < 0) { 758151787a6SHong Zhang SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: Cannot allocate required memory %d megabytes\n",-mumps->id.INFO(2)); 759151787a6SHong Zhang } else { 760151787a6SHong Zhang SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: Cannot allocate required memory %d bytes\n",mumps->id.INFO(2)); 761151787a6SHong Zhang } 762151787a6SHong Zhang } else SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: INFO(1)=%d, INFO(2)=%d\n",mumps->id.INFO(1),mumps->id.INFO(2)); 763397b6df1SKris Buschelman } 764a5e57a09SHong Zhang if (!mumps->myid && mumps->id.ICNTL(16) > 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB," mumps->id.ICNTL(16):=%d\n",mumps->id.INFOG(16)); 765397b6df1SKris Buschelman 766dcd589f8SShri Abhyankar (F)->assembled = PETSC_TRUE; 767a5e57a09SHong Zhang mumps->matstruc = SAME_NONZERO_PATTERN; 768a5e57a09SHong Zhang mumps->CleanUpMUMPS = PETSC_TRUE; 76967877ebaSShri Abhyankar 770a5e57a09SHong Zhang if (mumps->size > 1) { 77167877ebaSShri Abhyankar PetscInt lsol_loc; 77267877ebaSShri Abhyankar PetscScalar *sol_loc; 7732205254eSKarl Rupp 774c2093ab7SHong Zhang ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr); 775c2093ab7SHong Zhang if (isMPIAIJ) F_diag = ((Mat_MPIAIJ*)(F)->data)->A; 776c2093ab7SHong Zhang else F_diag = ((Mat_MPISBAIJ*)(F)->data)->A; 777c2093ab7SHong Zhang F_diag->assembled = PETSC_TRUE; 778c2093ab7SHong Zhang 779c2093ab7SHong Zhang /* distributed solution; Create x_seq=sol_loc for repeated use */ 780c2093ab7SHong Zhang if (mumps->x_seq) { 781c2093ab7SHong Zhang ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr); 782c2093ab7SHong Zhang ierr = PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);CHKERRQ(ierr); 783c2093ab7SHong Zhang ierr = VecDestroy(&mumps->x_seq);CHKERRQ(ierr); 784c2093ab7SHong Zhang } 785a5e57a09SHong Zhang lsol_loc = mumps->id.INFO(23); /* length of sol_loc */ 786dcca6d9dSJed Brown ierr = PetscMalloc2(lsol_loc,&sol_loc,lsol_loc,&mumps->id.isol_loc);CHKERRQ(ierr); 787a5e57a09SHong Zhang mumps->id.lsol_loc = lsol_loc; 78867877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 7892907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 790a5e57a09SHong Zhang mumps->id.sol_loc = (mumps_complex*)sol_loc; 7912907cef9SHong Zhang #else 792a5e57a09SHong Zhang mumps->id.sol_loc = (mumps_double_complex*)sol_loc; 7932907cef9SHong Zhang #endif 79467877ebaSShri Abhyankar #else 795a5e57a09SHong Zhang mumps->id.sol_loc = sol_loc; 79667877ebaSShri Abhyankar #endif 797a5e57a09SHong Zhang ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,lsol_loc,sol_loc,&mumps->x_seq);CHKERRQ(ierr); 79867877ebaSShri Abhyankar } 799397b6df1SKris Buschelman PetscFunctionReturn(0); 800397b6df1SKris Buschelman } 801397b6df1SKris Buschelman 8029a2535b5SHong Zhang /* Sets MUMPS options from the options database */ 803dcd589f8SShri Abhyankar #undef __FUNCT__ 8049a2535b5SHong Zhang #define __FUNCT__ "PetscSetMUMPSFromOptions" 8059a2535b5SHong Zhang PetscErrorCode PetscSetMUMPSFromOptions(Mat F, Mat A) 806dcd589f8SShri Abhyankar { 8079a2535b5SHong Zhang Mat_MUMPS *mumps = (Mat_MUMPS*)F->spptr; 808dcd589f8SShri Abhyankar PetscErrorCode ierr; 809dcd589f8SShri Abhyankar PetscInt icntl; 810ace3abfcSBarry Smith PetscBool flg; 811dcd589f8SShri Abhyankar 812dcd589f8SShri Abhyankar PetscFunctionBegin; 813ce94432eSBarry Smith ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MUMPS Options","Mat");CHKERRQ(ierr); 8149a2535b5SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_1","ICNTL(1): output stream for error messages","None",mumps->id.ICNTL(1),&icntl,&flg);CHKERRQ(ierr); 8159a2535b5SHong Zhang if (flg) mumps->id.ICNTL(1) = icntl; 8169a2535b5SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_2","ICNTL(2): output stream for diagnostic printing, statistics, and warning","None",mumps->id.ICNTL(2),&icntl,&flg);CHKERRQ(ierr); 8179a2535b5SHong Zhang if (flg) mumps->id.ICNTL(2) = icntl; 8189a2535b5SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_3","ICNTL(3): output stream for global information, collected on the host","None",mumps->id.ICNTL(3),&icntl,&flg);CHKERRQ(ierr); 8199a2535b5SHong Zhang if (flg) mumps->id.ICNTL(3) = icntl; 820dcd589f8SShri Abhyankar 8219a2535b5SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",mumps->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr); 8229a2535b5SHong Zhang if (flg) mumps->id.ICNTL(4) = icntl; 8239a2535b5SHong Zhang if (mumps->id.ICNTL(4) || PetscLogPrintInfo) mumps->id.ICNTL(3) = 6; /* resume MUMPS default id.ICNTL(3) = 6 */ 8249a2535b5SHong Zhang 825d341cd04SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_6","ICNTL(6): permutes to a zero-free diagonal and/or scale the matrix (0 to 7)","None",mumps->id.ICNTL(6),&icntl,&flg);CHKERRQ(ierr); 8269a2535b5SHong Zhang if (flg) mumps->id.ICNTL(6) = icntl; 8279a2535b5SHong Zhang 828d341cd04SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): computes a symmetric permutation in sequential analysis (0 to 7). 3=Scotch, 4=PORD, 5=Metis","None",mumps->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr); 829dcd589f8SShri Abhyankar if (flg) { 8302205254eSKarl Rupp if (icntl== 1 && mumps->size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"pivot order be set by the user in PERM_IN -- not supported by the PETSc/MUMPS interface\n"); 8312205254eSKarl Rupp else mumps->id.ICNTL(7) = icntl; 832dcd589f8SShri Abhyankar } 833e0b74bf9SHong Zhang 8340298fd71SBarry Smith ierr = PetscOptionsInt("-mat_mumps_icntl_8","ICNTL(8): scaling strategy (-2 to 8 or 77)","None",mumps->id.ICNTL(8),&mumps->id.ICNTL(8),NULL);CHKERRQ(ierr); 835d341cd04SHong Zhang /* ierr = PetscOptionsInt("-mat_mumps_icntl_9","ICNTL(9): computes the solution using A or A^T","None",mumps->id.ICNTL(9),&mumps->id.ICNTL(9),NULL);CHKERRQ(ierr); handled by MatSolveTranspose_MUMPS() */ 8360298fd71SBarry Smith ierr = PetscOptionsInt("-mat_mumps_icntl_10","ICNTL(10): max num of refinements","None",mumps->id.ICNTL(10),&mumps->id.ICNTL(10),NULL);CHKERRQ(ierr); 837d341cd04SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_11","ICNTL(11): statistics related to an error analysis (via -ksp_view)","None",mumps->id.ICNTL(11),&mumps->id.ICNTL(11),NULL);CHKERRQ(ierr); 838d341cd04SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): an ordering strategy for symmetric matrices (0 to 3)","None",mumps->id.ICNTL(12),&mumps->id.ICNTL(12),NULL);CHKERRQ(ierr); 839d341cd04SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): parallelism of the root node (enable ScaLAPACK) and its splitting","None",mumps->id.ICNTL(13),&mumps->id.ICNTL(13),NULL);CHKERRQ(ierr); 840d341cd04SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_14","ICNTL(14): percentage increase in the estimated working space","None",mumps->id.ICNTL(14),&mumps->id.ICNTL(14),NULL);CHKERRQ(ierr); 841d341cd04SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_19","ICNTL(19): computes the Schur complement","None",mumps->id.ICNTL(19),&mumps->id.ICNTL(19),NULL);CHKERRQ(ierr); 842*4e34a73bSHong Zhang /* ierr = PetscOptionsInt("-mat_mumps_icntl_20","ICNTL(20): the format (dense or sparse) of the right-hand sides","None",mumps->id.ICNTL(20),&mumps->id.ICNTL(20),NULL);CHKERRQ(ierr); -- sparse rhs is not supported in PETSc API */ 843d341cd04SHong Zhang /* ierr = PetscOptionsInt("-mat_mumps_icntl_21","ICNTL(21): the distribution (centralized or distributed) of the solution vectors","None",mumps->id.ICNTL(21),&mumps->id.ICNTL(21),NULL);CHKERRQ(ierr); we only use distributed solution vector */ 8449a2535b5SHong Zhang 845d341cd04SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_22","ICNTL(22): in-core/out-of-core factorization and solve (0 or 1)","None",mumps->id.ICNTL(22),&mumps->id.ICNTL(22),NULL);CHKERRQ(ierr); 8460298fd71SBarry Smith ierr = PetscOptionsInt("-mat_mumps_icntl_23","ICNTL(23): max size of the working memory (MB) that can allocate per processor","None",mumps->id.ICNTL(23),&mumps->id.ICNTL(23),NULL);CHKERRQ(ierr); 8470298fd71SBarry Smith ierr = PetscOptionsInt("-mat_mumps_icntl_24","ICNTL(24): detection of null pivot rows (0 or 1)","None",mumps->id.ICNTL(24),&mumps->id.ICNTL(24),NULL);CHKERRQ(ierr); 8489a2535b5SHong Zhang if (mumps->id.ICNTL(24)) { 8499a2535b5SHong Zhang mumps->id.ICNTL(13) = 1; /* turn-off ScaLAPACK to help with the correct detection of null pivots */ 850d7ebd59bSHong Zhang } 851d7ebd59bSHong Zhang 852d341cd04SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_25","ICNTL(25): compute a solution of a deficient matrix and a null space basis","None",mumps->id.ICNTL(25),&mumps->id.ICNTL(25),NULL);CHKERRQ(ierr); 853d341cd04SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_26","ICNTL(26): drives the solution phase if a Schur complement matrix","None",mumps->id.ICNTL(26),&mumps->id.ICNTL(26),NULL);CHKERRQ(ierr); 854*4e34a73bSHong Zhang /* ierr = PetscOptionsInt("-mat_mumps_icntl_27","ICNTL(27): the blocking size for multiple right-hand sides","None",mumps->id.ICNTL(27),&mumps->id.ICNTL(27),NULL);CHKERRQ(ierr); */ 8550298fd71SBarry Smith ierr = PetscOptionsInt("-mat_mumps_icntl_28","ICNTL(28): use 1 for sequential analysis and ictnl(7) ordering, or 2 for parallel analysis and ictnl(29) ordering","None",mumps->id.ICNTL(28),&mumps->id.ICNTL(28),NULL);CHKERRQ(ierr); 856d341cd04SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_29","ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis","None",mumps->id.ICNTL(29),&mumps->id.ICNTL(29),NULL);CHKERRQ(ierr); 8570298fd71SBarry Smith ierr = PetscOptionsInt("-mat_mumps_icntl_30","ICNTL(30): compute user-specified set of entries in inv(A)","None",mumps->id.ICNTL(30),&mumps->id.ICNTL(30),NULL);CHKERRQ(ierr); 858d341cd04SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_31","ICNTL(31): indicates which factors may be discarded during factorization","None",mumps->id.ICNTL(31),&mumps->id.ICNTL(31),NULL);CHKERRQ(ierr); 859*4e34a73bSHong Zhang /* ierr = PetscOptionsInt("-mat_mumps_icntl_32","ICNTL(32): performs the forward elemination of the right-hand sides during factorization","None",mumps->id.ICNTL(32),&mumps->id.ICNTL(32),NULL);CHKERRQ(ierr); -- not supported by PETSc API */ 8600298fd71SBarry Smith ierr = PetscOptionsInt("-mat_mumps_icntl_33","ICNTL(33): compute determinant","None",mumps->id.ICNTL(33),&mumps->id.ICNTL(33),NULL);CHKERRQ(ierr); 861dcd589f8SShri Abhyankar 8620298fd71SBarry Smith ierr = PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",mumps->id.CNTL(1),&mumps->id.CNTL(1),NULL);CHKERRQ(ierr); 8630298fd71SBarry Smith ierr = PetscOptionsReal("-mat_mumps_cntl_2","CNTL(2): stopping criterion of refinement","None",mumps->id.CNTL(2),&mumps->id.CNTL(2),NULL);CHKERRQ(ierr); 8640298fd71SBarry Smith ierr = PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",mumps->id.CNTL(3),&mumps->id.CNTL(3),NULL);CHKERRQ(ierr); 8650298fd71SBarry Smith ierr = PetscOptionsReal("-mat_mumps_cntl_4","CNTL(4): value for static pivoting","None",mumps->id.CNTL(4),&mumps->id.CNTL(4),NULL);CHKERRQ(ierr); 8660298fd71SBarry Smith ierr = PetscOptionsReal("-mat_mumps_cntl_5","CNTL(5): fixation for null pivots","None",mumps->id.CNTL(5),&mumps->id.CNTL(5),NULL);CHKERRQ(ierr); 867e5bb22a1SHong Zhang 8680298fd71SBarry Smith ierr = PetscOptionsString("-mat_mumps_ooc_tmpdir", "out of core directory", "None", mumps->id.ooc_tmpdir, mumps->id.ooc_tmpdir, 256, NULL); 869dcd589f8SShri Abhyankar PetscOptionsEnd(); 870dcd589f8SShri Abhyankar PetscFunctionReturn(0); 871dcd589f8SShri Abhyankar } 872dcd589f8SShri Abhyankar 873dcd589f8SShri Abhyankar #undef __FUNCT__ 874dcd589f8SShri Abhyankar #define __FUNCT__ "PetscInitializeMUMPS" 875f697e70eSHong Zhang PetscErrorCode PetscInitializeMUMPS(Mat A,Mat_MUMPS *mumps) 876dcd589f8SShri Abhyankar { 877dcd589f8SShri Abhyankar PetscErrorCode ierr; 878dcd589f8SShri Abhyankar 879dcd589f8SShri Abhyankar PetscFunctionBegin; 880ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)A), &mumps->myid); 881ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&mumps->size);CHKERRQ(ierr); 882ce94432eSBarry Smith ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)A),&(mumps->comm_mumps));CHKERRQ(ierr); 8832205254eSKarl Rupp 884f697e70eSHong Zhang mumps->id.comm_fortran = MPI_Comm_c2f(mumps->comm_mumps); 885f697e70eSHong Zhang 886f697e70eSHong Zhang mumps->id.job = JOB_INIT; 887f697e70eSHong Zhang mumps->id.par = 1; /* host participates factorizaton and solve */ 888f697e70eSHong Zhang mumps->id.sym = mumps->sym; 8892907cef9SHong Zhang PetscMUMPS_c(&mumps->id); 890f697e70eSHong Zhang 891f697e70eSHong Zhang mumps->CleanUpMUMPS = PETSC_FALSE; 8920298fd71SBarry Smith mumps->scat_rhs = NULL; 8930298fd71SBarry Smith mumps->scat_sol = NULL; 8949a2535b5SHong Zhang 89570544d5fSHong Zhang /* set PETSc-MUMPS default options - override MUMPS default */ 8969a2535b5SHong Zhang mumps->id.ICNTL(3) = 0; 8979a2535b5SHong Zhang mumps->id.ICNTL(4) = 0; 8989a2535b5SHong Zhang if (mumps->size == 1) { 8999a2535b5SHong Zhang mumps->id.ICNTL(18) = 0; /* centralized assembled matrix input */ 9009a2535b5SHong Zhang } else { 9019a2535b5SHong Zhang mumps->id.ICNTL(18) = 3; /* distributed assembled matrix input */ 902*4e34a73bSHong Zhang mumps->id.ICNTL(20) = 0; /* rhs is in dense format */ 90370544d5fSHong Zhang mumps->id.ICNTL(21) = 1; /* distributed solution */ 9049a2535b5SHong Zhang } 905dcd589f8SShri Abhyankar PetscFunctionReturn(0); 906dcd589f8SShri Abhyankar } 907dcd589f8SShri Abhyankar 908a5e57a09SHong Zhang /* Note Petsc r(=c) permutation is used when mumps->id.ICNTL(7)==1 with centralized assembled matrix input; otherwise r and c are ignored */ 909397b6df1SKris Buschelman #undef __FUNCT__ 910f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS" 9110481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 912b24902e0SBarry Smith { 913a5e57a09SHong Zhang Mat_MUMPS *mumps = (Mat_MUMPS*)F->spptr; 914dcd589f8SShri Abhyankar PetscErrorCode ierr; 91567877ebaSShri Abhyankar Vec b; 91667877ebaSShri Abhyankar IS is_iden; 91767877ebaSShri Abhyankar const PetscInt M = A->rmap->N; 918397b6df1SKris Buschelman 919397b6df1SKris Buschelman PetscFunctionBegin; 920a5e57a09SHong Zhang mumps->matstruc = DIFFERENT_NONZERO_PATTERN; 921dcd589f8SShri Abhyankar 9229a2535b5SHong Zhang /* Set MUMPS options from the options database */ 9239a2535b5SHong Zhang ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr); 924dcd589f8SShri Abhyankar 925a5e57a09SHong Zhang ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr); 926dcd589f8SShri Abhyankar 92767877ebaSShri Abhyankar /* analysis phase */ 92867877ebaSShri Abhyankar /*----------------*/ 929a5e57a09SHong Zhang mumps->id.job = JOB_FACTSYMBOLIC; 930a5e57a09SHong Zhang mumps->id.n = M; 931a5e57a09SHong Zhang switch (mumps->id.ICNTL(18)) { 93267877ebaSShri Abhyankar case 0: /* centralized assembled matrix input */ 933a5e57a09SHong Zhang if (!mumps->myid) { 934a5e57a09SHong Zhang mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn; 935a5e57a09SHong Zhang if (mumps->id.ICNTL(6)>1) { 93667877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 9372907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 938a5e57a09SHong Zhang mumps->id.a = (mumps_complex*)mumps->val; 9392907cef9SHong Zhang #else 940a5e57a09SHong Zhang mumps->id.a = (mumps_double_complex*)mumps->val; 9412907cef9SHong Zhang #endif 94267877ebaSShri Abhyankar #else 943a5e57a09SHong Zhang mumps->id.a = mumps->val; 94467877ebaSShri Abhyankar #endif 94567877ebaSShri Abhyankar } 946a5e57a09SHong Zhang if (mumps->id.ICNTL(7) == 1) { /* use user-provide matrix ordering - assuming r = c ordering */ 9475248a706SHong Zhang /* 9485248a706SHong Zhang PetscBool flag; 9495248a706SHong Zhang ierr = ISEqual(r,c,&flag);CHKERRQ(ierr); 9505248a706SHong Zhang if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"row_perm != col_perm"); 9515248a706SHong Zhang ierr = ISView(r,PETSC_VIEWER_STDOUT_SELF); 9525248a706SHong Zhang */ 953a5e57a09SHong Zhang if (!mumps->myid) { 954e0b74bf9SHong Zhang const PetscInt *idx; 955e0b74bf9SHong Zhang PetscInt i,*perm_in; 9562205254eSKarl Rupp 957785e854fSJed Brown ierr = PetscMalloc1(M,&perm_in);CHKERRQ(ierr); 958e0b74bf9SHong Zhang ierr = ISGetIndices(r,&idx);CHKERRQ(ierr); 9592205254eSKarl Rupp 960a5e57a09SHong Zhang mumps->id.perm_in = perm_in; 961e0b74bf9SHong Zhang for (i=0; i<M; i++) perm_in[i] = idx[i]+1; /* perm_in[]: start from 1, not 0! */ 962e0b74bf9SHong Zhang ierr = ISRestoreIndices(r,&idx);CHKERRQ(ierr); 963e0b74bf9SHong Zhang } 964e0b74bf9SHong Zhang } 96567877ebaSShri Abhyankar } 96667877ebaSShri Abhyankar break; 96767877ebaSShri Abhyankar case 3: /* distributed assembled matrix input (size>1) */ 968a5e57a09SHong Zhang mumps->id.nz_loc = mumps->nz; 969a5e57a09SHong Zhang mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn; 970a5e57a09SHong Zhang if (mumps->id.ICNTL(6)>1) { 97167877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 9722907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 973a5e57a09SHong Zhang mumps->id.a_loc = (mumps_complex*)mumps->val; 9742907cef9SHong Zhang #else 975a5e57a09SHong Zhang mumps->id.a_loc = (mumps_double_complex*)mumps->val; 9762907cef9SHong Zhang #endif 97767877ebaSShri Abhyankar #else 978a5e57a09SHong Zhang mumps->id.a_loc = mumps->val; 97967877ebaSShri Abhyankar #endif 98067877ebaSShri Abhyankar } 98167877ebaSShri Abhyankar /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */ 982a5e57a09SHong Zhang if (!mumps->myid) { 983a5e57a09SHong Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr); 98467877ebaSShri Abhyankar ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr); 98567877ebaSShri Abhyankar } else { 986a5e57a09SHong Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr); 98767877ebaSShri Abhyankar ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr); 98867877ebaSShri Abhyankar } 9892a7a6963SBarry Smith ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr); 990a5e57a09SHong Zhang ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr); 9916bf464f9SBarry Smith ierr = ISDestroy(&is_iden);CHKERRQ(ierr); 9926bf464f9SBarry Smith ierr = VecDestroy(&b);CHKERRQ(ierr); 99367877ebaSShri Abhyankar break; 99467877ebaSShri Abhyankar } 995a5e57a09SHong Zhang PetscMUMPS_c(&mumps->id); 996a5e57a09SHong Zhang if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",mumps->id.INFOG(1)); 99767877ebaSShri Abhyankar 998719d5645SBarry Smith F->ops->lufactornumeric = MatFactorNumeric_MUMPS; 999dcd589f8SShri Abhyankar F->ops->solve = MatSolve_MUMPS; 100051d5961aSHong Zhang F->ops->solvetranspose = MatSolveTranspose_MUMPS; 1001*4e34a73bSHong Zhang F->ops->matsolve = MatMatSolve_MUMPS; 1002b24902e0SBarry Smith PetscFunctionReturn(0); 1003b24902e0SBarry Smith } 1004b24902e0SBarry Smith 1005450b117fSShri Abhyankar /* Note the Petsc r and c permutations are ignored */ 1006450b117fSShri Abhyankar #undef __FUNCT__ 1007450b117fSShri Abhyankar #define __FUNCT__ "MatLUFactorSymbolic_BAIJMUMPS" 1008450b117fSShri Abhyankar PetscErrorCode MatLUFactorSymbolic_BAIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 1009450b117fSShri Abhyankar { 1010a5e57a09SHong Zhang Mat_MUMPS *mumps = (Mat_MUMPS*)F->spptr; 1011dcd589f8SShri Abhyankar PetscErrorCode ierr; 101267877ebaSShri Abhyankar Vec b; 101367877ebaSShri Abhyankar IS is_iden; 101467877ebaSShri Abhyankar const PetscInt M = A->rmap->N; 1015450b117fSShri Abhyankar 1016450b117fSShri Abhyankar PetscFunctionBegin; 1017a5e57a09SHong Zhang mumps->matstruc = DIFFERENT_NONZERO_PATTERN; 1018dcd589f8SShri Abhyankar 10199a2535b5SHong Zhang /* Set MUMPS options from the options database */ 10209a2535b5SHong Zhang ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr); 1021dcd589f8SShri Abhyankar 1022a5e57a09SHong Zhang ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr); 102367877ebaSShri Abhyankar 102467877ebaSShri Abhyankar /* analysis phase */ 102567877ebaSShri Abhyankar /*----------------*/ 1026a5e57a09SHong Zhang mumps->id.job = JOB_FACTSYMBOLIC; 1027a5e57a09SHong Zhang mumps->id.n = M; 1028a5e57a09SHong Zhang switch (mumps->id.ICNTL(18)) { 102967877ebaSShri Abhyankar case 0: /* centralized assembled matrix input */ 1030a5e57a09SHong Zhang if (!mumps->myid) { 1031a5e57a09SHong Zhang mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn; 1032a5e57a09SHong Zhang if (mumps->id.ICNTL(6)>1) { 103367877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 10342907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 1035a5e57a09SHong Zhang mumps->id.a = (mumps_complex*)mumps->val; 10362907cef9SHong Zhang #else 1037a5e57a09SHong Zhang mumps->id.a = (mumps_double_complex*)mumps->val; 10382907cef9SHong Zhang #endif 103967877ebaSShri Abhyankar #else 1040a5e57a09SHong Zhang mumps->id.a = mumps->val; 104167877ebaSShri Abhyankar #endif 104267877ebaSShri Abhyankar } 104367877ebaSShri Abhyankar } 104467877ebaSShri Abhyankar break; 104567877ebaSShri Abhyankar case 3: /* distributed assembled matrix input (size>1) */ 1046a5e57a09SHong Zhang mumps->id.nz_loc = mumps->nz; 1047a5e57a09SHong Zhang mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn; 1048a5e57a09SHong Zhang if (mumps->id.ICNTL(6)>1) { 104967877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 10502907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 1051a5e57a09SHong Zhang mumps->id.a_loc = (mumps_complex*)mumps->val; 10522907cef9SHong Zhang #else 1053a5e57a09SHong Zhang mumps->id.a_loc = (mumps_double_complex*)mumps->val; 10542907cef9SHong Zhang #endif 105567877ebaSShri Abhyankar #else 1056a5e57a09SHong Zhang mumps->id.a_loc = mumps->val; 105767877ebaSShri Abhyankar #endif 105867877ebaSShri Abhyankar } 105967877ebaSShri Abhyankar /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */ 1060a5e57a09SHong Zhang if (!mumps->myid) { 1061a5e57a09SHong Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr); 106267877ebaSShri Abhyankar ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr); 106367877ebaSShri Abhyankar } else { 1064a5e57a09SHong Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr); 106567877ebaSShri Abhyankar ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr); 106667877ebaSShri Abhyankar } 10672a7a6963SBarry Smith ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr); 1068a5e57a09SHong Zhang ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr); 10696bf464f9SBarry Smith ierr = ISDestroy(&is_iden);CHKERRQ(ierr); 10706bf464f9SBarry Smith ierr = VecDestroy(&b);CHKERRQ(ierr); 107167877ebaSShri Abhyankar break; 107267877ebaSShri Abhyankar } 1073a5e57a09SHong Zhang PetscMUMPS_c(&mumps->id); 1074a5e57a09SHong Zhang if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",mumps->id.INFOG(1)); 107567877ebaSShri Abhyankar 1076450b117fSShri Abhyankar F->ops->lufactornumeric = MatFactorNumeric_MUMPS; 1077dcd589f8SShri Abhyankar F->ops->solve = MatSolve_MUMPS; 107851d5961aSHong Zhang F->ops->solvetranspose = MatSolveTranspose_MUMPS; 1079450b117fSShri Abhyankar PetscFunctionReturn(0); 1080450b117fSShri Abhyankar } 1081b24902e0SBarry Smith 1082141f4205SHong Zhang /* Note the Petsc r permutation and factor info are ignored */ 1083397b6df1SKris Buschelman #undef __FUNCT__ 108467877ebaSShri Abhyankar #define __FUNCT__ "MatCholeskyFactorSymbolic_MUMPS" 108567877ebaSShri Abhyankar PetscErrorCode MatCholeskyFactorSymbolic_MUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info) 1086b24902e0SBarry Smith { 1087a5e57a09SHong Zhang Mat_MUMPS *mumps = (Mat_MUMPS*)F->spptr; 1088dcd589f8SShri Abhyankar PetscErrorCode ierr; 108967877ebaSShri Abhyankar Vec b; 109067877ebaSShri Abhyankar IS is_iden; 109167877ebaSShri Abhyankar const PetscInt M = A->rmap->N; 1092397b6df1SKris Buschelman 1093397b6df1SKris Buschelman PetscFunctionBegin; 1094a5e57a09SHong Zhang mumps->matstruc = DIFFERENT_NONZERO_PATTERN; 1095dcd589f8SShri Abhyankar 10969a2535b5SHong Zhang /* Set MUMPS options from the options database */ 10979a2535b5SHong Zhang ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr); 1098dcd589f8SShri Abhyankar 1099a5e57a09SHong Zhang ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr); 1100dcd589f8SShri Abhyankar 110167877ebaSShri Abhyankar /* analysis phase */ 110267877ebaSShri Abhyankar /*----------------*/ 1103a5e57a09SHong Zhang mumps->id.job = JOB_FACTSYMBOLIC; 1104a5e57a09SHong Zhang mumps->id.n = M; 1105a5e57a09SHong Zhang switch (mumps->id.ICNTL(18)) { 110667877ebaSShri Abhyankar case 0: /* centralized assembled matrix input */ 1107a5e57a09SHong Zhang if (!mumps->myid) { 1108a5e57a09SHong Zhang mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn; 1109a5e57a09SHong Zhang if (mumps->id.ICNTL(6)>1) { 111067877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 11112907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 1112a5e57a09SHong Zhang mumps->id.a = (mumps_complex*)mumps->val; 11132907cef9SHong Zhang #else 1114a5e57a09SHong Zhang mumps->id.a = (mumps_double_complex*)mumps->val; 11152907cef9SHong Zhang #endif 111667877ebaSShri Abhyankar #else 1117a5e57a09SHong Zhang mumps->id.a = mumps->val; 111867877ebaSShri Abhyankar #endif 111967877ebaSShri Abhyankar } 112067877ebaSShri Abhyankar } 112167877ebaSShri Abhyankar break; 112267877ebaSShri Abhyankar case 3: /* distributed assembled matrix input (size>1) */ 1123a5e57a09SHong Zhang mumps->id.nz_loc = mumps->nz; 1124a5e57a09SHong Zhang mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn; 1125a5e57a09SHong Zhang if (mumps->id.ICNTL(6)>1) { 112667877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 11272907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 1128a5e57a09SHong Zhang mumps->id.a_loc = (mumps_complex*)mumps->val; 11292907cef9SHong Zhang #else 1130a5e57a09SHong Zhang mumps->id.a_loc = (mumps_double_complex*)mumps->val; 11312907cef9SHong Zhang #endif 113267877ebaSShri Abhyankar #else 1133a5e57a09SHong Zhang mumps->id.a_loc = mumps->val; 113467877ebaSShri Abhyankar #endif 113567877ebaSShri Abhyankar } 113667877ebaSShri Abhyankar /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */ 1137a5e57a09SHong Zhang if (!mumps->myid) { 1138a5e57a09SHong Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr); 113967877ebaSShri Abhyankar ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr); 114067877ebaSShri Abhyankar } else { 1141a5e57a09SHong Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr); 114267877ebaSShri Abhyankar ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr); 114367877ebaSShri Abhyankar } 11442a7a6963SBarry Smith ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr); 1145a5e57a09SHong Zhang ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr); 11466bf464f9SBarry Smith ierr = ISDestroy(&is_iden);CHKERRQ(ierr); 11476bf464f9SBarry Smith ierr = VecDestroy(&b);CHKERRQ(ierr); 114867877ebaSShri Abhyankar break; 114967877ebaSShri Abhyankar } 1150a5e57a09SHong Zhang PetscMUMPS_c(&mumps->id); 1151a5e57a09SHong Zhang if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",mumps->id.INFOG(1)); 115267877ebaSShri Abhyankar 11532792810eSHong Zhang F->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS; 1154dcd589f8SShri Abhyankar F->ops->solve = MatSolve_MUMPS; 115551d5961aSHong Zhang F->ops->solvetranspose = MatSolve_MUMPS; 1156*4e34a73bSHong Zhang F->ops->matsolve = MatMatSolve_MUMPS; 1157*4e34a73bSHong Zhang #if defined(PETSC_USE_COMPLEX) 11580298fd71SBarry Smith F->ops->getinertia = NULL; 1159*4e34a73bSHong Zhang #else 1160*4e34a73bSHong Zhang F->ops->getinertia = MatGetInertia_SBAIJMUMPS; 1161db4efbfdSBarry Smith #endif 1162b24902e0SBarry Smith PetscFunctionReturn(0); 1163b24902e0SBarry Smith } 1164b24902e0SBarry Smith 1165*4e34a73bSHong Zhang //update!!! 1166397b6df1SKris Buschelman #undef __FUNCT__ 116764e6c443SBarry Smith #define __FUNCT__ "MatView_MUMPS" 116864e6c443SBarry Smith PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer) 116974ed9c26SBarry Smith { 1170f6c57405SHong Zhang PetscErrorCode ierr; 117164e6c443SBarry Smith PetscBool iascii; 117264e6c443SBarry Smith PetscViewerFormat format; 1173a5e57a09SHong Zhang Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 1174f6c57405SHong Zhang 1175f6c57405SHong Zhang PetscFunctionBegin; 117664e6c443SBarry Smith /* check if matrix is mumps type */ 117764e6c443SBarry Smith if (A->ops->solve != MatSolve_MUMPS) PetscFunctionReturn(0); 117864e6c443SBarry Smith 1179251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 118064e6c443SBarry Smith if (iascii) { 118164e6c443SBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 118264e6c443SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO) { 118364e6c443SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr); 1184a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," SYM (matrix type): %d \n",mumps->id.sym);CHKERRQ(ierr); 1185a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," PAR (host participation): %d \n",mumps->id.par);CHKERRQ(ierr); 1186a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(1) (output for error): %d \n",mumps->id.ICNTL(1));CHKERRQ(ierr); 1187a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(2) (output of diagnostic msg): %d \n",mumps->id.ICNTL(2));CHKERRQ(ierr); 1188a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(3) (output for global info): %d \n",mumps->id.ICNTL(3));CHKERRQ(ierr); 1189a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(4) (level of printing): %d \n",mumps->id.ICNTL(4));CHKERRQ(ierr); 1190a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(5) (input mat struct): %d \n",mumps->id.ICNTL(5));CHKERRQ(ierr); 1191a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(6) (matrix prescaling): %d \n",mumps->id.ICNTL(6));CHKERRQ(ierr); 1192a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(7) (sequentia matrix ordering):%d \n",mumps->id.ICNTL(7));CHKERRQ(ierr); 1193a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(8) (scalling strategy): %d \n",mumps->id.ICNTL(8));CHKERRQ(ierr); 1194a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(10) (max num of refinements): %d \n",mumps->id.ICNTL(10));CHKERRQ(ierr); 1195a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(11) (error analysis): %d \n",mumps->id.ICNTL(11));CHKERRQ(ierr); 1196a5e57a09SHong Zhang if (mumps->id.ICNTL(11)>0) { 1197a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(4) (inf norm of input mat): %g\n",mumps->id.RINFOG(4));CHKERRQ(ierr); 1198a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(5) (inf norm of solution): %g\n",mumps->id.RINFOG(5));CHKERRQ(ierr); 1199a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(6) (inf norm of residual): %g\n",mumps->id.RINFOG(6));CHKERRQ(ierr); 1200a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",mumps->id.RINFOG(7),mumps->id.RINFOG(8));CHKERRQ(ierr); 1201a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(9) (error estimate): %g \n",mumps->id.RINFOG(9));CHKERRQ(ierr); 1202a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",mumps->id.RINFOG(10),mumps->id.RINFOG(11));CHKERRQ(ierr); 1203f6c57405SHong Zhang } 1204a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(12) (efficiency control): %d \n",mumps->id.ICNTL(12));CHKERRQ(ierr); 1205a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(13) (efficiency control): %d \n",mumps->id.ICNTL(13));CHKERRQ(ierr); 1206a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(14) (percentage of estimated workspace increase): %d \n",mumps->id.ICNTL(14));CHKERRQ(ierr); 1207f6c57405SHong Zhang /* ICNTL(15-17) not used */ 1208a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(18) (input mat struct): %d \n",mumps->id.ICNTL(18));CHKERRQ(ierr); 1209a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(19) (Shur complement info): %d \n",mumps->id.ICNTL(19));CHKERRQ(ierr); 1210a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(20) (rhs sparse pattern): %d \n",mumps->id.ICNTL(20));CHKERRQ(ierr); 1211a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(21) (somumpstion struct): %d \n",mumps->id.ICNTL(21));CHKERRQ(ierr); 1212a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(22) (in-core/out-of-core facility): %d \n",mumps->id.ICNTL(22));CHKERRQ(ierr); 1213a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(23) (max size of memory can be allocated locally):%d \n",mumps->id.ICNTL(23));CHKERRQ(ierr); 1214c0165424SHong Zhang 1215a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(24) (detection of null pivot rows): %d \n",mumps->id.ICNTL(24));CHKERRQ(ierr); 1216a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(25) (computation of a null space basis): %d \n",mumps->id.ICNTL(25));CHKERRQ(ierr); 1217a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(26) (Schur options for rhs or solution): %d \n",mumps->id.ICNTL(26));CHKERRQ(ierr); 1218a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(27) (experimental parameter): %d \n",mumps->id.ICNTL(27));CHKERRQ(ierr); 1219a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(28) (use parallel or sequential ordering): %d \n",mumps->id.ICNTL(28));CHKERRQ(ierr); 1220a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(29) (parallel ordering): %d \n",mumps->id.ICNTL(29));CHKERRQ(ierr); 122142179a6aSHong Zhang 1222a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(30) (user-specified set of entries in inv(A)): %d \n",mumps->id.ICNTL(30));CHKERRQ(ierr); 1223a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(31) (factors is discarded in the solve phase): %d \n",mumps->id.ICNTL(31));CHKERRQ(ierr); 1224a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(33) (compute determinant): %d \n",mumps->id.ICNTL(33));CHKERRQ(ierr); 1225f6c57405SHong Zhang 1226a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(1) (relative pivoting threshold): %g \n",mumps->id.CNTL(1));CHKERRQ(ierr); 1227a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(2) (stopping criterion of refinement): %g \n",mumps->id.CNTL(2));CHKERRQ(ierr); 1228a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(3) (absomumpste pivoting threshold): %g \n",mumps->id.CNTL(3));CHKERRQ(ierr); 1229a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(4) (vamumpse of static pivoting): %g \n",mumps->id.CNTL(4));CHKERRQ(ierr); 1230a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(5) (fixation for null pivots): %g \n",mumps->id.CNTL(5));CHKERRQ(ierr); 1231f6c57405SHong Zhang 1232f6c57405SHong Zhang /* infomation local to each processor */ 123334ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer, " RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr); 12347b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 1235a5e57a09SHong Zhang ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %g \n",mumps->myid,mumps->id.RINFO(1));CHKERRQ(ierr); 123634ed7027SBarry Smith ierr = PetscViewerFlush(viewer); 123734ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer, " RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr); 1238a5e57a09SHong Zhang ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %g \n",mumps->myid,mumps->id.RINFO(2));CHKERRQ(ierr); 123934ed7027SBarry Smith ierr = PetscViewerFlush(viewer); 124034ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer, " RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr); 1241a5e57a09SHong Zhang ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %g \n",mumps->myid,mumps->id.RINFO(3));CHKERRQ(ierr); 124234ed7027SBarry Smith ierr = PetscViewerFlush(viewer); 1243f6c57405SHong Zhang 124434ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer, " INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization): \n");CHKERRQ(ierr); 1245a5e57a09SHong Zhang ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %d \n",mumps->myid,mumps->id.INFO(15));CHKERRQ(ierr); 124634ed7027SBarry Smith ierr = PetscViewerFlush(viewer); 1247f6c57405SHong Zhang 124834ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer, " INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization): \n");CHKERRQ(ierr); 1249a5e57a09SHong Zhang ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %d \n",mumps->myid,mumps->id.INFO(16));CHKERRQ(ierr); 125034ed7027SBarry Smith ierr = PetscViewerFlush(viewer); 1251f6c57405SHong Zhang 125234ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer, " INFO(23) (num of pivots eliminated on this processor after factorization): \n");CHKERRQ(ierr); 1253a5e57a09SHong Zhang ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %d \n",mumps->myid,mumps->id.INFO(23));CHKERRQ(ierr); 125434ed7027SBarry Smith ierr = PetscViewerFlush(viewer); 12557b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr); 1256f6c57405SHong Zhang 1257a5e57a09SHong Zhang if (!mumps->myid) { /* information from the host */ 1258a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",mumps->id.RINFOG(1));CHKERRQ(ierr); 1259a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",mumps->id.RINFOG(2));CHKERRQ(ierr); 1260a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",mumps->id.RINFOG(3));CHKERRQ(ierr); 1261a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," (RINFOG(12) RINFOG(13))*2^INFOG(34) (determinant): (%g,%g)*(2^%d)\n",mumps->id.RINFOG(12),mumps->id.RINFOG(13),mumps->id.INFOG(34));CHKERRQ(ierr); 1262f6c57405SHong Zhang 1263a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(3));CHKERRQ(ierr); 1264a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(4));CHKERRQ(ierr); 1265a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(5) (estimated maximum front size in the complete tree): %d \n",mumps->id.INFOG(5));CHKERRQ(ierr); 1266a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(6) (number of nodes in the complete tree): %d \n",mumps->id.INFOG(6));CHKERRQ(ierr); 1267a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(7) (ordering option effectively use after analysis): %d \n",mumps->id.INFOG(7));CHKERRQ(ierr); 1268a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",mumps->id.INFOG(8));CHKERRQ(ierr); 1269a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(9) (total real/complex workspace to store the matrix factors after factorization): %d \n",mumps->id.INFOG(9));CHKERRQ(ierr); 1270a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(10) (total integer space store the matrix factors after factorization): %d \n",mumps->id.INFOG(10));CHKERRQ(ierr); 1271a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(11) (order of largest frontal matrix after factorization): %d \n",mumps->id.INFOG(11));CHKERRQ(ierr); 1272a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(12) (number of off-diagonal pivots): %d \n",mumps->id.INFOG(12));CHKERRQ(ierr); 1273a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(13) (number of delayed pivots after factorization): %d \n",mumps->id.INFOG(13));CHKERRQ(ierr); 1274a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(14) (number of memory compress after factorization): %d \n",mumps->id.INFOG(14));CHKERRQ(ierr); 1275a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(15) (number of steps of iterative refinement after solution): %d \n",mumps->id.INFOG(15));CHKERRQ(ierr); 1276a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(16) (estimated size (in MB) of all MUMPS internal data for factorization after analysis: value on the most memory consuming processor): %d \n",mumps->id.INFOG(16));CHKERRQ(ierr); 1277a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(17) (estimated size of all MUMPS internal data for factorization after analysis: sum over all processors): %d \n",mumps->id.INFOG(17));CHKERRQ(ierr); 1278a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(18) (size of all MUMPS internal data allocated during factorization: value on the most memory consuming processor): %d \n",mumps->id.INFOG(18));CHKERRQ(ierr); 1279a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(19) (size of all MUMPS internal data allocated during factorization: sum over all processors): %d \n",mumps->id.INFOG(19));CHKERRQ(ierr); 1280a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(20) (estimated number of entries in the factors): %d \n",mumps->id.INFOG(20));CHKERRQ(ierr); 1281a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(21) (size in MB of memory effectively used during factorization - value on the most memory consuming processor): %d \n",mumps->id.INFOG(21));CHKERRQ(ierr); 1282a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(22) (size in MB of memory effectively used during factorization - sum over all processors): %d \n",mumps->id.INFOG(22));CHKERRQ(ierr); 1283a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",mumps->id.INFOG(23));CHKERRQ(ierr); 1284a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",mumps->id.INFOG(24));CHKERRQ(ierr); 1285a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",mumps->id.INFOG(25));CHKERRQ(ierr); 128640d435e3SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(28) (after factorization: number of null pivots encountered): %d\n",mumps->id.INFOG(28));CHKERRQ(ierr); 128740d435e3SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(29) (after factorization: effective number of entries in the factors (sum over all processors)): %d\n",mumps->id.INFOG(29));CHKERRQ(ierr); 128840d435e3SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(30, 31) (after solution: size in Mbytes of memory used during solution phase): %d, %d\n",mumps->id.INFOG(30),mumps->id.INFOG(31));CHKERRQ(ierr); 128940d435e3SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(32) (after analysis: type of analysis done): %d\n",mumps->id.INFOG(32));CHKERRQ(ierr); 129040d435e3SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(33) (value used for ICNTL(8)): %d\n",mumps->id.INFOG(33));CHKERRQ(ierr); 129140d435e3SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(34) (exponent of the determinant if determinant is requested): %d\n",mumps->id.INFOG(34));CHKERRQ(ierr); 1292f6c57405SHong Zhang } 1293f6c57405SHong Zhang } 1294cb828f0fSHong Zhang } 1295f6c57405SHong Zhang PetscFunctionReturn(0); 1296f6c57405SHong Zhang } 1297f6c57405SHong Zhang 129835bd34faSBarry Smith #undef __FUNCT__ 129935bd34faSBarry Smith #define __FUNCT__ "MatGetInfo_MUMPS" 130035bd34faSBarry Smith PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info) 130135bd34faSBarry Smith { 1302cb828f0fSHong Zhang Mat_MUMPS *mumps =(Mat_MUMPS*)A->spptr; 130335bd34faSBarry Smith 130435bd34faSBarry Smith PetscFunctionBegin; 130535bd34faSBarry Smith info->block_size = 1.0; 1306cb828f0fSHong Zhang info->nz_allocated = mumps->id.INFOG(20); 1307cb828f0fSHong Zhang info->nz_used = mumps->id.INFOG(20); 130835bd34faSBarry Smith info->nz_unneeded = 0.0; 130935bd34faSBarry Smith info->assemblies = 0.0; 131035bd34faSBarry Smith info->mallocs = 0.0; 131135bd34faSBarry Smith info->memory = 0.0; 131235bd34faSBarry Smith info->fill_ratio_given = 0; 131335bd34faSBarry Smith info->fill_ratio_needed = 0; 131435bd34faSBarry Smith info->factor_mallocs = 0; 131535bd34faSBarry Smith PetscFunctionReturn(0); 131635bd34faSBarry Smith } 131735bd34faSBarry Smith 13185ccb76cbSHong Zhang /* -------------------------------------------------------------------------------------------*/ 13195ccb76cbSHong Zhang #undef __FUNCT__ 13205ccb76cbSHong Zhang #define __FUNCT__ "MatMumpsSetIcntl_MUMPS" 13215ccb76cbSHong Zhang PetscErrorCode MatMumpsSetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt ival) 13225ccb76cbSHong Zhang { 1323a5e57a09SHong Zhang Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 13245ccb76cbSHong Zhang 13255ccb76cbSHong Zhang PetscFunctionBegin; 1326a5e57a09SHong Zhang mumps->id.ICNTL(icntl) = ival; 13275ccb76cbSHong Zhang PetscFunctionReturn(0); 13285ccb76cbSHong Zhang } 13295ccb76cbSHong Zhang 13305ccb76cbSHong Zhang #undef __FUNCT__ 1331bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetIcntl_MUMPS" 1332bc6112feSHong Zhang PetscErrorCode MatMumpsGetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt *ival) 1333bc6112feSHong Zhang { 1334bc6112feSHong Zhang Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 1335bc6112feSHong Zhang 1336bc6112feSHong Zhang PetscFunctionBegin; 1337bc6112feSHong Zhang *ival = mumps->id.ICNTL(icntl); 1338bc6112feSHong Zhang PetscFunctionReturn(0); 1339bc6112feSHong Zhang } 1340bc6112feSHong Zhang 1341bc6112feSHong Zhang #undef __FUNCT__ 13425ccb76cbSHong Zhang #define __FUNCT__ "MatMumpsSetIcntl" 13435ccb76cbSHong Zhang /*@ 13445ccb76cbSHong Zhang MatMumpsSetIcntl - Set MUMPS parameter ICNTL() 13455ccb76cbSHong Zhang 13465ccb76cbSHong Zhang Logically Collective on Mat 13475ccb76cbSHong Zhang 13485ccb76cbSHong Zhang Input Parameters: 13495ccb76cbSHong Zhang + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 13505ccb76cbSHong Zhang . icntl - index of MUMPS parameter array ICNTL() 13515ccb76cbSHong Zhang - ival - value of MUMPS ICNTL(icntl) 13525ccb76cbSHong Zhang 13535ccb76cbSHong Zhang Options Database: 13545ccb76cbSHong Zhang . -mat_mumps_icntl_<icntl> <ival> 13555ccb76cbSHong Zhang 13565ccb76cbSHong Zhang Level: beginner 13575ccb76cbSHong Zhang 13585ccb76cbSHong Zhang References: MUMPS Users' Guide 13595ccb76cbSHong Zhang 13605ccb76cbSHong Zhang .seealso: MatGetFactor() 13615ccb76cbSHong Zhang @*/ 13625ccb76cbSHong Zhang PetscErrorCode MatMumpsSetIcntl(Mat F,PetscInt icntl,PetscInt ival) 13635ccb76cbSHong Zhang { 13645ccb76cbSHong Zhang PetscErrorCode ierr; 13655ccb76cbSHong Zhang 13665ccb76cbSHong Zhang PetscFunctionBegin; 13675ccb76cbSHong Zhang PetscValidLogicalCollectiveInt(F,icntl,2); 13685ccb76cbSHong Zhang PetscValidLogicalCollectiveInt(F,ival,3); 13695ccb76cbSHong Zhang ierr = PetscTryMethod(F,"MatMumpsSetIcntl_C",(Mat,PetscInt,PetscInt),(F,icntl,ival));CHKERRQ(ierr); 13705ccb76cbSHong Zhang PetscFunctionReturn(0); 13715ccb76cbSHong Zhang } 13725ccb76cbSHong Zhang 1373bc6112feSHong Zhang #undef __FUNCT__ 1374bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetIcntl" 1375a21f80fcSHong Zhang /*@ 1376a21f80fcSHong Zhang MatMumpsGetIcntl - Get MUMPS parameter ICNTL() 1377a21f80fcSHong Zhang 1378a21f80fcSHong Zhang Logically Collective on Mat 1379a21f80fcSHong Zhang 1380a21f80fcSHong Zhang Input Parameters: 1381a21f80fcSHong Zhang + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 1382a21f80fcSHong Zhang - icntl - index of MUMPS parameter array ICNTL() 1383a21f80fcSHong Zhang 1384a21f80fcSHong Zhang Output Parameter: 1385a21f80fcSHong Zhang . ival - value of MUMPS ICNTL(icntl) 1386a21f80fcSHong Zhang 1387a21f80fcSHong Zhang Level: beginner 1388a21f80fcSHong Zhang 1389a21f80fcSHong Zhang References: MUMPS Users' Guide 1390a21f80fcSHong Zhang 1391a21f80fcSHong Zhang .seealso: MatGetFactor() 1392a21f80fcSHong Zhang @*/ 1393bc6112feSHong Zhang PetscErrorCode MatMumpsGetIcntl(Mat F,PetscInt icntl,PetscInt *ival) 1394bc6112feSHong Zhang { 1395bc6112feSHong Zhang PetscErrorCode ierr; 1396bc6112feSHong Zhang 1397bc6112feSHong Zhang PetscFunctionBegin; 1398bc6112feSHong Zhang PetscValidLogicalCollectiveInt(F,icntl,2); 1399bc6112feSHong Zhang PetscValidIntPointer(ival,3); 1400bc6112feSHong Zhang ierr = PetscTryMethod(F,"MatMumpsGetIcntl_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr); 1401bc6112feSHong Zhang PetscFunctionReturn(0); 1402bc6112feSHong Zhang } 1403bc6112feSHong Zhang 14048928b65cSHong Zhang /* -------------------------------------------------------------------------------------------*/ 14058928b65cSHong Zhang #undef __FUNCT__ 14068928b65cSHong Zhang #define __FUNCT__ "MatMumpsSetCntl_MUMPS" 14078928b65cSHong Zhang PetscErrorCode MatMumpsSetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal val) 14088928b65cSHong Zhang { 14098928b65cSHong Zhang Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 14108928b65cSHong Zhang 14118928b65cSHong Zhang PetscFunctionBegin; 14128928b65cSHong Zhang mumps->id.CNTL(icntl) = val; 14138928b65cSHong Zhang PetscFunctionReturn(0); 14148928b65cSHong Zhang } 14158928b65cSHong Zhang 14168928b65cSHong Zhang #undef __FUNCT__ 1417bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetCntl_MUMPS" 1418bc6112feSHong Zhang PetscErrorCode MatMumpsGetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal *val) 1419bc6112feSHong Zhang { 1420bc6112feSHong Zhang Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 1421bc6112feSHong Zhang 1422bc6112feSHong Zhang PetscFunctionBegin; 1423bc6112feSHong Zhang *val = mumps->id.CNTL(icntl); 1424bc6112feSHong Zhang PetscFunctionReturn(0); 1425bc6112feSHong Zhang } 1426bc6112feSHong Zhang 1427bc6112feSHong Zhang #undef __FUNCT__ 14288928b65cSHong Zhang #define __FUNCT__ "MatMumpsSetCntl" 14298928b65cSHong Zhang /*@ 14308928b65cSHong Zhang MatMumpsSetCntl - Set MUMPS parameter CNTL() 14318928b65cSHong Zhang 14328928b65cSHong Zhang Logically Collective on Mat 14338928b65cSHong Zhang 14348928b65cSHong Zhang Input Parameters: 14358928b65cSHong Zhang + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 14368928b65cSHong Zhang . icntl - index of MUMPS parameter array CNTL() 14378928b65cSHong Zhang - val - value of MUMPS CNTL(icntl) 14388928b65cSHong Zhang 14398928b65cSHong Zhang Options Database: 14408928b65cSHong Zhang . -mat_mumps_cntl_<icntl> <val> 14418928b65cSHong Zhang 14428928b65cSHong Zhang Level: beginner 14438928b65cSHong Zhang 14448928b65cSHong Zhang References: MUMPS Users' Guide 14458928b65cSHong Zhang 14468928b65cSHong Zhang .seealso: MatGetFactor() 14478928b65cSHong Zhang @*/ 14488928b65cSHong Zhang PetscErrorCode MatMumpsSetCntl(Mat F,PetscInt icntl,PetscReal val) 14498928b65cSHong Zhang { 14508928b65cSHong Zhang PetscErrorCode ierr; 14518928b65cSHong Zhang 14528928b65cSHong Zhang PetscFunctionBegin; 14538928b65cSHong Zhang PetscValidLogicalCollectiveInt(F,icntl,2); 1454bc6112feSHong Zhang PetscValidLogicalCollectiveReal(F,val,3); 14558928b65cSHong Zhang ierr = PetscTryMethod(F,"MatMumpsSetCntl_C",(Mat,PetscInt,PetscReal),(F,icntl,val));CHKERRQ(ierr); 14568928b65cSHong Zhang PetscFunctionReturn(0); 14578928b65cSHong Zhang } 14588928b65cSHong Zhang 1459bc6112feSHong Zhang #undef __FUNCT__ 1460bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetCntl" 1461a21f80fcSHong Zhang /*@ 1462a21f80fcSHong Zhang MatMumpsGetCntl - Get MUMPS parameter CNTL() 1463a21f80fcSHong Zhang 1464a21f80fcSHong Zhang Logically Collective on Mat 1465a21f80fcSHong Zhang 1466a21f80fcSHong Zhang Input Parameters: 1467a21f80fcSHong Zhang + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 1468a21f80fcSHong Zhang - icntl - index of MUMPS parameter array CNTL() 1469a21f80fcSHong Zhang 1470a21f80fcSHong Zhang Output Parameter: 1471a21f80fcSHong Zhang . val - value of MUMPS CNTL(icntl) 1472a21f80fcSHong Zhang 1473a21f80fcSHong Zhang Level: beginner 1474a21f80fcSHong Zhang 1475a21f80fcSHong Zhang References: MUMPS Users' Guide 1476a21f80fcSHong Zhang 1477a21f80fcSHong Zhang .seealso: MatGetFactor() 1478a21f80fcSHong Zhang @*/ 1479bc6112feSHong Zhang PetscErrorCode MatMumpsGetCntl(Mat F,PetscInt icntl,PetscReal *val) 1480bc6112feSHong Zhang { 1481bc6112feSHong Zhang PetscErrorCode ierr; 1482bc6112feSHong Zhang 1483bc6112feSHong Zhang PetscFunctionBegin; 1484bc6112feSHong Zhang PetscValidLogicalCollectiveInt(F,icntl,2); 1485bc6112feSHong Zhang PetscValidRealPointer(val,3); 1486bc6112feSHong Zhang ierr = PetscTryMethod(F,"MatMumpsGetCntl_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr); 1487bc6112feSHong Zhang PetscFunctionReturn(0); 1488bc6112feSHong Zhang } 1489bc6112feSHong Zhang 1490bc6112feSHong Zhang #undef __FUNCT__ 1491ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfo_MUMPS" 1492ca810319SHong Zhang PetscErrorCode MatMumpsGetInfo_MUMPS(Mat F,PetscInt icntl,PetscInt *info) 1493bc6112feSHong Zhang { 1494bc6112feSHong Zhang Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 1495bc6112feSHong Zhang 1496bc6112feSHong Zhang PetscFunctionBegin; 1497bc6112feSHong Zhang *info = mumps->id.INFO(icntl); 1498bc6112feSHong Zhang PetscFunctionReturn(0); 1499bc6112feSHong Zhang } 1500bc6112feSHong Zhang 1501bc6112feSHong Zhang #undef __FUNCT__ 1502ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfog_MUMPS" 1503ca810319SHong Zhang PetscErrorCode MatMumpsGetInfog_MUMPS(Mat F,PetscInt icntl,PetscInt *infog) 1504bc6112feSHong Zhang { 1505bc6112feSHong Zhang Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 1506bc6112feSHong Zhang 1507bc6112feSHong Zhang PetscFunctionBegin; 1508bc6112feSHong Zhang *infog = mumps->id.INFOG(icntl); 1509bc6112feSHong Zhang PetscFunctionReturn(0); 1510bc6112feSHong Zhang } 1511bc6112feSHong Zhang 1512bc6112feSHong Zhang #undef __FUNCT__ 1513ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfo_MUMPS" 1514ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfo_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfo) 1515bc6112feSHong Zhang { 1516bc6112feSHong Zhang Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 1517bc6112feSHong Zhang 1518bc6112feSHong Zhang PetscFunctionBegin; 1519bc6112feSHong Zhang *rinfo = mumps->id.RINFO(icntl); 1520bc6112feSHong Zhang PetscFunctionReturn(0); 1521bc6112feSHong Zhang } 1522bc6112feSHong Zhang 1523bc6112feSHong Zhang #undef __FUNCT__ 1524ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfog_MUMPS" 1525ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfog_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfog) 1526bc6112feSHong Zhang { 1527bc6112feSHong Zhang Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 1528bc6112feSHong Zhang 1529bc6112feSHong Zhang PetscFunctionBegin; 1530bc6112feSHong Zhang *rinfog = mumps->id.RINFOG(icntl); 1531bc6112feSHong Zhang PetscFunctionReturn(0); 1532bc6112feSHong Zhang } 1533bc6112feSHong Zhang 1534bc6112feSHong Zhang #undef __FUNCT__ 1535ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfo" 1536a21f80fcSHong Zhang /*@ 1537a21f80fcSHong Zhang MatMumpsGetInfo - Get MUMPS parameter INFO() 1538a21f80fcSHong Zhang 1539a21f80fcSHong Zhang Logically Collective on Mat 1540a21f80fcSHong Zhang 1541a21f80fcSHong Zhang Input Parameters: 1542a21f80fcSHong Zhang + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 1543a21f80fcSHong Zhang - icntl - index of MUMPS parameter array INFO() 1544a21f80fcSHong Zhang 1545a21f80fcSHong Zhang Output Parameter: 1546a21f80fcSHong Zhang . ival - value of MUMPS INFO(icntl) 1547a21f80fcSHong Zhang 1548a21f80fcSHong Zhang Level: beginner 1549a21f80fcSHong Zhang 1550a21f80fcSHong Zhang References: MUMPS Users' Guide 1551a21f80fcSHong Zhang 1552a21f80fcSHong Zhang .seealso: MatGetFactor() 1553a21f80fcSHong Zhang @*/ 1554ca810319SHong Zhang PetscErrorCode MatMumpsGetInfo(Mat F,PetscInt icntl,PetscInt *ival) 1555bc6112feSHong Zhang { 1556bc6112feSHong Zhang PetscErrorCode ierr; 1557bc6112feSHong Zhang 1558bc6112feSHong Zhang PetscFunctionBegin; 1559ca810319SHong Zhang PetscValidIntPointer(ival,3); 1560ca810319SHong Zhang ierr = PetscTryMethod(F,"MatMumpsGetInfo_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr); 1561bc6112feSHong Zhang PetscFunctionReturn(0); 1562bc6112feSHong Zhang } 1563bc6112feSHong Zhang 1564bc6112feSHong Zhang #undef __FUNCT__ 1565ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfog" 1566a21f80fcSHong Zhang /*@ 1567a21f80fcSHong Zhang MatMumpsGetInfog - Get MUMPS parameter INFOG() 1568a21f80fcSHong Zhang 1569a21f80fcSHong Zhang Logically Collective on Mat 1570a21f80fcSHong Zhang 1571a21f80fcSHong Zhang Input Parameters: 1572a21f80fcSHong Zhang + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 1573a21f80fcSHong Zhang - icntl - index of MUMPS parameter array INFOG() 1574a21f80fcSHong Zhang 1575a21f80fcSHong Zhang Output Parameter: 1576a21f80fcSHong Zhang . ival - value of MUMPS INFOG(icntl) 1577a21f80fcSHong Zhang 1578a21f80fcSHong Zhang Level: beginner 1579a21f80fcSHong Zhang 1580a21f80fcSHong Zhang References: MUMPS Users' Guide 1581a21f80fcSHong Zhang 1582a21f80fcSHong Zhang .seealso: MatGetFactor() 1583a21f80fcSHong Zhang @*/ 1584ca810319SHong Zhang PetscErrorCode MatMumpsGetInfog(Mat F,PetscInt icntl,PetscInt *ival) 1585bc6112feSHong Zhang { 1586bc6112feSHong Zhang PetscErrorCode ierr; 1587bc6112feSHong Zhang 1588bc6112feSHong Zhang PetscFunctionBegin; 1589ca810319SHong Zhang PetscValidIntPointer(ival,3); 1590ca810319SHong Zhang ierr = PetscTryMethod(F,"MatMumpsGetInfog_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr); 1591bc6112feSHong Zhang PetscFunctionReturn(0); 1592bc6112feSHong Zhang } 1593bc6112feSHong Zhang 1594bc6112feSHong Zhang #undef __FUNCT__ 1595ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfo" 1596a21f80fcSHong Zhang /*@ 1597a21f80fcSHong Zhang MatMumpsGetRinfo - Get MUMPS parameter RINFO() 1598a21f80fcSHong Zhang 1599a21f80fcSHong Zhang Logically Collective on Mat 1600a21f80fcSHong Zhang 1601a21f80fcSHong Zhang Input Parameters: 1602a21f80fcSHong Zhang + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 1603a21f80fcSHong Zhang - icntl - index of MUMPS parameter array RINFO() 1604a21f80fcSHong Zhang 1605a21f80fcSHong Zhang Output Parameter: 1606a21f80fcSHong Zhang . val - value of MUMPS RINFO(icntl) 1607a21f80fcSHong Zhang 1608a21f80fcSHong Zhang Level: beginner 1609a21f80fcSHong Zhang 1610a21f80fcSHong Zhang References: MUMPS Users' Guide 1611a21f80fcSHong Zhang 1612a21f80fcSHong Zhang .seealso: MatGetFactor() 1613a21f80fcSHong Zhang @*/ 1614ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfo(Mat F,PetscInt icntl,PetscReal *val) 1615bc6112feSHong Zhang { 1616bc6112feSHong Zhang PetscErrorCode ierr; 1617bc6112feSHong Zhang 1618bc6112feSHong Zhang PetscFunctionBegin; 1619bc6112feSHong Zhang PetscValidRealPointer(val,3); 1620ca810319SHong Zhang ierr = PetscTryMethod(F,"MatMumpsGetRinfo_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr); 1621bc6112feSHong Zhang PetscFunctionReturn(0); 1622bc6112feSHong Zhang } 1623bc6112feSHong Zhang 1624bc6112feSHong Zhang #undef __FUNCT__ 1625ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfog" 1626a21f80fcSHong Zhang /*@ 1627a21f80fcSHong Zhang MatMumpsGetRinfog - Get MUMPS parameter RINFOG() 1628a21f80fcSHong Zhang 1629a21f80fcSHong Zhang Logically Collective on Mat 1630a21f80fcSHong Zhang 1631a21f80fcSHong Zhang Input Parameters: 1632a21f80fcSHong Zhang + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 1633a21f80fcSHong Zhang - icntl - index of MUMPS parameter array RINFOG() 1634a21f80fcSHong Zhang 1635a21f80fcSHong Zhang Output Parameter: 1636a21f80fcSHong Zhang . val - value of MUMPS RINFOG(icntl) 1637a21f80fcSHong Zhang 1638a21f80fcSHong Zhang Level: beginner 1639a21f80fcSHong Zhang 1640a21f80fcSHong Zhang References: MUMPS Users' Guide 1641a21f80fcSHong Zhang 1642a21f80fcSHong Zhang .seealso: MatGetFactor() 1643a21f80fcSHong Zhang @*/ 1644ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfog(Mat F,PetscInt icntl,PetscReal *val) 1645bc6112feSHong Zhang { 1646bc6112feSHong Zhang PetscErrorCode ierr; 1647bc6112feSHong Zhang 1648bc6112feSHong Zhang PetscFunctionBegin; 1649bc6112feSHong Zhang PetscValidRealPointer(val,3); 1650ca810319SHong Zhang ierr = PetscTryMethod(F,"MatMumpsGetRinfog_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr); 1651bc6112feSHong Zhang PetscFunctionReturn(0); 1652bc6112feSHong Zhang } 1653bc6112feSHong Zhang 165424b6179bSKris Buschelman /*MC 16552692d6eeSBarry Smith MATSOLVERMUMPS - A matrix type providing direct solvers (LU and Cholesky) for 165624b6179bSKris Buschelman distributed and sequential matrices via the external package MUMPS. 165724b6179bSKris Buschelman 165841c8de11SBarry Smith Works with MATAIJ and MATSBAIJ matrices 165924b6179bSKris Buschelman 166024b6179bSKris Buschelman Options Database Keys: 1661*4e34a73bSHong Zhang + -mat_mumps_icntl_1 <6>: ICNTL(1): output stream for error messages (None) 1662*4e34a73bSHong Zhang . -mat_mumps_icntl_2 <0>: ICNTL(2): output stream for diagnostic printing, statistics, and warning (None) 1663*4e34a73bSHong Zhang . -mat_mumps_icntl_3 <0>: ICNTL(3): output stream for global information, collected on the host (None) 1664*4e34a73bSHong Zhang . -mat_mumps_icntl_4 <0>: ICNTL(4): level of printing (0 to 4) (None) 1665*4e34a73bSHong Zhang . -mat_mumps_icntl_6 <7>: ICNTL(6): permutes to a zero-free diagonal and/or scale the matrix (0 to 7) (None) 1666*4e34a73bSHong Zhang . -mat_mumps_icntl_7 <7>: ICNTL(7): computes a symmetric permutation in sequential analysis (0 to 7). 3=Scotch, 4=PORD, 5=Metis (None) 1667*4e34a73bSHong Zhang . -mat_mumps_icntl_8 <77>: ICNTL(8): scaling strategy (-2 to 8 or 77) (None) 1668*4e34a73bSHong Zhang . -mat_mumps_icntl_10 <0>: ICNTL(10): max num of refinements (None) 1669*4e34a73bSHong Zhang . -mat_mumps_icntl_11 <0>: ICNTL(11): statistics related to an error analysis (via -ksp_view) (None) 1670*4e34a73bSHong Zhang . -mat_mumps_icntl_12 <1>: ICNTL(12): an ordering strategy for symmetric matrices (0 to 3) (None) 1671*4e34a73bSHong Zhang . -mat_mumps_icntl_13 <0>: ICNTL(13): parallelism of the root node (enable ScaLAPACK) and its splitting (None) 1672*4e34a73bSHong Zhang . -mat_mumps_icntl_14 <20>: ICNTL(14): percentage increase in the estimated working space (None) 1673*4e34a73bSHong Zhang . -mat_mumps_icntl_19 <0>: ICNTL(19): computes the Schur complement (None) 1674*4e34a73bSHong Zhang . -mat_mumps_icntl_22 <0>: ICNTL(22): in-core/out-of-core factorization and solve (0 or 1) (None) 1675*4e34a73bSHong Zhang . -mat_mumps_icntl_23 <0>: ICNTL(23): max size of the working memory (MB) that can allocate per processor (None) 1676*4e34a73bSHong Zhang . -mat_mumps_icntl_24 <0>: ICNTL(24): detection of null pivot rows (0 or 1) (None) 1677*4e34a73bSHong Zhang . -mat_mumps_icntl_25 <0>: ICNTL(25): compute a solution of a deficient matrix and a null space basis (None) 1678*4e34a73bSHong Zhang . -mat_mumps_icntl_26 <0>: ICNTL(26): drives the solution phase if a Schur complement matrix (None) 1679*4e34a73bSHong Zhang . -mat_mumps_icntl_28 <1>: ICNTL(28): use 1 for sequential analysis and ictnl(7) ordering, or 2 for parallel analysis and ictnl(29) ordering (None) 1680*4e34a73bSHong Zhang . -mat_mumps_icntl_29 <0>: ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis (None) 1681*4e34a73bSHong Zhang . -mat_mumps_icntl_30 <0>: ICNTL(30): compute user-specified set of entries in inv(A) (None) 1682*4e34a73bSHong Zhang . -mat_mumps_icntl_31 <0>: ICNTL(31): indicates which factors may be discarded during factorization (None) 1683*4e34a73bSHong Zhang . -mat_mumps_icntl_33 <0>: ICNTL(33): compute determinant (None) 1684*4e34a73bSHong Zhang . -mat_mumps_cntl_1 <0.01>: CNTL(1): relative pivoting threshold (None) 1685*4e34a73bSHong Zhang . -mat_mumps_cntl_2 <1.49012e-08>: CNTL(2): stopping criterion of refinement (None) 1686*4e34a73bSHong Zhang . -mat_mumps_cntl_3 <0>: CNTL(3): absolute pivoting threshold (None) 1687*4e34a73bSHong Zhang . -mat_mumps_cntl_4 <-1>: CNTL(4): value for static pivoting (None) 1688*4e34a73bSHong Zhang - -mat_mumps_cntl_5 <0>: CNTL(5): fixation for null pivots (None) 168924b6179bSKris Buschelman 169024b6179bSKris Buschelman Level: beginner 169124b6179bSKris Buschelman 169241c8de11SBarry Smith .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage 169341c8de11SBarry Smith 169424b6179bSKris Buschelman M*/ 169524b6179bSKris Buschelman 169635bd34faSBarry Smith #undef __FUNCT__ 169735bd34faSBarry Smith #define __FUNCT__ "MatFactorGetSolverPackage_mumps" 1698f7a08781SBarry Smith static PetscErrorCode MatFactorGetSolverPackage_mumps(Mat A,const MatSolverPackage *type) 169935bd34faSBarry Smith { 170035bd34faSBarry Smith PetscFunctionBegin; 17012692d6eeSBarry Smith *type = MATSOLVERMUMPS; 170235bd34faSBarry Smith PetscFunctionReturn(0); 170335bd34faSBarry Smith } 170435bd34faSBarry Smith 1705bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI AIJ matrices */ 17062877fffaSHong Zhang #undef __FUNCT__ 1707bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_aij_mumps" 17088cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_aij_mumps(Mat A,MatFactorType ftype,Mat *F) 17092877fffaSHong Zhang { 17102877fffaSHong Zhang Mat B; 17112877fffaSHong Zhang PetscErrorCode ierr; 17122877fffaSHong Zhang Mat_MUMPS *mumps; 1713ace3abfcSBarry Smith PetscBool isSeqAIJ; 17142877fffaSHong Zhang 17152877fffaSHong Zhang PetscFunctionBegin; 17162877fffaSHong Zhang /* Create the factorization matrix */ 1717251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr); 1718ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 17192877fffaSHong Zhang ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 17202877fffaSHong Zhang ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 1721bccb9932SShri Abhyankar if (isSeqAIJ) { 17220298fd71SBarry Smith ierr = MatSeqAIJSetPreallocation(B,0,NULL);CHKERRQ(ierr); 1723bccb9932SShri Abhyankar } else { 17240298fd71SBarry Smith ierr = MatMPIAIJSetPreallocation(B,0,NULL,0,NULL);CHKERRQ(ierr); 1725bccb9932SShri Abhyankar } 17262877fffaSHong Zhang 1727b00a9115SJed Brown ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr); 17282205254eSKarl Rupp 17292877fffaSHong Zhang B->ops->view = MatView_MUMPS; 173035bd34faSBarry Smith B->ops->getinfo = MatGetInfo_MUMPS; 173120be8e61SHong Zhang B->ops->getdiagonal = MatGetDiagonal_MUMPS; 17322205254eSKarl Rupp 1733bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr); 1734bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr); 1735bc6112feSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr); 1736bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr); 1737bc6112feSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr); 1738bc6112feSHong Zhang 1739ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr); 1740ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr); 1741ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr); 1742ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr); 1743450b117fSShri Abhyankar if (ftype == MAT_FACTOR_LU) { 1744450b117fSShri Abhyankar B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS; 1745d5f3da31SBarry Smith B->factortype = MAT_FACTOR_LU; 1746bccb9932SShri Abhyankar if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqaij; 1747bccb9932SShri Abhyankar else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpiaij; 1748746480a1SHong Zhang mumps->sym = 0; 1749dcd589f8SShri Abhyankar } else { 175067877ebaSShri Abhyankar B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS; 1751450b117fSShri Abhyankar B->factortype = MAT_FACTOR_CHOLESKY; 1752bccb9932SShri Abhyankar if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqsbaij; 1753bccb9932SShri Abhyankar else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpisbaij; 17546fdc2a6dSBarry Smith if (A->spd_set && A->spd) mumps->sym = 1; 17556fdc2a6dSBarry Smith else mumps->sym = 2; 1756450b117fSShri Abhyankar } 17572877fffaSHong Zhang 17582877fffaSHong Zhang mumps->isAIJ = PETSC_TRUE; 1759bf0cc555SLisandro Dalcin mumps->Destroy = B->ops->destroy; 17602877fffaSHong Zhang B->ops->destroy = MatDestroy_MUMPS; 17612877fffaSHong Zhang B->spptr = (void*)mumps; 17622205254eSKarl Rupp 1763f697e70eSHong Zhang ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr); 1764746480a1SHong Zhang 17652877fffaSHong Zhang *F = B; 17662877fffaSHong Zhang PetscFunctionReturn(0); 17672877fffaSHong Zhang } 17682877fffaSHong Zhang 1769bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI SBAIJ matrices */ 17702877fffaSHong Zhang #undef __FUNCT__ 1771bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_sbaij_mumps" 17728cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_sbaij_mumps(Mat A,MatFactorType ftype,Mat *F) 17732877fffaSHong Zhang { 17742877fffaSHong Zhang Mat B; 17752877fffaSHong Zhang PetscErrorCode ierr; 17762877fffaSHong Zhang Mat_MUMPS *mumps; 1777ace3abfcSBarry Smith PetscBool isSeqSBAIJ; 17782877fffaSHong Zhang 17792877fffaSHong Zhang PetscFunctionBegin; 1780ce94432eSBarry Smith if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix"); 1781ce94432eSBarry Smith if (A->rmap->bs > 1) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with block size > 1 with MUMPS Cholesky, use AIJ matrix instead"); 1782251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr); 17832877fffaSHong Zhang /* Create the factorization matrix */ 1784ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 17852877fffaSHong Zhang ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 17862877fffaSHong Zhang ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 1787b00a9115SJed Brown ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr); 1788bccb9932SShri Abhyankar if (isSeqSBAIJ) { 17890298fd71SBarry Smith ierr = MatSeqSBAIJSetPreallocation(B,1,0,NULL);CHKERRQ(ierr); 17902205254eSKarl Rupp 179116ebf90aSShri Abhyankar mumps->ConvertToTriples = MatConvertToTriples_seqsbaij_seqsbaij; 1792dcd589f8SShri Abhyankar } else { 17930298fd71SBarry Smith ierr = MatMPISBAIJSetPreallocation(B,1,0,NULL,0,NULL);CHKERRQ(ierr); 17942205254eSKarl Rupp 1795bccb9932SShri Abhyankar mumps->ConvertToTriples = MatConvertToTriples_mpisbaij_mpisbaij; 1796bccb9932SShri Abhyankar } 1797bccb9932SShri Abhyankar 179867877ebaSShri Abhyankar B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS; 1799bccb9932SShri Abhyankar B->ops->view = MatView_MUMPS; 180020be8e61SHong Zhang B->ops->getdiagonal = MatGetDiagonal_MUMPS; 18012205254eSKarl Rupp 1802bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr); 1803b13644aeSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr); 1804b13644aeSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr); 1805b13644aeSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr); 1806b13644aeSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr); 1807bc6112feSHong Zhang 1808ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr); 1809ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr); 1810ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr); 1811ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr); 18122205254eSKarl Rupp 1813f4762488SHong Zhang B->factortype = MAT_FACTOR_CHOLESKY; 18146fdc2a6dSBarry Smith if (A->spd_set && A->spd) mumps->sym = 1; 18156fdc2a6dSBarry Smith else mumps->sym = 2; 1816a214ac2aSShri Abhyankar 1817bccb9932SShri Abhyankar mumps->isAIJ = PETSC_FALSE; 1818bf0cc555SLisandro Dalcin mumps->Destroy = B->ops->destroy; 1819f3c0ef26SHong Zhang B->ops->destroy = MatDestroy_MUMPS; 18202877fffaSHong Zhang B->spptr = (void*)mumps; 18212205254eSKarl Rupp 1822f697e70eSHong Zhang ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr); 1823746480a1SHong Zhang 18242877fffaSHong Zhang *F = B; 18252877fffaSHong Zhang PetscFunctionReturn(0); 18262877fffaSHong Zhang } 182797969023SHong Zhang 1828450b117fSShri Abhyankar #undef __FUNCT__ 1829bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_baij_mumps" 18308cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_baij_mumps(Mat A,MatFactorType ftype,Mat *F) 183167877ebaSShri Abhyankar { 183267877ebaSShri Abhyankar Mat B; 183367877ebaSShri Abhyankar PetscErrorCode ierr; 183467877ebaSShri Abhyankar Mat_MUMPS *mumps; 1835ace3abfcSBarry Smith PetscBool isSeqBAIJ; 183667877ebaSShri Abhyankar 183767877ebaSShri Abhyankar PetscFunctionBegin; 183867877ebaSShri Abhyankar /* Create the factorization matrix */ 1839251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);CHKERRQ(ierr); 1840ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 184167877ebaSShri Abhyankar ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 184267877ebaSShri Abhyankar ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 1843bccb9932SShri Abhyankar if (isSeqBAIJ) { 18440298fd71SBarry Smith ierr = MatSeqBAIJSetPreallocation(B,A->rmap->bs,0,NULL);CHKERRQ(ierr); 1845bccb9932SShri Abhyankar } else { 18460298fd71SBarry Smith ierr = MatMPIBAIJSetPreallocation(B,A->rmap->bs,0,NULL,0,NULL);CHKERRQ(ierr); 1847bccb9932SShri Abhyankar } 1848450b117fSShri Abhyankar 1849b00a9115SJed Brown ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr); 1850450b117fSShri Abhyankar if (ftype == MAT_FACTOR_LU) { 1851450b117fSShri Abhyankar B->ops->lufactorsymbolic = MatLUFactorSymbolic_BAIJMUMPS; 1852450b117fSShri Abhyankar B->factortype = MAT_FACTOR_LU; 1853bccb9932SShri Abhyankar if (isSeqBAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqbaij_seqaij; 1854bccb9932SShri Abhyankar else mumps->ConvertToTriples = MatConvertToTriples_mpibaij_mpiaij; 1855746480a1SHong Zhang mumps->sym = 0; 1856f23aa3ddSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc BAIJ matrices with MUMPS Cholesky, use SBAIJ or AIJ matrix instead\n"); 1857bccb9932SShri Abhyankar 1858450b117fSShri Abhyankar B->ops->view = MatView_MUMPS; 185920be8e61SHong Zhang B->ops->getdiagonal = MatGetDiagonal_MUMPS; 18602205254eSKarl Rupp 1861bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr); 1862bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr); 1863bc6112feSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr); 1864bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr); 1865bc6112feSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr); 1866bc6112feSHong Zhang 1867ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr); 1868ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr); 1869ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr); 1870ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr); 1871450b117fSShri Abhyankar 1872450b117fSShri Abhyankar mumps->isAIJ = PETSC_TRUE; 1873bf0cc555SLisandro Dalcin mumps->Destroy = B->ops->destroy; 1874450b117fSShri Abhyankar B->ops->destroy = MatDestroy_MUMPS; 1875450b117fSShri Abhyankar B->spptr = (void*)mumps; 18762205254eSKarl Rupp 1877f697e70eSHong Zhang ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr); 1878746480a1SHong Zhang 1879450b117fSShri Abhyankar *F = B; 1880450b117fSShri Abhyankar PetscFunctionReturn(0); 1881450b117fSShri Abhyankar } 188242c9c57cSBarry Smith 188342c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor_aij_mumps(Mat,MatFactorType,Mat*); 188442c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor_baij_mumps(Mat,MatFactorType,Mat*); 188542c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor_sbaij_mumps(Mat,MatFactorType,Mat*); 188642c9c57cSBarry Smith 188742c9c57cSBarry Smith #undef __FUNCT__ 188842c9c57cSBarry Smith #define __FUNCT__ "MatSolverPackageRegister_MUMPS" 188929b38603SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_MUMPS(void) 189042c9c57cSBarry Smith { 189142c9c57cSBarry Smith PetscErrorCode ierr; 189242c9c57cSBarry Smith 189342c9c57cSBarry Smith PetscFunctionBegin; 189442c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIAIJ, MAT_FACTOR_LU,MatGetFactor_aij_mumps);CHKERRQ(ierr); 189542c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIAIJ, MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);CHKERRQ(ierr); 189642c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIBAIJ, MAT_FACTOR_LU,MatGetFactor_baij_mumps);CHKERRQ(ierr); 189742c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIBAIJ, MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);CHKERRQ(ierr); 189842c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPISBAIJ, MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);CHKERRQ(ierr); 189942c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQAIJ, MAT_FACTOR_LU,MatGetFactor_aij_mumps);CHKERRQ(ierr); 190042c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQAIJ, MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);CHKERRQ(ierr); 190142c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQBAIJ, MAT_FACTOR_LU,MatGetFactor_baij_mumps);CHKERRQ(ierr); 190242c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQBAIJ, MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);CHKERRQ(ierr); 190342c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQSBAIJ, MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);CHKERRQ(ierr); 190442c9c57cSBarry Smith PetscFunctionReturn(0); 190542c9c57cSBarry Smith } 190642c9c57cSBarry Smith 1907