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) 11c6db04a5SJed Brown #include <zmumps_c.h> 12397b6df1SKris Buschelman #else 13c6db04a5SJed Brown #include <dmumps_c.h> 14397b6df1SKris Buschelman #endif 15397b6df1SKris Buschelman EXTERN_C_END 16397b6df1SKris Buschelman #define JOB_INIT -1 173d472b54SHong Zhang #define JOB_FACTSYMBOLIC 1 183d472b54SHong Zhang #define JOB_FACTNUMERIC 2 193d472b54SHong Zhang #define JOB_SOLVE 3 20397b6df1SKris Buschelman #define JOB_END -2 213d472b54SHong Zhang 223d472b54SHong Zhang 23397b6df1SKris Buschelman /* macros s.t. indices match MUMPS documentation */ 24397b6df1SKris Buschelman #define ICNTL(I) icntl[(I)-1] 25397b6df1SKris Buschelman #define CNTL(I) cntl[(I)-1] 26397b6df1SKris Buschelman #define INFOG(I) infog[(I)-1] 27a7aca84bSHong Zhang #define INFO(I) info[(I)-1] 28397b6df1SKris Buschelman #define RINFOG(I) rinfog[(I)-1] 29adc1d99fSHong Zhang #define RINFO(I) rinfo[(I)-1] 30397b6df1SKris Buschelman 31397b6df1SKris Buschelman typedef struct { 32397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 33397b6df1SKris Buschelman ZMUMPS_STRUC_C id; 34397b6df1SKris Buschelman #else 35397b6df1SKris Buschelman DMUMPS_STRUC_C id; 36397b6df1SKris Buschelman #endif 37397b6df1SKris Buschelman MatStructure matstruc; 38c1490034SHong Zhang PetscMPIInt myid,size; 3916ebf90aSShri Abhyankar PetscInt *irn,*jcn,nz,sym,nSolve; 40397b6df1SKris Buschelman PetscScalar *val; 41397b6df1SKris Buschelman MPI_Comm comm_mumps; 42329ec9b3SHong Zhang VecScatter scat_rhs, scat_sol; 4364e6c443SBarry Smith PetscBool isAIJ,CleanUpMUMPS; 44329ec9b3SHong Zhang Vec b_seq,x_seq; 45bf0cc555SLisandro Dalcin PetscErrorCode (*Destroy)(Mat); 46bccb9932SShri Abhyankar PetscErrorCode (*ConvertToTriples)(Mat, int, MatReuse, int*, int**, int**, PetscScalar**); 47f0c56d0fSKris Buschelman } Mat_MUMPS; 48f0c56d0fSKris Buschelman 4909573ac7SBarry Smith extern PetscErrorCode MatDuplicate_MUMPS(Mat,MatDuplicateOption,Mat*); 50b24902e0SBarry Smith 5167877ebaSShri Abhyankar 5267877ebaSShri Abhyankar /* MatConvertToTriples_A_B */ 5367877ebaSShri Abhyankar /*convert Petsc matrix to triples: row[nz], col[nz], val[nz] */ 54397b6df1SKris Buschelman /* 55397b6df1SKris Buschelman input: 5667877ebaSShri Abhyankar A - matrix in aij,baij or sbaij (bs=1) format 57397b6df1SKris Buschelman shift - 0: C style output triple; 1: Fortran style output triple. 58bccb9932SShri Abhyankar reuse - MAT_INITIAL_MATRIX: spaces are allocated and values are set for the triple 59bccb9932SShri Abhyankar MAT_REUSE_MATRIX: only the values in v array are updated 60397b6df1SKris Buschelman output: 61397b6df1SKris Buschelman nnz - dim of r, c, and v (number of local nonzero entries of A) 62397b6df1SKris Buschelman r, c, v - row and col index, matrix values (matrix triples) 63397b6df1SKris Buschelman */ 6416ebf90aSShri Abhyankar 6516ebf90aSShri Abhyankar #undef __FUNCT__ 6616ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqaij_seqaij" 67bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 68b24902e0SBarry Smith { 69185f6596SHong Zhang const PetscInt *ai,*aj,*ajj,M=A->rmap->n; 7067877ebaSShri Abhyankar PetscInt nz,rnz,i,j; 71dfbe8321SBarry Smith PetscErrorCode ierr; 72c1490034SHong Zhang PetscInt *row,*col; 7316ebf90aSShri Abhyankar Mat_SeqAIJ *aa=(Mat_SeqAIJ*)A->data; 74397b6df1SKris Buschelman 75397b6df1SKris Buschelman PetscFunctionBegin; 7616ebf90aSShri Abhyankar *v=aa->a; 77bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX){ 7816ebf90aSShri Abhyankar nz = aa->nz; ai = aa->i; aj = aa->j; 7916ebf90aSShri Abhyankar *nnz = nz; 80185f6596SHong Zhang ierr = PetscMalloc(2*nz*sizeof(PetscInt), &row);CHKERRQ(ierr); 81185f6596SHong Zhang col = row + nz; 82185f6596SHong Zhang 8316ebf90aSShri Abhyankar nz = 0; 8416ebf90aSShri Abhyankar for(i=0; i<M; i++) { 8516ebf90aSShri Abhyankar rnz = ai[i+1] - ai[i]; 8667877ebaSShri Abhyankar ajj = aj + ai[i]; 8767877ebaSShri Abhyankar for(j=0; j<rnz; j++) { 8867877ebaSShri Abhyankar row[nz] = i+shift; col[nz++] = ajj[j] + shift; 8916ebf90aSShri Abhyankar } 9016ebf90aSShri Abhyankar } 9116ebf90aSShri Abhyankar *r = row; *c = col; 9216ebf90aSShri Abhyankar } 9316ebf90aSShri Abhyankar PetscFunctionReturn(0); 9416ebf90aSShri Abhyankar } 95397b6df1SKris Buschelman 9616ebf90aSShri Abhyankar #undef __FUNCT__ 9767877ebaSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqbaij_seqaij" 98bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqbaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 9967877ebaSShri Abhyankar { 10067877ebaSShri Abhyankar Mat_SeqBAIJ *aa=(Mat_SeqBAIJ*)A->data; 10167877ebaSShri Abhyankar const PetscInt *ai,*aj,*ajj,bs=A->rmap->bs,bs2=aa->bs2,M=A->rmap->N/bs; 1020ad0caddSJed Brown PetscInt nz,idx=0,rnz,i,j,k,m; 10367877ebaSShri Abhyankar PetscErrorCode ierr; 10467877ebaSShri Abhyankar PetscInt *row,*col; 10567877ebaSShri Abhyankar 10667877ebaSShri Abhyankar PetscFunctionBegin; 107cf3759fdSShri Abhyankar *v = aa->a; 108bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX){ 109cf3759fdSShri Abhyankar ai = aa->i; aj = aa->j; 11067877ebaSShri Abhyankar nz = bs2*aa->nz; 11167877ebaSShri Abhyankar *nnz = nz; 112185f6596SHong Zhang ierr = PetscMalloc(2*nz*sizeof(PetscInt), &row);CHKERRQ(ierr); 113185f6596SHong Zhang col = row + nz; 114185f6596SHong Zhang 11567877ebaSShri Abhyankar for(i=0; i<M; i++) { 11667877ebaSShri Abhyankar ajj = aj + ai[i]; 11767877ebaSShri Abhyankar rnz = ai[i+1] - ai[i]; 11867877ebaSShri Abhyankar for(k=0; k<rnz; k++) { 11967877ebaSShri Abhyankar for(j=0; j<bs; j++) { 12067877ebaSShri Abhyankar for(m=0; m<bs; m++) { 12167877ebaSShri Abhyankar row[idx] = i*bs + m + shift; 122cf3759fdSShri Abhyankar col[idx++] = bs*(ajj[k]) + j + shift; 12367877ebaSShri Abhyankar } 12467877ebaSShri Abhyankar } 12567877ebaSShri Abhyankar } 12667877ebaSShri Abhyankar } 127cf3759fdSShri Abhyankar *r = row; *c = col; 12867877ebaSShri Abhyankar } 12967877ebaSShri Abhyankar PetscFunctionReturn(0); 13067877ebaSShri Abhyankar } 13167877ebaSShri Abhyankar 13267877ebaSShri Abhyankar #undef __FUNCT__ 13316ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqsbaij_seqsbaij" 134bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqsbaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 13516ebf90aSShri Abhyankar { 13667877ebaSShri Abhyankar const PetscInt *ai, *aj,*ajj,M=A->rmap->n; 13767877ebaSShri Abhyankar PetscInt nz,rnz,i,j; 13816ebf90aSShri Abhyankar PetscErrorCode ierr; 13916ebf90aSShri Abhyankar PetscInt *row,*col; 14016ebf90aSShri Abhyankar Mat_SeqSBAIJ *aa=(Mat_SeqSBAIJ*)A->data; 14116ebf90aSShri Abhyankar 14216ebf90aSShri Abhyankar PetscFunctionBegin; 143882afa5aSHong Zhang *v = aa->a; 144bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX){ 14516ebf90aSShri Abhyankar nz = aa->nz;ai=aa->i; aj=aa->j;*v=aa->a; 14616ebf90aSShri Abhyankar *nnz = nz; 147185f6596SHong Zhang ierr = PetscMalloc(2*nz*sizeof(PetscInt), &row);CHKERRQ(ierr); 148185f6596SHong Zhang col = row + nz; 149185f6596SHong Zhang 15016ebf90aSShri Abhyankar nz = 0; 15116ebf90aSShri Abhyankar for(i=0; i<M; i++) { 15216ebf90aSShri Abhyankar rnz = ai[i+1] - ai[i]; 15367877ebaSShri Abhyankar ajj = aj + ai[i]; 15467877ebaSShri Abhyankar for(j=0; j<rnz; j++) { 15567877ebaSShri Abhyankar row[nz] = i+shift; col[nz++] = ajj[j] + shift; 15616ebf90aSShri Abhyankar } 15716ebf90aSShri Abhyankar } 15816ebf90aSShri Abhyankar *r = row; *c = col; 15916ebf90aSShri Abhyankar } 16016ebf90aSShri Abhyankar PetscFunctionReturn(0); 16116ebf90aSShri Abhyankar } 16216ebf90aSShri Abhyankar 16316ebf90aSShri Abhyankar #undef __FUNCT__ 16416ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqaij_seqsbaij" 165bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 16616ebf90aSShri Abhyankar { 16767877ebaSShri Abhyankar const PetscInt *ai,*aj,*ajj,*adiag,M=A->rmap->n; 16867877ebaSShri Abhyankar PetscInt nz,rnz,i,j; 16967877ebaSShri Abhyankar const PetscScalar *av,*v1; 17016ebf90aSShri Abhyankar PetscScalar *val; 17116ebf90aSShri Abhyankar PetscErrorCode ierr; 17216ebf90aSShri Abhyankar PetscInt *row,*col; 17316ebf90aSShri Abhyankar Mat_SeqSBAIJ *aa=(Mat_SeqSBAIJ*)A->data; 17416ebf90aSShri Abhyankar 17516ebf90aSShri Abhyankar PetscFunctionBegin; 17616ebf90aSShri Abhyankar ai=aa->i; aj=aa->j;av=aa->a; 17716ebf90aSShri Abhyankar adiag=aa->diag; 178bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX){ 17916ebf90aSShri Abhyankar nz = M + (aa->nz-M)/2; 18016ebf90aSShri Abhyankar *nnz = nz; 181185f6596SHong Zhang ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr); 182185f6596SHong Zhang col = row + nz; 183185f6596SHong Zhang val = (PetscScalar*)(col + nz); 184185f6596SHong Zhang 18516ebf90aSShri Abhyankar nz = 0; 18616ebf90aSShri Abhyankar for(i=0; i<M; i++) { 18716ebf90aSShri Abhyankar rnz = ai[i+1] - adiag[i]; 18867877ebaSShri Abhyankar ajj = aj + adiag[i]; 189cf3759fdSShri Abhyankar v1 = av + adiag[i]; 19067877ebaSShri Abhyankar for(j=0; j<rnz; j++) { 19167877ebaSShri Abhyankar row[nz] = i+shift; col[nz] = ajj[j] + shift; val[nz++] = v1[j]; 19216ebf90aSShri Abhyankar } 19316ebf90aSShri Abhyankar } 19416ebf90aSShri Abhyankar *r = row; *c = col; *v = val; 195397b6df1SKris Buschelman } else { 19616ebf90aSShri Abhyankar nz = 0; val = *v; 19716ebf90aSShri Abhyankar for(i=0; i <M; i++) { 19816ebf90aSShri Abhyankar rnz = ai[i+1] - adiag[i]; 19967877ebaSShri Abhyankar ajj = aj + adiag[i]; 20067877ebaSShri Abhyankar v1 = av + adiag[i]; 20167877ebaSShri Abhyankar for(j=0; j<rnz; j++) { 20267877ebaSShri Abhyankar val[nz++] = v1[j]; 20316ebf90aSShri Abhyankar } 20416ebf90aSShri Abhyankar } 20516ebf90aSShri Abhyankar } 20616ebf90aSShri Abhyankar PetscFunctionReturn(0); 20716ebf90aSShri Abhyankar } 20816ebf90aSShri Abhyankar 20916ebf90aSShri Abhyankar #undef __FUNCT__ 21016ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpisbaij_mpisbaij" 211bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpisbaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 21216ebf90aSShri Abhyankar { 21316ebf90aSShri Abhyankar const PetscInt *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj; 21416ebf90aSShri Abhyankar PetscErrorCode ierr; 21516ebf90aSShri Abhyankar PetscInt rstart,nz,i,j,jj,irow,countA,countB; 21616ebf90aSShri Abhyankar PetscInt *row,*col; 21716ebf90aSShri Abhyankar const PetscScalar *av, *bv,*v1,*v2; 21816ebf90aSShri Abhyankar PetscScalar *val; 219397b6df1SKris Buschelman Mat_MPISBAIJ *mat = (Mat_MPISBAIJ*)A->data; 220397b6df1SKris Buschelman Mat_SeqSBAIJ *aa=(Mat_SeqSBAIJ*)(mat->A)->data; 221397b6df1SKris Buschelman Mat_SeqBAIJ *bb=(Mat_SeqBAIJ*)(mat->B)->data; 22216ebf90aSShri Abhyankar 22316ebf90aSShri Abhyankar PetscFunctionBegin; 224d0f46423SBarry Smith ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart; 225397b6df1SKris Buschelman garray = mat->garray; 226397b6df1SKris Buschelman av=aa->a; bv=bb->a; 227397b6df1SKris Buschelman 228bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX){ 22916ebf90aSShri Abhyankar nz = aa->nz + bb->nz; 23016ebf90aSShri Abhyankar *nnz = nz; 231185f6596SHong Zhang ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr); 232185f6596SHong Zhang col = row + nz; 233185f6596SHong Zhang val = (PetscScalar*)(col + nz); 234185f6596SHong Zhang 235397b6df1SKris Buschelman *r = row; *c = col; *v = val; 236397b6df1SKris Buschelman } else { 237397b6df1SKris Buschelman row = *r; col = *c; val = *v; 238397b6df1SKris Buschelman } 239397b6df1SKris Buschelman 240028e57e8SHong Zhang jj = 0; irow = rstart; 241397b6df1SKris Buschelman for ( i=0; i<m; i++ ) { 242397b6df1SKris Buschelman ajj = aj + ai[i]; /* ptr to the beginning of this row */ 243397b6df1SKris Buschelman countA = ai[i+1] - ai[i]; 244397b6df1SKris Buschelman countB = bi[i+1] - bi[i]; 245397b6df1SKris Buschelman bjj = bj + bi[i]; 24616ebf90aSShri Abhyankar v1 = av + ai[i]; 24716ebf90aSShri Abhyankar v2 = bv + bi[i]; 248397b6df1SKris Buschelman 249397b6df1SKris Buschelman /* A-part */ 250397b6df1SKris Buschelman for (j=0; j<countA; j++){ 251bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 252397b6df1SKris Buschelman row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift; 253397b6df1SKris Buschelman } 25416ebf90aSShri Abhyankar val[jj++] = v1[j]; 255397b6df1SKris Buschelman } 25616ebf90aSShri Abhyankar 25716ebf90aSShri Abhyankar /* B-part */ 25816ebf90aSShri Abhyankar for(j=0; j < countB; j++){ 259bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 260397b6df1SKris Buschelman row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift; 261397b6df1SKris Buschelman } 26216ebf90aSShri Abhyankar val[jj++] = v2[j]; 26316ebf90aSShri Abhyankar } 26416ebf90aSShri Abhyankar irow++; 26516ebf90aSShri Abhyankar } 26616ebf90aSShri Abhyankar PetscFunctionReturn(0); 26716ebf90aSShri Abhyankar } 26816ebf90aSShri Abhyankar 26916ebf90aSShri Abhyankar #undef __FUNCT__ 27016ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpiaij_mpiaij" 271bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpiaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 27216ebf90aSShri Abhyankar { 27316ebf90aSShri Abhyankar const PetscInt *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj; 27416ebf90aSShri Abhyankar PetscErrorCode ierr; 27516ebf90aSShri Abhyankar PetscInt rstart,nz,i,j,jj,irow,countA,countB; 27616ebf90aSShri Abhyankar PetscInt *row,*col; 27716ebf90aSShri Abhyankar const PetscScalar *av, *bv,*v1,*v2; 27816ebf90aSShri Abhyankar PetscScalar *val; 27916ebf90aSShri Abhyankar Mat_MPIAIJ *mat = (Mat_MPIAIJ*)A->data; 28016ebf90aSShri Abhyankar Mat_SeqAIJ *aa=(Mat_SeqAIJ*)(mat->A)->data; 28116ebf90aSShri Abhyankar Mat_SeqAIJ *bb=(Mat_SeqAIJ*)(mat->B)->data; 28216ebf90aSShri Abhyankar 28316ebf90aSShri Abhyankar PetscFunctionBegin; 28416ebf90aSShri Abhyankar ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart; 28516ebf90aSShri Abhyankar garray = mat->garray; 28616ebf90aSShri Abhyankar av=aa->a; bv=bb->a; 28716ebf90aSShri Abhyankar 288bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX){ 28916ebf90aSShri Abhyankar nz = aa->nz + bb->nz; 29016ebf90aSShri Abhyankar *nnz = nz; 291185f6596SHong Zhang ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr); 292185f6596SHong Zhang col = row + nz; 293185f6596SHong Zhang val = (PetscScalar*)(col + nz); 294185f6596SHong Zhang 29516ebf90aSShri Abhyankar *r = row; *c = col; *v = val; 29616ebf90aSShri Abhyankar } else { 29716ebf90aSShri Abhyankar row = *r; col = *c; val = *v; 29816ebf90aSShri Abhyankar } 29916ebf90aSShri Abhyankar 30016ebf90aSShri Abhyankar jj = 0; irow = rstart; 30116ebf90aSShri Abhyankar for ( i=0; i<m; i++ ) { 30216ebf90aSShri Abhyankar ajj = aj + ai[i]; /* ptr to the beginning of this row */ 30316ebf90aSShri Abhyankar countA = ai[i+1] - ai[i]; 30416ebf90aSShri Abhyankar countB = bi[i+1] - bi[i]; 30516ebf90aSShri Abhyankar bjj = bj + bi[i]; 30616ebf90aSShri Abhyankar v1 = av + ai[i]; 30716ebf90aSShri Abhyankar v2 = bv + bi[i]; 30816ebf90aSShri Abhyankar 30916ebf90aSShri Abhyankar /* A-part */ 31016ebf90aSShri Abhyankar for (j=0; j<countA; j++){ 311bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX){ 31216ebf90aSShri Abhyankar row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift; 31316ebf90aSShri Abhyankar } 31416ebf90aSShri Abhyankar val[jj++] = v1[j]; 31516ebf90aSShri Abhyankar } 31616ebf90aSShri Abhyankar 31716ebf90aSShri Abhyankar /* B-part */ 31816ebf90aSShri Abhyankar for(j=0; j < countB; j++){ 319bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX){ 32016ebf90aSShri Abhyankar row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift; 32116ebf90aSShri Abhyankar } 32216ebf90aSShri Abhyankar val[jj++] = v2[j]; 32316ebf90aSShri Abhyankar } 32416ebf90aSShri Abhyankar irow++; 32516ebf90aSShri Abhyankar } 32616ebf90aSShri Abhyankar PetscFunctionReturn(0); 32716ebf90aSShri Abhyankar } 32816ebf90aSShri Abhyankar 32916ebf90aSShri Abhyankar #undef __FUNCT__ 33067877ebaSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpibaij_mpiaij" 331bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpibaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 33267877ebaSShri Abhyankar { 33367877ebaSShri Abhyankar Mat_MPIBAIJ *mat = (Mat_MPIBAIJ*)A->data; 33467877ebaSShri Abhyankar Mat_SeqBAIJ *aa=(Mat_SeqBAIJ*)(mat->A)->data; 33567877ebaSShri Abhyankar Mat_SeqBAIJ *bb=(Mat_SeqBAIJ*)(mat->B)->data; 33667877ebaSShri Abhyankar const PetscInt *ai = aa->i, *bi = bb->i, *aj = aa->j, *bj = bb->j,*ajj, *bjj; 337d985c460SShri Abhyankar const PetscInt *garray = mat->garray,mbs=mat->mbs,rstart=A->rmap->rstart; 33867877ebaSShri Abhyankar const PetscInt bs = A->rmap->bs,bs2=mat->bs2; 33967877ebaSShri Abhyankar PetscErrorCode ierr; 34067877ebaSShri Abhyankar PetscInt nz,i,j,k,n,jj,irow,countA,countB,idx; 34167877ebaSShri Abhyankar PetscInt *row,*col; 34267877ebaSShri Abhyankar const PetscScalar *av=aa->a, *bv=bb->a,*v1,*v2; 34367877ebaSShri Abhyankar PetscScalar *val; 34467877ebaSShri Abhyankar 34567877ebaSShri Abhyankar PetscFunctionBegin; 34667877ebaSShri Abhyankar 347bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 34867877ebaSShri Abhyankar nz = bs2*(aa->nz + bb->nz); 34967877ebaSShri Abhyankar *nnz = nz; 350185f6596SHong Zhang ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr); 351185f6596SHong Zhang col = row + nz; 352185f6596SHong Zhang val = (PetscScalar*)(col + nz); 353185f6596SHong Zhang 35467877ebaSShri Abhyankar *r = row; *c = col; *v = val; 35567877ebaSShri Abhyankar } else { 35667877ebaSShri Abhyankar row = *r; col = *c; val = *v; 35767877ebaSShri Abhyankar } 35867877ebaSShri Abhyankar 359d985c460SShri Abhyankar jj = 0; irow = rstart; 36067877ebaSShri Abhyankar for ( i=0; i<mbs; i++ ) { 36167877ebaSShri Abhyankar countA = ai[i+1] - ai[i]; 36267877ebaSShri Abhyankar countB = bi[i+1] - bi[i]; 36367877ebaSShri Abhyankar ajj = aj + ai[i]; 36467877ebaSShri Abhyankar bjj = bj + bi[i]; 36567877ebaSShri Abhyankar v1 = av + bs2*ai[i]; 36667877ebaSShri Abhyankar v2 = bv + bs2*bi[i]; 36767877ebaSShri Abhyankar 36867877ebaSShri Abhyankar idx = 0; 36967877ebaSShri Abhyankar /* A-part */ 37067877ebaSShri Abhyankar for (k=0; k<countA; k++){ 37167877ebaSShri Abhyankar for (j=0; j<bs; j++) { 37267877ebaSShri Abhyankar for (n=0; n<bs; n++) { 373bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX){ 374d985c460SShri Abhyankar row[jj] = irow + n + shift; 375d985c460SShri Abhyankar col[jj] = rstart + bs*ajj[k] + j + shift; 37667877ebaSShri Abhyankar } 37767877ebaSShri Abhyankar val[jj++] = v1[idx++]; 37867877ebaSShri Abhyankar } 37967877ebaSShri Abhyankar } 38067877ebaSShri Abhyankar } 38167877ebaSShri Abhyankar 38267877ebaSShri Abhyankar idx = 0; 38367877ebaSShri Abhyankar /* B-part */ 38467877ebaSShri Abhyankar for(k=0; k<countB; k++){ 38567877ebaSShri Abhyankar for (j=0; j<bs; j++) { 38667877ebaSShri Abhyankar for (n=0; n<bs; n++) { 387bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX){ 388d985c460SShri Abhyankar row[jj] = irow + n + shift; 389d985c460SShri Abhyankar col[jj] = bs*garray[bjj[k]] + j + shift; 39067877ebaSShri Abhyankar } 391d985c460SShri Abhyankar val[jj++] = v2[idx++]; 39267877ebaSShri Abhyankar } 39367877ebaSShri Abhyankar } 39467877ebaSShri Abhyankar } 395d985c460SShri Abhyankar irow += bs; 39667877ebaSShri Abhyankar } 39767877ebaSShri Abhyankar PetscFunctionReturn(0); 39867877ebaSShri Abhyankar } 39967877ebaSShri Abhyankar 40067877ebaSShri Abhyankar #undef __FUNCT__ 40116ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpiaij_mpisbaij" 402bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpiaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 40316ebf90aSShri Abhyankar { 40416ebf90aSShri Abhyankar const PetscInt *ai, *aj,*adiag, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj; 40516ebf90aSShri Abhyankar PetscErrorCode ierr; 406e0bace9bSHong Zhang PetscInt rstart,nz,nza,nzb,i,j,jj,irow,countA,countB; 40716ebf90aSShri Abhyankar PetscInt *row,*col; 40816ebf90aSShri Abhyankar const PetscScalar *av, *bv,*v1,*v2; 40916ebf90aSShri Abhyankar PetscScalar *val; 41016ebf90aSShri Abhyankar Mat_MPIAIJ *mat = (Mat_MPIAIJ*)A->data; 41116ebf90aSShri Abhyankar Mat_SeqAIJ *aa=(Mat_SeqAIJ*)(mat->A)->data; 41216ebf90aSShri Abhyankar Mat_SeqAIJ *bb=(Mat_SeqAIJ*)(mat->B)->data; 41316ebf90aSShri Abhyankar 41416ebf90aSShri Abhyankar PetscFunctionBegin; 41516ebf90aSShri Abhyankar ai=aa->i; aj=aa->j; adiag=aa->diag; 41616ebf90aSShri Abhyankar bi=bb->i; bj=bb->j; garray = mat->garray; 41716ebf90aSShri Abhyankar av=aa->a; bv=bb->a; 41816ebf90aSShri Abhyankar rstart = A->rmap->rstart; 41916ebf90aSShri Abhyankar 420bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 421e0bace9bSHong Zhang nza = 0; /* num of upper triangular entries in mat->A, including diagonals */ 422e0bace9bSHong Zhang nzb = 0; /* num of upper triangular entries in mat->B */ 42316ebf90aSShri Abhyankar for(i=0; i<m; i++){ 424e0bace9bSHong Zhang nza += (ai[i+1] - adiag[i]); 42516ebf90aSShri Abhyankar countB = bi[i+1] - bi[i]; 42616ebf90aSShri Abhyankar bjj = bj + bi[i]; 427e0bace9bSHong Zhang for (j=0; j<countB; j++){ 428e0bace9bSHong Zhang if (garray[bjj[j]] > rstart) nzb++; 429e0bace9bSHong Zhang } 430e0bace9bSHong Zhang } 43116ebf90aSShri Abhyankar 432e0bace9bSHong Zhang nz = nza + nzb; /* total nz of upper triangular part of mat */ 43316ebf90aSShri Abhyankar *nnz = nz; 434185f6596SHong Zhang ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr); 435185f6596SHong Zhang col = row + nz; 436185f6596SHong Zhang val = (PetscScalar*)(col + nz); 437185f6596SHong Zhang 43816ebf90aSShri Abhyankar *r = row; *c = col; *v = val; 43916ebf90aSShri Abhyankar } else { 44016ebf90aSShri Abhyankar row = *r; col = *c; val = *v; 44116ebf90aSShri Abhyankar } 44216ebf90aSShri Abhyankar 44316ebf90aSShri Abhyankar jj = 0; irow = rstart; 44416ebf90aSShri Abhyankar for ( i=0; i<m; i++ ) { 44516ebf90aSShri Abhyankar ajj = aj + adiag[i]; /* ptr to the beginning of the diagonal of this row */ 44616ebf90aSShri Abhyankar v1 = av + adiag[i]; 44716ebf90aSShri Abhyankar countA = ai[i+1] - adiag[i]; 44816ebf90aSShri Abhyankar countB = bi[i+1] - bi[i]; 44916ebf90aSShri Abhyankar bjj = bj + bi[i]; 45016ebf90aSShri Abhyankar v2 = bv + bi[i]; 45116ebf90aSShri Abhyankar 45216ebf90aSShri Abhyankar /* A-part */ 45316ebf90aSShri Abhyankar for (j=0; j<countA; j++){ 454bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 45516ebf90aSShri Abhyankar row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift; 45616ebf90aSShri Abhyankar } 45716ebf90aSShri Abhyankar val[jj++] = v1[j]; 45816ebf90aSShri Abhyankar } 45916ebf90aSShri Abhyankar 46016ebf90aSShri Abhyankar /* B-part */ 46116ebf90aSShri Abhyankar for(j=0; j < countB; j++){ 46216ebf90aSShri Abhyankar if (garray[bjj[j]] > rstart) { 463bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 46416ebf90aSShri Abhyankar row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift; 46516ebf90aSShri Abhyankar } 46616ebf90aSShri Abhyankar val[jj++] = v2[j]; 46716ebf90aSShri Abhyankar } 468397b6df1SKris Buschelman } 469397b6df1SKris Buschelman irow++; 470397b6df1SKris Buschelman } 471397b6df1SKris Buschelman PetscFunctionReturn(0); 472397b6df1SKris Buschelman } 473397b6df1SKris Buschelman 474397b6df1SKris Buschelman #undef __FUNCT__ 4753924e44cSKris Buschelman #define __FUNCT__ "MatDestroy_MUMPS" 476dfbe8321SBarry Smith PetscErrorCode MatDestroy_MUMPS(Mat A) 477dfbe8321SBarry Smith { 478f0c56d0fSKris Buschelman Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 479dfbe8321SBarry Smith PetscErrorCode ierr; 480b24902e0SBarry Smith 481397b6df1SKris Buschelman PetscFunctionBegin; 482bf0cc555SLisandro Dalcin if (lu && lu->CleanUpMUMPS) { 483397b6df1SKris Buschelman /* Terminate instance, deallocate memories */ 4846bf464f9SBarry Smith ierr = PetscFree2(lu->id.sol_loc,lu->id.isol_loc);CHKERRQ(ierr); 4856bf464f9SBarry Smith ierr = VecScatterDestroy(&lu->scat_rhs);CHKERRQ(ierr); 4866bf464f9SBarry Smith ierr = VecDestroy(&lu->b_seq);CHKERRQ(ierr); 487bf0cc555SLisandro Dalcin ierr = VecScatterDestroy(&lu->scat_sol);CHKERRQ(ierr); 488bf0cc555SLisandro Dalcin ierr = VecDestroy(&lu->x_seq);CHKERRQ(ierr); 4896bf464f9SBarry Smith ierr=PetscFree(lu->id.perm_in);CHKERRQ(ierr); 490185f6596SHong Zhang ierr = PetscFree(lu->irn);CHKERRQ(ierr); 491397b6df1SKris Buschelman lu->id.job=JOB_END; 492397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 493397b6df1SKris Buschelman zmumps_c(&lu->id); 494397b6df1SKris Buschelman #else 495397b6df1SKris Buschelman dmumps_c(&lu->id); 496397b6df1SKris Buschelman #endif 497397b6df1SKris Buschelman ierr = MPI_Comm_free(&(lu->comm_mumps));CHKERRQ(ierr); 498397b6df1SKris Buschelman } 499bf0cc555SLisandro Dalcin if (lu && lu->Destroy) { 500bf0cc555SLisandro Dalcin ierr = (lu->Destroy)(A);CHKERRQ(ierr); 501bf0cc555SLisandro Dalcin } 502bf0cc555SLisandro Dalcin ierr = PetscFree(A->spptr);CHKERRQ(ierr); 503bf0cc555SLisandro Dalcin 50497969023SHong Zhang /* clear composed functions */ 50597969023SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatFactorGetSolverPackage_C","",PETSC_NULL);CHKERRQ(ierr); 506f250808bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMumpsSetIcntl_C","",PETSC_NULL);CHKERRQ(ierr); 507397b6df1SKris Buschelman PetscFunctionReturn(0); 508397b6df1SKris Buschelman } 509397b6df1SKris Buschelman 510397b6df1SKris Buschelman #undef __FUNCT__ 511f6c57405SHong Zhang #define __FUNCT__ "MatSolve_MUMPS" 512b24902e0SBarry Smith PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x) 513b24902e0SBarry Smith { 514f0c56d0fSKris Buschelman Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 515d54de34fSKris Buschelman PetscScalar *array; 51667877ebaSShri Abhyankar Vec b_seq; 517329ec9b3SHong Zhang IS is_iden,is_petsc; 518dfbe8321SBarry Smith PetscErrorCode ierr; 519329ec9b3SHong Zhang PetscInt i; 520397b6df1SKris Buschelman 521397b6df1SKris Buschelman PetscFunctionBegin; 522329ec9b3SHong Zhang lu->id.nrhs = 1; 52367877ebaSShri Abhyankar b_seq = lu->b_seq; 524397b6df1SKris Buschelman if (lu->size > 1){ 525329ec9b3SHong Zhang /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */ 52667877ebaSShri Abhyankar ierr = VecScatterBegin(lu->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 52767877ebaSShri Abhyankar ierr = VecScatterEnd(lu->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 52867877ebaSShri Abhyankar if (!lu->myid) {ierr = VecGetArray(b_seq,&array);CHKERRQ(ierr);} 529397b6df1SKris Buschelman } else { /* size == 1 */ 530397b6df1SKris Buschelman ierr = VecCopy(b,x);CHKERRQ(ierr); 531397b6df1SKris Buschelman ierr = VecGetArray(x,&array);CHKERRQ(ierr); 532397b6df1SKris Buschelman } 533397b6df1SKris Buschelman if (!lu->myid) { /* define rhs on the host */ 5348278f211SHong Zhang lu->id.nrhs = 1; 535397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 536397b6df1SKris Buschelman lu->id.rhs = (mumps_double_complex*)array; 537397b6df1SKris Buschelman #else 538397b6df1SKris Buschelman lu->id.rhs = array; 539397b6df1SKris Buschelman #endif 540397b6df1SKris Buschelman } 541397b6df1SKris Buschelman 542397b6df1SKris Buschelman /* solve phase */ 543329ec9b3SHong Zhang /*-------------*/ 5443d472b54SHong Zhang lu->id.job = JOB_SOLVE; 545397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 546397b6df1SKris Buschelman zmumps_c(&lu->id); 547397b6df1SKris Buschelman #else 548397b6df1SKris Buschelman dmumps_c(&lu->id); 549397b6df1SKris Buschelman #endif 55065e19b50SBarry Smith if (lu->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",lu->id.INFOG(1)); 551397b6df1SKris Buschelman 552329ec9b3SHong Zhang if (lu->size > 1) { /* convert mumps distributed solution to petsc mpi x */ 553329ec9b3SHong Zhang if (!lu->nSolve){ /* create scatter scat_sol */ 554329ec9b3SHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,lu->id.lsol_loc,0,1,&is_iden);CHKERRQ(ierr); /* from */ 555329ec9b3SHong Zhang for (i=0; i<lu->id.lsol_loc; i++){ 556329ec9b3SHong Zhang lu->id.isol_loc[i] -= 1; /* change Fortran style to C style */ 557397b6df1SKris Buschelman } 55870b3c8c7SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,lu->id.lsol_loc,lu->id.isol_loc,PETSC_COPY_VALUES,&is_petsc);CHKERRQ(ierr); /* to */ 559329ec9b3SHong Zhang ierr = VecScatterCreate(lu->x_seq,is_iden,x,is_petsc,&lu->scat_sol);CHKERRQ(ierr); 5606bf464f9SBarry Smith ierr = ISDestroy(&is_iden);CHKERRQ(ierr); 5616bf464f9SBarry Smith ierr = ISDestroy(&is_petsc);CHKERRQ(ierr); 562397b6df1SKris Buschelman } 563ca9f406cSSatish Balay ierr = VecScatterBegin(lu->scat_sol,lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 564ca9f406cSSatish Balay ierr = VecScatterEnd(lu->scat_sol,lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 565329ec9b3SHong Zhang } 566329ec9b3SHong Zhang lu->nSolve++; 567397b6df1SKris Buschelman PetscFunctionReturn(0); 568397b6df1SKris Buschelman } 569397b6df1SKris Buschelman 57051d5961aSHong Zhang #undef __FUNCT__ 57151d5961aSHong Zhang #define __FUNCT__ "MatSolveTranspose_MUMPS" 57251d5961aSHong Zhang PetscErrorCode MatSolveTranspose_MUMPS(Mat A,Vec b,Vec x) 57351d5961aSHong Zhang { 57451d5961aSHong Zhang Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 57551d5961aSHong Zhang PetscErrorCode ierr; 57651d5961aSHong Zhang 57751d5961aSHong Zhang PetscFunctionBegin; 57851d5961aSHong Zhang lu->id.ICNTL(9) = 0; 5790ad0caddSJed Brown ierr = MatSolve_MUMPS(A,b,x);CHKERRQ(ierr); 58051d5961aSHong Zhang lu->id.ICNTL(9) = 1; 58151d5961aSHong Zhang PetscFunctionReturn(0); 58251d5961aSHong Zhang } 58351d5961aSHong Zhang 584e0b74bf9SHong Zhang #undef __FUNCT__ 585e0b74bf9SHong Zhang #define __FUNCT__ "MatMatSolve_MUMPS" 586e0b74bf9SHong Zhang PetscErrorCode MatMatSolve_MUMPS(Mat A,Mat B,Mat X) 587e0b74bf9SHong Zhang { 588bda8bf91SBarry Smith PetscErrorCode ierr; 589bda8bf91SBarry Smith PetscBool flg; 590bda8bf91SBarry Smith 591e0b74bf9SHong Zhang PetscFunctionBegin; 592251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,PETSC_NULL);CHKERRQ(ierr); 593bda8bf91SBarry Smith if (!flg) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix"); 594251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,PETSC_NULL);CHKERRQ(ierr); 595bda8bf91SBarry Smith if (!flg) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix"); SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatSolve_MUMPS() is not implemented yet"); 596e0b74bf9SHong Zhang PetscFunctionReturn(0); 597e0b74bf9SHong Zhang } 598e0b74bf9SHong Zhang 599ace3df97SHong Zhang #if !defined(PETSC_USE_COMPLEX) 600a58c3f20SHong Zhang /* 601a58c3f20SHong Zhang input: 602a58c3f20SHong Zhang F: numeric factor 603a58c3f20SHong Zhang output: 604a58c3f20SHong Zhang nneg: total number of negative pivots 605a58c3f20SHong Zhang nzero: 0 606a58c3f20SHong Zhang npos: (global dimension of F) - nneg 607a58c3f20SHong Zhang */ 608a58c3f20SHong Zhang 609a58c3f20SHong Zhang #undef __FUNCT__ 610a58c3f20SHong Zhang #define __FUNCT__ "MatGetInertia_SBAIJMUMPS" 611dfbe8321SBarry Smith PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos) 612a58c3f20SHong Zhang { 613a58c3f20SHong Zhang Mat_MUMPS *lu =(Mat_MUMPS*)F->spptr; 614dfbe8321SBarry Smith PetscErrorCode ierr; 615c1490034SHong Zhang PetscMPIInt size; 616a58c3f20SHong Zhang 617a58c3f20SHong Zhang PetscFunctionBegin; 6187adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)F)->comm,&size);CHKERRQ(ierr); 619bcb30aebSHong 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 */ 62065e19b50SBarry Smith if (size > 1 && lu->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",lu->id.INFOG(13)); 621a58c3f20SHong Zhang if (nneg){ 622a58c3f20SHong Zhang if (!lu->myid){ 623a58c3f20SHong Zhang *nneg = lu->id.INFOG(12); 624a58c3f20SHong Zhang } 625bcb30aebSHong Zhang ierr = MPI_Bcast(nneg,1,MPI_INT,0,lu->comm_mumps);CHKERRQ(ierr); 626a58c3f20SHong Zhang } 627a58c3f20SHong Zhang if (nzero) *nzero = 0; 628d0f46423SBarry Smith if (npos) *npos = F->rmap->N - (*nneg); 629a58c3f20SHong Zhang PetscFunctionReturn(0); 630a58c3f20SHong Zhang } 631ace3df97SHong Zhang #endif /* !defined(PETSC_USE_COMPLEX) */ 632a58c3f20SHong Zhang 633397b6df1SKris Buschelman #undef __FUNCT__ 634f6c57405SHong Zhang #define __FUNCT__ "MatFactorNumeric_MUMPS" 6350481f469SBarry Smith PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info) 636af281ebdSHong Zhang { 637dcd589f8SShri Abhyankar Mat_MUMPS *lu =(Mat_MUMPS*)(F)->spptr; 6386849ba73SBarry Smith PetscErrorCode ierr; 639e09efc27SHong Zhang Mat F_diag; 640ace3abfcSBarry Smith PetscBool isMPIAIJ; 641397b6df1SKris Buschelman 642397b6df1SKris Buschelman PetscFunctionBegin; 643eb85809bSHong Zhang ierr = (*lu->ConvertToTriples)(A, 1, MAT_REUSE_MATRIX, &lu->nz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr); 644397b6df1SKris Buschelman 645397b6df1SKris Buschelman /* numerical factorization phase */ 646329ec9b3SHong Zhang /*-------------------------------*/ 6473d472b54SHong Zhang lu->id.job = JOB_FACTNUMERIC; 648958c9bccSBarry Smith if(!lu->id.ICNTL(18)) { 649a7aca84bSHong Zhang if (!lu->myid) { 650397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 651397b6df1SKris Buschelman lu->id.a = (mumps_double_complex*)lu->val; 652397b6df1SKris Buschelman #else 653397b6df1SKris Buschelman lu->id.a = lu->val; 654397b6df1SKris Buschelman #endif 655397b6df1SKris Buschelman } 656397b6df1SKris Buschelman } else { 657397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 658397b6df1SKris Buschelman lu->id.a_loc = (mumps_double_complex*)lu->val; 659397b6df1SKris Buschelman #else 660397b6df1SKris Buschelman lu->id.a_loc = lu->val; 661397b6df1SKris Buschelman #endif 662397b6df1SKris Buschelman } 663397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 664397b6df1SKris Buschelman zmumps_c(&lu->id); 665397b6df1SKris Buschelman #else 666397b6df1SKris Buschelman dmumps_c(&lu->id); 667397b6df1SKris Buschelman #endif 668397b6df1SKris Buschelman if (lu->id.INFOG(1) < 0) { 66965e19b50SBarry Smith if (lu->id.INFO(1) == -13) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: Cannot allocate required memory %d megabytes\n",lu->id.INFO(2)); 67065e19b50SBarry Smith else SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: INFO(1)=%d, INFO(2)=%d\n",lu->id.INFO(1),lu->id.INFO(2)); 671397b6df1SKris Buschelman } 67265e19b50SBarry Smith if (!lu->myid && lu->id.ICNTL(16) > 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB," lu->id.ICNTL(16):=%d\n",lu->id.INFOG(16)); 673397b6df1SKris Buschelman 6748ada1bb4SHong Zhang if (lu->size > 1){ 675251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr); 676a214ac2aSShri Abhyankar if(isMPIAIJ) { 677dcd589f8SShri Abhyankar F_diag = ((Mat_MPIAIJ *)(F)->data)->A; 678e09efc27SHong Zhang } else { 679dcd589f8SShri Abhyankar F_diag = ((Mat_MPISBAIJ *)(F)->data)->A; 680e09efc27SHong Zhang } 681e09efc27SHong Zhang F_diag->assembled = PETSC_TRUE; 682329ec9b3SHong Zhang if (lu->nSolve){ 6836bf464f9SBarry Smith ierr = VecScatterDestroy(&lu->scat_sol);CHKERRQ(ierr); 6840e83c824SBarry Smith ierr = PetscFree2(lu->id.sol_loc,lu->id.isol_loc);CHKERRQ(ierr); 6856bf464f9SBarry Smith ierr = VecDestroy(&lu->x_seq);CHKERRQ(ierr); 686329ec9b3SHong Zhang } 6878ada1bb4SHong Zhang } 688dcd589f8SShri Abhyankar (F)->assembled = PETSC_TRUE; 689397b6df1SKris Buschelman lu->matstruc = SAME_NONZERO_PATTERN; 690ace87b0dSHong Zhang lu->CleanUpMUMPS = PETSC_TRUE; 691329ec9b3SHong Zhang lu->nSolve = 0; 69267877ebaSShri Abhyankar 69367877ebaSShri Abhyankar if (lu->size > 1){ 69467877ebaSShri Abhyankar /* distributed solution */ 69567877ebaSShri Abhyankar if (!lu->nSolve){ 69667877ebaSShri Abhyankar /* Create x_seq=sol_loc for repeated use */ 69767877ebaSShri Abhyankar PetscInt lsol_loc; 69867877ebaSShri Abhyankar PetscScalar *sol_loc; 69967877ebaSShri Abhyankar lsol_loc = lu->id.INFO(23); /* length of sol_loc */ 70067877ebaSShri Abhyankar ierr = PetscMalloc2(lsol_loc,PetscScalar,&sol_loc,lsol_loc,PetscInt,&lu->id.isol_loc);CHKERRQ(ierr); 70167877ebaSShri Abhyankar lu->id.lsol_loc = lsol_loc; 70267877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 70367877ebaSShri Abhyankar lu->id.sol_loc = (mumps_double_complex*)sol_loc; 70467877ebaSShri Abhyankar #else 70567877ebaSShri Abhyankar lu->id.sol_loc = sol_loc; 70667877ebaSShri Abhyankar #endif 707778a2246SBarry Smith ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,lsol_loc,sol_loc,&lu->x_seq);CHKERRQ(ierr); 70867877ebaSShri Abhyankar } 70967877ebaSShri Abhyankar } 710397b6df1SKris Buschelman PetscFunctionReturn(0); 711397b6df1SKris Buschelman } 712397b6df1SKris Buschelman 7139a2535b5SHong Zhang /* Sets MUMPS options from the options database */ 714dcd589f8SShri Abhyankar #undef __FUNCT__ 7159a2535b5SHong Zhang #define __FUNCT__ "PetscSetMUMPSFromOptions" 7169a2535b5SHong Zhang PetscErrorCode PetscSetMUMPSFromOptions(Mat F, Mat A) 717dcd589f8SShri Abhyankar { 7189a2535b5SHong Zhang Mat_MUMPS *mumps = (Mat_MUMPS*)F->spptr; 719dcd589f8SShri Abhyankar PetscErrorCode ierr; 720dcd589f8SShri Abhyankar PetscInt icntl; 721ace3abfcSBarry Smith PetscBool flg; 722dcd589f8SShri Abhyankar 723dcd589f8SShri Abhyankar PetscFunctionBegin; 724dcd589f8SShri Abhyankar ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"MUMPS Options","Mat");CHKERRQ(ierr); 7259a2535b5SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_1","ICNTL(1): output stream for error messages","None",mumps->id.ICNTL(1),&icntl,&flg);CHKERRQ(ierr); 7269a2535b5SHong Zhang if (flg) mumps->id.ICNTL(1) = icntl; 7279a2535b5SHong 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); 7289a2535b5SHong Zhang if (flg) mumps->id.ICNTL(2) = icntl; 7299a2535b5SHong 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); 7309a2535b5SHong Zhang if (flg) mumps->id.ICNTL(3) = icntl; 731dcd589f8SShri Abhyankar 7329a2535b5SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",mumps->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr); 7339a2535b5SHong Zhang if (flg) mumps->id.ICNTL(4) = icntl; 7349a2535b5SHong Zhang if (mumps->id.ICNTL(4) || PetscLogPrintInfo ) mumps->id.ICNTL(3) = 6; /* resume MUMPS default id.ICNTL(3) = 6 */ 7359a2535b5SHong Zhang 7369a2535b5SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_6","ICNTL(6): permuting and/or scaling the matrix (0 to 7)","None",mumps->id.ICNTL(6),&icntl,&flg);CHKERRQ(ierr); 7379a2535b5SHong Zhang if (flg) mumps->id.ICNTL(6) = icntl; 7389a2535b5SHong Zhang 7399a2535b5SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): matrix ordering (0 to 7). 3=Scotch, 4=PORD, 5=Metis","None",mumps->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr); 740dcd589f8SShri Abhyankar if (flg) { 7419a2535b5SHong Zhang if (icntl== 1 && mumps->size > 1){ 742e32f2f54SBarry Smith 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"); 743dcd589f8SShri Abhyankar } else { 7449a2535b5SHong Zhang mumps->id.ICNTL(7) = icntl; 745dcd589f8SShri Abhyankar } 746dcd589f8SShri Abhyankar } 747e0b74bf9SHong Zhang 74870544d5fSHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_8","ICNTL(8): scaling strategy (-2 to 8 or 77)","None",mumps->id.ICNTL(8),&mumps->id.ICNTL(8),PETSC_NULL);CHKERRQ(ierr); 7499a2535b5SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_10","ICNTL(10): max num of refinements","None",mumps->id.ICNTL(10),&mumps->id.ICNTL(10),PETSC_NULL);CHKERRQ(ierr); 7509a2535b5SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_11","ICNTL(11): statistics related to the linear system solved (via -ksp_view)","None",mumps->id.ICNTL(11),&mumps->id.ICNTL(11),PETSC_NULL);CHKERRQ(ierr); 75170544d5fSHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): efficiency control: defines the ordering strategy with scaling constraints (0 to 3)","None",mumps->id.ICNTL(12),&mumps->id.ICNTL(12),PETSC_NULL);CHKERRQ(ierr); 7529a2535b5SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): efficiency control: with or without ScaLAPACK","None",mumps->id.ICNTL(13),&mumps->id.ICNTL(13),PETSC_NULL);CHKERRQ(ierr); 7539a2535b5SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_14","ICNTL(14): percentage of estimated workspace increase","None",mumps->id.ICNTL(14),&mumps->id.ICNTL(14),PETSC_NULL);CHKERRQ(ierr); 7549a2535b5SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_19","ICNTL(19): Schur complement","None",mumps->id.ICNTL(19),&mumps->id.ICNTL(19),PETSC_NULL);CHKERRQ(ierr); 7559a2535b5SHong Zhang 7569a2535b5SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_22","ICNTL(22): in-core/out-of-core facility (0 or 1)","None",mumps->id.ICNTL(22),&mumps->id.ICNTL(22),PETSC_NULL);CHKERRQ(ierr); 7579a2535b5SHong Zhang 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),PETSC_NULL);CHKERRQ(ierr); 7589a2535b5SHong Zhang 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),PETSC_NULL);CHKERRQ(ierr); 7599a2535b5SHong Zhang if (mumps->id.ICNTL(24)){ 7609a2535b5SHong Zhang mumps->id.ICNTL(13) = 1; /* turn-off ScaLAPACK to help with the correct detection of null pivots */ 761d7ebd59bSHong Zhang } 762d7ebd59bSHong Zhang 7639a2535b5SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_25","ICNTL(25): computation of a null space basis","None",mumps->id.ICNTL(25),&mumps->id.ICNTL(25),PETSC_NULL);CHKERRQ(ierr); 7649a2535b5SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_26","ICNTL(26): Schur options for right-hand side or solution vector","None",mumps->id.ICNTL(26),&mumps->id.ICNTL(26),PETSC_NULL);CHKERRQ(ierr); 7659a2535b5SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_27","ICNTL(27): experimental parameter","None",mumps->id.ICNTL(27),&mumps->id.ICNTL(27),PETSC_NULL);CHKERRQ(ierr); 7669a2535b5SHong Zhang 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),PETSC_NULL);CHKERRQ(ierr); 7679a2535b5SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_29","ICNTL(29): parallel ordering 1 = ptscotch 2 = parmetis","None",mumps->id.ICNTL(29),&mumps->id.ICNTL(29),PETSC_NULL);CHKERRQ(ierr); 7689a2535b5SHong Zhang 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),PETSC_NULL);CHKERRQ(ierr); 7699a2535b5SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_31","ICNTL(31): factors can be discarded in the solve phase","None",mumps->id.ICNTL(31),&mumps->id.ICNTL(31),PETSC_NULL);CHKERRQ(ierr); 7709a2535b5SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_33","ICNTL(33): compute determinant","None",mumps->id.ICNTL(33),&mumps->id.ICNTL(33),PETSC_NULL);CHKERRQ(ierr); 771dcd589f8SShri Abhyankar 7729a2535b5SHong Zhang ierr = PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",mumps->id.CNTL(1),&mumps->id.CNTL(1),PETSC_NULL);CHKERRQ(ierr); 7739a2535b5SHong Zhang ierr = PetscOptionsReal("-mat_mumps_cntl_2","CNTL(2): stopping criterion of refinement","None",mumps->id.CNTL(2),&mumps->id.CNTL(2),PETSC_NULL);CHKERRQ(ierr); 7749a2535b5SHong Zhang ierr = PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",mumps->id.CNTL(3),&mumps->id.CNTL(3),PETSC_NULL);CHKERRQ(ierr); 7759a2535b5SHong Zhang ierr = PetscOptionsReal("-mat_mumps_cntl_4","CNTL(4): value for static pivoting","None",mumps->id.CNTL(4),&mumps->id.CNTL(4),PETSC_NULL);CHKERRQ(ierr); 7769a2535b5SHong Zhang ierr = PetscOptionsReal("-mat_mumps_cntl_5","CNTL(5): fixation for null pivots","None",mumps->id.CNTL(5),&mumps->id.CNTL(5),PETSC_NULL);CHKERRQ(ierr); 777*e5bb22a1SHong Zhang 778*e5bb22a1SHong Zhang ierr = PetscOptionsString("-mat_mumps_ooc_tmpdir", "out of core directory", "None", mumps->id.ooc_tmpdir, mumps->id.ooc_tmpdir, 256, PETSC_NULL); 779dcd589f8SShri Abhyankar PetscOptionsEnd(); 780dcd589f8SShri Abhyankar PetscFunctionReturn(0); 781dcd589f8SShri Abhyankar } 782dcd589f8SShri Abhyankar 783dcd589f8SShri Abhyankar #undef __FUNCT__ 784dcd589f8SShri Abhyankar #define __FUNCT__ "PetscInitializeMUMPS" 785f697e70eSHong Zhang PetscErrorCode PetscInitializeMUMPS(Mat A,Mat_MUMPS* mumps) 786dcd589f8SShri Abhyankar { 787dcd589f8SShri Abhyankar PetscErrorCode ierr; 788dcd589f8SShri Abhyankar 789dcd589f8SShri Abhyankar PetscFunctionBegin; 790f697e70eSHong Zhang ierr = MPI_Comm_rank(((PetscObject)A)->comm, &mumps->myid); 791f697e70eSHong Zhang ierr = MPI_Comm_size(((PetscObject)A)->comm,&mumps->size);CHKERRQ(ierr); 792f697e70eSHong Zhang ierr = MPI_Comm_dup(((PetscObject)A)->comm,&(mumps->comm_mumps));CHKERRQ(ierr); 793f697e70eSHong Zhang mumps->id.comm_fortran = MPI_Comm_c2f(mumps->comm_mumps); 794f697e70eSHong Zhang 795f697e70eSHong Zhang mumps->id.job = JOB_INIT; 796f697e70eSHong Zhang mumps->id.par = 1; /* host participates factorizaton and solve */ 797f697e70eSHong Zhang mumps->id.sym = mumps->sym; 798dcd589f8SShri Abhyankar #if defined(PETSC_USE_COMPLEX) 799f697e70eSHong Zhang zmumps_c(&mumps->id); 800dcd589f8SShri Abhyankar #else 801f697e70eSHong Zhang dmumps_c(&mumps->id); 802dcd589f8SShri Abhyankar #endif 803f697e70eSHong Zhang 804f697e70eSHong Zhang mumps->CleanUpMUMPS = PETSC_FALSE; 805f697e70eSHong Zhang mumps->scat_rhs = PETSC_NULL; 806f697e70eSHong Zhang mumps->scat_sol = PETSC_NULL; 807f697e70eSHong Zhang mumps->nSolve = 0; 8089a2535b5SHong Zhang 80970544d5fSHong Zhang /* set PETSc-MUMPS default options - override MUMPS default */ 8109a2535b5SHong Zhang mumps->id.ICNTL(3) = 0; 8119a2535b5SHong Zhang mumps->id.ICNTL(4) = 0; 8129a2535b5SHong Zhang if (mumps->size == 1){ 8139a2535b5SHong Zhang mumps->id.ICNTL(18) = 0; /* centralized assembled matrix input */ 8149a2535b5SHong Zhang } else { 8159a2535b5SHong Zhang mumps->id.ICNTL(18) = 3; /* distributed assembled matrix input */ 81670544d5fSHong Zhang mumps->id.ICNTL(21) = 1; /* distributed solution */ 8179a2535b5SHong Zhang } 818dcd589f8SShri Abhyankar PetscFunctionReturn(0); 819dcd589f8SShri Abhyankar } 820dcd589f8SShri Abhyankar 821397b6df1SKris Buschelman /* Note the Petsc r and c permutations are ignored */ 822397b6df1SKris Buschelman #undef __FUNCT__ 823f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS" 8240481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 825b24902e0SBarry Smith { 826719d5645SBarry Smith Mat_MUMPS *lu = (Mat_MUMPS*)F->spptr; 827dcd589f8SShri Abhyankar PetscErrorCode ierr; 82867877ebaSShri Abhyankar Vec b; 82967877ebaSShri Abhyankar IS is_iden; 83067877ebaSShri Abhyankar const PetscInt M = A->rmap->N; 831397b6df1SKris Buschelman 832397b6df1SKris Buschelman PetscFunctionBegin; 833b24902e0SBarry Smith lu->matstruc = DIFFERENT_NONZERO_PATTERN; 834dcd589f8SShri Abhyankar 8359a2535b5SHong Zhang /* Set MUMPS options from the options database */ 8369a2535b5SHong Zhang ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr); 837dcd589f8SShri Abhyankar 838eb85809bSHong Zhang ierr = (*lu->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &lu->nz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr); 839dcd589f8SShri Abhyankar 84067877ebaSShri Abhyankar /* analysis phase */ 84167877ebaSShri Abhyankar /*----------------*/ 8423d472b54SHong Zhang lu->id.job = JOB_FACTSYMBOLIC; 84367877ebaSShri Abhyankar lu->id.n = M; 84467877ebaSShri Abhyankar switch (lu->id.ICNTL(18)){ 84567877ebaSShri Abhyankar case 0: /* centralized assembled matrix input */ 84667877ebaSShri Abhyankar if (!lu->myid) { 84767877ebaSShri Abhyankar lu->id.nz =lu->nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn; 84867877ebaSShri Abhyankar if (lu->id.ICNTL(6)>1){ 84967877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 85067877ebaSShri Abhyankar lu->id.a = (mumps_double_complex*)lu->val; 85167877ebaSShri Abhyankar #else 85267877ebaSShri Abhyankar lu->id.a = lu->val; 85367877ebaSShri Abhyankar #endif 85467877ebaSShri Abhyankar } 855e0b74bf9SHong Zhang if (lu->id.ICNTL(7) == 1){ /* use user-provide matrix ordering */ 856e0b74bf9SHong Zhang if (!lu->myid) { 857e0b74bf9SHong Zhang const PetscInt *idx; 858e0b74bf9SHong Zhang PetscInt i,*perm_in; 859e0b74bf9SHong Zhang ierr = PetscMalloc(M*sizeof(PetscInt),&perm_in);CHKERRQ(ierr); 860e0b74bf9SHong Zhang ierr = ISGetIndices(r,&idx);CHKERRQ(ierr); 861e0b74bf9SHong Zhang lu->id.perm_in = perm_in; 862e0b74bf9SHong Zhang for (i=0; i<M; i++) perm_in[i] = idx[i]+1; /* perm_in[]: start from 1, not 0! */ 863e0b74bf9SHong Zhang ierr = ISRestoreIndices(r,&idx);CHKERRQ(ierr); 864e0b74bf9SHong Zhang } 865e0b74bf9SHong Zhang } 86667877ebaSShri Abhyankar } 86767877ebaSShri Abhyankar break; 86867877ebaSShri Abhyankar case 3: /* distributed assembled matrix input (size>1) */ 86967877ebaSShri Abhyankar lu->id.nz_loc = lu->nz; 87067877ebaSShri Abhyankar lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn; 87167877ebaSShri Abhyankar if (lu->id.ICNTL(6)>1) { 87267877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 87367877ebaSShri Abhyankar lu->id.a_loc = (mumps_double_complex*)lu->val; 87467877ebaSShri Abhyankar #else 87567877ebaSShri Abhyankar lu->id.a_loc = lu->val; 87667877ebaSShri Abhyankar #endif 87767877ebaSShri Abhyankar } 87867877ebaSShri Abhyankar /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */ 87967877ebaSShri Abhyankar if (!lu->myid){ 88067877ebaSShri Abhyankar ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&lu->b_seq);CHKERRQ(ierr); 88167877ebaSShri Abhyankar ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr); 88267877ebaSShri Abhyankar } else { 88367877ebaSShri Abhyankar ierr = VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);CHKERRQ(ierr); 88467877ebaSShri Abhyankar ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr); 88567877ebaSShri Abhyankar } 88667877ebaSShri Abhyankar ierr = VecCreate(((PetscObject)A)->comm,&b);CHKERRQ(ierr); 88767877ebaSShri Abhyankar ierr = VecSetSizes(b,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr); 88867877ebaSShri Abhyankar ierr = VecSetFromOptions(b);CHKERRQ(ierr); 88967877ebaSShri Abhyankar 89067877ebaSShri Abhyankar ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr); 8916bf464f9SBarry Smith ierr = ISDestroy(&is_iden);CHKERRQ(ierr); 8926bf464f9SBarry Smith ierr = VecDestroy(&b);CHKERRQ(ierr); 89367877ebaSShri Abhyankar break; 89467877ebaSShri Abhyankar } 89567877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 89667877ebaSShri Abhyankar zmumps_c(&lu->id); 89767877ebaSShri Abhyankar #else 89867877ebaSShri Abhyankar dmumps_c(&lu->id); 89967877ebaSShri Abhyankar #endif 90067877ebaSShri Abhyankar if (lu->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1)); 90167877ebaSShri Abhyankar 902719d5645SBarry Smith F->ops->lufactornumeric = MatFactorNumeric_MUMPS; 903dcd589f8SShri Abhyankar F->ops->solve = MatSolve_MUMPS; 90451d5961aSHong Zhang F->ops->solvetranspose = MatSolveTranspose_MUMPS; 905e0b74bf9SHong Zhang F->ops->matsolve = MatMatSolve_MUMPS; 906b24902e0SBarry Smith PetscFunctionReturn(0); 907b24902e0SBarry Smith } 908b24902e0SBarry Smith 909450b117fSShri Abhyankar /* Note the Petsc r and c permutations are ignored */ 910450b117fSShri Abhyankar #undef __FUNCT__ 911450b117fSShri Abhyankar #define __FUNCT__ "MatLUFactorSymbolic_BAIJMUMPS" 912450b117fSShri Abhyankar PetscErrorCode MatLUFactorSymbolic_BAIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 913450b117fSShri Abhyankar { 914dcd589f8SShri Abhyankar 915450b117fSShri Abhyankar Mat_MUMPS *lu = (Mat_MUMPS*)F->spptr; 916dcd589f8SShri Abhyankar PetscErrorCode ierr; 91767877ebaSShri Abhyankar Vec b; 91867877ebaSShri Abhyankar IS is_iden; 91967877ebaSShri Abhyankar const PetscInt M = A->rmap->N; 920450b117fSShri Abhyankar 921450b117fSShri Abhyankar PetscFunctionBegin; 922450b117fSShri Abhyankar lu->matstruc = DIFFERENT_NONZERO_PATTERN; 923dcd589f8SShri Abhyankar 9249a2535b5SHong Zhang /* Set MUMPS options from the options database */ 9259a2535b5SHong Zhang ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr); 926dcd589f8SShri Abhyankar 927eb85809bSHong Zhang ierr = (*lu->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &lu->nz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr); 92867877ebaSShri Abhyankar 92967877ebaSShri Abhyankar /* analysis phase */ 93067877ebaSShri Abhyankar /*----------------*/ 9313d472b54SHong Zhang lu->id.job = JOB_FACTSYMBOLIC; 93267877ebaSShri Abhyankar lu->id.n = M; 93367877ebaSShri Abhyankar switch (lu->id.ICNTL(18)){ 93467877ebaSShri Abhyankar case 0: /* centralized assembled matrix input */ 93567877ebaSShri Abhyankar if (!lu->myid) { 93667877ebaSShri Abhyankar lu->id.nz =lu->nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn; 93767877ebaSShri Abhyankar if (lu->id.ICNTL(6)>1){ 93867877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 93967877ebaSShri Abhyankar lu->id.a = (mumps_double_complex*)lu->val; 94067877ebaSShri Abhyankar #else 94167877ebaSShri Abhyankar lu->id.a = lu->val; 94267877ebaSShri Abhyankar #endif 94367877ebaSShri Abhyankar } 94467877ebaSShri Abhyankar } 94567877ebaSShri Abhyankar break; 94667877ebaSShri Abhyankar case 3: /* distributed assembled matrix input (size>1) */ 94767877ebaSShri Abhyankar lu->id.nz_loc = lu->nz; 94867877ebaSShri Abhyankar lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn; 94967877ebaSShri Abhyankar if (lu->id.ICNTL(6)>1) { 95067877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 95167877ebaSShri Abhyankar lu->id.a_loc = (mumps_double_complex*)lu->val; 95267877ebaSShri Abhyankar #else 95367877ebaSShri Abhyankar lu->id.a_loc = lu->val; 95467877ebaSShri Abhyankar #endif 95567877ebaSShri Abhyankar } 95667877ebaSShri Abhyankar /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */ 95767877ebaSShri Abhyankar if (!lu->myid){ 95867877ebaSShri Abhyankar ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&lu->b_seq);CHKERRQ(ierr); 95967877ebaSShri Abhyankar ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr); 96067877ebaSShri Abhyankar } else { 96167877ebaSShri Abhyankar ierr = VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);CHKERRQ(ierr); 96267877ebaSShri Abhyankar ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr); 96367877ebaSShri Abhyankar } 96467877ebaSShri Abhyankar ierr = VecCreate(((PetscObject)A)->comm,&b);CHKERRQ(ierr); 96567877ebaSShri Abhyankar ierr = VecSetSizes(b,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr); 96667877ebaSShri Abhyankar ierr = VecSetFromOptions(b);CHKERRQ(ierr); 96767877ebaSShri Abhyankar 96867877ebaSShri Abhyankar ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr); 9696bf464f9SBarry Smith ierr = ISDestroy(&is_iden);CHKERRQ(ierr); 9706bf464f9SBarry Smith ierr = VecDestroy(&b);CHKERRQ(ierr); 97167877ebaSShri Abhyankar break; 97267877ebaSShri Abhyankar } 97367877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 97467877ebaSShri Abhyankar zmumps_c(&lu->id); 97567877ebaSShri Abhyankar #else 97667877ebaSShri Abhyankar dmumps_c(&lu->id); 97767877ebaSShri Abhyankar #endif 97867877ebaSShri Abhyankar if (lu->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1)); 97967877ebaSShri Abhyankar 980450b117fSShri Abhyankar F->ops->lufactornumeric = MatFactorNumeric_MUMPS; 981dcd589f8SShri Abhyankar F->ops->solve = MatSolve_MUMPS; 98251d5961aSHong Zhang F->ops->solvetranspose = MatSolveTranspose_MUMPS; 983450b117fSShri Abhyankar PetscFunctionReturn(0); 984450b117fSShri Abhyankar } 985b24902e0SBarry Smith 986141f4205SHong Zhang /* Note the Petsc r permutation and factor info are ignored */ 987397b6df1SKris Buschelman #undef __FUNCT__ 98867877ebaSShri Abhyankar #define __FUNCT__ "MatCholeskyFactorSymbolic_MUMPS" 98967877ebaSShri Abhyankar PetscErrorCode MatCholeskyFactorSymbolic_MUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info) 990b24902e0SBarry Smith { 9912792810eSHong Zhang Mat_MUMPS *lu = (Mat_MUMPS*)F->spptr; 992dcd589f8SShri Abhyankar PetscErrorCode ierr; 99367877ebaSShri Abhyankar Vec b; 99467877ebaSShri Abhyankar IS is_iden; 99567877ebaSShri Abhyankar const PetscInt M = A->rmap->N; 996397b6df1SKris Buschelman 997397b6df1SKris Buschelman PetscFunctionBegin; 998b24902e0SBarry Smith lu->matstruc = DIFFERENT_NONZERO_PATTERN; 999dcd589f8SShri Abhyankar 10009a2535b5SHong Zhang /* Set MUMPS options from the options database */ 10019a2535b5SHong Zhang ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr); 1002dcd589f8SShri Abhyankar 1003eb85809bSHong Zhang ierr = (*lu->ConvertToTriples)(A, 1 , MAT_INITIAL_MATRIX, &lu->nz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr); 1004dcd589f8SShri Abhyankar 100567877ebaSShri Abhyankar /* analysis phase */ 100667877ebaSShri Abhyankar /*----------------*/ 10073d472b54SHong Zhang lu->id.job = JOB_FACTSYMBOLIC; 100867877ebaSShri Abhyankar lu->id.n = M; 100967877ebaSShri Abhyankar switch (lu->id.ICNTL(18)){ 101067877ebaSShri Abhyankar case 0: /* centralized assembled matrix input */ 101167877ebaSShri Abhyankar if (!lu->myid) { 101267877ebaSShri Abhyankar lu->id.nz =lu->nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn; 101367877ebaSShri Abhyankar if (lu->id.ICNTL(6)>1){ 101467877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 101567877ebaSShri Abhyankar lu->id.a = (mumps_double_complex*)lu->val; 101667877ebaSShri Abhyankar #else 101767877ebaSShri Abhyankar lu->id.a = lu->val; 101867877ebaSShri Abhyankar #endif 101967877ebaSShri Abhyankar } 102067877ebaSShri Abhyankar } 102167877ebaSShri Abhyankar break; 102267877ebaSShri Abhyankar case 3: /* distributed assembled matrix input (size>1) */ 102367877ebaSShri Abhyankar lu->id.nz_loc = lu->nz; 102467877ebaSShri Abhyankar lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn; 102567877ebaSShri Abhyankar if (lu->id.ICNTL(6)>1) { 102667877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 102767877ebaSShri Abhyankar lu->id.a_loc = (mumps_double_complex*)lu->val; 102867877ebaSShri Abhyankar #else 102967877ebaSShri Abhyankar lu->id.a_loc = lu->val; 103067877ebaSShri Abhyankar #endif 103167877ebaSShri Abhyankar } 103267877ebaSShri Abhyankar /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */ 103367877ebaSShri Abhyankar if (!lu->myid){ 103467877ebaSShri Abhyankar ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&lu->b_seq);CHKERRQ(ierr); 103567877ebaSShri Abhyankar ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr); 103667877ebaSShri Abhyankar } else { 103767877ebaSShri Abhyankar ierr = VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);CHKERRQ(ierr); 103867877ebaSShri Abhyankar ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr); 103967877ebaSShri Abhyankar } 104067877ebaSShri Abhyankar ierr = VecCreate(((PetscObject)A)->comm,&b);CHKERRQ(ierr); 104167877ebaSShri Abhyankar ierr = VecSetSizes(b,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr); 104267877ebaSShri Abhyankar ierr = VecSetFromOptions(b);CHKERRQ(ierr); 104367877ebaSShri Abhyankar 104467877ebaSShri Abhyankar ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr); 10456bf464f9SBarry Smith ierr = ISDestroy(&is_iden);CHKERRQ(ierr); 10466bf464f9SBarry Smith ierr = VecDestroy(&b);CHKERRQ(ierr); 104767877ebaSShri Abhyankar break; 104867877ebaSShri Abhyankar } 104967877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 105067877ebaSShri Abhyankar zmumps_c(&lu->id); 105167877ebaSShri Abhyankar #else 105267877ebaSShri Abhyankar dmumps_c(&lu->id); 105367877ebaSShri Abhyankar #endif 105467877ebaSShri Abhyankar if (lu->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1)); 105567877ebaSShri Abhyankar 10562792810eSHong Zhang F->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS; 1057dcd589f8SShri Abhyankar F->ops->solve = MatSolve_MUMPS; 105851d5961aSHong Zhang F->ops->solvetranspose = MatSolve_MUMPS; 1059db4efbfdSBarry Smith #if !defined(PETSC_USE_COMPLEX) 106005aa0992SJose Roman F->ops->getinertia = MatGetInertia_SBAIJMUMPS; 106105aa0992SJose Roman #else 106205aa0992SJose Roman F->ops->getinertia = PETSC_NULL; 1063db4efbfdSBarry Smith #endif 1064b24902e0SBarry Smith PetscFunctionReturn(0); 1065b24902e0SBarry Smith } 1066b24902e0SBarry Smith 1067397b6df1SKris Buschelman #undef __FUNCT__ 106864e6c443SBarry Smith #define __FUNCT__ "MatView_MUMPS" 106964e6c443SBarry Smith PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer) 107074ed9c26SBarry Smith { 1071f6c57405SHong Zhang PetscErrorCode ierr; 107264e6c443SBarry Smith PetscBool iascii; 107364e6c443SBarry Smith PetscViewerFormat format; 107464e6c443SBarry Smith Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 1075f6c57405SHong Zhang 1076f6c57405SHong Zhang PetscFunctionBegin; 107764e6c443SBarry Smith /* check if matrix is mumps type */ 107864e6c443SBarry Smith if (A->ops->solve != MatSolve_MUMPS) PetscFunctionReturn(0); 107964e6c443SBarry Smith 1080251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 108164e6c443SBarry Smith if (iascii) { 108264e6c443SBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 108364e6c443SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO){ 108464e6c443SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr); 108564e6c443SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," SYM (matrix type): %d \n",lu->id.sym);CHKERRQ(ierr); 108664e6c443SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," PAR (host participation): %d \n",lu->id.par);CHKERRQ(ierr); 1087f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(1) (output for error): %d \n",lu->id.ICNTL(1));CHKERRQ(ierr); 1088f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(2) (output of diagnostic msg): %d \n",lu->id.ICNTL(2));CHKERRQ(ierr); 1089f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(3) (output for global info): %d \n",lu->id.ICNTL(3));CHKERRQ(ierr); 1090f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(4) (level of printing): %d \n",lu->id.ICNTL(4));CHKERRQ(ierr); 1091f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(5) (input mat struct): %d \n",lu->id.ICNTL(5));CHKERRQ(ierr); 1092f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(6) (matrix prescaling): %d \n",lu->id.ICNTL(6));CHKERRQ(ierr); 1093d06efebcSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ICNTL(7) (sequentia matrix ordering):%d \n",lu->id.ICNTL(7));CHKERRQ(ierr); 1094f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(8) (scalling strategy): %d \n",lu->id.ICNTL(8));CHKERRQ(ierr); 1095f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(10) (max num of refinements): %d \n",lu->id.ICNTL(10));CHKERRQ(ierr); 1096f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(11) (error analysis): %d \n",lu->id.ICNTL(11));CHKERRQ(ierr); 109734ed7027SBarry Smith if (lu->id.ICNTL(11)>0) { 109834ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," RINFOG(4) (inf norm of input mat): %g\n",lu->id.RINFOG(4));CHKERRQ(ierr); 109934ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," RINFOG(5) (inf norm of solution): %g\n",lu->id.RINFOG(5));CHKERRQ(ierr); 110034ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," RINFOG(6) (inf norm of residual): %g\n",lu->id.RINFOG(6));CHKERRQ(ierr); 110134ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",lu->id.RINFOG(7),lu->id.RINFOG(8));CHKERRQ(ierr); 110234ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," RINFOG(9) (error estimate): %g \n",lu->id.RINFOG(9));CHKERRQ(ierr); 110334ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",lu->id.RINFOG(10),lu->id.RINFOG(11));CHKERRQ(ierr); 1104f6c57405SHong Zhang } 1105f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(12) (efficiency control): %d \n",lu->id.ICNTL(12));CHKERRQ(ierr); 1106f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(13) (efficiency control): %d \n",lu->id.ICNTL(13));CHKERRQ(ierr); 1107f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(14) (percentage of estimated workspace increase): %d \n",lu->id.ICNTL(14));CHKERRQ(ierr); 1108f6c57405SHong Zhang /* ICNTL(15-17) not used */ 1109f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(18) (input mat struct): %d \n",lu->id.ICNTL(18));CHKERRQ(ierr); 1110f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(19) (Shur complement info): %d \n",lu->id.ICNTL(19));CHKERRQ(ierr); 1111f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(20) (rhs sparse pattern): %d \n",lu->id.ICNTL(20));CHKERRQ(ierr); 1112f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(21) (solution struct): %d \n",lu->id.ICNTL(21));CHKERRQ(ierr); 1113c0165424SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(22) (in-core/out-of-core facility): %d \n",lu->id.ICNTL(22));CHKERRQ(ierr); 1114c0165424SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(23) (max size of memory can be allocated locally):%d \n",lu->id.ICNTL(23));CHKERRQ(ierr); 1115c0165424SHong Zhang 1116c0165424SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(24) (detection of null pivot rows): %d \n",lu->id.ICNTL(24));CHKERRQ(ierr); 1117c0165424SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(25) (computation of a null space basis): %d \n",lu->id.ICNTL(25));CHKERRQ(ierr); 1118c0165424SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(26) (Schur options for rhs or solution): %d \n",lu->id.ICNTL(26));CHKERRQ(ierr); 1119c0165424SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(27) (experimental parameter): %d \n",lu->id.ICNTL(27));CHKERRQ(ierr); 112042179a6aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(28) (use parallel or sequential ordering): %d \n",lu->id.ICNTL(28));CHKERRQ(ierr); 112142179a6aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(29) (parallel ordering): %d \n",lu->id.ICNTL(29));CHKERRQ(ierr); 112242179a6aSHong Zhang 112342179a6aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(30) (user-specified set of entries in inv(A)): %d \n",lu->id.ICNTL(30));CHKERRQ(ierr); 112442179a6aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(31) (factors is discarded in the solve phase): %d \n",lu->id.ICNTL(31));CHKERRQ(ierr); 112542179a6aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(33) (compute determinant): %d \n",lu->id.ICNTL(33));CHKERRQ(ierr); 1126f6c57405SHong Zhang 1127f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(1) (relative pivoting threshold): %g \n",lu->id.CNTL(1));CHKERRQ(ierr); 1128f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(2) (stopping criterion of refinement): %g \n",lu->id.CNTL(2));CHKERRQ(ierr); 1129f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(3) (absolute pivoting threshold): %g \n",lu->id.CNTL(3));CHKERRQ(ierr); 1130f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(4) (value of static pivoting): %g \n",lu->id.CNTL(4));CHKERRQ(ierr); 1131c0165424SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(5) (fixation for null pivots): %g \n",lu->id.CNTL(5));CHKERRQ(ierr); 1132f6c57405SHong Zhang 1133f6c57405SHong Zhang /* infomation local to each processor */ 113434ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer, " RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr); 11357b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 113634ed7027SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %g \n",lu->myid,lu->id.RINFO(1));CHKERRQ(ierr); 113734ed7027SBarry Smith ierr = PetscViewerFlush(viewer); 113834ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer, " RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr); 113934ed7027SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %g \n",lu->myid,lu->id.RINFO(2));CHKERRQ(ierr); 114034ed7027SBarry Smith ierr = PetscViewerFlush(viewer); 114134ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer, " RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr); 114234ed7027SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %g \n",lu->myid,lu->id.RINFO(3));CHKERRQ(ierr); 114334ed7027SBarry Smith ierr = PetscViewerFlush(viewer); 1144f6c57405SHong Zhang 114534ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer, " INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization): \n");CHKERRQ(ierr); 114634ed7027SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %d \n",lu->myid,lu->id.INFO(15));CHKERRQ(ierr); 114734ed7027SBarry Smith ierr = PetscViewerFlush(viewer); 1148f6c57405SHong Zhang 114934ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer, " INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization): \n");CHKERRQ(ierr); 115034ed7027SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %d \n",lu->myid,lu->id.INFO(16));CHKERRQ(ierr); 115134ed7027SBarry Smith ierr = PetscViewerFlush(viewer); 1152f6c57405SHong Zhang 115334ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer, " INFO(23) (num of pivots eliminated on this processor after factorization): \n");CHKERRQ(ierr); 115434ed7027SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %d \n",lu->myid,lu->id.INFO(23));CHKERRQ(ierr); 115534ed7027SBarry Smith ierr = PetscViewerFlush(viewer); 11567b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr); 1157f6c57405SHong Zhang 1158f6c57405SHong Zhang if (!lu->myid){ /* information from the host */ 1159f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",lu->id.RINFOG(1));CHKERRQ(ierr); 1160f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",lu->id.RINFOG(2));CHKERRQ(ierr); 1161f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",lu->id.RINFOG(3));CHKERRQ(ierr); 116242179a6aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," (RINFOG(12) RINFOG(13))*2^INFOG(34) (determinant): (%g,%g)*(2^%d)\n",lu->id.RINFOG(12),lu->id.RINFOG(13),lu->id.INFOG(34));CHKERRQ(ierr); 1163f6c57405SHong Zhang 1164f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(3));CHKERRQ(ierr); 1165f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(4));CHKERRQ(ierr); 1166f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(5) (estimated maximum front size in the complete tree): %d \n",lu->id.INFOG(5));CHKERRQ(ierr); 1167f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(6) (number of nodes in the complete tree): %d \n",lu->id.INFOG(6));CHKERRQ(ierr); 11682bd8dccdSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," INFOG(7) (ordering option effectively use after analysis): %d \n",lu->id.INFOG(7));CHKERRQ(ierr); 1169f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",lu->id.INFOG(8));CHKERRQ(ierr); 1170f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(9) (total real/complex workspace to store the matrix factors after factorization): %d \n",lu->id.INFOG(9));CHKERRQ(ierr); 1171f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(10) (total integer space store the matrix factors after factorization): %d \n",lu->id.INFOG(10));CHKERRQ(ierr); 1172f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(11) (order of largest frontal matrix after factorization): %d \n",lu->id.INFOG(11));CHKERRQ(ierr); 1173f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(12) (number of off-diagonal pivots): %d \n",lu->id.INFOG(12));CHKERRQ(ierr); 1174f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(13) (number of delayed pivots after factorization): %d \n",lu->id.INFOG(13));CHKERRQ(ierr); 1175f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(14) (number of memory compress after factorization): %d \n",lu->id.INFOG(14));CHKERRQ(ierr); 1176f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(15) (number of steps of iterative refinement after solution): %d \n",lu->id.INFOG(15));CHKERRQ(ierr); 1177f6c57405SHong 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",lu->id.INFOG(16));CHKERRQ(ierr); 1178f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(17) (estimated size of all MUMPS internal data for factorization after analysis: sum over all processors): %d \n",lu->id.INFOG(17));CHKERRQ(ierr); 1179f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(18) (size of all MUMPS internal data allocated during factorization: value on the most memory consuming processor): %d \n",lu->id.INFOG(18));CHKERRQ(ierr); 1180f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(19) (size of all MUMPS internal data allocated during factorization: sum over all processors): %d \n",lu->id.INFOG(19));CHKERRQ(ierr); 1181f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(20) (estimated number of entries in the factors): %d \n",lu->id.INFOG(20));CHKERRQ(ierr); 1182f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(21) (size in MB of memory effectively used during factorization - value on the most memory consuming processor): %d \n",lu->id.INFOG(21));CHKERRQ(ierr); 1183f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(22) (size in MB of memory effectively used during factorization - sum over all processors): %d \n",lu->id.INFOG(22));CHKERRQ(ierr); 1184f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",lu->id.INFOG(23));CHKERRQ(ierr); 1185f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",lu->id.INFOG(24));CHKERRQ(ierr); 1186f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",lu->id.INFOG(25));CHKERRQ(ierr); 1187f6c57405SHong Zhang } 1188f6c57405SHong Zhang } 1189cb828f0fSHong Zhang } 1190f6c57405SHong Zhang PetscFunctionReturn(0); 1191f6c57405SHong Zhang } 1192f6c57405SHong Zhang 119335bd34faSBarry Smith #undef __FUNCT__ 119435bd34faSBarry Smith #define __FUNCT__ "MatGetInfo_MUMPS" 119535bd34faSBarry Smith PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info) 119635bd34faSBarry Smith { 1197cb828f0fSHong Zhang Mat_MUMPS *mumps =(Mat_MUMPS*)A->spptr; 119835bd34faSBarry Smith 119935bd34faSBarry Smith PetscFunctionBegin; 120035bd34faSBarry Smith info->block_size = 1.0; 1201cb828f0fSHong Zhang info->nz_allocated = mumps->id.INFOG(20); 1202cb828f0fSHong Zhang info->nz_used = mumps->id.INFOG(20); 120335bd34faSBarry Smith info->nz_unneeded = 0.0; 120435bd34faSBarry Smith info->assemblies = 0.0; 120535bd34faSBarry Smith info->mallocs = 0.0; 120635bd34faSBarry Smith info->memory = 0.0; 120735bd34faSBarry Smith info->fill_ratio_given = 0; 120835bd34faSBarry Smith info->fill_ratio_needed = 0; 120935bd34faSBarry Smith info->factor_mallocs = 0; 121035bd34faSBarry Smith PetscFunctionReturn(0); 121135bd34faSBarry Smith } 121235bd34faSBarry Smith 12135ccb76cbSHong Zhang /* -------------------------------------------------------------------------------------------*/ 12145ccb76cbSHong Zhang #undef __FUNCT__ 12155ccb76cbSHong Zhang #define __FUNCT__ "MatMumpsSetIcntl_MUMPS" 12165ccb76cbSHong Zhang PetscErrorCode MatMumpsSetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt ival) 12175ccb76cbSHong Zhang { 12185ccb76cbSHong Zhang Mat_MUMPS *lu =(Mat_MUMPS*)F->spptr; 12195ccb76cbSHong Zhang 12205ccb76cbSHong Zhang PetscFunctionBegin; 12215ccb76cbSHong Zhang lu->id.ICNTL(icntl) = ival; 12225ccb76cbSHong Zhang PetscFunctionReturn(0); 12235ccb76cbSHong Zhang } 12245ccb76cbSHong Zhang 12255ccb76cbSHong Zhang #undef __FUNCT__ 12265ccb76cbSHong Zhang #define __FUNCT__ "MatMumpsSetIcntl" 12275ccb76cbSHong Zhang /*@ 12285ccb76cbSHong Zhang MatMumpsSetIcntl - Set MUMPS parameter ICNTL() 12295ccb76cbSHong Zhang 12305ccb76cbSHong Zhang Logically Collective on Mat 12315ccb76cbSHong Zhang 12325ccb76cbSHong Zhang Input Parameters: 12335ccb76cbSHong Zhang + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 12345ccb76cbSHong Zhang . icntl - index of MUMPS parameter array ICNTL() 12355ccb76cbSHong Zhang - ival - value of MUMPS ICNTL(icntl) 12365ccb76cbSHong Zhang 12375ccb76cbSHong Zhang Options Database: 12385ccb76cbSHong Zhang . -mat_mumps_icntl_<icntl> <ival> 12395ccb76cbSHong Zhang 12405ccb76cbSHong Zhang Level: beginner 12415ccb76cbSHong Zhang 12425ccb76cbSHong Zhang References: MUMPS Users' Guide 12435ccb76cbSHong Zhang 12445ccb76cbSHong Zhang .seealso: MatGetFactor() 12455ccb76cbSHong Zhang @*/ 12465ccb76cbSHong Zhang PetscErrorCode MatMumpsSetIcntl(Mat F,PetscInt icntl,PetscInt ival) 12475ccb76cbSHong Zhang { 12485ccb76cbSHong Zhang PetscErrorCode ierr; 12495ccb76cbSHong Zhang 12505ccb76cbSHong Zhang PetscFunctionBegin; 12515ccb76cbSHong Zhang PetscValidLogicalCollectiveInt(F,icntl,2); 12525ccb76cbSHong Zhang PetscValidLogicalCollectiveInt(F,ival,3); 12535ccb76cbSHong Zhang ierr = PetscTryMethod(F,"MatMumpsSetIcntl_C",(Mat,PetscInt,PetscInt),(F,icntl,ival));CHKERRQ(ierr); 12545ccb76cbSHong Zhang PetscFunctionReturn(0); 12555ccb76cbSHong Zhang } 12565ccb76cbSHong Zhang 125724b6179bSKris Buschelman /*MC 12582692d6eeSBarry Smith MATSOLVERMUMPS - A matrix type providing direct solvers (LU and Cholesky) for 125924b6179bSKris Buschelman distributed and sequential matrices via the external package MUMPS. 126024b6179bSKris Buschelman 126141c8de11SBarry Smith Works with MATAIJ and MATSBAIJ matrices 126224b6179bSKris Buschelman 126324b6179bSKris Buschelman Options Database Keys: 1264fb8376fbSHong Zhang + -mat_mumps_icntl_4 <0,...,4> - print level 126524b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide) 126664e6c443SBarry Smith . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guidec) 126724b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T 126824b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements 126994b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view 127024b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide) 127124b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide) 127224b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide) 127324b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide) 127424b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold 127524b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement 127624b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold 127724b6179bSKris Buschelman 127824b6179bSKris Buschelman Level: beginner 127924b6179bSKris Buschelman 128041c8de11SBarry Smith .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage 128141c8de11SBarry Smith 128224b6179bSKris Buschelman M*/ 128324b6179bSKris Buschelman 12842877fffaSHong Zhang EXTERN_C_BEGIN 128535bd34faSBarry Smith #undef __FUNCT__ 128635bd34faSBarry Smith #define __FUNCT__ "MatFactorGetSolverPackage_mumps" 128735bd34faSBarry Smith PetscErrorCode MatFactorGetSolverPackage_mumps(Mat A,const MatSolverPackage *type) 128835bd34faSBarry Smith { 128935bd34faSBarry Smith PetscFunctionBegin; 12902692d6eeSBarry Smith *type = MATSOLVERMUMPS; 129135bd34faSBarry Smith PetscFunctionReturn(0); 129235bd34faSBarry Smith } 129335bd34faSBarry Smith EXTERN_C_END 129435bd34faSBarry Smith 129535bd34faSBarry Smith EXTERN_C_BEGIN 1296bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI AIJ matrices */ 12972877fffaSHong Zhang #undef __FUNCT__ 1298bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_aij_mumps" 1299bccb9932SShri Abhyankar PetscErrorCode MatGetFactor_aij_mumps(Mat A,MatFactorType ftype,Mat *F) 13002877fffaSHong Zhang { 13012877fffaSHong Zhang Mat B; 13022877fffaSHong Zhang PetscErrorCode ierr; 13032877fffaSHong Zhang Mat_MUMPS *mumps; 1304ace3abfcSBarry Smith PetscBool isSeqAIJ; 13052877fffaSHong Zhang 13062877fffaSHong Zhang PetscFunctionBegin; 13072877fffaSHong Zhang /* Create the factorization matrix */ 1308251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr); 13092877fffaSHong Zhang ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 13102877fffaSHong Zhang ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 13112877fffaSHong Zhang ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 1312bccb9932SShri Abhyankar if (isSeqAIJ) { 13132877fffaSHong Zhang ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 1314bccb9932SShri Abhyankar } else { 1315bccb9932SShri Abhyankar ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 1316bccb9932SShri Abhyankar } 13172877fffaSHong Zhang 131816ebf90aSShri Abhyankar ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr); 13192877fffaSHong Zhang B->ops->view = MatView_MUMPS; 132035bd34faSBarry Smith B->ops->getinfo = MatGetInfo_MUMPS; 132135bd34faSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr); 13225ccb76cbSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMumpsSetIcntl_C","MatMumpsSetIcntl_MUMPS",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr); 1323450b117fSShri Abhyankar if (ftype == MAT_FACTOR_LU) { 1324450b117fSShri Abhyankar B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS; 1325d5f3da31SBarry Smith B->factortype = MAT_FACTOR_LU; 1326bccb9932SShri Abhyankar if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqaij; 1327bccb9932SShri Abhyankar else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpiaij; 1328746480a1SHong Zhang mumps->sym = 0; 1329dcd589f8SShri Abhyankar } else { 133067877ebaSShri Abhyankar B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS; 1331450b117fSShri Abhyankar B->factortype = MAT_FACTOR_CHOLESKY; 1332bccb9932SShri Abhyankar if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqsbaij; 1333bccb9932SShri Abhyankar else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpisbaij; 13346fdc2a6dSBarry Smith if (A->spd_set && A->spd) mumps->sym = 1; 13356fdc2a6dSBarry Smith else mumps->sym = 2; 1336450b117fSShri Abhyankar } 13372877fffaSHong Zhang 13382877fffaSHong Zhang mumps->isAIJ = PETSC_TRUE; 1339bf0cc555SLisandro Dalcin mumps->Destroy = B->ops->destroy; 13402877fffaSHong Zhang B->ops->destroy = MatDestroy_MUMPS; 13412877fffaSHong Zhang B->spptr = (void*)mumps; 1342f697e70eSHong Zhang ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr); 1343746480a1SHong Zhang 13442877fffaSHong Zhang *F = B; 13452877fffaSHong Zhang PetscFunctionReturn(0); 13462877fffaSHong Zhang } 13472877fffaSHong Zhang EXTERN_C_END 13482877fffaSHong Zhang 1349bccb9932SShri Abhyankar 13502877fffaSHong Zhang EXTERN_C_BEGIN 1351bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI SBAIJ matrices */ 13522877fffaSHong Zhang #undef __FUNCT__ 1353bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_sbaij_mumps" 1354bccb9932SShri Abhyankar PetscErrorCode MatGetFactor_sbaij_mumps(Mat A,MatFactorType ftype,Mat *F) 13552877fffaSHong Zhang { 13562877fffaSHong Zhang Mat B; 13572877fffaSHong Zhang PetscErrorCode ierr; 13582877fffaSHong Zhang Mat_MUMPS *mumps; 1359ace3abfcSBarry Smith PetscBool isSeqSBAIJ; 13602877fffaSHong Zhang 13612877fffaSHong Zhang PetscFunctionBegin; 1362e7e72b3dSBarry Smith if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix"); 1363bccb9932SShri Abhyankar if(A->rmap->bs > 1) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with block size > 1 with MUMPS Cholesky, use AIJ matrix instead"); 1364251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr); 13652877fffaSHong Zhang /* Create the factorization matrix */ 13662877fffaSHong Zhang ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 13672877fffaSHong Zhang ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 13682877fffaSHong Zhang ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 136916ebf90aSShri Abhyankar ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr); 1370bccb9932SShri Abhyankar if (isSeqSBAIJ) { 1371bccb9932SShri Abhyankar ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr); 137216ebf90aSShri Abhyankar mumps->ConvertToTriples = MatConvertToTriples_seqsbaij_seqsbaij; 1373dcd589f8SShri Abhyankar } else { 1374bccb9932SShri Abhyankar ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 1375bccb9932SShri Abhyankar mumps->ConvertToTriples = MatConvertToTriples_mpisbaij_mpisbaij; 1376bccb9932SShri Abhyankar } 1377bccb9932SShri Abhyankar 137867877ebaSShri Abhyankar B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS; 1379bccb9932SShri Abhyankar B->ops->view = MatView_MUMPS; 1380bccb9932SShri Abhyankar ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr); 1381f250808bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMumpsSetIcntl_C","MatMumpsSetIcntl",MatMumpsSetIcntl);CHKERRQ(ierr); 1382f4762488SHong Zhang B->factortype = MAT_FACTOR_CHOLESKY; 13836fdc2a6dSBarry Smith if (A->spd_set && A->spd) mumps->sym = 1; 13846fdc2a6dSBarry Smith else mumps->sym = 2; 1385a214ac2aSShri Abhyankar 1386bccb9932SShri Abhyankar mumps->isAIJ = PETSC_FALSE; 1387bf0cc555SLisandro Dalcin mumps->Destroy = B->ops->destroy; 1388f3c0ef26SHong Zhang B->ops->destroy = MatDestroy_MUMPS; 13892877fffaSHong Zhang B->spptr = (void*)mumps; 1390f697e70eSHong Zhang ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr); 1391746480a1SHong Zhang 13922877fffaSHong Zhang *F = B; 13932877fffaSHong Zhang PetscFunctionReturn(0); 13942877fffaSHong Zhang } 13952877fffaSHong Zhang EXTERN_C_END 139697969023SHong Zhang 1397450b117fSShri Abhyankar EXTERN_C_BEGIN 1398450b117fSShri Abhyankar #undef __FUNCT__ 1399bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_baij_mumps" 1400bccb9932SShri Abhyankar PetscErrorCode MatGetFactor_baij_mumps(Mat A,MatFactorType ftype,Mat *F) 140167877ebaSShri Abhyankar { 140267877ebaSShri Abhyankar Mat B; 140367877ebaSShri Abhyankar PetscErrorCode ierr; 140467877ebaSShri Abhyankar Mat_MUMPS *mumps; 1405ace3abfcSBarry Smith PetscBool isSeqBAIJ; 140667877ebaSShri Abhyankar 140767877ebaSShri Abhyankar PetscFunctionBegin; 140867877ebaSShri Abhyankar /* Create the factorization matrix */ 1409251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);CHKERRQ(ierr); 141067877ebaSShri Abhyankar ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 141167877ebaSShri Abhyankar ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 141267877ebaSShri Abhyankar ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 1413bccb9932SShri Abhyankar if (isSeqBAIJ) { 141467877ebaSShri Abhyankar ierr = MatSeqBAIJSetPreallocation(B,A->rmap->bs,0,PETSC_NULL);CHKERRQ(ierr); 1415bccb9932SShri Abhyankar } else { 141667877ebaSShri Abhyankar ierr = MatMPIBAIJSetPreallocation(B,A->rmap->bs,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 1417bccb9932SShri Abhyankar } 1418450b117fSShri Abhyankar 141967877ebaSShri Abhyankar ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr); 1420450b117fSShri Abhyankar if (ftype == MAT_FACTOR_LU) { 1421450b117fSShri Abhyankar B->ops->lufactorsymbolic = MatLUFactorSymbolic_BAIJMUMPS; 1422450b117fSShri Abhyankar B->factortype = MAT_FACTOR_LU; 1423bccb9932SShri Abhyankar if (isSeqBAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqbaij_seqaij; 1424bccb9932SShri Abhyankar else mumps->ConvertToTriples = MatConvertToTriples_mpibaij_mpiaij; 1425746480a1SHong Zhang mumps->sym = 0; 1426746480a1SHong Zhang } else { 1427746480a1SHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc BAIJ matrices with MUMPS Cholesky, use SBAIJ or AIJ matrix instead\n"); 1428450b117fSShri Abhyankar } 1429bccb9932SShri Abhyankar 1430450b117fSShri Abhyankar B->ops->view = MatView_MUMPS; 1431450b117fSShri Abhyankar ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr); 14325ccb76cbSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMumpsSetIcntl_C","MatMumpsSetIcntl_MUMPS",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr); 1433450b117fSShri Abhyankar 1434450b117fSShri Abhyankar mumps->isAIJ = PETSC_TRUE; 1435bf0cc555SLisandro Dalcin mumps->Destroy = B->ops->destroy; 1436450b117fSShri Abhyankar B->ops->destroy = MatDestroy_MUMPS; 1437450b117fSShri Abhyankar B->spptr = (void*)mumps; 1438f697e70eSHong Zhang ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr); 1439746480a1SHong Zhang 1440450b117fSShri Abhyankar *F = B; 1441450b117fSShri Abhyankar PetscFunctionReturn(0); 1442450b117fSShri Abhyankar } 1443450b117fSShri Abhyankar EXTERN_C_END 1444a214ac2aSShri Abhyankar 1445