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; 45*bf0cc555SLisandro 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; 10267877ebaSShri Abhyankar PetscInt nz,idx=0,rnz,i,j,k,m,ii; 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 ii = 0; 11767877ebaSShri Abhyankar ajj = aj + ai[i]; 11867877ebaSShri Abhyankar rnz = ai[i+1] - ai[i]; 11967877ebaSShri Abhyankar for(k=0; k<rnz; k++) { 12067877ebaSShri Abhyankar for(j=0; j<bs; j++) { 12167877ebaSShri Abhyankar for(m=0; m<bs; m++) { 12267877ebaSShri Abhyankar row[idx] = i*bs + m + shift; 123cf3759fdSShri Abhyankar col[idx++] = bs*(ajj[k]) + j + shift; 12467877ebaSShri Abhyankar } 12567877ebaSShri Abhyankar } 12667877ebaSShri Abhyankar } 12767877ebaSShri Abhyankar } 128cf3759fdSShri Abhyankar *r = row; *c = col; 12967877ebaSShri Abhyankar } 13067877ebaSShri Abhyankar PetscFunctionReturn(0); 13167877ebaSShri Abhyankar } 13267877ebaSShri Abhyankar 13367877ebaSShri Abhyankar #undef __FUNCT__ 13416ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqsbaij_seqsbaij" 135bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqsbaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 13616ebf90aSShri Abhyankar { 13767877ebaSShri Abhyankar const PetscInt *ai, *aj,*ajj,M=A->rmap->n; 13867877ebaSShri Abhyankar PetscInt nz,rnz,i,j; 13916ebf90aSShri Abhyankar PetscErrorCode ierr; 14016ebf90aSShri Abhyankar PetscInt *row,*col; 14116ebf90aSShri Abhyankar Mat_SeqSBAIJ *aa=(Mat_SeqSBAIJ*)A->data; 14216ebf90aSShri Abhyankar 14316ebf90aSShri Abhyankar PetscFunctionBegin; 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; 40616ebf90aSShri Abhyankar PetscInt rstart,nz,nza,nzb_low,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) { 42116ebf90aSShri Abhyankar nza = 0;nzb_low = 0; 42216ebf90aSShri Abhyankar for(i=0; i<m; i++){ 42316ebf90aSShri Abhyankar nza = nza + (ai[i+1] - adiag[i]); 42416ebf90aSShri Abhyankar countB = bi[i+1] - bi[i]; 42516ebf90aSShri Abhyankar bjj = bj + bi[i]; 42616ebf90aSShri Abhyankar 42716ebf90aSShri Abhyankar j = 0; 42816ebf90aSShri Abhyankar while(garray[bjj[j]] < rstart) { 42916ebf90aSShri Abhyankar if(j == countB) break; 43016ebf90aSShri Abhyankar j++;nzb_low++; 43116ebf90aSShri Abhyankar } 43216ebf90aSShri Abhyankar } 43316ebf90aSShri Abhyankar /* Total nz = nz for the upper triangular A part + nz for the 2nd B part */ 43416ebf90aSShri Abhyankar nz = nza + (bb->nz - nzb_low); 43516ebf90aSShri Abhyankar *nnz = nz; 436185f6596SHong Zhang ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr); 437185f6596SHong Zhang col = row + nz; 438185f6596SHong Zhang val = (PetscScalar*)(col + nz); 439185f6596SHong Zhang 44016ebf90aSShri Abhyankar *r = row; *c = col; *v = val; 44116ebf90aSShri Abhyankar } else { 44216ebf90aSShri Abhyankar row = *r; col = *c; val = *v; 44316ebf90aSShri Abhyankar } 44416ebf90aSShri Abhyankar 44516ebf90aSShri Abhyankar jj = 0; irow = rstart; 44616ebf90aSShri Abhyankar for ( i=0; i<m; i++ ) { 44716ebf90aSShri Abhyankar ajj = aj + adiag[i]; /* ptr to the beginning of the diagonal of this row */ 44816ebf90aSShri Abhyankar v1 = av + adiag[i]; 44916ebf90aSShri Abhyankar countA = ai[i+1] - adiag[i]; 45016ebf90aSShri Abhyankar countB = bi[i+1] - bi[i]; 45116ebf90aSShri Abhyankar bjj = bj + bi[i]; 45216ebf90aSShri Abhyankar v2 = bv + bi[i]; 45316ebf90aSShri Abhyankar 45416ebf90aSShri Abhyankar /* A-part */ 45516ebf90aSShri Abhyankar for (j=0; j<countA; j++){ 456bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 45716ebf90aSShri Abhyankar row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift; 45816ebf90aSShri Abhyankar } 45916ebf90aSShri Abhyankar val[jj++] = v1[j]; 46016ebf90aSShri Abhyankar } 46116ebf90aSShri Abhyankar 46216ebf90aSShri Abhyankar /* B-part */ 46316ebf90aSShri Abhyankar for(j=0; j < countB; j++){ 46416ebf90aSShri Abhyankar if (garray[bjj[j]] > rstart) { 465bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 46616ebf90aSShri Abhyankar row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift; 46716ebf90aSShri Abhyankar } 46816ebf90aSShri Abhyankar val[jj++] = v2[j]; 46916ebf90aSShri Abhyankar } 470397b6df1SKris Buschelman } 471397b6df1SKris Buschelman irow++; 472397b6df1SKris Buschelman } 473397b6df1SKris Buschelman PetscFunctionReturn(0); 474397b6df1SKris Buschelman } 475397b6df1SKris Buschelman 476397b6df1SKris Buschelman #undef __FUNCT__ 4773924e44cSKris Buschelman #define __FUNCT__ "MatDestroy_MUMPS" 478dfbe8321SBarry Smith PetscErrorCode MatDestroy_MUMPS(Mat A) 479dfbe8321SBarry Smith { 480f0c56d0fSKris Buschelman Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 481dfbe8321SBarry Smith PetscErrorCode ierr; 482b24902e0SBarry Smith 483397b6df1SKris Buschelman PetscFunctionBegin; 484*bf0cc555SLisandro Dalcin if (lu && lu->CleanUpMUMPS) { 485397b6df1SKris Buschelman /* Terminate instance, deallocate memories */ 4866bf464f9SBarry Smith ierr = PetscFree2(lu->id.sol_loc,lu->id.isol_loc);CHKERRQ(ierr); 4876bf464f9SBarry Smith ierr = VecScatterDestroy(&lu->scat_rhs);CHKERRQ(ierr); 4886bf464f9SBarry Smith ierr = VecDestroy(&lu->b_seq);CHKERRQ(ierr); 489*bf0cc555SLisandro Dalcin ierr = VecScatterDestroy(&lu->scat_sol);CHKERRQ(ierr); 490*bf0cc555SLisandro Dalcin ierr = VecDestroy(&lu->x_seq);CHKERRQ(ierr); 4916bf464f9SBarry Smith ierr=PetscFree(lu->id.perm_in);CHKERRQ(ierr); 492185f6596SHong Zhang ierr = PetscFree(lu->irn);CHKERRQ(ierr); 493397b6df1SKris Buschelman lu->id.job=JOB_END; 494397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 495397b6df1SKris Buschelman zmumps_c(&lu->id); 496397b6df1SKris Buschelman #else 497397b6df1SKris Buschelman dmumps_c(&lu->id); 498397b6df1SKris Buschelman #endif 499397b6df1SKris Buschelman ierr = MPI_Comm_free(&(lu->comm_mumps));CHKERRQ(ierr); 500397b6df1SKris Buschelman } 501*bf0cc555SLisandro Dalcin if (lu && lu->Destroy) { 502*bf0cc555SLisandro Dalcin ierr = (lu->Destroy)(A);CHKERRQ(ierr); 503*bf0cc555SLisandro Dalcin } 504*bf0cc555SLisandro Dalcin ierr = PetscFree(A->spptr);CHKERRQ(ierr); 505*bf0cc555SLisandro Dalcin 50697969023SHong Zhang /* clear composed functions */ 50797969023SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatFactorGetSolverPackage_C","",PETSC_NULL);CHKERRQ(ierr); 508f250808bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMumpsSetIcntl_C","",PETSC_NULL);CHKERRQ(ierr); 509397b6df1SKris Buschelman PetscFunctionReturn(0); 510397b6df1SKris Buschelman } 511397b6df1SKris Buschelman 512397b6df1SKris Buschelman #undef __FUNCT__ 513f6c57405SHong Zhang #define __FUNCT__ "MatSolve_MUMPS" 514b24902e0SBarry Smith PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x) 515b24902e0SBarry Smith { 516f0c56d0fSKris Buschelman Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 517d54de34fSKris Buschelman PetscScalar *array; 51867877ebaSShri Abhyankar Vec b_seq; 519329ec9b3SHong Zhang IS is_iden,is_petsc; 520dfbe8321SBarry Smith PetscErrorCode ierr; 521329ec9b3SHong Zhang PetscInt i; 522397b6df1SKris Buschelman 523397b6df1SKris Buschelman PetscFunctionBegin; 524329ec9b3SHong Zhang lu->id.nrhs = 1; 52567877ebaSShri Abhyankar b_seq = lu->b_seq; 526397b6df1SKris Buschelman if (lu->size > 1){ 527329ec9b3SHong Zhang /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */ 52867877ebaSShri Abhyankar ierr = VecScatterBegin(lu->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 52967877ebaSShri Abhyankar ierr = VecScatterEnd(lu->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 53067877ebaSShri Abhyankar if (!lu->myid) {ierr = VecGetArray(b_seq,&array);CHKERRQ(ierr);} 531397b6df1SKris Buschelman } else { /* size == 1 */ 532397b6df1SKris Buschelman ierr = VecCopy(b,x);CHKERRQ(ierr); 533397b6df1SKris Buschelman ierr = VecGetArray(x,&array);CHKERRQ(ierr); 534397b6df1SKris Buschelman } 535397b6df1SKris Buschelman if (!lu->myid) { /* define rhs on the host */ 5368278f211SHong Zhang lu->id.nrhs = 1; 537397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 538397b6df1SKris Buschelman lu->id.rhs = (mumps_double_complex*)array; 539397b6df1SKris Buschelman #else 540397b6df1SKris Buschelman lu->id.rhs = array; 541397b6df1SKris Buschelman #endif 542397b6df1SKris Buschelman } 543397b6df1SKris Buschelman 544397b6df1SKris Buschelman /* solve phase */ 545329ec9b3SHong Zhang /*-------------*/ 5463d472b54SHong Zhang lu->id.job = JOB_SOLVE; 547397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 548397b6df1SKris Buschelman zmumps_c(&lu->id); 549397b6df1SKris Buschelman #else 550397b6df1SKris Buschelman dmumps_c(&lu->id); 551397b6df1SKris Buschelman #endif 55265e19b50SBarry 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)); 553397b6df1SKris Buschelman 554329ec9b3SHong Zhang if (lu->size > 1) { /* convert mumps distributed solution to petsc mpi x */ 555329ec9b3SHong Zhang if (!lu->nSolve){ /* create scatter scat_sol */ 556329ec9b3SHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,lu->id.lsol_loc,0,1,&is_iden);CHKERRQ(ierr); /* from */ 557329ec9b3SHong Zhang for (i=0; i<lu->id.lsol_loc; i++){ 558329ec9b3SHong Zhang lu->id.isol_loc[i] -= 1; /* change Fortran style to C style */ 559397b6df1SKris Buschelman } 56070b3c8c7SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,lu->id.lsol_loc,lu->id.isol_loc,PETSC_COPY_VALUES,&is_petsc);CHKERRQ(ierr); /* to */ 561329ec9b3SHong Zhang ierr = VecScatterCreate(lu->x_seq,is_iden,x,is_petsc,&lu->scat_sol);CHKERRQ(ierr); 5626bf464f9SBarry Smith ierr = ISDestroy(&is_iden);CHKERRQ(ierr); 5636bf464f9SBarry Smith ierr = ISDestroy(&is_petsc);CHKERRQ(ierr); 564397b6df1SKris Buschelman } 565ca9f406cSSatish Balay ierr = VecScatterBegin(lu->scat_sol,lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 566ca9f406cSSatish Balay ierr = VecScatterEnd(lu->scat_sol,lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 567329ec9b3SHong Zhang } 568329ec9b3SHong Zhang lu->nSolve++; 569397b6df1SKris Buschelman PetscFunctionReturn(0); 570397b6df1SKris Buschelman } 571397b6df1SKris Buschelman 57251d5961aSHong Zhang #undef __FUNCT__ 57351d5961aSHong Zhang #define __FUNCT__ "MatSolveTranspose_MUMPS" 57451d5961aSHong Zhang PetscErrorCode MatSolveTranspose_MUMPS(Mat A,Vec b,Vec x) 57551d5961aSHong Zhang { 57651d5961aSHong Zhang Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 57751d5961aSHong Zhang PetscErrorCode ierr; 57851d5961aSHong Zhang 57951d5961aSHong Zhang PetscFunctionBegin; 58051d5961aSHong Zhang lu->id.ICNTL(9) = 0; 58151d5961aSHong Zhang ierr = MatSolve_MUMPS(A,b,x); 58251d5961aSHong Zhang lu->id.ICNTL(9) = 1; 58351d5961aSHong Zhang PetscFunctionReturn(0); 58451d5961aSHong Zhang } 58551d5961aSHong Zhang 586e0b74bf9SHong Zhang #undef __FUNCT__ 587e0b74bf9SHong Zhang #define __FUNCT__ "MatMatSolve_MUMPS" 588e0b74bf9SHong Zhang PetscErrorCode MatMatSolve_MUMPS(Mat A,Mat B,Mat X) 589e0b74bf9SHong Zhang { 590e0b74bf9SHong Zhang PetscFunctionBegin; 591e0b74bf9SHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatSolve_MUMPS() is not implemented yet"); 592e0b74bf9SHong Zhang PetscFunctionReturn(0); 593e0b74bf9SHong Zhang } 594e0b74bf9SHong Zhang 595ace3df97SHong Zhang #if !defined(PETSC_USE_COMPLEX) 596a58c3f20SHong Zhang /* 597a58c3f20SHong Zhang input: 598a58c3f20SHong Zhang F: numeric factor 599a58c3f20SHong Zhang output: 600a58c3f20SHong Zhang nneg: total number of negative pivots 601a58c3f20SHong Zhang nzero: 0 602a58c3f20SHong Zhang npos: (global dimension of F) - nneg 603a58c3f20SHong Zhang */ 604a58c3f20SHong Zhang 605a58c3f20SHong Zhang #undef __FUNCT__ 606a58c3f20SHong Zhang #define __FUNCT__ "MatGetInertia_SBAIJMUMPS" 607dfbe8321SBarry Smith PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos) 608a58c3f20SHong Zhang { 609a58c3f20SHong Zhang Mat_MUMPS *lu =(Mat_MUMPS*)F->spptr; 610dfbe8321SBarry Smith PetscErrorCode ierr; 611c1490034SHong Zhang PetscMPIInt size; 612a58c3f20SHong Zhang 613a58c3f20SHong Zhang PetscFunctionBegin; 6147adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)F)->comm,&size);CHKERRQ(ierr); 615bcb30aebSHong 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 */ 61665e19b50SBarry 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)); 617a58c3f20SHong Zhang if (nneg){ 618a58c3f20SHong Zhang if (!lu->myid){ 619a58c3f20SHong Zhang *nneg = lu->id.INFOG(12); 620a58c3f20SHong Zhang } 621bcb30aebSHong Zhang ierr = MPI_Bcast(nneg,1,MPI_INT,0,lu->comm_mumps);CHKERRQ(ierr); 622a58c3f20SHong Zhang } 623a58c3f20SHong Zhang if (nzero) *nzero = 0; 624d0f46423SBarry Smith if (npos) *npos = F->rmap->N - (*nneg); 625a58c3f20SHong Zhang PetscFunctionReturn(0); 626a58c3f20SHong Zhang } 627ace3df97SHong Zhang #endif /* !defined(PETSC_USE_COMPLEX) */ 628a58c3f20SHong Zhang 629397b6df1SKris Buschelman #undef __FUNCT__ 630f6c57405SHong Zhang #define __FUNCT__ "MatFactorNumeric_MUMPS" 6310481f469SBarry Smith PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info) 632af281ebdSHong Zhang { 633dcd589f8SShri Abhyankar Mat_MUMPS *lu =(Mat_MUMPS*)(F)->spptr; 6346849ba73SBarry Smith PetscErrorCode ierr; 635bccb9932SShri Abhyankar MatReuse reuse; 636e09efc27SHong Zhang Mat F_diag; 637ace3abfcSBarry Smith PetscBool isMPIAIJ; 638397b6df1SKris Buschelman 639397b6df1SKris Buschelman PetscFunctionBegin; 640bccb9932SShri Abhyankar reuse = MAT_REUSE_MATRIX; 641bccb9932SShri Abhyankar ierr = (*lu->ConvertToTriples)(A, 1, reuse, &lu->nz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr); 642397b6df1SKris Buschelman 643397b6df1SKris Buschelman /* numerical factorization phase */ 644329ec9b3SHong Zhang /*-------------------------------*/ 6453d472b54SHong Zhang lu->id.job = JOB_FACTNUMERIC; 646958c9bccSBarry Smith if(!lu->id.ICNTL(18)) { 647a7aca84bSHong Zhang if (!lu->myid) { 648397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 649397b6df1SKris Buschelman lu->id.a = (mumps_double_complex*)lu->val; 650397b6df1SKris Buschelman #else 651397b6df1SKris Buschelman lu->id.a = lu->val; 652397b6df1SKris Buschelman #endif 653397b6df1SKris Buschelman } 654397b6df1SKris Buschelman } else { 655397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 656397b6df1SKris Buschelman lu->id.a_loc = (mumps_double_complex*)lu->val; 657397b6df1SKris Buschelman #else 658397b6df1SKris Buschelman lu->id.a_loc = lu->val; 659397b6df1SKris Buschelman #endif 660397b6df1SKris Buschelman } 661397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 662397b6df1SKris Buschelman zmumps_c(&lu->id); 663397b6df1SKris Buschelman #else 664397b6df1SKris Buschelman dmumps_c(&lu->id); 665397b6df1SKris Buschelman #endif 666397b6df1SKris Buschelman if (lu->id.INFOG(1) < 0) { 66765e19b50SBarry 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)); 66865e19b50SBarry 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)); 669397b6df1SKris Buschelman } 67065e19b50SBarry 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)); 671397b6df1SKris Buschelman 6728ada1bb4SHong Zhang if (lu->size > 1){ 67316ebf90aSShri Abhyankar ierr = PetscTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr); 674a214ac2aSShri Abhyankar if(isMPIAIJ) { 675dcd589f8SShri Abhyankar F_diag = ((Mat_MPIAIJ *)(F)->data)->A; 676e09efc27SHong Zhang } else { 677dcd589f8SShri Abhyankar F_diag = ((Mat_MPISBAIJ *)(F)->data)->A; 678e09efc27SHong Zhang } 679e09efc27SHong Zhang F_diag->assembled = PETSC_TRUE; 680329ec9b3SHong Zhang if (lu->nSolve){ 6816bf464f9SBarry Smith ierr = VecScatterDestroy(&lu->scat_sol);CHKERRQ(ierr); 6820e83c824SBarry Smith ierr = PetscFree2(lu->id.sol_loc,lu->id.isol_loc);CHKERRQ(ierr); 6836bf464f9SBarry Smith ierr = VecDestroy(&lu->x_seq);CHKERRQ(ierr); 684329ec9b3SHong Zhang } 6858ada1bb4SHong Zhang } 686dcd589f8SShri Abhyankar (F)->assembled = PETSC_TRUE; 687397b6df1SKris Buschelman lu->matstruc = SAME_NONZERO_PATTERN; 688ace87b0dSHong Zhang lu->CleanUpMUMPS = PETSC_TRUE; 689329ec9b3SHong Zhang lu->nSolve = 0; 69067877ebaSShri Abhyankar 69167877ebaSShri Abhyankar if (lu->size > 1){ 69267877ebaSShri Abhyankar /* distributed solution */ 69367877ebaSShri Abhyankar lu->id.ICNTL(21) = 1; 69467877ebaSShri Abhyankar if (!lu->nSolve){ 69567877ebaSShri Abhyankar /* Create x_seq=sol_loc for repeated use */ 69667877ebaSShri Abhyankar PetscInt lsol_loc; 69767877ebaSShri Abhyankar PetscScalar *sol_loc; 69867877ebaSShri Abhyankar lsol_loc = lu->id.INFO(23); /* length of sol_loc */ 69967877ebaSShri Abhyankar ierr = PetscMalloc2(lsol_loc,PetscScalar,&sol_loc,lsol_loc,PetscInt,&lu->id.isol_loc);CHKERRQ(ierr); 70067877ebaSShri Abhyankar lu->id.lsol_loc = lsol_loc; 70167877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 70267877ebaSShri Abhyankar lu->id.sol_loc = (mumps_double_complex*)sol_loc; 70367877ebaSShri Abhyankar #else 70467877ebaSShri Abhyankar lu->id.sol_loc = sol_loc; 70567877ebaSShri Abhyankar #endif 70667877ebaSShri Abhyankar ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,lsol_loc,sol_loc,&lu->x_seq);CHKERRQ(ierr); 70767877ebaSShri Abhyankar } 70867877ebaSShri Abhyankar } 709397b6df1SKris Buschelman PetscFunctionReturn(0); 710397b6df1SKris Buschelman } 711397b6df1SKris Buschelman 712dcd589f8SShri Abhyankar #undef __FUNCT__ 713dcd589f8SShri Abhyankar #define __FUNCT__ "PetscSetMUMPSOptions" 714dcd589f8SShri Abhyankar PetscErrorCode PetscSetMUMPSOptions(Mat F, Mat A) 715dcd589f8SShri Abhyankar { 716dcd589f8SShri Abhyankar Mat_MUMPS *lu = (Mat_MUMPS*)F->spptr; 717dcd589f8SShri Abhyankar PetscErrorCode ierr; 718dcd589f8SShri Abhyankar PetscInt icntl; 719ace3abfcSBarry Smith PetscBool flg; 720dcd589f8SShri Abhyankar 721dcd589f8SShri Abhyankar PetscFunctionBegin; 722dcd589f8SShri Abhyankar ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"MUMPS Options","Mat");CHKERRQ(ierr); 723dcd589f8SShri Abhyankar if (lu->size == 1){ 724dcd589f8SShri Abhyankar lu->id.ICNTL(18) = 0; /* centralized assembled matrix input */ 725dcd589f8SShri Abhyankar } else { 726dcd589f8SShri Abhyankar lu->id.ICNTL(18) = 3; /* distributed assembled matrix input */ 727dcd589f8SShri Abhyankar } 728dcd589f8SShri Abhyankar 729dcd589f8SShri Abhyankar icntl=-1; 730dcd589f8SShri Abhyankar lu->id.ICNTL(4) = 0; /* level of printing; overwrite mumps default ICNTL(4)=2 */ 731dcd589f8SShri Abhyankar ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",lu->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr); 732dcd589f8SShri Abhyankar if ((flg && icntl > 0) || PetscLogPrintInfo) { 733dcd589f8SShri Abhyankar lu->id.ICNTL(4)=icntl; /* and use mumps default icntl(i), i=1,2,3 */ 734dcd589f8SShri Abhyankar } else { /* no output */ 735dcd589f8SShri Abhyankar lu->id.ICNTL(1) = 0; /* error message, default= 6 */ 736dcd589f8SShri Abhyankar lu->id.ICNTL(2) = 0; /* output stream for diagnostic printing, statistics, and warning. default=0 */ 737dcd589f8SShri Abhyankar lu->id.ICNTL(3) = 0; /* output stream for global information, default=6 */ 738dcd589f8SShri Abhyankar } 739dcd589f8SShri Abhyankar ierr = PetscOptionsInt("-mat_mumps_icntl_6","ICNTL(6): column permutation and/or scaling to get a zero-free diagonal (0 to 7)","None",lu->id.ICNTL(6),&lu->id.ICNTL(6),PETSC_NULL);CHKERRQ(ierr); 740dcd589f8SShri Abhyankar icntl=-1; 741292fb18eSBarry Smith ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): sequential matrix ordering (0 to 7) 3 = Scotch, 5 = Metis","None",lu->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr); 742dcd589f8SShri Abhyankar if (flg) { 743e0b74bf9SHong Zhang if (icntl== 1 && lu->size > 1){ 744e32f2f54SBarry 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"); 745dcd589f8SShri Abhyankar } else { 746dcd589f8SShri Abhyankar lu->id.ICNTL(7) = icntl; 747dcd589f8SShri Abhyankar } 748dcd589f8SShri Abhyankar } 749e0b74bf9SHong Zhang 750dcd589f8SShri Abhyankar ierr = PetscOptionsInt("-mat_mumps_icntl_8","ICNTL(8): scaling strategy (-2 to 7 or 77)","None",lu->id.ICNTL(8),&lu->id.ICNTL(8),PETSC_NULL);CHKERRQ(ierr); 751dcd589f8SShri Abhyankar ierr = PetscOptionsInt("-mat_mumps_icntl_10","ICNTL(10): max num of refinements","None",lu->id.ICNTL(10),&lu->id.ICNTL(10),PETSC_NULL);CHKERRQ(ierr); 752dcd589f8SShri Abhyankar ierr = PetscOptionsInt("-mat_mumps_icntl_11","ICNTL(11): statistics related to the linear system solved (via -ksp_view)","None",lu->id.ICNTL(11),&lu->id.ICNTL(11),PETSC_NULL);CHKERRQ(ierr); 753dcd589f8SShri Abhyankar ierr = PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): efficiency control: defines the ordering strategy with scaling constraints (0 to 3","None",lu->id.ICNTL(12),&lu->id.ICNTL(12),PETSC_NULL);CHKERRQ(ierr); 754dcd589f8SShri Abhyankar ierr = PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): efficiency control: with or without ScaLAPACK","None",lu->id.ICNTL(13),&lu->id.ICNTL(13),PETSC_NULL);CHKERRQ(ierr); 755dcd589f8SShri Abhyankar ierr = PetscOptionsInt("-mat_mumps_icntl_14","ICNTL(14): percentage of estimated workspace increase","None",lu->id.ICNTL(14),&lu->id.ICNTL(14),PETSC_NULL);CHKERRQ(ierr); 756dcd589f8SShri Abhyankar ierr = PetscOptionsInt("-mat_mumps_icntl_19","ICNTL(19): Schur complement","None",lu->id.ICNTL(19),&lu->id.ICNTL(19),PETSC_NULL);CHKERRQ(ierr); 757dcd589f8SShri Abhyankar 758dcd589f8SShri Abhyankar ierr = PetscOptionsInt("-mat_mumps_icntl_22","ICNTL(22): in-core/out-of-core facility (0 or 1)","None",lu->id.ICNTL(22),&lu->id.ICNTL(22),PETSC_NULL);CHKERRQ(ierr); 759dcd589f8SShri Abhyankar ierr = PetscOptionsInt("-mat_mumps_icntl_23","ICNTL(23): max size of the working memory (MB) that can allocate per processor","None",lu->id.ICNTL(23),&lu->id.ICNTL(23),PETSC_NULL);CHKERRQ(ierr); 760dcd589f8SShri Abhyankar ierr = PetscOptionsInt("-mat_mumps_icntl_24","ICNTL(24): detection of null pivot rows (0 or 1)","None",lu->id.ICNTL(24),&lu->id.ICNTL(24),PETSC_NULL);CHKERRQ(ierr); 761d7ebd59bSHong Zhang if (lu->id.ICNTL(24)){ 762d7ebd59bSHong Zhang lu->id.ICNTL(13) = 1; /* turn-off ScaLAPACK to help with the correct detection of null pivots */ 763d7ebd59bSHong Zhang } 764d7ebd59bSHong Zhang 765dcd589f8SShri Abhyankar ierr = PetscOptionsInt("-mat_mumps_icntl_25","ICNTL(25): computation of a null space basis","None",lu->id.ICNTL(25),&lu->id.ICNTL(25),PETSC_NULL);CHKERRQ(ierr); 766dcd589f8SShri Abhyankar ierr = PetscOptionsInt("-mat_mumps_icntl_26","ICNTL(26): Schur options for right-hand side or solution vector","None",lu->id.ICNTL(26),&lu->id.ICNTL(26),PETSC_NULL);CHKERRQ(ierr); 767dcd589f8SShri Abhyankar ierr = PetscOptionsInt("-mat_mumps_icntl_27","ICNTL(27): experimental parameter","None",lu->id.ICNTL(27),&lu->id.ICNTL(27),PETSC_NULL);CHKERRQ(ierr); 768d7ebd59bSHong 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",lu->id.ICNTL(28),&lu->id.ICNTL(28),PETSC_NULL);CHKERRQ(ierr); 769292fb18eSBarry Smith ierr = PetscOptionsInt("-mat_mumps_icntl_29","ICNTL(29): parallel ordering 1 = ptscotch 2 = parmetis","None",lu->id.ICNTL(29),&lu->id.ICNTL(29),PETSC_NULL);CHKERRQ(ierr); 770dcd589f8SShri Abhyankar 771dcd589f8SShri Abhyankar ierr = PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",lu->id.CNTL(1),&lu->id.CNTL(1),PETSC_NULL);CHKERRQ(ierr); 772dcd589f8SShri Abhyankar ierr = PetscOptionsReal("-mat_mumps_cntl_2","CNTL(2): stopping criterion of refinement","None",lu->id.CNTL(2),&lu->id.CNTL(2),PETSC_NULL);CHKERRQ(ierr); 773dcd589f8SShri Abhyankar ierr = PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",lu->id.CNTL(3),&lu->id.CNTL(3),PETSC_NULL);CHKERRQ(ierr); 774dcd589f8SShri Abhyankar ierr = PetscOptionsReal("-mat_mumps_cntl_4","CNTL(4): value for static pivoting","None",lu->id.CNTL(4),&lu->id.CNTL(4),PETSC_NULL);CHKERRQ(ierr); 775dcd589f8SShri Abhyankar ierr = PetscOptionsReal("-mat_mumps_cntl_5","CNTL(5): fixation for null pivots","None",lu->id.CNTL(5),&lu->id.CNTL(5),PETSC_NULL);CHKERRQ(ierr); 776dcd589f8SShri Abhyankar PetscOptionsEnd(); 777dcd589f8SShri Abhyankar PetscFunctionReturn(0); 778dcd589f8SShri Abhyankar } 779dcd589f8SShri Abhyankar 780dcd589f8SShri Abhyankar #undef __FUNCT__ 781dcd589f8SShri Abhyankar #define __FUNCT__ "PetscInitializeMUMPS" 782f697e70eSHong Zhang PetscErrorCode PetscInitializeMUMPS(Mat A,Mat_MUMPS* mumps) 783dcd589f8SShri Abhyankar { 784dcd589f8SShri Abhyankar PetscErrorCode ierr; 785dcd589f8SShri Abhyankar 786dcd589f8SShri Abhyankar PetscFunctionBegin; 787f697e70eSHong Zhang ierr = MPI_Comm_rank(((PetscObject)A)->comm, &mumps->myid); 788f697e70eSHong Zhang ierr = MPI_Comm_size(((PetscObject)A)->comm,&mumps->size);CHKERRQ(ierr); 789f697e70eSHong Zhang ierr = MPI_Comm_dup(((PetscObject)A)->comm,&(mumps->comm_mumps));CHKERRQ(ierr); 790f697e70eSHong Zhang mumps->id.comm_fortran = MPI_Comm_c2f(mumps->comm_mumps); 791f697e70eSHong Zhang 792f697e70eSHong Zhang mumps->id.job = JOB_INIT; 793f697e70eSHong Zhang mumps->id.par = 1; /* host participates factorizaton and solve */ 794f697e70eSHong Zhang mumps->id.sym = mumps->sym; 795dcd589f8SShri Abhyankar #if defined(PETSC_USE_COMPLEX) 796f697e70eSHong Zhang zmumps_c(&mumps->id); 797dcd589f8SShri Abhyankar #else 798f697e70eSHong Zhang dmumps_c(&mumps->id); 799dcd589f8SShri Abhyankar #endif 800f697e70eSHong Zhang 801f697e70eSHong Zhang mumps->CleanUpMUMPS = PETSC_FALSE; 802f697e70eSHong Zhang mumps->scat_rhs = PETSC_NULL; 803f697e70eSHong Zhang mumps->scat_sol = PETSC_NULL; 804f697e70eSHong Zhang mumps->nSolve = 0; 805dcd589f8SShri Abhyankar PetscFunctionReturn(0); 806dcd589f8SShri Abhyankar } 807dcd589f8SShri Abhyankar 808397b6df1SKris Buschelman /* Note the Petsc r and c permutations are ignored */ 809397b6df1SKris Buschelman #undef __FUNCT__ 810f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS" 8110481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 812b24902e0SBarry Smith { 813719d5645SBarry Smith Mat_MUMPS *lu = (Mat_MUMPS*)F->spptr; 814dcd589f8SShri Abhyankar PetscErrorCode ierr; 815bccb9932SShri Abhyankar MatReuse reuse; 81667877ebaSShri Abhyankar Vec b; 81767877ebaSShri Abhyankar IS is_iden; 81867877ebaSShri Abhyankar const PetscInt M = A->rmap->N; 819397b6df1SKris Buschelman 820397b6df1SKris Buschelman PetscFunctionBegin; 821b24902e0SBarry Smith lu->matstruc = DIFFERENT_NONZERO_PATTERN; 822dcd589f8SShri Abhyankar 823dcd589f8SShri Abhyankar /* Set MUMPS options */ 824dcd589f8SShri Abhyankar ierr = PetscSetMUMPSOptions(F,A);CHKERRQ(ierr); 825dcd589f8SShri Abhyankar 826bccb9932SShri Abhyankar reuse = MAT_INITIAL_MATRIX; 827bccb9932SShri Abhyankar ierr = (*lu->ConvertToTriples)(A, 1, reuse, &lu->nz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr); 828dcd589f8SShri Abhyankar 82967877ebaSShri Abhyankar /* analysis phase */ 83067877ebaSShri Abhyankar /*----------------*/ 8313d472b54SHong Zhang lu->id.job = JOB_FACTSYMBOLIC; 83267877ebaSShri Abhyankar lu->id.n = M; 83367877ebaSShri Abhyankar switch (lu->id.ICNTL(18)){ 83467877ebaSShri Abhyankar case 0: /* centralized assembled matrix input */ 83567877ebaSShri Abhyankar if (!lu->myid) { 83667877ebaSShri Abhyankar lu->id.nz =lu->nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn; 83767877ebaSShri Abhyankar if (lu->id.ICNTL(6)>1){ 83867877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 83967877ebaSShri Abhyankar lu->id.a = (mumps_double_complex*)lu->val; 84067877ebaSShri Abhyankar #else 84167877ebaSShri Abhyankar lu->id.a = lu->val; 84267877ebaSShri Abhyankar #endif 84367877ebaSShri Abhyankar } 844e0b74bf9SHong Zhang if (lu->id.ICNTL(7) == 1){ /* use user-provide matrix ordering */ 845e0b74bf9SHong Zhang if (!lu->myid) { 846e0b74bf9SHong Zhang const PetscInt *idx; 847e0b74bf9SHong Zhang PetscInt i,*perm_in; 848e0b74bf9SHong Zhang ierr = PetscMalloc(M*sizeof(PetscInt),&perm_in);CHKERRQ(ierr); 849e0b74bf9SHong Zhang ierr = ISGetIndices(r,&idx);CHKERRQ(ierr); 850e0b74bf9SHong Zhang lu->id.perm_in = perm_in; 851e0b74bf9SHong Zhang for (i=0; i<M; i++) perm_in[i] = idx[i]+1; /* perm_in[]: start from 1, not 0! */ 852e0b74bf9SHong Zhang ierr = ISRestoreIndices(r,&idx);CHKERRQ(ierr); 853e0b74bf9SHong Zhang } 854e0b74bf9SHong Zhang } 85567877ebaSShri Abhyankar } 85667877ebaSShri Abhyankar break; 85767877ebaSShri Abhyankar case 3: /* distributed assembled matrix input (size>1) */ 85867877ebaSShri Abhyankar lu->id.nz_loc = lu->nz; 85967877ebaSShri Abhyankar lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn; 86067877ebaSShri Abhyankar if (lu->id.ICNTL(6)>1) { 86167877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 86267877ebaSShri Abhyankar lu->id.a_loc = (mumps_double_complex*)lu->val; 86367877ebaSShri Abhyankar #else 86467877ebaSShri Abhyankar lu->id.a_loc = lu->val; 86567877ebaSShri Abhyankar #endif 86667877ebaSShri Abhyankar } 86767877ebaSShri Abhyankar /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */ 86867877ebaSShri Abhyankar if (!lu->myid){ 86967877ebaSShri Abhyankar ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&lu->b_seq);CHKERRQ(ierr); 87067877ebaSShri Abhyankar ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr); 87167877ebaSShri Abhyankar } else { 87267877ebaSShri Abhyankar ierr = VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);CHKERRQ(ierr); 87367877ebaSShri Abhyankar ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr); 87467877ebaSShri Abhyankar } 87567877ebaSShri Abhyankar ierr = VecCreate(((PetscObject)A)->comm,&b);CHKERRQ(ierr); 87667877ebaSShri Abhyankar ierr = VecSetSizes(b,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr); 87767877ebaSShri Abhyankar ierr = VecSetFromOptions(b);CHKERRQ(ierr); 87867877ebaSShri Abhyankar 87967877ebaSShri Abhyankar ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr); 8806bf464f9SBarry Smith ierr = ISDestroy(&is_iden);CHKERRQ(ierr); 8816bf464f9SBarry Smith ierr = VecDestroy(&b);CHKERRQ(ierr); 88267877ebaSShri Abhyankar break; 88367877ebaSShri Abhyankar } 88467877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 88567877ebaSShri Abhyankar zmumps_c(&lu->id); 88667877ebaSShri Abhyankar #else 88767877ebaSShri Abhyankar dmumps_c(&lu->id); 88867877ebaSShri Abhyankar #endif 88967877ebaSShri 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)); 89067877ebaSShri Abhyankar 891719d5645SBarry Smith F->ops->lufactornumeric = MatFactorNumeric_MUMPS; 892dcd589f8SShri Abhyankar F->ops->solve = MatSolve_MUMPS; 89351d5961aSHong Zhang F->ops->solvetranspose = MatSolveTranspose_MUMPS; 894e0b74bf9SHong Zhang F->ops->matsolve = MatMatSolve_MUMPS; 895b24902e0SBarry Smith PetscFunctionReturn(0); 896b24902e0SBarry Smith } 897b24902e0SBarry Smith 898450b117fSShri Abhyankar /* Note the Petsc r and c permutations are ignored */ 899450b117fSShri Abhyankar #undef __FUNCT__ 900450b117fSShri Abhyankar #define __FUNCT__ "MatLUFactorSymbolic_BAIJMUMPS" 901450b117fSShri Abhyankar PetscErrorCode MatLUFactorSymbolic_BAIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 902450b117fSShri Abhyankar { 903dcd589f8SShri Abhyankar 904450b117fSShri Abhyankar Mat_MUMPS *lu = (Mat_MUMPS*)F->spptr; 905dcd589f8SShri Abhyankar PetscErrorCode ierr; 906bccb9932SShri Abhyankar MatReuse reuse; 90767877ebaSShri Abhyankar Vec b; 90867877ebaSShri Abhyankar IS is_iden; 90967877ebaSShri Abhyankar const PetscInt M = A->rmap->N; 910450b117fSShri Abhyankar 911450b117fSShri Abhyankar PetscFunctionBegin; 912450b117fSShri Abhyankar lu->matstruc = DIFFERENT_NONZERO_PATTERN; 913dcd589f8SShri Abhyankar 914dcd589f8SShri Abhyankar /* Set MUMPS options */ 915dcd589f8SShri Abhyankar ierr = PetscSetMUMPSOptions(F,A);CHKERRQ(ierr); 916dcd589f8SShri Abhyankar 917bccb9932SShri Abhyankar reuse = MAT_INITIAL_MATRIX; 918bccb9932SShri Abhyankar ierr = (*lu->ConvertToTriples)(A, 1, reuse, &lu->nz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr); 91967877ebaSShri Abhyankar 92067877ebaSShri Abhyankar /* analysis phase */ 92167877ebaSShri Abhyankar /*----------------*/ 9223d472b54SHong Zhang lu->id.job = JOB_FACTSYMBOLIC; 92367877ebaSShri Abhyankar lu->id.n = M; 92467877ebaSShri Abhyankar switch (lu->id.ICNTL(18)){ 92567877ebaSShri Abhyankar case 0: /* centralized assembled matrix input */ 92667877ebaSShri Abhyankar if (!lu->myid) { 92767877ebaSShri Abhyankar lu->id.nz =lu->nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn; 92867877ebaSShri Abhyankar if (lu->id.ICNTL(6)>1){ 92967877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 93067877ebaSShri Abhyankar lu->id.a = (mumps_double_complex*)lu->val; 93167877ebaSShri Abhyankar #else 93267877ebaSShri Abhyankar lu->id.a = lu->val; 93367877ebaSShri Abhyankar #endif 93467877ebaSShri Abhyankar } 93567877ebaSShri Abhyankar } 93667877ebaSShri Abhyankar break; 93767877ebaSShri Abhyankar case 3: /* distributed assembled matrix input (size>1) */ 93867877ebaSShri Abhyankar lu->id.nz_loc = lu->nz; 93967877ebaSShri Abhyankar lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn; 94067877ebaSShri Abhyankar if (lu->id.ICNTL(6)>1) { 94167877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 94267877ebaSShri Abhyankar lu->id.a_loc = (mumps_double_complex*)lu->val; 94367877ebaSShri Abhyankar #else 94467877ebaSShri Abhyankar lu->id.a_loc = lu->val; 94567877ebaSShri Abhyankar #endif 94667877ebaSShri Abhyankar } 94767877ebaSShri Abhyankar /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */ 94867877ebaSShri Abhyankar if (!lu->myid){ 94967877ebaSShri Abhyankar ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&lu->b_seq);CHKERRQ(ierr); 95067877ebaSShri Abhyankar ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr); 95167877ebaSShri Abhyankar } else { 95267877ebaSShri Abhyankar ierr = VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);CHKERRQ(ierr); 95367877ebaSShri Abhyankar ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr); 95467877ebaSShri Abhyankar } 95567877ebaSShri Abhyankar ierr = VecCreate(((PetscObject)A)->comm,&b);CHKERRQ(ierr); 95667877ebaSShri Abhyankar ierr = VecSetSizes(b,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr); 95767877ebaSShri Abhyankar ierr = VecSetFromOptions(b);CHKERRQ(ierr); 95867877ebaSShri Abhyankar 95967877ebaSShri Abhyankar ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr); 9606bf464f9SBarry Smith ierr = ISDestroy(&is_iden);CHKERRQ(ierr); 9616bf464f9SBarry Smith ierr = VecDestroy(&b);CHKERRQ(ierr); 96267877ebaSShri Abhyankar break; 96367877ebaSShri Abhyankar } 96467877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 96567877ebaSShri Abhyankar zmumps_c(&lu->id); 96667877ebaSShri Abhyankar #else 96767877ebaSShri Abhyankar dmumps_c(&lu->id); 96867877ebaSShri Abhyankar #endif 96967877ebaSShri 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)); 97067877ebaSShri Abhyankar 971450b117fSShri Abhyankar F->ops->lufactornumeric = MatFactorNumeric_MUMPS; 972dcd589f8SShri Abhyankar F->ops->solve = MatSolve_MUMPS; 97351d5961aSHong Zhang F->ops->solvetranspose = MatSolveTranspose_MUMPS; 974450b117fSShri Abhyankar PetscFunctionReturn(0); 975450b117fSShri Abhyankar } 976b24902e0SBarry Smith 977141f4205SHong Zhang /* Note the Petsc r permutation and factor info are ignored */ 978397b6df1SKris Buschelman #undef __FUNCT__ 97967877ebaSShri Abhyankar #define __FUNCT__ "MatCholeskyFactorSymbolic_MUMPS" 98067877ebaSShri Abhyankar PetscErrorCode MatCholeskyFactorSymbolic_MUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info) 981b24902e0SBarry Smith { 9822792810eSHong Zhang Mat_MUMPS *lu = (Mat_MUMPS*)F->spptr; 983dcd589f8SShri Abhyankar PetscErrorCode ierr; 984bccb9932SShri Abhyankar MatReuse reuse; 98567877ebaSShri Abhyankar Vec b; 98667877ebaSShri Abhyankar IS is_iden; 98767877ebaSShri Abhyankar const PetscInt M = A->rmap->N; 988397b6df1SKris Buschelman 989397b6df1SKris Buschelman PetscFunctionBegin; 990b24902e0SBarry Smith lu->matstruc = DIFFERENT_NONZERO_PATTERN; 991dcd589f8SShri Abhyankar 992dcd589f8SShri Abhyankar /* Set MUMPS options */ 993dcd589f8SShri Abhyankar ierr = PetscSetMUMPSOptions(F,A);CHKERRQ(ierr); 994dcd589f8SShri Abhyankar 995bccb9932SShri Abhyankar reuse = MAT_INITIAL_MATRIX; 996bccb9932SShri Abhyankar ierr = (*lu->ConvertToTriples)(A, 1 , reuse, &lu->nz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr); 997dcd589f8SShri Abhyankar 99867877ebaSShri Abhyankar /* analysis phase */ 99967877ebaSShri Abhyankar /*----------------*/ 10003d472b54SHong Zhang lu->id.job = JOB_FACTSYMBOLIC; 100167877ebaSShri Abhyankar lu->id.n = M; 100267877ebaSShri Abhyankar switch (lu->id.ICNTL(18)){ 100367877ebaSShri Abhyankar case 0: /* centralized assembled matrix input */ 100467877ebaSShri Abhyankar if (!lu->myid) { 100567877ebaSShri Abhyankar lu->id.nz =lu->nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn; 100667877ebaSShri Abhyankar if (lu->id.ICNTL(6)>1){ 100767877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 100867877ebaSShri Abhyankar lu->id.a = (mumps_double_complex*)lu->val; 100967877ebaSShri Abhyankar #else 101067877ebaSShri Abhyankar lu->id.a = lu->val; 101167877ebaSShri Abhyankar #endif 101267877ebaSShri Abhyankar } 101367877ebaSShri Abhyankar } 101467877ebaSShri Abhyankar break; 101567877ebaSShri Abhyankar case 3: /* distributed assembled matrix input (size>1) */ 101667877ebaSShri Abhyankar lu->id.nz_loc = lu->nz; 101767877ebaSShri Abhyankar lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn; 101867877ebaSShri Abhyankar if (lu->id.ICNTL(6)>1) { 101967877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 102067877ebaSShri Abhyankar lu->id.a_loc = (mumps_double_complex*)lu->val; 102167877ebaSShri Abhyankar #else 102267877ebaSShri Abhyankar lu->id.a_loc = lu->val; 102367877ebaSShri Abhyankar #endif 102467877ebaSShri Abhyankar } 102567877ebaSShri Abhyankar /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */ 102667877ebaSShri Abhyankar if (!lu->myid){ 102767877ebaSShri Abhyankar ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&lu->b_seq);CHKERRQ(ierr); 102867877ebaSShri Abhyankar ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr); 102967877ebaSShri Abhyankar } else { 103067877ebaSShri Abhyankar ierr = VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);CHKERRQ(ierr); 103167877ebaSShri Abhyankar ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr); 103267877ebaSShri Abhyankar } 103367877ebaSShri Abhyankar ierr = VecCreate(((PetscObject)A)->comm,&b);CHKERRQ(ierr); 103467877ebaSShri Abhyankar ierr = VecSetSizes(b,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr); 103567877ebaSShri Abhyankar ierr = VecSetFromOptions(b);CHKERRQ(ierr); 103667877ebaSShri Abhyankar 103767877ebaSShri Abhyankar ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr); 10386bf464f9SBarry Smith ierr = ISDestroy(&is_iden);CHKERRQ(ierr); 10396bf464f9SBarry Smith ierr = VecDestroy(&b);CHKERRQ(ierr); 104067877ebaSShri Abhyankar break; 104167877ebaSShri Abhyankar } 104267877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 104367877ebaSShri Abhyankar zmumps_c(&lu->id); 104467877ebaSShri Abhyankar #else 104567877ebaSShri Abhyankar dmumps_c(&lu->id); 104667877ebaSShri Abhyankar #endif 104767877ebaSShri 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)); 104867877ebaSShri Abhyankar 10492792810eSHong Zhang F->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS; 1050dcd589f8SShri Abhyankar F->ops->solve = MatSolve_MUMPS; 105151d5961aSHong Zhang F->ops->solvetranspose = MatSolve_MUMPS; 1052db4efbfdSBarry Smith #if !defined(PETSC_USE_COMPLEX) 1053dcd589f8SShri Abhyankar (F)->ops->getinertia = MatGetInertia_SBAIJMUMPS; 1054db4efbfdSBarry Smith #endif 1055b24902e0SBarry Smith PetscFunctionReturn(0); 1056b24902e0SBarry Smith } 1057b24902e0SBarry Smith 1058397b6df1SKris Buschelman #undef __FUNCT__ 105964e6c443SBarry Smith #define __FUNCT__ "MatView_MUMPS" 106064e6c443SBarry Smith PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer) 106174ed9c26SBarry Smith { 1062f6c57405SHong Zhang PetscErrorCode ierr; 106364e6c443SBarry Smith PetscBool iascii; 106464e6c443SBarry Smith PetscViewerFormat format; 106564e6c443SBarry Smith Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 1066f6c57405SHong Zhang 1067f6c57405SHong Zhang PetscFunctionBegin; 106864e6c443SBarry Smith /* check if matrix is mumps type */ 106964e6c443SBarry Smith if (A->ops->solve != MatSolve_MUMPS) PetscFunctionReturn(0); 107064e6c443SBarry Smith 107164e6c443SBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 107264e6c443SBarry Smith if (iascii) { 107364e6c443SBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 107464e6c443SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO){ 107564e6c443SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr); 107664e6c443SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," SYM (matrix type): %d \n",lu->id.sym);CHKERRQ(ierr); 107764e6c443SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," PAR (host participation): %d \n",lu->id.par);CHKERRQ(ierr); 1078f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(1) (output for error): %d \n",lu->id.ICNTL(1));CHKERRQ(ierr); 1079f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(2) (output of diagnostic msg): %d \n",lu->id.ICNTL(2));CHKERRQ(ierr); 1080f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(3) (output for global info): %d \n",lu->id.ICNTL(3));CHKERRQ(ierr); 1081f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(4) (level of printing): %d \n",lu->id.ICNTL(4));CHKERRQ(ierr); 1082f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(5) (input mat struct): %d \n",lu->id.ICNTL(5));CHKERRQ(ierr); 1083f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(6) (matrix prescaling): %d \n",lu->id.ICNTL(6));CHKERRQ(ierr); 1084d06efebcSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ICNTL(7) (sequentia matrix ordering):%d \n",lu->id.ICNTL(7));CHKERRQ(ierr); 1085f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(8) (scalling strategy): %d \n",lu->id.ICNTL(8));CHKERRQ(ierr); 1086f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(10) (max num of refinements): %d \n",lu->id.ICNTL(10));CHKERRQ(ierr); 1087f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(11) (error analysis): %d \n",lu->id.ICNTL(11));CHKERRQ(ierr); 108834ed7027SBarry Smith if (lu->id.ICNTL(11)>0) { 108934ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," RINFOG(4) (inf norm of input mat): %g\n",lu->id.RINFOG(4));CHKERRQ(ierr); 109034ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," RINFOG(5) (inf norm of solution): %g\n",lu->id.RINFOG(5));CHKERRQ(ierr); 109134ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," RINFOG(6) (inf norm of residual): %g\n",lu->id.RINFOG(6));CHKERRQ(ierr); 109234ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",lu->id.RINFOG(7),lu->id.RINFOG(8));CHKERRQ(ierr); 109334ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," RINFOG(9) (error estimate): %g \n",lu->id.RINFOG(9));CHKERRQ(ierr); 109434ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",lu->id.RINFOG(10),lu->id.RINFOG(11));CHKERRQ(ierr); 1095f6c57405SHong Zhang } 1096f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(12) (efficiency control): %d \n",lu->id.ICNTL(12));CHKERRQ(ierr); 1097f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(13) (efficiency control): %d \n",lu->id.ICNTL(13));CHKERRQ(ierr); 1098f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(14) (percentage of estimated workspace increase): %d \n",lu->id.ICNTL(14));CHKERRQ(ierr); 1099f6c57405SHong Zhang /* ICNTL(15-17) not used */ 1100f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(18) (input mat struct): %d \n",lu->id.ICNTL(18));CHKERRQ(ierr); 1101f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(19) (Shur complement info): %d \n",lu->id.ICNTL(19));CHKERRQ(ierr); 1102f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(20) (rhs sparse pattern): %d \n",lu->id.ICNTL(20));CHKERRQ(ierr); 1103f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(21) (solution struct): %d \n",lu->id.ICNTL(21));CHKERRQ(ierr); 1104c0165424SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(22) (in-core/out-of-core facility): %d \n",lu->id.ICNTL(22));CHKERRQ(ierr); 1105c0165424SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(23) (max size of memory can be allocated locally):%d \n",lu->id.ICNTL(23));CHKERRQ(ierr); 1106c0165424SHong Zhang 1107c0165424SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(24) (detection of null pivot rows): %d \n",lu->id.ICNTL(24));CHKERRQ(ierr); 1108c0165424SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(25) (computation of a null space basis): %d \n",lu->id.ICNTL(25));CHKERRQ(ierr); 1109c0165424SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(26) (Schur options for rhs or solution): %d \n",lu->id.ICNTL(26));CHKERRQ(ierr); 1110c0165424SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(27) (experimental parameter): %d \n",lu->id.ICNTL(27));CHKERRQ(ierr); 1111d06efebcSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ICNTL(28) (Use parallel or sequential ordering): %d \n",lu->id.ICNTL(28));CHKERRQ(ierr); 1112d06efebcSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ICNTL(29) (Parallel ordering): %d \n",lu->id.ICNTL(29));CHKERRQ(ierr); 1113f6c57405SHong Zhang 1114f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(1) (relative pivoting threshold): %g \n",lu->id.CNTL(1));CHKERRQ(ierr); 1115f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(2) (stopping criterion of refinement): %g \n",lu->id.CNTL(2));CHKERRQ(ierr); 1116f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(3) (absolute pivoting threshold): %g \n",lu->id.CNTL(3));CHKERRQ(ierr); 1117f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(4) (value of static pivoting): %g \n",lu->id.CNTL(4));CHKERRQ(ierr); 1118c0165424SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(5) (fixation for null pivots): %g \n",lu->id.CNTL(5));CHKERRQ(ierr); 1119f6c57405SHong Zhang 1120f6c57405SHong Zhang /* infomation local to each processor */ 112134ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer, " RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr); 11227b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 112334ed7027SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %g \n",lu->myid,lu->id.RINFO(1));CHKERRQ(ierr); 112434ed7027SBarry Smith ierr = PetscViewerFlush(viewer); 112534ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer, " RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr); 112634ed7027SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %g \n",lu->myid,lu->id.RINFO(2));CHKERRQ(ierr); 112734ed7027SBarry Smith ierr = PetscViewerFlush(viewer); 112834ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer, " RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr); 112934ed7027SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %g \n",lu->myid,lu->id.RINFO(3));CHKERRQ(ierr); 113034ed7027SBarry Smith ierr = PetscViewerFlush(viewer); 1131f6c57405SHong Zhang 113234ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer, " INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization): \n");CHKERRQ(ierr); 113334ed7027SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %d \n",lu->myid,lu->id.INFO(15));CHKERRQ(ierr); 113434ed7027SBarry Smith ierr = PetscViewerFlush(viewer); 1135f6c57405SHong Zhang 113634ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer, " INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization): \n");CHKERRQ(ierr); 113734ed7027SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %d \n",lu->myid,lu->id.INFO(16));CHKERRQ(ierr); 113834ed7027SBarry Smith ierr = PetscViewerFlush(viewer); 1139f6c57405SHong Zhang 114034ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer, " INFO(23) (num of pivots eliminated on this processor after factorization): \n");CHKERRQ(ierr); 114134ed7027SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %d \n",lu->myid,lu->id.INFO(23));CHKERRQ(ierr); 114234ed7027SBarry Smith ierr = PetscViewerFlush(viewer); 11437b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr); 1144f6c57405SHong Zhang 1145f6c57405SHong Zhang if (!lu->myid){ /* information from the host */ 1146f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",lu->id.RINFOG(1));CHKERRQ(ierr); 1147f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",lu->id.RINFOG(2));CHKERRQ(ierr); 1148f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",lu->id.RINFOG(3));CHKERRQ(ierr); 1149f6c57405SHong Zhang 1150f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(3));CHKERRQ(ierr); 1151f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(4));CHKERRQ(ierr); 1152f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(5) (estimated maximum front size in the complete tree): %d \n",lu->id.INFOG(5));CHKERRQ(ierr); 1153f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(6) (number of nodes in the complete tree): %d \n",lu->id.INFOG(6));CHKERRQ(ierr); 1154f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(7) (ordering option effectively uese after analysis): %d \n",lu->id.INFOG(7));CHKERRQ(ierr); 1155f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",lu->id.INFOG(8));CHKERRQ(ierr); 1156f6c57405SHong 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); 1157f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(10) (total integer space store the matrix factors after factorization): %d \n",lu->id.INFOG(10));CHKERRQ(ierr); 1158f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(11) (order of largest frontal matrix after factorization): %d \n",lu->id.INFOG(11));CHKERRQ(ierr); 1159f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(12) (number of off-diagonal pivots): %d \n",lu->id.INFOG(12));CHKERRQ(ierr); 1160f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(13) (number of delayed pivots after factorization): %d \n",lu->id.INFOG(13));CHKERRQ(ierr); 1161f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(14) (number of memory compress after factorization): %d \n",lu->id.INFOG(14));CHKERRQ(ierr); 1162f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(15) (number of steps of iterative refinement after solution): %d \n",lu->id.INFOG(15));CHKERRQ(ierr); 1163f6c57405SHong 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); 1164f6c57405SHong 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); 1165f6c57405SHong 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); 1166f6c57405SHong 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); 1167f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(20) (estimated number of entries in the factors): %d \n",lu->id.INFOG(20));CHKERRQ(ierr); 1168f6c57405SHong 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); 1169f6c57405SHong 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); 1170f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",lu->id.INFOG(23));CHKERRQ(ierr); 1171f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",lu->id.INFOG(24));CHKERRQ(ierr); 1172f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",lu->id.INFOG(25));CHKERRQ(ierr); 1173f6c57405SHong Zhang } 1174f6c57405SHong Zhang } 1175cb828f0fSHong Zhang } 1176f6c57405SHong Zhang PetscFunctionReturn(0); 1177f6c57405SHong Zhang } 1178f6c57405SHong Zhang 117935bd34faSBarry Smith #undef __FUNCT__ 118035bd34faSBarry Smith #define __FUNCT__ "MatGetInfo_MUMPS" 118135bd34faSBarry Smith PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info) 118235bd34faSBarry Smith { 1183cb828f0fSHong Zhang Mat_MUMPS *mumps =(Mat_MUMPS*)A->spptr; 118435bd34faSBarry Smith 118535bd34faSBarry Smith PetscFunctionBegin; 118635bd34faSBarry Smith info->block_size = 1.0; 1187cb828f0fSHong Zhang info->nz_allocated = mumps->id.INFOG(20); 1188cb828f0fSHong Zhang info->nz_used = mumps->id.INFOG(20); 118935bd34faSBarry Smith info->nz_unneeded = 0.0; 119035bd34faSBarry Smith info->assemblies = 0.0; 119135bd34faSBarry Smith info->mallocs = 0.0; 119235bd34faSBarry Smith info->memory = 0.0; 119335bd34faSBarry Smith info->fill_ratio_given = 0; 119435bd34faSBarry Smith info->fill_ratio_needed = 0; 119535bd34faSBarry Smith info->factor_mallocs = 0; 119635bd34faSBarry Smith PetscFunctionReturn(0); 119735bd34faSBarry Smith } 119835bd34faSBarry Smith 11995ccb76cbSHong Zhang /* -------------------------------------------------------------------------------------------*/ 12005ccb76cbSHong Zhang #undef __FUNCT__ 12015ccb76cbSHong Zhang #define __FUNCT__ "MatMumpsSetIcntl_MUMPS" 12025ccb76cbSHong Zhang PetscErrorCode MatMumpsSetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt ival) 12035ccb76cbSHong Zhang { 12045ccb76cbSHong Zhang Mat_MUMPS *lu =(Mat_MUMPS*)F->spptr; 12055ccb76cbSHong Zhang 12065ccb76cbSHong Zhang PetscFunctionBegin; 12075ccb76cbSHong Zhang lu->id.ICNTL(icntl) = ival; 12085ccb76cbSHong Zhang PetscFunctionReturn(0); 12095ccb76cbSHong Zhang } 12105ccb76cbSHong Zhang 12115ccb76cbSHong Zhang #undef __FUNCT__ 12125ccb76cbSHong Zhang #define __FUNCT__ "MatMumpsSetIcntl" 12135ccb76cbSHong Zhang /*@ 12145ccb76cbSHong Zhang MatMumpsSetIcntl - Set MUMPS parameter ICNTL() 12155ccb76cbSHong Zhang 12165ccb76cbSHong Zhang Logically Collective on Mat 12175ccb76cbSHong Zhang 12185ccb76cbSHong Zhang Input Parameters: 12195ccb76cbSHong Zhang + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 12205ccb76cbSHong Zhang . icntl - index of MUMPS parameter array ICNTL() 12215ccb76cbSHong Zhang - ival - value of MUMPS ICNTL(icntl) 12225ccb76cbSHong Zhang 12235ccb76cbSHong Zhang Options Database: 12245ccb76cbSHong Zhang . -mat_mumps_icntl_<icntl> <ival> 12255ccb76cbSHong Zhang 12265ccb76cbSHong Zhang Level: beginner 12275ccb76cbSHong Zhang 12285ccb76cbSHong Zhang References: MUMPS Users' Guide 12295ccb76cbSHong Zhang 12305ccb76cbSHong Zhang .seealso: MatGetFactor() 12315ccb76cbSHong Zhang @*/ 12325ccb76cbSHong Zhang PetscErrorCode MatMumpsSetIcntl(Mat F,PetscInt icntl,PetscInt ival) 12335ccb76cbSHong Zhang { 12345ccb76cbSHong Zhang PetscErrorCode ierr; 12355ccb76cbSHong Zhang 12365ccb76cbSHong Zhang PetscFunctionBegin; 12375ccb76cbSHong Zhang PetscValidLogicalCollectiveInt(F,icntl,2); 12385ccb76cbSHong Zhang PetscValidLogicalCollectiveInt(F,ival,3); 12395ccb76cbSHong Zhang ierr = PetscTryMethod(F,"MatMumpsSetIcntl_C",(Mat,PetscInt,PetscInt),(F,icntl,ival));CHKERRQ(ierr); 12405ccb76cbSHong Zhang PetscFunctionReturn(0); 12415ccb76cbSHong Zhang } 12425ccb76cbSHong Zhang 124324b6179bSKris Buschelman /*MC 12442692d6eeSBarry Smith MATSOLVERMUMPS - A matrix type providing direct solvers (LU and Cholesky) for 124524b6179bSKris Buschelman distributed and sequential matrices via the external package MUMPS. 124624b6179bSKris Buschelman 124741c8de11SBarry Smith Works with MATAIJ and MATSBAIJ matrices 124824b6179bSKris Buschelman 124924b6179bSKris Buschelman Options Database Keys: 1250fb8376fbSHong Zhang + -mat_mumps_icntl_4 <0,...,4> - print level 125124b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide) 125264e6c443SBarry Smith . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guidec) 125324b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T 125424b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements 125594b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view 125624b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide) 125724b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide) 125824b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide) 125924b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide) 126024b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold 126124b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement 126224b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold 126324b6179bSKris Buschelman 126424b6179bSKris Buschelman Level: beginner 126524b6179bSKris Buschelman 126641c8de11SBarry Smith .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage 126741c8de11SBarry Smith 126824b6179bSKris Buschelman M*/ 126924b6179bSKris Buschelman 12702877fffaSHong Zhang EXTERN_C_BEGIN 127135bd34faSBarry Smith #undef __FUNCT__ 127235bd34faSBarry Smith #define __FUNCT__ "MatFactorGetSolverPackage_mumps" 127335bd34faSBarry Smith PetscErrorCode MatFactorGetSolverPackage_mumps(Mat A,const MatSolverPackage *type) 127435bd34faSBarry Smith { 127535bd34faSBarry Smith PetscFunctionBegin; 12762692d6eeSBarry Smith *type = MATSOLVERMUMPS; 127735bd34faSBarry Smith PetscFunctionReturn(0); 127835bd34faSBarry Smith } 127935bd34faSBarry Smith EXTERN_C_END 128035bd34faSBarry Smith 128135bd34faSBarry Smith EXTERN_C_BEGIN 1282bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI AIJ matrices */ 12832877fffaSHong Zhang #undef __FUNCT__ 1284bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_aij_mumps" 1285bccb9932SShri Abhyankar PetscErrorCode MatGetFactor_aij_mumps(Mat A,MatFactorType ftype,Mat *F) 12862877fffaSHong Zhang { 12872877fffaSHong Zhang Mat B; 12882877fffaSHong Zhang PetscErrorCode ierr; 12892877fffaSHong Zhang Mat_MUMPS *mumps; 1290ace3abfcSBarry Smith PetscBool isSeqAIJ; 12912877fffaSHong Zhang 12922877fffaSHong Zhang PetscFunctionBegin; 12932877fffaSHong Zhang /* Create the factorization matrix */ 1294bccb9932SShri Abhyankar ierr = PetscTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr); 12952877fffaSHong Zhang ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 12962877fffaSHong Zhang ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 12972877fffaSHong Zhang ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 1298bccb9932SShri Abhyankar if (isSeqAIJ) { 12992877fffaSHong Zhang ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 1300bccb9932SShri Abhyankar } else { 1301bccb9932SShri Abhyankar ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 1302bccb9932SShri Abhyankar } 13032877fffaSHong Zhang 130416ebf90aSShri Abhyankar ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr); 13052877fffaSHong Zhang B->ops->view = MatView_MUMPS; 130635bd34faSBarry Smith B->ops->getinfo = MatGetInfo_MUMPS; 130735bd34faSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr); 13085ccb76cbSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMumpsSetIcntl_C","MatMumpsSetIcntl_MUMPS",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr); 1309450b117fSShri Abhyankar if (ftype == MAT_FACTOR_LU) { 1310450b117fSShri Abhyankar B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS; 1311d5f3da31SBarry Smith B->factortype = MAT_FACTOR_LU; 1312bccb9932SShri Abhyankar if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqaij; 1313bccb9932SShri Abhyankar else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpiaij; 1314746480a1SHong Zhang mumps->sym = 0; 1315dcd589f8SShri Abhyankar } else { 131667877ebaSShri Abhyankar B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS; 1317450b117fSShri Abhyankar B->factortype = MAT_FACTOR_CHOLESKY; 1318bccb9932SShri Abhyankar if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqsbaij; 1319bccb9932SShri Abhyankar else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpisbaij; 13206fdc2a6dSBarry Smith if (A->spd_set && A->spd) mumps->sym = 1; 13216fdc2a6dSBarry Smith else mumps->sym = 2; 1322450b117fSShri Abhyankar } 13232877fffaSHong Zhang 13242877fffaSHong Zhang mumps->isAIJ = PETSC_TRUE; 1325*bf0cc555SLisandro Dalcin mumps->Destroy = B->ops->destroy; 13262877fffaSHong Zhang B->ops->destroy = MatDestroy_MUMPS; 13272877fffaSHong Zhang B->spptr = (void*)mumps; 1328f697e70eSHong Zhang ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr); 1329746480a1SHong Zhang 13302877fffaSHong Zhang *F = B; 13312877fffaSHong Zhang PetscFunctionReturn(0); 13322877fffaSHong Zhang } 13332877fffaSHong Zhang EXTERN_C_END 13342877fffaSHong Zhang 1335bccb9932SShri Abhyankar 13362877fffaSHong Zhang EXTERN_C_BEGIN 1337bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI SBAIJ matrices */ 13382877fffaSHong Zhang #undef __FUNCT__ 1339bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_sbaij_mumps" 1340bccb9932SShri Abhyankar PetscErrorCode MatGetFactor_sbaij_mumps(Mat A,MatFactorType ftype,Mat *F) 13412877fffaSHong Zhang { 13422877fffaSHong Zhang Mat B; 13432877fffaSHong Zhang PetscErrorCode ierr; 13442877fffaSHong Zhang Mat_MUMPS *mumps; 1345ace3abfcSBarry Smith PetscBool isSeqSBAIJ; 13462877fffaSHong Zhang 13472877fffaSHong Zhang PetscFunctionBegin; 1348e7e72b3dSBarry Smith if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix"); 1349bccb9932SShri 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"); 1350bccb9932SShri Abhyankar ierr = PetscTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr); 13512877fffaSHong Zhang /* Create the factorization matrix */ 13522877fffaSHong Zhang ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 13532877fffaSHong Zhang ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 13542877fffaSHong Zhang ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 135516ebf90aSShri Abhyankar ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr); 1356bccb9932SShri Abhyankar if (isSeqSBAIJ) { 1357bccb9932SShri Abhyankar ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr); 135816ebf90aSShri Abhyankar mumps->ConvertToTriples = MatConvertToTriples_seqsbaij_seqsbaij; 1359dcd589f8SShri Abhyankar } else { 1360bccb9932SShri Abhyankar ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 1361bccb9932SShri Abhyankar mumps->ConvertToTriples = MatConvertToTriples_mpisbaij_mpisbaij; 1362bccb9932SShri Abhyankar } 1363bccb9932SShri Abhyankar 136467877ebaSShri Abhyankar B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS; 1365bccb9932SShri Abhyankar B->ops->view = MatView_MUMPS; 1366bccb9932SShri Abhyankar ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr); 1367f250808bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMumpsSetIcntl_C","MatMumpsSetIcntl",MatMumpsSetIcntl);CHKERRQ(ierr); 1368f4762488SHong Zhang B->factortype = MAT_FACTOR_CHOLESKY; 13696fdc2a6dSBarry Smith if (A->spd_set && A->spd) mumps->sym = 1; 13706fdc2a6dSBarry Smith else mumps->sym = 2; 1371a214ac2aSShri Abhyankar 1372bccb9932SShri Abhyankar mumps->isAIJ = PETSC_FALSE; 1373*bf0cc555SLisandro Dalcin mumps->Destroy = B->ops->destroy; 1374f3c0ef26SHong Zhang B->ops->destroy = MatDestroy_MUMPS; 13752877fffaSHong Zhang B->spptr = (void*)mumps; 1376f697e70eSHong Zhang ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr); 1377746480a1SHong Zhang 13782877fffaSHong Zhang *F = B; 13792877fffaSHong Zhang PetscFunctionReturn(0); 13802877fffaSHong Zhang } 13812877fffaSHong Zhang EXTERN_C_END 138297969023SHong Zhang 1383450b117fSShri Abhyankar EXTERN_C_BEGIN 1384450b117fSShri Abhyankar #undef __FUNCT__ 1385bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_baij_mumps" 1386bccb9932SShri Abhyankar PetscErrorCode MatGetFactor_baij_mumps(Mat A,MatFactorType ftype,Mat *F) 138767877ebaSShri Abhyankar { 138867877ebaSShri Abhyankar Mat B; 138967877ebaSShri Abhyankar PetscErrorCode ierr; 139067877ebaSShri Abhyankar Mat_MUMPS *mumps; 1391ace3abfcSBarry Smith PetscBool isSeqBAIJ; 139267877ebaSShri Abhyankar 139367877ebaSShri Abhyankar PetscFunctionBegin; 139467877ebaSShri Abhyankar /* Create the factorization matrix */ 1395bccb9932SShri Abhyankar ierr = PetscTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);CHKERRQ(ierr); 139667877ebaSShri Abhyankar ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 139767877ebaSShri Abhyankar ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 139867877ebaSShri Abhyankar ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 1399bccb9932SShri Abhyankar if (isSeqBAIJ) { 140067877ebaSShri Abhyankar ierr = MatSeqBAIJSetPreallocation(B,A->rmap->bs,0,PETSC_NULL);CHKERRQ(ierr); 1401bccb9932SShri Abhyankar } else { 140267877ebaSShri Abhyankar ierr = MatMPIBAIJSetPreallocation(B,A->rmap->bs,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 1403bccb9932SShri Abhyankar } 1404450b117fSShri Abhyankar 140567877ebaSShri Abhyankar ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr); 1406450b117fSShri Abhyankar if (ftype == MAT_FACTOR_LU) { 1407450b117fSShri Abhyankar B->ops->lufactorsymbolic = MatLUFactorSymbolic_BAIJMUMPS; 1408450b117fSShri Abhyankar B->factortype = MAT_FACTOR_LU; 1409bccb9932SShri Abhyankar if (isSeqBAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqbaij_seqaij; 1410bccb9932SShri Abhyankar else mumps->ConvertToTriples = MatConvertToTriples_mpibaij_mpiaij; 1411746480a1SHong Zhang mumps->sym = 0; 1412746480a1SHong Zhang } else { 1413746480a1SHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc BAIJ matrices with MUMPS Cholesky, use SBAIJ or AIJ matrix instead\n"); 1414450b117fSShri Abhyankar } 1415bccb9932SShri Abhyankar 1416450b117fSShri Abhyankar B->ops->view = MatView_MUMPS; 1417450b117fSShri Abhyankar ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr); 14185ccb76cbSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMumpsSetIcntl_C","MatMumpsSetIcntl_MUMPS",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr); 1419450b117fSShri Abhyankar 1420450b117fSShri Abhyankar mumps->isAIJ = PETSC_TRUE; 1421*bf0cc555SLisandro Dalcin mumps->Destroy = B->ops->destroy; 1422450b117fSShri Abhyankar B->ops->destroy = MatDestroy_MUMPS; 1423450b117fSShri Abhyankar B->spptr = (void*)mumps; 1424f697e70eSHong Zhang ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr); 1425746480a1SHong Zhang 1426450b117fSShri Abhyankar *F = B; 1427450b117fSShri Abhyankar PetscFunctionReturn(0); 1428450b117fSShri Abhyankar } 1429450b117fSShri Abhyankar EXTERN_C_END 1430a214ac2aSShri Abhyankar 1431