1be1d678aSKris Buschelman #define PETSCMAT_DLL 21c2a3de1SBarry Smith 3397b6df1SKris Buschelman /* 4c2b5dc30SHong Zhang Provides an interface to the MUMPS sparse solver 5397b6df1SKris Buschelman */ 64e1dd20eSHong Zhang #include "../src/mat/impls/aij/seq/aij.h" /*I "petscmat.h" I*/ 77c4f633dSBarry Smith #include "../src/mat/impls/aij/mpi/mpiaij.h" 87c4f633dSBarry Smith #include "../src/mat/impls/sbaij/seq/sbaij.h" 97c4f633dSBarry Smith #include "../src/mat/impls/sbaij/mpi/mpisbaij.h" 10450b117fSShri Abhyankar #include "../src/mat/impls/baij/seq/baij.h" 11450b117fSShri Abhyankar #include "../src/mat/impls/baij/mpi/mpibaij.h" 12397b6df1SKris Buschelman 13397b6df1SKris Buschelman EXTERN_C_BEGIN 14397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 15397b6df1SKris Buschelman #include "zmumps_c.h" 16397b6df1SKris Buschelman #else 17397b6df1SKris Buschelman #include "dmumps_c.h" 18397b6df1SKris Buschelman #endif 19397b6df1SKris Buschelman EXTERN_C_END 20397b6df1SKris Buschelman #define JOB_INIT -1 21397b6df1SKris Buschelman #define JOB_END -2 22397b6df1SKris Buschelman /* macros s.t. indices match MUMPS documentation */ 23397b6df1SKris Buschelman #define ICNTL(I) icntl[(I)-1] 24397b6df1SKris Buschelman #define CNTL(I) cntl[(I)-1] 25397b6df1SKris Buschelman #define INFOG(I) infog[(I)-1] 26a7aca84bSHong Zhang #define INFO(I) info[(I)-1] 27397b6df1SKris Buschelman #define RINFOG(I) rinfog[(I)-1] 28adc1d99fSHong Zhang #define RINFO(I) rinfo[(I)-1] 29397b6df1SKris Buschelman 30397b6df1SKris Buschelman typedef struct { 31397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 32397b6df1SKris Buschelman ZMUMPS_STRUC_C id; 33397b6df1SKris Buschelman #else 34397b6df1SKris Buschelman DMUMPS_STRUC_C id; 35397b6df1SKris Buschelman #endif 36397b6df1SKris Buschelman MatStructure matstruc; 37c1490034SHong Zhang PetscMPIInt myid,size; 3816ebf90aSShri Abhyankar PetscInt *irn,*jcn,nz,sym,nSolve; 39397b6df1SKris Buschelman PetscScalar *val; 40397b6df1SKris Buschelman MPI_Comm comm_mumps; 41329ec9b3SHong Zhang VecScatter scat_rhs, scat_sol; 42cb828f0fSHong Zhang PetscTruth isAIJ,CleanUpMUMPS,mumpsview; 43329ec9b3SHong Zhang Vec b_seq,x_seq; 4467334b25SHong Zhang PetscErrorCode (*MatDestroy)(Mat); 45bccb9932SShri Abhyankar PetscErrorCode (*ConvertToTriples)(Mat, int, MatReuse, int*, int**, int**, PetscScalar**); 46f0c56d0fSKris Buschelman } Mat_MUMPS; 47f0c56d0fSKris Buschelman 48dfbe8321SBarry Smith EXTERN PetscErrorCode MatDuplicate_MUMPS(Mat,MatDuplicateOption,Mat*); 49b24902e0SBarry Smith 5067877ebaSShri Abhyankar 5167877ebaSShri Abhyankar /* MatConvertToTriples_A_B */ 5267877ebaSShri Abhyankar /*convert Petsc matrix to triples: row[nz], col[nz], val[nz] */ 53397b6df1SKris Buschelman /* 54397b6df1SKris Buschelman input: 5567877ebaSShri Abhyankar A - matrix in aij,baij or sbaij (bs=1) format 56397b6df1SKris Buschelman shift - 0: C style output triple; 1: Fortran style output triple. 57bccb9932SShri Abhyankar reuse - MAT_INITIAL_MATRIX: spaces are allocated and values are set for the triple 58bccb9932SShri Abhyankar MAT_REUSE_MATRIX: only the values in v array are updated 59397b6df1SKris Buschelman output: 60397b6df1SKris Buschelman nnz - dim of r, c, and v (number of local nonzero entries of A) 61397b6df1SKris Buschelman r, c, v - row and col index, matrix values (matrix triples) 62397b6df1SKris Buschelman */ 6316ebf90aSShri Abhyankar 6416ebf90aSShri Abhyankar #undef __FUNCT__ 6516ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqaij_seqaij" 66bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 67b24902e0SBarry Smith { 6867877ebaSShri Abhyankar const PetscInt *ai,*aj,*ajj,M=A->rmap->n;; 6967877ebaSShri Abhyankar PetscInt nz,rnz,i,j; 70dfbe8321SBarry Smith PetscErrorCode ierr; 71c1490034SHong Zhang PetscInt *row,*col; 7216ebf90aSShri Abhyankar Mat_SeqAIJ *aa=(Mat_SeqAIJ*)A->data; 73397b6df1SKris Buschelman 74397b6df1SKris Buschelman PetscFunctionBegin; 7516ebf90aSShri Abhyankar *v=aa->a; 76bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX){ 7716ebf90aSShri Abhyankar nz = aa->nz; ai = aa->i; aj = aa->j; 7816ebf90aSShri Abhyankar *nnz = nz; 79*cf3759fdSShri Abhyankar ierr = PetscMalloc2(nz,PetscInt, &row,nz,PetscInt,&col);CHKERRQ(ierr); 8016ebf90aSShri Abhyankar nz = 0; 8116ebf90aSShri Abhyankar for(i=0; i<M; i++) { 8216ebf90aSShri Abhyankar rnz = ai[i+1] - ai[i]; 8367877ebaSShri Abhyankar ajj = aj + ai[i]; 8467877ebaSShri Abhyankar for(j=0; j<rnz; j++) { 8567877ebaSShri Abhyankar row[nz] = i+shift; col[nz++] = ajj[j] + shift; 8616ebf90aSShri Abhyankar } 8716ebf90aSShri Abhyankar } 8816ebf90aSShri Abhyankar *r = row; *c = col; 8916ebf90aSShri Abhyankar } 9016ebf90aSShri Abhyankar PetscFunctionReturn(0); 9116ebf90aSShri Abhyankar } 92397b6df1SKris Buschelman 9316ebf90aSShri Abhyankar #undef __FUNCT__ 9467877ebaSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqbaij_seqaij" 95bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqbaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 9667877ebaSShri Abhyankar { 9767877ebaSShri Abhyankar Mat_SeqBAIJ *aa=(Mat_SeqBAIJ*)A->data; 9867877ebaSShri Abhyankar const PetscInt *ai,*aj,*ajj,bs=A->rmap->bs,bs2=aa->bs2,M=A->rmap->N/bs; 9967877ebaSShri Abhyankar PetscInt nz,idx=0,rnz,i,j,k,m,ii; 10067877ebaSShri Abhyankar PetscErrorCode ierr; 10167877ebaSShri Abhyankar PetscInt *row,*col; 10267877ebaSShri Abhyankar 10367877ebaSShri Abhyankar PetscFunctionBegin; 104*cf3759fdSShri Abhyankar *v = aa->a; 105bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX){ 106*cf3759fdSShri Abhyankar ai = aa->i; aj = aa->j; 10767877ebaSShri Abhyankar nz = bs2*aa->nz; 10867877ebaSShri Abhyankar *nnz = nz; 109*cf3759fdSShri Abhyankar ierr = PetscMalloc2(nz,PetscInt, &row,nz,PetscInt,&col);CHKERRQ(ierr); 11067877ebaSShri Abhyankar for(i=0; i<M; i++) { 11167877ebaSShri Abhyankar ii = 0; 11267877ebaSShri Abhyankar ajj = aj + ai[i]; 11367877ebaSShri Abhyankar rnz = ai[i+1] - ai[i]; 11467877ebaSShri Abhyankar for(k=0; k<rnz; k++) { 11567877ebaSShri Abhyankar for(j=0; j<bs; j++) { 11667877ebaSShri Abhyankar for(m=0; m<bs; m++) { 11767877ebaSShri Abhyankar row[idx] = i*bs + m + shift; 118*cf3759fdSShri Abhyankar col[idx++] = bs*(ajj[k]) + j + shift; 11967877ebaSShri Abhyankar } 12067877ebaSShri Abhyankar } 12167877ebaSShri Abhyankar } 12267877ebaSShri Abhyankar } 123*cf3759fdSShri Abhyankar *r = row; *c = col; 12467877ebaSShri Abhyankar } 12567877ebaSShri Abhyankar PetscFunctionReturn(0); 12667877ebaSShri Abhyankar } 12767877ebaSShri Abhyankar 12867877ebaSShri Abhyankar #undef __FUNCT__ 12916ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqsbaij_seqsbaij" 130bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqsbaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 13116ebf90aSShri Abhyankar { 13267877ebaSShri Abhyankar const PetscInt *ai, *aj,*ajj,M=A->rmap->n; 13367877ebaSShri Abhyankar PetscInt nz,rnz,i,j; 13416ebf90aSShri Abhyankar PetscErrorCode ierr; 13516ebf90aSShri Abhyankar PetscInt *row,*col; 13616ebf90aSShri Abhyankar Mat_SeqSBAIJ *aa=(Mat_SeqSBAIJ*)A->data; 13716ebf90aSShri Abhyankar 13816ebf90aSShri Abhyankar PetscFunctionBegin; 139bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX){ 14016ebf90aSShri Abhyankar nz = aa->nz;ai=aa->i; aj=aa->j;*v=aa->a; 14116ebf90aSShri Abhyankar *nnz = nz; 142*cf3759fdSShri Abhyankar ierr = PetscMalloc2(nz,PetscInt, &row,nz,PetscInt,&col);CHKERRQ(ierr); 14316ebf90aSShri Abhyankar nz = 0; 14416ebf90aSShri Abhyankar for(i=0; i<M; i++) { 14516ebf90aSShri Abhyankar rnz = ai[i+1] - ai[i]; 14667877ebaSShri Abhyankar ajj = aj + ai[i]; 14767877ebaSShri Abhyankar for(j=0; j<rnz; j++) { 14867877ebaSShri Abhyankar row[nz] = i+shift; col[nz++] = ajj[j] + shift; 14916ebf90aSShri Abhyankar } 15016ebf90aSShri Abhyankar } 15116ebf90aSShri Abhyankar *r = row; *c = col; 15216ebf90aSShri Abhyankar } 15316ebf90aSShri Abhyankar PetscFunctionReturn(0); 15416ebf90aSShri Abhyankar } 15516ebf90aSShri Abhyankar 15616ebf90aSShri Abhyankar #undef __FUNCT__ 15716ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqaij_seqsbaij" 158bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 15916ebf90aSShri Abhyankar { 16067877ebaSShri Abhyankar const PetscInt *ai,*aj,*ajj,*adiag,M=A->rmap->n; 16167877ebaSShri Abhyankar PetscInt nz,rnz,i,j; 16267877ebaSShri Abhyankar const PetscScalar *av,*v1; 16316ebf90aSShri Abhyankar PetscScalar *val; 16416ebf90aSShri Abhyankar PetscErrorCode ierr; 16516ebf90aSShri Abhyankar PetscInt *row,*col; 16616ebf90aSShri Abhyankar Mat_SeqSBAIJ *aa=(Mat_SeqSBAIJ*)A->data; 16716ebf90aSShri Abhyankar 16816ebf90aSShri Abhyankar PetscFunctionBegin; 16916ebf90aSShri Abhyankar ai=aa->i; aj=aa->j;av=aa->a; 17016ebf90aSShri Abhyankar adiag=aa->diag; 171bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX){ 17216ebf90aSShri Abhyankar nz = M + (aa->nz-M)/2; 17316ebf90aSShri Abhyankar *nnz = nz; 174*cf3759fdSShri Abhyankar ierr = PetscMalloc3(nz,PetscInt, &row,nz,PetscInt,&col,nz,PetscScalar,&val);CHKERRQ(ierr); 17516ebf90aSShri Abhyankar nz = 0; 17616ebf90aSShri Abhyankar for(i=0; i<M; i++) { 17716ebf90aSShri Abhyankar rnz = ai[i+1] - adiag[i]; 17867877ebaSShri Abhyankar ajj = aj + adiag[i]; 179*cf3759fdSShri Abhyankar v1 = av + adiag[i]; 18067877ebaSShri Abhyankar for(j=0; j<rnz; j++) { 18167877ebaSShri Abhyankar row[nz] = i+shift; col[nz] = ajj[j] + shift; val[nz++] = v1[j]; 18216ebf90aSShri Abhyankar } 18316ebf90aSShri Abhyankar } 18416ebf90aSShri Abhyankar *r = row; *c = col; *v = val; 185397b6df1SKris Buschelman } else { 18616ebf90aSShri Abhyankar nz = 0; val = *v; 18716ebf90aSShri Abhyankar for(i=0; i <M; i++) { 18816ebf90aSShri Abhyankar rnz = ai[i+1] - adiag[i]; 18967877ebaSShri Abhyankar ajj = aj + adiag[i]; 19067877ebaSShri Abhyankar v1 = av + adiag[i]; 19167877ebaSShri Abhyankar for(j=0; j<rnz; j++) { 19267877ebaSShri Abhyankar val[nz++] = v1[j]; 19316ebf90aSShri Abhyankar } 19416ebf90aSShri Abhyankar } 19516ebf90aSShri Abhyankar } 19616ebf90aSShri Abhyankar PetscFunctionReturn(0); 19716ebf90aSShri Abhyankar } 19816ebf90aSShri Abhyankar 19916ebf90aSShri Abhyankar #undef __FUNCT__ 20016ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpisbaij_mpisbaij" 201bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpisbaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 20216ebf90aSShri Abhyankar { 20316ebf90aSShri Abhyankar const PetscInt *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj; 20416ebf90aSShri Abhyankar PetscErrorCode ierr; 20516ebf90aSShri Abhyankar PetscInt rstart,nz,i,j,jj,irow,countA,countB; 20616ebf90aSShri Abhyankar PetscInt *row,*col; 20716ebf90aSShri Abhyankar const PetscScalar *av, *bv,*v1,*v2; 20816ebf90aSShri Abhyankar PetscScalar *val; 209397b6df1SKris Buschelman Mat_MPISBAIJ *mat = (Mat_MPISBAIJ*)A->data; 210397b6df1SKris Buschelman Mat_SeqSBAIJ *aa=(Mat_SeqSBAIJ*)(mat->A)->data; 211397b6df1SKris Buschelman Mat_SeqBAIJ *bb=(Mat_SeqBAIJ*)(mat->B)->data; 21216ebf90aSShri Abhyankar 21316ebf90aSShri Abhyankar PetscFunctionBegin; 214d0f46423SBarry Smith ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart; 215397b6df1SKris Buschelman garray = mat->garray; 216397b6df1SKris Buschelman av=aa->a; bv=bb->a; 217397b6df1SKris Buschelman 218bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX){ 21916ebf90aSShri Abhyankar nz = aa->nz + bb->nz; 22016ebf90aSShri Abhyankar *nnz = nz; 221*cf3759fdSShri Abhyankar ierr = PetscMalloc3(nz,PetscInt, &row,nz,PetscInt,&col,nz,PetscScalar,&val);CHKERRQ(ierr); 222397b6df1SKris Buschelman *r = row; *c = col; *v = val; 223397b6df1SKris Buschelman } else { 224397b6df1SKris Buschelman row = *r; col = *c; val = *v; 225397b6df1SKris Buschelman } 226397b6df1SKris Buschelman 227028e57e8SHong Zhang jj = 0; irow = rstart; 228397b6df1SKris Buschelman for ( i=0; i<m; i++ ) { 229397b6df1SKris Buschelman ajj = aj + ai[i]; /* ptr to the beginning of this row */ 230397b6df1SKris Buschelman countA = ai[i+1] - ai[i]; 231397b6df1SKris Buschelman countB = bi[i+1] - bi[i]; 232397b6df1SKris Buschelman bjj = bj + bi[i]; 23316ebf90aSShri Abhyankar v1 = av + ai[i]; 23416ebf90aSShri Abhyankar v2 = bv + bi[i]; 235397b6df1SKris Buschelman 236397b6df1SKris Buschelman /* A-part */ 237397b6df1SKris Buschelman for (j=0; j<countA; j++){ 238bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 239397b6df1SKris Buschelman row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift; 240397b6df1SKris Buschelman } 24116ebf90aSShri Abhyankar val[jj++] = v1[j]; 242397b6df1SKris Buschelman } 24316ebf90aSShri Abhyankar 24416ebf90aSShri Abhyankar /* B-part */ 24516ebf90aSShri Abhyankar for(j=0; j < countB; j++){ 246bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 247397b6df1SKris Buschelman row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift; 248397b6df1SKris Buschelman } 24916ebf90aSShri Abhyankar val[jj++] = v2[j]; 25016ebf90aSShri Abhyankar } 25116ebf90aSShri Abhyankar irow++; 25216ebf90aSShri Abhyankar } 25316ebf90aSShri Abhyankar PetscFunctionReturn(0); 25416ebf90aSShri Abhyankar } 25516ebf90aSShri Abhyankar 25616ebf90aSShri Abhyankar #undef __FUNCT__ 25716ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpiaij_mpiaij" 258bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpiaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 25916ebf90aSShri Abhyankar { 26016ebf90aSShri Abhyankar const PetscInt *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj; 26116ebf90aSShri Abhyankar PetscErrorCode ierr; 26216ebf90aSShri Abhyankar PetscInt rstart,nz,i,j,jj,irow,countA,countB; 26316ebf90aSShri Abhyankar PetscInt *row,*col; 26416ebf90aSShri Abhyankar const PetscScalar *av, *bv,*v1,*v2; 26516ebf90aSShri Abhyankar PetscScalar *val; 26616ebf90aSShri Abhyankar Mat_MPIAIJ *mat = (Mat_MPIAIJ*)A->data; 26716ebf90aSShri Abhyankar Mat_SeqAIJ *aa=(Mat_SeqAIJ*)(mat->A)->data; 26816ebf90aSShri Abhyankar Mat_SeqAIJ *bb=(Mat_SeqAIJ*)(mat->B)->data; 26916ebf90aSShri Abhyankar 27016ebf90aSShri Abhyankar PetscFunctionBegin; 27116ebf90aSShri Abhyankar ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart; 27216ebf90aSShri Abhyankar garray = mat->garray; 27316ebf90aSShri Abhyankar av=aa->a; bv=bb->a; 27416ebf90aSShri Abhyankar 275bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX){ 27616ebf90aSShri Abhyankar nz = aa->nz + bb->nz; 27716ebf90aSShri Abhyankar *nnz = nz; 278*cf3759fdSShri Abhyankar ierr = PetscMalloc3(nz,PetscInt, &row,nz,PetscInt,&col,nz,PetscScalar,&val);CHKERRQ(ierr); 27916ebf90aSShri Abhyankar *r = row; *c = col; *v = val; 28016ebf90aSShri Abhyankar } else { 28116ebf90aSShri Abhyankar row = *r; col = *c; val = *v; 28216ebf90aSShri Abhyankar } 28316ebf90aSShri Abhyankar 28416ebf90aSShri Abhyankar jj = 0; irow = rstart; 28516ebf90aSShri Abhyankar for ( i=0; i<m; i++ ) { 28616ebf90aSShri Abhyankar ajj = aj + ai[i]; /* ptr to the beginning of this row */ 28716ebf90aSShri Abhyankar countA = ai[i+1] - ai[i]; 28816ebf90aSShri Abhyankar countB = bi[i+1] - bi[i]; 28916ebf90aSShri Abhyankar bjj = bj + bi[i]; 29016ebf90aSShri Abhyankar v1 = av + ai[i]; 29116ebf90aSShri Abhyankar v2 = bv + bi[i]; 29216ebf90aSShri Abhyankar 29316ebf90aSShri Abhyankar /* A-part */ 29416ebf90aSShri Abhyankar for (j=0; j<countA; j++){ 295bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX){ 29616ebf90aSShri Abhyankar row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift; 29716ebf90aSShri Abhyankar } 29816ebf90aSShri Abhyankar val[jj++] = v1[j]; 29916ebf90aSShri Abhyankar } 30016ebf90aSShri Abhyankar 30116ebf90aSShri Abhyankar /* B-part */ 30216ebf90aSShri Abhyankar for(j=0; j < countB; j++){ 303bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX){ 30416ebf90aSShri Abhyankar row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift; 30516ebf90aSShri Abhyankar } 30616ebf90aSShri Abhyankar val[jj++] = v2[j]; 30716ebf90aSShri Abhyankar } 30816ebf90aSShri Abhyankar irow++; 30916ebf90aSShri Abhyankar } 31016ebf90aSShri Abhyankar PetscFunctionReturn(0); 31116ebf90aSShri Abhyankar } 31216ebf90aSShri Abhyankar 31316ebf90aSShri Abhyankar #undef __FUNCT__ 31467877ebaSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpibaij_mpiaij" 315bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpibaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 31667877ebaSShri Abhyankar { 31767877ebaSShri Abhyankar Mat_MPIBAIJ *mat = (Mat_MPIBAIJ*)A->data; 31867877ebaSShri Abhyankar Mat_SeqBAIJ *aa=(Mat_SeqBAIJ*)(mat->A)->data; 31967877ebaSShri Abhyankar Mat_SeqBAIJ *bb=(Mat_SeqBAIJ*)(mat->B)->data; 32067877ebaSShri Abhyankar const PetscInt *ai = aa->i, *bi = bb->i, *aj = aa->j, *bj = bb->j,*ajj, *bjj; 32167877ebaSShri Abhyankar const PetscInt *garray = mat->garray,mbs=mat->mbs,rstartbs=mat->rstartbs; 32267877ebaSShri Abhyankar const PetscInt bs = A->rmap->bs,bs2=mat->bs2; 32367877ebaSShri Abhyankar PetscErrorCode ierr; 32467877ebaSShri Abhyankar PetscInt nz,i,j,k,n,jj,irow,countA,countB,idx; 32567877ebaSShri Abhyankar PetscInt *row,*col; 32667877ebaSShri Abhyankar const PetscScalar *av=aa->a, *bv=bb->a,*v1,*v2; 32767877ebaSShri Abhyankar PetscScalar *val; 32867877ebaSShri Abhyankar 32967877ebaSShri Abhyankar PetscFunctionBegin; 33067877ebaSShri Abhyankar 331bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 33267877ebaSShri Abhyankar nz = bs2*(aa->nz + bb->nz); 33367877ebaSShri Abhyankar *nnz = nz; 334*cf3759fdSShri Abhyankar ierr = PetscMalloc3(nz,PetscInt, &row,nz,PetscInt,&col,nz,PetscScalar,&val);CHKERRQ(ierr); 33567877ebaSShri Abhyankar *r = row; *c = col; *v = val; 33667877ebaSShri Abhyankar } else { 33767877ebaSShri Abhyankar row = *r; col = *c; val = *v; 33867877ebaSShri Abhyankar } 33967877ebaSShri Abhyankar 34067877ebaSShri Abhyankar jj = 0; irow = rstartbs; 34167877ebaSShri Abhyankar for ( i=0; i<mbs; i++ ) { 34267877ebaSShri Abhyankar countA = ai[i+1] - ai[i]; 34367877ebaSShri Abhyankar countB = bi[i+1] - bi[i]; 34467877ebaSShri Abhyankar ajj = aj + ai[i]; 34567877ebaSShri Abhyankar bjj = bj + bi[i]; 34667877ebaSShri Abhyankar v1 = av + bs2*ai[i]; 34767877ebaSShri Abhyankar v2 = bv + bs2*bi[i]; 34867877ebaSShri Abhyankar 34967877ebaSShri Abhyankar idx = 0; 35067877ebaSShri Abhyankar /* A-part */ 35167877ebaSShri Abhyankar for (k=0; k<countA; k++){ 35267877ebaSShri Abhyankar for (j=0; j<bs; j++) { 35367877ebaSShri Abhyankar for (n=0; n<bs; n++) { 354bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX){ 35567877ebaSShri Abhyankar row[jj] = bs*irow + n + shift; 35667877ebaSShri Abhyankar col[jj] = bs*(rstartbs + ajj[k]) + j + shift; 35767877ebaSShri Abhyankar } 35867877ebaSShri Abhyankar val[jj++] = v1[idx++]; 35967877ebaSShri Abhyankar } 36067877ebaSShri Abhyankar } 36167877ebaSShri Abhyankar } 36267877ebaSShri Abhyankar 36367877ebaSShri Abhyankar idx = 0; 36467877ebaSShri Abhyankar /* B-part */ 36567877ebaSShri Abhyankar for(k=0; k<countB; k++){ 36667877ebaSShri Abhyankar for (j=0; j<bs; j++) { 36767877ebaSShri Abhyankar for (n=0; n<bs; n++) { 368bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX){ 36967877ebaSShri Abhyankar row[jj] = bs*irow + n + shift; 37067877ebaSShri Abhyankar col[jj] = bs*(garray[bjj[k]]) + j + shift; 37167877ebaSShri Abhyankar } 37267877ebaSShri Abhyankar val[jj++] = bv[idx++]; 37367877ebaSShri Abhyankar } 37467877ebaSShri Abhyankar } 37567877ebaSShri Abhyankar } 37667877ebaSShri Abhyankar irow++; 37767877ebaSShri Abhyankar } 37867877ebaSShri Abhyankar PetscFunctionReturn(0); 37967877ebaSShri Abhyankar } 38067877ebaSShri Abhyankar 38167877ebaSShri Abhyankar #undef __FUNCT__ 38216ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpiaij_mpisbaij" 383bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpiaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 38416ebf90aSShri Abhyankar { 38516ebf90aSShri Abhyankar const PetscInt *ai, *aj,*adiag, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj; 38616ebf90aSShri Abhyankar PetscErrorCode ierr; 38716ebf90aSShri Abhyankar PetscInt rstart,nz,nza,nzb_low,i,j,jj,irow,countA,countB; 38816ebf90aSShri Abhyankar PetscInt *row,*col; 38916ebf90aSShri Abhyankar const PetscScalar *av, *bv,*v1,*v2; 39016ebf90aSShri Abhyankar PetscScalar *val; 39116ebf90aSShri Abhyankar Mat_MPIAIJ *mat = (Mat_MPIAIJ*)A->data; 39216ebf90aSShri Abhyankar Mat_SeqAIJ *aa=(Mat_SeqAIJ*)(mat->A)->data; 39316ebf90aSShri Abhyankar Mat_SeqAIJ *bb=(Mat_SeqAIJ*)(mat->B)->data; 39416ebf90aSShri Abhyankar 39516ebf90aSShri Abhyankar PetscFunctionBegin; 39616ebf90aSShri Abhyankar ai=aa->i; aj=aa->j; adiag=aa->diag; 39716ebf90aSShri Abhyankar bi=bb->i; bj=bb->j; garray = mat->garray; 39816ebf90aSShri Abhyankar av=aa->a; bv=bb->a; 39916ebf90aSShri Abhyankar rstart = A->rmap->rstart; 40016ebf90aSShri Abhyankar 401bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 40216ebf90aSShri Abhyankar nza = 0;nzb_low = 0; 40316ebf90aSShri Abhyankar for(i=0; i<m; i++){ 40416ebf90aSShri Abhyankar nza = nza + (ai[i+1] - adiag[i]); 40516ebf90aSShri Abhyankar countB = bi[i+1] - bi[i]; 40616ebf90aSShri Abhyankar bjj = bj + bi[i]; 40716ebf90aSShri Abhyankar 40816ebf90aSShri Abhyankar j = 0; 40916ebf90aSShri Abhyankar while(garray[bjj[j]] < rstart) { 41016ebf90aSShri Abhyankar if(j == countB) break; 41116ebf90aSShri Abhyankar j++;nzb_low++; 41216ebf90aSShri Abhyankar } 41316ebf90aSShri Abhyankar } 41416ebf90aSShri Abhyankar /* Total nz = nz for the upper triangular A part + nz for the 2nd B part */ 41516ebf90aSShri Abhyankar nz = nza + (bb->nz - nzb_low); 41616ebf90aSShri Abhyankar *nnz = nz; 417*cf3759fdSShri Abhyankar ierr = PetscMalloc3(nz,PetscInt, &row,nz,PetscInt,&col,nz,PetscScalar,&val);CHKERRQ(ierr); 41816ebf90aSShri Abhyankar *r = row; *c = col; *v = val; 41916ebf90aSShri Abhyankar } else { 42016ebf90aSShri Abhyankar row = *r; col = *c; val = *v; 42116ebf90aSShri Abhyankar } 42216ebf90aSShri Abhyankar 42316ebf90aSShri Abhyankar jj = 0; irow = rstart; 42416ebf90aSShri Abhyankar for ( i=0; i<m; i++ ) { 42516ebf90aSShri Abhyankar ajj = aj + adiag[i]; /* ptr to the beginning of the diagonal of this row */ 42616ebf90aSShri Abhyankar v1 = av + adiag[i]; 42716ebf90aSShri Abhyankar countA = ai[i+1] - adiag[i]; 42816ebf90aSShri Abhyankar countB = bi[i+1] - bi[i]; 42916ebf90aSShri Abhyankar bjj = bj + bi[i]; 43016ebf90aSShri Abhyankar v2 = bv + bi[i]; 43116ebf90aSShri Abhyankar 43216ebf90aSShri Abhyankar /* A-part */ 43316ebf90aSShri Abhyankar for (j=0; j<countA; j++){ 434bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 43516ebf90aSShri Abhyankar row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift; 43616ebf90aSShri Abhyankar } 43716ebf90aSShri Abhyankar val[jj++] = v1[j]; 43816ebf90aSShri Abhyankar } 43916ebf90aSShri Abhyankar 44016ebf90aSShri Abhyankar /* B-part */ 44116ebf90aSShri Abhyankar for(j=0; j < countB; j++){ 44216ebf90aSShri Abhyankar if (garray[bjj[j]] > rstart) { 443bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 44416ebf90aSShri Abhyankar row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift; 44516ebf90aSShri Abhyankar } 44616ebf90aSShri Abhyankar val[jj++] = v2[j]; 44716ebf90aSShri Abhyankar } 448397b6df1SKris Buschelman } 449397b6df1SKris Buschelman irow++; 450397b6df1SKris Buschelman } 451397b6df1SKris Buschelman PetscFunctionReturn(0); 452397b6df1SKris Buschelman } 453397b6df1SKris Buschelman 454397b6df1SKris Buschelman #undef __FUNCT__ 4553924e44cSKris Buschelman #define __FUNCT__ "MatDestroy_MUMPS" 456dfbe8321SBarry Smith PetscErrorCode MatDestroy_MUMPS(Mat A) 457dfbe8321SBarry Smith { 458f0c56d0fSKris Buschelman Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 459dfbe8321SBarry Smith PetscErrorCode ierr; 460c1490034SHong Zhang PetscMPIInt size=lu->size; 461b24902e0SBarry Smith 462397b6df1SKris Buschelman PetscFunctionBegin; 463397b6df1SKris Buschelman if (lu->CleanUpMUMPS) { 464397b6df1SKris Buschelman /* Terminate instance, deallocate memories */ 465329ec9b3SHong Zhang if (size > 1){ 46668653410SBarry Smith ierr = PetscFree2(lu->id.sol_loc,lu->id.isol_loc);CHKERRQ(ierr); 467329ec9b3SHong Zhang ierr = VecScatterDestroy(lu->scat_rhs);CHKERRQ(ierr); 468329ec9b3SHong Zhang ierr = VecDestroy(lu->b_seq);CHKERRQ(ierr); 4692750af12SHong Zhang if (lu->nSolve && lu->scat_sol){ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr);} 4702750af12SHong Zhang if (lu->nSolve && lu->x_seq){ierr = VecDestroy(lu->x_seq);CHKERRQ(ierr);} 471329ec9b3SHong Zhang ierr = PetscFree(lu->val);CHKERRQ(ierr); 472329ec9b3SHong Zhang } 47316ebf90aSShri Abhyankar if( size == 1 && A->factortype == MAT_FACTOR_CHOLESKY && lu->isAIJ) { 474450b117fSShri Abhyankar ierr = PetscFree(lu->val);CHKERRQ(ierr); 475450b117fSShri Abhyankar } 476397b6df1SKris Buschelman lu->id.job=JOB_END; 477397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 478397b6df1SKris Buschelman zmumps_c(&lu->id); 479397b6df1SKris Buschelman #else 480397b6df1SKris Buschelman dmumps_c(&lu->id); 481397b6df1SKris Buschelman #endif 482c338a77dSKris Buschelman ierr = PetscFree(lu->irn);CHKERRQ(ierr); 483c338a77dSKris Buschelman ierr = PetscFree(lu->jcn);CHKERRQ(ierr); 484397b6df1SKris Buschelman ierr = MPI_Comm_free(&(lu->comm_mumps));CHKERRQ(ierr); 485397b6df1SKris Buschelman } 48697969023SHong Zhang /* clear composed functions */ 48797969023SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatFactorGetSolverPackage_C","",PETSC_NULL);CHKERRQ(ierr); 48897969023SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMumpsSetIcntl_C","",PETSC_NULL);CHKERRQ(ierr); 48967334b25SHong Zhang ierr = (lu->MatDestroy)(A);CHKERRQ(ierr); 490397b6df1SKris Buschelman PetscFunctionReturn(0); 491397b6df1SKris Buschelman } 492397b6df1SKris Buschelman 493397b6df1SKris Buschelman #undef __FUNCT__ 494f6c57405SHong Zhang #define __FUNCT__ "MatSolve_MUMPS" 495b24902e0SBarry Smith PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x) 496b24902e0SBarry Smith { 497f0c56d0fSKris Buschelman Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 498d54de34fSKris Buschelman PetscScalar *array; 49967877ebaSShri Abhyankar Vec b_seq; 500329ec9b3SHong Zhang IS is_iden,is_petsc; 501dfbe8321SBarry Smith PetscErrorCode ierr; 502329ec9b3SHong Zhang PetscInt i; 503397b6df1SKris Buschelman 504397b6df1SKris Buschelman PetscFunctionBegin; 505329ec9b3SHong Zhang lu->id.nrhs = 1; 50667877ebaSShri Abhyankar b_seq = lu->b_seq; 507397b6df1SKris Buschelman if (lu->size > 1){ 508329ec9b3SHong Zhang /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */ 50967877ebaSShri Abhyankar ierr = VecScatterBegin(lu->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 51067877ebaSShri Abhyankar ierr = VecScatterEnd(lu->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 51167877ebaSShri Abhyankar if (!lu->myid) {ierr = VecGetArray(b_seq,&array);CHKERRQ(ierr);} 512397b6df1SKris Buschelman } else { /* size == 1 */ 513397b6df1SKris Buschelman ierr = VecCopy(b,x);CHKERRQ(ierr); 514397b6df1SKris Buschelman ierr = VecGetArray(x,&array);CHKERRQ(ierr); 515397b6df1SKris Buschelman } 516397b6df1SKris Buschelman if (!lu->myid) { /* define rhs on the host */ 5178278f211SHong Zhang lu->id.nrhs = 1; 518397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 519397b6df1SKris Buschelman lu->id.rhs = (mumps_double_complex*)array; 520397b6df1SKris Buschelman #else 521397b6df1SKris Buschelman lu->id.rhs = array; 522397b6df1SKris Buschelman #endif 523397b6df1SKris Buschelman } 524397b6df1SKris Buschelman 525397b6df1SKris Buschelman /* solve phase */ 526329ec9b3SHong Zhang /*-------------*/ 527397b6df1SKris Buschelman lu->id.job = 3; 528397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 529397b6df1SKris Buschelman zmumps_c(&lu->id); 530397b6df1SKris Buschelman #else 531397b6df1SKris Buschelman dmumps_c(&lu->id); 532397b6df1SKris Buschelman #endif 53365e19b50SBarry 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)); 534397b6df1SKris Buschelman 535329ec9b3SHong Zhang if (lu->size > 1) { /* convert mumps distributed solution to petsc mpi x */ 536329ec9b3SHong Zhang if (!lu->nSolve){ /* create scatter scat_sol */ 537329ec9b3SHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,lu->id.lsol_loc,0,1,&is_iden);CHKERRQ(ierr); /* from */ 538329ec9b3SHong Zhang for (i=0; i<lu->id.lsol_loc; i++){ 539329ec9b3SHong Zhang lu->id.isol_loc[i] -= 1; /* change Fortran style to C style */ 540397b6df1SKris Buschelman } 541329ec9b3SHong Zhang ierr = ISCreateGeneral(PETSC_COMM_SELF,lu->id.lsol_loc,lu->id.isol_loc,&is_petsc);CHKERRQ(ierr); /* to */ 542329ec9b3SHong Zhang ierr = VecScatterCreate(lu->x_seq,is_iden,x,is_petsc,&lu->scat_sol);CHKERRQ(ierr); 543329ec9b3SHong Zhang ierr = ISDestroy(is_iden);CHKERRQ(ierr); 544329ec9b3SHong Zhang ierr = ISDestroy(is_petsc);CHKERRQ(ierr); 545397b6df1SKris Buschelman } 546ca9f406cSSatish Balay ierr = VecScatterBegin(lu->scat_sol,lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 547ca9f406cSSatish Balay ierr = VecScatterEnd(lu->scat_sol,lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 548329ec9b3SHong Zhang } 549329ec9b3SHong Zhang lu->nSolve++; 550397b6df1SKris Buschelman PetscFunctionReturn(0); 551397b6df1SKris Buschelman } 552397b6df1SKris Buschelman 553ace3df97SHong Zhang #if !defined(PETSC_USE_COMPLEX) 554a58c3f20SHong Zhang /* 555a58c3f20SHong Zhang input: 556a58c3f20SHong Zhang F: numeric factor 557a58c3f20SHong Zhang output: 558a58c3f20SHong Zhang nneg: total number of negative pivots 559a58c3f20SHong Zhang nzero: 0 560a58c3f20SHong Zhang npos: (global dimension of F) - nneg 561a58c3f20SHong Zhang */ 562a58c3f20SHong Zhang 563a58c3f20SHong Zhang #undef __FUNCT__ 564a58c3f20SHong Zhang #define __FUNCT__ "MatGetInertia_SBAIJMUMPS" 565dfbe8321SBarry Smith PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos) 566a58c3f20SHong Zhang { 567a58c3f20SHong Zhang Mat_MUMPS *lu =(Mat_MUMPS*)F->spptr; 568dfbe8321SBarry Smith PetscErrorCode ierr; 569c1490034SHong Zhang PetscMPIInt size; 570a58c3f20SHong Zhang 571a58c3f20SHong Zhang PetscFunctionBegin; 5727adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)F)->comm,&size);CHKERRQ(ierr); 573bcb30aebSHong 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 */ 57465e19b50SBarry 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)); 575a58c3f20SHong Zhang if (nneg){ 576a58c3f20SHong Zhang if (!lu->myid){ 577a58c3f20SHong Zhang *nneg = lu->id.INFOG(12); 578a58c3f20SHong Zhang } 579bcb30aebSHong Zhang ierr = MPI_Bcast(nneg,1,MPI_INT,0,lu->comm_mumps);CHKERRQ(ierr); 580a58c3f20SHong Zhang } 581a58c3f20SHong Zhang if (nzero) *nzero = 0; 582d0f46423SBarry Smith if (npos) *npos = F->rmap->N - (*nneg); 583a58c3f20SHong Zhang PetscFunctionReturn(0); 584a58c3f20SHong Zhang } 585ace3df97SHong Zhang #endif /* !defined(PETSC_USE_COMPLEX) */ 586a58c3f20SHong Zhang 587397b6df1SKris Buschelman #undef __FUNCT__ 588f6c57405SHong Zhang #define __FUNCT__ "MatFactorNumeric_MUMPS" 5890481f469SBarry Smith PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info) 590af281ebdSHong Zhang { 591dcd589f8SShri Abhyankar Mat_MUMPS *lu =(Mat_MUMPS*)(F)->spptr; 5926849ba73SBarry Smith PetscErrorCode ierr; 593bccb9932SShri Abhyankar MatReuse reuse; 594e09efc27SHong Zhang Mat F_diag; 59516ebf90aSShri Abhyankar PetscTruth isMPIAIJ; 596397b6df1SKris Buschelman 597397b6df1SKris Buschelman PetscFunctionBegin; 598bccb9932SShri Abhyankar reuse = MAT_REUSE_MATRIX; 599bccb9932SShri Abhyankar ierr = (*lu->ConvertToTriples)(A, 1, reuse, &lu->nz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr); 600397b6df1SKris Buschelman 601397b6df1SKris Buschelman /* numerical factorization phase */ 602329ec9b3SHong Zhang /*-------------------------------*/ 603329ec9b3SHong Zhang lu->id.job = 2; 604958c9bccSBarry Smith if(!lu->id.ICNTL(18)) { 605a7aca84bSHong Zhang if (!lu->myid) { 606397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 607397b6df1SKris Buschelman lu->id.a = (mumps_double_complex*)lu->val; 608397b6df1SKris Buschelman #else 609397b6df1SKris Buschelman lu->id.a = lu->val; 610397b6df1SKris Buschelman #endif 611397b6df1SKris Buschelman } 612397b6df1SKris Buschelman } else { 613397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 614397b6df1SKris Buschelman lu->id.a_loc = (mumps_double_complex*)lu->val; 615397b6df1SKris Buschelman #else 616397b6df1SKris Buschelman lu->id.a_loc = lu->val; 617397b6df1SKris Buschelman #endif 618397b6df1SKris Buschelman } 619397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 620397b6df1SKris Buschelman zmumps_c(&lu->id); 621397b6df1SKris Buschelman #else 622397b6df1SKris Buschelman dmumps_c(&lu->id); 623397b6df1SKris Buschelman #endif 624397b6df1SKris Buschelman if (lu->id.INFOG(1) < 0) { 62565e19b50SBarry 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)); 62665e19b50SBarry 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)); 627397b6df1SKris Buschelman } 62865e19b50SBarry 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)); 629397b6df1SKris Buschelman 6308ada1bb4SHong Zhang if (lu->size > 1){ 63116ebf90aSShri Abhyankar ierr = PetscTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr); 632a214ac2aSShri Abhyankar if(isMPIAIJ) { 633dcd589f8SShri Abhyankar F_diag = ((Mat_MPIAIJ *)(F)->data)->A; 634e09efc27SHong Zhang } else { 635dcd589f8SShri Abhyankar F_diag = ((Mat_MPISBAIJ *)(F)->data)->A; 636e09efc27SHong Zhang } 637e09efc27SHong Zhang F_diag->assembled = PETSC_TRUE; 638329ec9b3SHong Zhang if (lu->nSolve){ 639329ec9b3SHong Zhang ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr); 6400e83c824SBarry Smith ierr = PetscFree2(lu->id.sol_loc,lu->id.isol_loc);CHKERRQ(ierr); 641329ec9b3SHong Zhang ierr = VecDestroy(lu->x_seq);CHKERRQ(ierr); 642329ec9b3SHong Zhang } 6438ada1bb4SHong Zhang } 644dcd589f8SShri Abhyankar (F)->assembled = PETSC_TRUE; 645397b6df1SKris Buschelman lu->matstruc = SAME_NONZERO_PATTERN; 646ace87b0dSHong Zhang lu->CleanUpMUMPS = PETSC_TRUE; 647329ec9b3SHong Zhang lu->nSolve = 0; 64867877ebaSShri Abhyankar 64967877ebaSShri Abhyankar if (lu->size > 1){ 65067877ebaSShri Abhyankar /* distributed solution */ 65167877ebaSShri Abhyankar lu->id.ICNTL(21) = 1; 65267877ebaSShri Abhyankar if (!lu->nSolve){ 65367877ebaSShri Abhyankar /* Create x_seq=sol_loc for repeated use */ 65467877ebaSShri Abhyankar PetscInt lsol_loc; 65567877ebaSShri Abhyankar PetscScalar *sol_loc; 65667877ebaSShri Abhyankar lsol_loc = lu->id.INFO(23); /* length of sol_loc */ 65767877ebaSShri Abhyankar ierr = PetscMalloc2(lsol_loc,PetscScalar,&sol_loc,lsol_loc,PetscInt,&lu->id.isol_loc);CHKERRQ(ierr); 65867877ebaSShri Abhyankar lu->id.lsol_loc = lsol_loc; 65967877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 66067877ebaSShri Abhyankar lu->id.sol_loc = (mumps_double_complex*)sol_loc; 66167877ebaSShri Abhyankar #else 66267877ebaSShri Abhyankar lu->id.sol_loc = sol_loc; 66367877ebaSShri Abhyankar #endif 66467877ebaSShri Abhyankar ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,lsol_loc,sol_loc,&lu->x_seq);CHKERRQ(ierr); 66567877ebaSShri Abhyankar } 66667877ebaSShri Abhyankar } 667397b6df1SKris Buschelman PetscFunctionReturn(0); 668397b6df1SKris Buschelman } 669397b6df1SKris Buschelman 670dcd589f8SShri Abhyankar #undef __FUNCT__ 671dcd589f8SShri Abhyankar #define __FUNCT__ "PetscSetMUMPSOptions" 672dcd589f8SShri Abhyankar PetscErrorCode PetscSetMUMPSOptions(Mat F, Mat A) 673dcd589f8SShri Abhyankar { 674dcd589f8SShri Abhyankar Mat_MUMPS *lu = (Mat_MUMPS*)F->spptr; 675dcd589f8SShri Abhyankar PetscErrorCode ierr; 676dcd589f8SShri Abhyankar PetscInt icntl; 677dcd589f8SShri Abhyankar PetscTruth flg; 678dcd589f8SShri Abhyankar 679dcd589f8SShri Abhyankar PetscFunctionBegin; 680dcd589f8SShri Abhyankar ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"MUMPS Options","Mat");CHKERRQ(ierr); 681cb828f0fSHong Zhang ierr = PetscOptionsTruth("-mat_mumps_view","View MUMPS parameters","None",lu->mumpsview,&lu->mumpsview,PETSC_NULL);CHKERRQ(ierr); 682dcd589f8SShri Abhyankar if (lu->size == 1){ 683dcd589f8SShri Abhyankar lu->id.ICNTL(18) = 0; /* centralized assembled matrix input */ 684dcd589f8SShri Abhyankar } else { 685dcd589f8SShri Abhyankar lu->id.ICNTL(18) = 3; /* distributed assembled matrix input */ 686dcd589f8SShri Abhyankar } 687dcd589f8SShri Abhyankar 688dcd589f8SShri Abhyankar icntl=-1; 689dcd589f8SShri Abhyankar lu->id.ICNTL(4) = 0; /* level of printing; overwrite mumps default ICNTL(4)=2 */ 690dcd589f8SShri Abhyankar ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",lu->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr); 691dcd589f8SShri Abhyankar if ((flg && icntl > 0) || PetscLogPrintInfo) { 692dcd589f8SShri Abhyankar lu->id.ICNTL(4)=icntl; /* and use mumps default icntl(i), i=1,2,3 */ 693dcd589f8SShri Abhyankar } else { /* no output */ 694dcd589f8SShri Abhyankar lu->id.ICNTL(1) = 0; /* error message, default= 6 */ 695dcd589f8SShri Abhyankar lu->id.ICNTL(2) = 0; /* output stream for diagnostic printing, statistics, and warning. default=0 */ 696dcd589f8SShri Abhyankar lu->id.ICNTL(3) = 0; /* output stream for global information, default=6 */ 697dcd589f8SShri Abhyankar } 698dcd589f8SShri 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); 699dcd589f8SShri Abhyankar icntl=-1; 700dcd589f8SShri Abhyankar ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): matrix ordering (0 to 7)","None",lu->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr); 701dcd589f8SShri Abhyankar if (flg) { 702dcd589f8SShri Abhyankar if (icntl== 1){ 703e32f2f54SBarry 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"); 704dcd589f8SShri Abhyankar } else { 705dcd589f8SShri Abhyankar lu->id.ICNTL(7) = icntl; 706dcd589f8SShri Abhyankar } 707dcd589f8SShri Abhyankar } 708dcd589f8SShri 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); 709dcd589f8SShri Abhyankar ierr = PetscOptionsInt("-mat_mumps_icntl_9","ICNTL(9): A or A^T x=b to be solved. 1: A; otherwise: A^T","None",lu->id.ICNTL(9),&lu->id.ICNTL(9),PETSC_NULL);CHKERRQ(ierr); 710dcd589f8SShri 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); 711dcd589f8SShri 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); 712dcd589f8SShri 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); 713dcd589f8SShri 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); 714dcd589f8SShri 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); 715dcd589f8SShri Abhyankar ierr = PetscOptionsInt("-mat_mumps_icntl_19","ICNTL(19): Schur complement","None",lu->id.ICNTL(19),&lu->id.ICNTL(19),PETSC_NULL);CHKERRQ(ierr); 716dcd589f8SShri Abhyankar 717dcd589f8SShri 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); 718dcd589f8SShri 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); 719dcd589f8SShri 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); 720dcd589f8SShri 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); 721dcd589f8SShri 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); 722dcd589f8SShri Abhyankar ierr = PetscOptionsInt("-mat_mumps_icntl_27","ICNTL(27): experimental parameter","None",lu->id.ICNTL(27),&lu->id.ICNTL(27),PETSC_NULL);CHKERRQ(ierr); 723dcd589f8SShri Abhyankar 724dcd589f8SShri 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); 725dcd589f8SShri 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); 726dcd589f8SShri 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); 727dcd589f8SShri 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); 728dcd589f8SShri 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); 729dcd589f8SShri Abhyankar PetscOptionsEnd(); 730dcd589f8SShri Abhyankar PetscFunctionReturn(0); 731dcd589f8SShri Abhyankar } 732dcd589f8SShri Abhyankar 733dcd589f8SShri Abhyankar #undef __FUNCT__ 734dcd589f8SShri Abhyankar #define __FUNCT__ "PetscInitializeMUMPS" 735dcd589f8SShri Abhyankar PetscErrorCode PetscInitializeMUMPS(Mat F) 736dcd589f8SShri Abhyankar { 737dcd589f8SShri Abhyankar Mat_MUMPS *lu = (Mat_MUMPS*)F->spptr; 738dcd589f8SShri Abhyankar PetscErrorCode ierr; 739dcd589f8SShri Abhyankar PetscInt icntl; 740dcd589f8SShri Abhyankar PetscTruth flg; 741dcd589f8SShri Abhyankar 742dcd589f8SShri Abhyankar PetscFunctionBegin; 743dcd589f8SShri Abhyankar lu->id.job = JOB_INIT; 744dcd589f8SShri Abhyankar lu->id.par=1; /* host participates factorizaton and solve */ 745dcd589f8SShri Abhyankar lu->id.sym=lu->sym; 746dcd589f8SShri Abhyankar if (lu->sym == 2){ 747dcd589f8SShri Abhyankar ierr = PetscOptionsInt("-mat_mumps_sym","SYM: (1,2)","None",lu->id.sym,&icntl,&flg);CHKERRQ(ierr); 748dcd589f8SShri Abhyankar if (flg && icntl == 1) lu->id.sym=icntl; /* matrix is spd */ 749dcd589f8SShri Abhyankar } 750dcd589f8SShri Abhyankar #if defined(PETSC_USE_COMPLEX) 751dcd589f8SShri Abhyankar zmumps_c(&lu->id); 752dcd589f8SShri Abhyankar #else 753dcd589f8SShri Abhyankar dmumps_c(&lu->id); 754dcd589f8SShri Abhyankar #endif 755dcd589f8SShri Abhyankar PetscFunctionReturn(0); 756dcd589f8SShri Abhyankar } 757dcd589f8SShri Abhyankar 758397b6df1SKris Buschelman /* Note the Petsc r and c permutations are ignored */ 759397b6df1SKris Buschelman #undef __FUNCT__ 760f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS" 7610481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 762b24902e0SBarry Smith { 763719d5645SBarry Smith Mat_MUMPS *lu = (Mat_MUMPS*)F->spptr; 764dcd589f8SShri Abhyankar PetscErrorCode ierr; 765bccb9932SShri Abhyankar MatReuse reuse; 76667877ebaSShri Abhyankar Vec b; 76767877ebaSShri Abhyankar IS is_iden; 76867877ebaSShri Abhyankar const PetscInt M = A->rmap->N; 769397b6df1SKris Buschelman 770397b6df1SKris Buschelman PetscFunctionBegin; 771b24902e0SBarry Smith lu->sym = 0; 772b24902e0SBarry Smith lu->matstruc = DIFFERENT_NONZERO_PATTERN; 773dcd589f8SShri Abhyankar 774dcd589f8SShri Abhyankar ierr = MPI_Comm_rank(((PetscObject)A)->comm, &lu->myid); 775dcd589f8SShri Abhyankar ierr = MPI_Comm_size(((PetscObject)A)->comm,&lu->size);CHKERRQ(ierr); 776dcd589f8SShri Abhyankar ierr = MPI_Comm_dup(((PetscObject)A)->comm,&(lu->comm_mumps));CHKERRQ(ierr); 777dcd589f8SShri Abhyankar lu->id.comm_fortran = MPI_Comm_c2f(lu->comm_mumps); 778dcd589f8SShri Abhyankar 779dcd589f8SShri Abhyankar /* Initialize a MUMPS instance */ 780dcd589f8SShri Abhyankar ierr = PetscInitializeMUMPS(F);CHKERRQ(ierr); 781dcd589f8SShri Abhyankar /* Set MUMPS options */ 782dcd589f8SShri Abhyankar ierr = PetscSetMUMPSOptions(F,A);CHKERRQ(ierr); 783dcd589f8SShri Abhyankar 784bccb9932SShri Abhyankar reuse = MAT_INITIAL_MATRIX; 785bccb9932SShri Abhyankar ierr = (*lu->ConvertToTriples)(A, 1, reuse, &lu->nz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr); 786dcd589f8SShri Abhyankar 78767877ebaSShri Abhyankar /* analysis phase */ 78867877ebaSShri Abhyankar /*----------------*/ 78967877ebaSShri Abhyankar lu->id.job = 1; 79067877ebaSShri Abhyankar lu->id.n = M; 79167877ebaSShri Abhyankar switch (lu->id.ICNTL(18)){ 79267877ebaSShri Abhyankar case 0: /* centralized assembled matrix input */ 79367877ebaSShri Abhyankar if (!lu->myid) { 79467877ebaSShri Abhyankar lu->id.nz =lu->nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn; 79567877ebaSShri Abhyankar if (lu->id.ICNTL(6)>1){ 79667877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 79767877ebaSShri Abhyankar lu->id.a = (mumps_double_complex*)lu->val; 79867877ebaSShri Abhyankar #else 79967877ebaSShri Abhyankar lu->id.a = lu->val; 80067877ebaSShri Abhyankar #endif 80167877ebaSShri Abhyankar } 80267877ebaSShri Abhyankar } 80367877ebaSShri Abhyankar break; 80467877ebaSShri Abhyankar case 3: /* distributed assembled matrix input (size>1) */ 80567877ebaSShri Abhyankar lu->id.nz_loc = lu->nz; 80667877ebaSShri Abhyankar lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn; 80767877ebaSShri Abhyankar if (lu->id.ICNTL(6)>1) { 80867877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 80967877ebaSShri Abhyankar lu->id.a_loc = (mumps_double_complex*)lu->val; 81067877ebaSShri Abhyankar #else 81167877ebaSShri Abhyankar lu->id.a_loc = lu->val; 81267877ebaSShri Abhyankar #endif 81367877ebaSShri Abhyankar } 81467877ebaSShri Abhyankar /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */ 81567877ebaSShri Abhyankar if (!lu->myid){ 81667877ebaSShri Abhyankar ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&lu->b_seq);CHKERRQ(ierr); 81767877ebaSShri Abhyankar ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr); 81867877ebaSShri Abhyankar } else { 81967877ebaSShri Abhyankar ierr = VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);CHKERRQ(ierr); 82067877ebaSShri Abhyankar ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr); 82167877ebaSShri Abhyankar } 82267877ebaSShri Abhyankar ierr = VecCreate(((PetscObject)A)->comm,&b);CHKERRQ(ierr); 82367877ebaSShri Abhyankar ierr = VecSetSizes(b,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr); 82467877ebaSShri Abhyankar ierr = VecSetFromOptions(b);CHKERRQ(ierr); 82567877ebaSShri Abhyankar 82667877ebaSShri Abhyankar ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr); 82767877ebaSShri Abhyankar ierr = ISDestroy(is_iden);CHKERRQ(ierr); 82867877ebaSShri Abhyankar ierr = VecDestroy(b);CHKERRQ(ierr); 82967877ebaSShri Abhyankar break; 83067877ebaSShri Abhyankar } 83167877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 83267877ebaSShri Abhyankar zmumps_c(&lu->id); 83367877ebaSShri Abhyankar #else 83467877ebaSShri Abhyankar dmumps_c(&lu->id); 83567877ebaSShri Abhyankar #endif 83667877ebaSShri 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)); 83767877ebaSShri Abhyankar 838719d5645SBarry Smith F->ops->lufactornumeric = MatFactorNumeric_MUMPS; 839dcd589f8SShri Abhyankar F->ops->solve = MatSolve_MUMPS; 840b24902e0SBarry Smith PetscFunctionReturn(0); 841b24902e0SBarry Smith } 842b24902e0SBarry Smith 843450b117fSShri Abhyankar /* Note the Petsc r and c permutations are ignored */ 844450b117fSShri Abhyankar #undef __FUNCT__ 845450b117fSShri Abhyankar #define __FUNCT__ "MatLUFactorSymbolic_BAIJMUMPS" 846450b117fSShri Abhyankar PetscErrorCode MatLUFactorSymbolic_BAIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 847450b117fSShri Abhyankar { 848dcd589f8SShri Abhyankar 849450b117fSShri Abhyankar Mat_MUMPS *lu = (Mat_MUMPS*)F->spptr; 850dcd589f8SShri Abhyankar PetscErrorCode ierr; 851bccb9932SShri Abhyankar MatReuse reuse; 85267877ebaSShri Abhyankar Vec b; 85367877ebaSShri Abhyankar IS is_iden; 85467877ebaSShri Abhyankar const PetscInt M = A->rmap->N; 855450b117fSShri Abhyankar 856450b117fSShri Abhyankar PetscFunctionBegin; 857450b117fSShri Abhyankar lu->sym = 0; 858450b117fSShri Abhyankar lu->matstruc = DIFFERENT_NONZERO_PATTERN; 859dcd589f8SShri Abhyankar ierr = MPI_Comm_rank(((PetscObject)A)->comm, &lu->myid); 860dcd589f8SShri Abhyankar ierr = MPI_Comm_size(((PetscObject)A)->comm,&lu->size);CHKERRQ(ierr); 861dcd589f8SShri Abhyankar ierr = MPI_Comm_dup(((PetscObject)A)->comm,&(lu->comm_mumps));CHKERRQ(ierr); 862dcd589f8SShri Abhyankar lu->id.comm_fortran = MPI_Comm_c2f(lu->comm_mumps); 863dcd589f8SShri Abhyankar 864dcd589f8SShri Abhyankar /* Initialize a MUMPS instance */ 865dcd589f8SShri Abhyankar ierr = PetscInitializeMUMPS(F);CHKERRQ(ierr); 866dcd589f8SShri Abhyankar /* Set MUMPS options */ 867dcd589f8SShri Abhyankar ierr = PetscSetMUMPSOptions(F,A);CHKERRQ(ierr); 868dcd589f8SShri Abhyankar 869bccb9932SShri Abhyankar reuse = MAT_INITIAL_MATRIX; 870bccb9932SShri Abhyankar ierr = (*lu->ConvertToTriples)(A, 1, reuse, &lu->nz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr); 87167877ebaSShri Abhyankar 87267877ebaSShri Abhyankar /* analysis phase */ 87367877ebaSShri Abhyankar /*----------------*/ 87467877ebaSShri Abhyankar lu->id.job = 1; 87567877ebaSShri Abhyankar lu->id.n = M; 87667877ebaSShri Abhyankar switch (lu->id.ICNTL(18)){ 87767877ebaSShri Abhyankar case 0: /* centralized assembled matrix input */ 87867877ebaSShri Abhyankar if (!lu->myid) { 87967877ebaSShri Abhyankar lu->id.nz =lu->nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn; 88067877ebaSShri Abhyankar if (lu->id.ICNTL(6)>1){ 88167877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 88267877ebaSShri Abhyankar lu->id.a = (mumps_double_complex*)lu->val; 88367877ebaSShri Abhyankar #else 88467877ebaSShri Abhyankar lu->id.a = lu->val; 88567877ebaSShri Abhyankar #endif 88667877ebaSShri Abhyankar } 88767877ebaSShri Abhyankar } 88867877ebaSShri Abhyankar break; 88967877ebaSShri Abhyankar case 3: /* distributed assembled matrix input (size>1) */ 89067877ebaSShri Abhyankar lu->id.nz_loc = lu->nz; 89167877ebaSShri Abhyankar lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn; 89267877ebaSShri Abhyankar if (lu->id.ICNTL(6)>1) { 89367877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 89467877ebaSShri Abhyankar lu->id.a_loc = (mumps_double_complex*)lu->val; 89567877ebaSShri Abhyankar #else 89667877ebaSShri Abhyankar lu->id.a_loc = lu->val; 89767877ebaSShri Abhyankar #endif 89867877ebaSShri Abhyankar } 89967877ebaSShri Abhyankar /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */ 90067877ebaSShri Abhyankar if (!lu->myid){ 90167877ebaSShri Abhyankar ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&lu->b_seq);CHKERRQ(ierr); 90267877ebaSShri Abhyankar ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr); 90367877ebaSShri Abhyankar } else { 90467877ebaSShri Abhyankar ierr = VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);CHKERRQ(ierr); 90567877ebaSShri Abhyankar ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr); 90667877ebaSShri Abhyankar } 90767877ebaSShri Abhyankar ierr = VecCreate(((PetscObject)A)->comm,&b);CHKERRQ(ierr); 90867877ebaSShri Abhyankar ierr = VecSetSizes(b,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr); 90967877ebaSShri Abhyankar ierr = VecSetFromOptions(b);CHKERRQ(ierr); 91067877ebaSShri Abhyankar 91167877ebaSShri Abhyankar ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr); 91267877ebaSShri Abhyankar ierr = ISDestroy(is_iden);CHKERRQ(ierr); 91367877ebaSShri Abhyankar ierr = VecDestroy(b);CHKERRQ(ierr); 91467877ebaSShri Abhyankar break; 91567877ebaSShri Abhyankar } 91667877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 91767877ebaSShri Abhyankar zmumps_c(&lu->id); 91867877ebaSShri Abhyankar #else 91967877ebaSShri Abhyankar dmumps_c(&lu->id); 92067877ebaSShri Abhyankar #endif 92167877ebaSShri 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)); 92267877ebaSShri Abhyankar 923450b117fSShri Abhyankar F->ops->lufactornumeric = MatFactorNumeric_MUMPS; 924dcd589f8SShri Abhyankar F->ops->solve = MatSolve_MUMPS; 925450b117fSShri Abhyankar PetscFunctionReturn(0); 926450b117fSShri Abhyankar } 927b24902e0SBarry Smith 928397b6df1SKris Buschelman /* Note the Petsc r permutation is ignored */ 929397b6df1SKris Buschelman #undef __FUNCT__ 93067877ebaSShri Abhyankar #define __FUNCT__ "MatCholeskyFactorSymbolic_MUMPS" 93167877ebaSShri Abhyankar PetscErrorCode MatCholeskyFactorSymbolic_MUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info) 932b24902e0SBarry Smith { 9332792810eSHong Zhang Mat_MUMPS *lu = (Mat_MUMPS*)F->spptr; 934dcd589f8SShri Abhyankar PetscErrorCode ierr; 935bccb9932SShri Abhyankar MatReuse reuse; 93667877ebaSShri Abhyankar Vec b; 93767877ebaSShri Abhyankar IS is_iden; 93867877ebaSShri Abhyankar const PetscInt M = A->rmap->N; 939397b6df1SKris Buschelman 940397b6df1SKris Buschelman PetscFunctionBegin; 941b24902e0SBarry Smith lu->sym = 2; 942b24902e0SBarry Smith lu->matstruc = DIFFERENT_NONZERO_PATTERN; 943dcd589f8SShri Abhyankar ierr = MPI_Comm_rank(((PetscObject)A)->comm, &lu->myid); 944dcd589f8SShri Abhyankar ierr = MPI_Comm_size(((PetscObject)A)->comm,&lu->size);CHKERRQ(ierr); 945dcd589f8SShri Abhyankar ierr = MPI_Comm_dup(((PetscObject)A)->comm,&(lu->comm_mumps));CHKERRQ(ierr); 946dcd589f8SShri Abhyankar lu->id.comm_fortran = MPI_Comm_c2f(lu->comm_mumps); 947dcd589f8SShri Abhyankar 948dcd589f8SShri Abhyankar /* Initialize a MUMPS instance */ 949dcd589f8SShri Abhyankar ierr = PetscInitializeMUMPS(F);CHKERRQ(ierr); 950dcd589f8SShri Abhyankar /* Set MUMPS options */ 951dcd589f8SShri Abhyankar ierr = PetscSetMUMPSOptions(F,A);CHKERRQ(ierr); 952dcd589f8SShri Abhyankar 953bccb9932SShri Abhyankar reuse = MAT_INITIAL_MATRIX; 954bccb9932SShri Abhyankar ierr = (*lu->ConvertToTriples)(A, 1 , reuse, &lu->nz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr); 955dcd589f8SShri Abhyankar 95667877ebaSShri Abhyankar /* analysis phase */ 95767877ebaSShri Abhyankar /*----------------*/ 95867877ebaSShri Abhyankar lu->id.job = 1; 95967877ebaSShri Abhyankar lu->id.n = M; 96067877ebaSShri Abhyankar switch (lu->id.ICNTL(18)){ 96167877ebaSShri Abhyankar case 0: /* centralized assembled matrix input */ 96267877ebaSShri Abhyankar if (!lu->myid) { 96367877ebaSShri Abhyankar lu->id.nz =lu->nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn; 96467877ebaSShri Abhyankar if (lu->id.ICNTL(6)>1){ 96567877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 96667877ebaSShri Abhyankar lu->id.a = (mumps_double_complex*)lu->val; 96767877ebaSShri Abhyankar #else 96867877ebaSShri Abhyankar lu->id.a = lu->val; 96967877ebaSShri Abhyankar #endif 97067877ebaSShri Abhyankar } 97167877ebaSShri Abhyankar } 97267877ebaSShri Abhyankar break; 97367877ebaSShri Abhyankar case 3: /* distributed assembled matrix input (size>1) */ 97467877ebaSShri Abhyankar lu->id.nz_loc = lu->nz; 97567877ebaSShri Abhyankar lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn; 97667877ebaSShri Abhyankar if (lu->id.ICNTL(6)>1) { 97767877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 97867877ebaSShri Abhyankar lu->id.a_loc = (mumps_double_complex*)lu->val; 97967877ebaSShri Abhyankar #else 98067877ebaSShri Abhyankar lu->id.a_loc = lu->val; 98167877ebaSShri Abhyankar #endif 98267877ebaSShri Abhyankar } 98367877ebaSShri Abhyankar /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */ 98467877ebaSShri Abhyankar if (!lu->myid){ 98567877ebaSShri Abhyankar ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&lu->b_seq);CHKERRQ(ierr); 98667877ebaSShri Abhyankar ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr); 98767877ebaSShri Abhyankar } else { 98867877ebaSShri Abhyankar ierr = VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);CHKERRQ(ierr); 98967877ebaSShri Abhyankar ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr); 99067877ebaSShri Abhyankar } 99167877ebaSShri Abhyankar ierr = VecCreate(((PetscObject)A)->comm,&b);CHKERRQ(ierr); 99267877ebaSShri Abhyankar ierr = VecSetSizes(b,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr); 99367877ebaSShri Abhyankar ierr = VecSetFromOptions(b);CHKERRQ(ierr); 99467877ebaSShri Abhyankar 99567877ebaSShri Abhyankar ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr); 99667877ebaSShri Abhyankar ierr = ISDestroy(is_iden);CHKERRQ(ierr); 99767877ebaSShri Abhyankar ierr = VecDestroy(b);CHKERRQ(ierr); 99867877ebaSShri Abhyankar break; 99967877ebaSShri Abhyankar } 100067877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 100167877ebaSShri Abhyankar zmumps_c(&lu->id); 100267877ebaSShri Abhyankar #else 100367877ebaSShri Abhyankar dmumps_c(&lu->id); 100467877ebaSShri Abhyankar #endif 100567877ebaSShri 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)); 100667877ebaSShri Abhyankar 10072792810eSHong Zhang F->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS; 1008dcd589f8SShri Abhyankar F->ops->solve = MatSolve_MUMPS; 1009db4efbfdSBarry Smith #if !defined(PETSC_USE_COMPLEX) 1010dcd589f8SShri Abhyankar (F)->ops->getinertia = MatGetInertia_SBAIJMUMPS; 1011db4efbfdSBarry Smith #endif 1012b24902e0SBarry Smith PetscFunctionReturn(0); 1013b24902e0SBarry Smith } 1014b24902e0SBarry Smith 1015397b6df1SKris Buschelman #undef __FUNCT__ 1016f6c57405SHong Zhang #define __FUNCT__ "MatFactorInfo_MUMPS" 101774ed9c26SBarry Smith PetscErrorCode MatFactorInfo_MUMPS(Mat A,PetscViewer viewer) 101874ed9c26SBarry Smith { 1019f6c57405SHong Zhang Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 1020f6c57405SHong Zhang PetscErrorCode ierr; 1021f6c57405SHong Zhang 1022f6c57405SHong Zhang PetscFunctionBegin; 1023f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(1) (output for error): %d \n",lu->id.ICNTL(1));CHKERRQ(ierr); 1024f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(2) (output of diagnostic msg):%d \n",lu->id.ICNTL(2));CHKERRQ(ierr); 1025f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(3) (output for global info): %d \n",lu->id.ICNTL(3));CHKERRQ(ierr); 1026f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(4) (level of printing): %d \n",lu->id.ICNTL(4));CHKERRQ(ierr); 1027f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(5) (input mat struct): %d \n",lu->id.ICNTL(5));CHKERRQ(ierr); 1028f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(6) (matrix prescaling): %d \n",lu->id.ICNTL(6));CHKERRQ(ierr); 1029f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(7) (matrix ordering): %d \n",lu->id.ICNTL(7));CHKERRQ(ierr); 1030f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(8) (scalling strategy): %d \n",lu->id.ICNTL(8));CHKERRQ(ierr); 1031f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(9) (A/A^T x=b is solved): %d \n",lu->id.ICNTL(9));CHKERRQ(ierr); 1032f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(10) (max num of refinements): %d \n",lu->id.ICNTL(10));CHKERRQ(ierr); 1033f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(11) (error analysis): %d \n",lu->id.ICNTL(11));CHKERRQ(ierr); 1034f6c57405SHong Zhang if (!lu->myid && lu->id.ICNTL(11)>0) { 1035f6c57405SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(4) (inf norm of input mat): %g\n",lu->id.RINFOG(4));CHKERRQ(ierr); 1036f6c57405SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(5) (inf norm of solution): %g\n",lu->id.RINFOG(5));CHKERRQ(ierr); 1037f6c57405SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(6) (inf norm of residual): %g\n",lu->id.RINFOG(6));CHKERRQ(ierr); 1038f6c57405SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",lu->id.RINFOG(7),lu->id.RINFOG(8));CHKERRQ(ierr); 1039f6c57405SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(9) (error estimate): %g \n",lu->id.RINFOG(9));CHKERRQ(ierr); 1040f6c57405SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",lu->id.RINFOG(10),lu->id.RINFOG(11));CHKERRQ(ierr); 1041f6c57405SHong Zhang 1042f6c57405SHong Zhang } 1043f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(12) (efficiency control): %d \n",lu->id.ICNTL(12));CHKERRQ(ierr); 1044f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(13) (efficiency control): %d \n",lu->id.ICNTL(13));CHKERRQ(ierr); 1045f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(14) (percentage of estimated workspace increase): %d \n",lu->id.ICNTL(14));CHKERRQ(ierr); 1046f6c57405SHong Zhang /* ICNTL(15-17) not used */ 1047f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(18) (input mat struct): %d \n",lu->id.ICNTL(18));CHKERRQ(ierr); 1048f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(19) (Shur complement info): %d \n",lu->id.ICNTL(19));CHKERRQ(ierr); 1049f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(20) (rhs sparse pattern): %d \n",lu->id.ICNTL(20));CHKERRQ(ierr); 1050f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(21) (solution struct): %d \n",lu->id.ICNTL(21));CHKERRQ(ierr); 1051c0165424SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(22) (in-core/out-of-core facility): %d \n",lu->id.ICNTL(22));CHKERRQ(ierr); 1052c0165424SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(23) (max size of memory can be allocated locally):%d \n",lu->id.ICNTL(23));CHKERRQ(ierr); 1053c0165424SHong Zhang 1054c0165424SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(24) (detection of null pivot rows): %d \n",lu->id.ICNTL(24));CHKERRQ(ierr); 1055c0165424SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(25) (computation of a null space basis): %d \n",lu->id.ICNTL(25));CHKERRQ(ierr); 1056c0165424SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(26) (Schur options for rhs or solution): %d \n",lu->id.ICNTL(26));CHKERRQ(ierr); 1057c0165424SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(27) (experimental parameter): %d \n",lu->id.ICNTL(27));CHKERRQ(ierr); 1058f6c57405SHong Zhang 1059f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(1) (relative pivoting threshold): %g \n",lu->id.CNTL(1));CHKERRQ(ierr); 1060f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(2) (stopping criterion of refinement): %g \n",lu->id.CNTL(2));CHKERRQ(ierr); 1061f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(3) (absolute pivoting threshold): %g \n",lu->id.CNTL(3));CHKERRQ(ierr); 1062f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(4) (value of static pivoting): %g \n",lu->id.CNTL(4));CHKERRQ(ierr); 1063c0165424SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(5) (fixation for null pivots): %g \n",lu->id.CNTL(5));CHKERRQ(ierr); 1064f6c57405SHong Zhang 1065f6c57405SHong Zhang /* infomation local to each processor */ 1066f6c57405SHong Zhang if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, " RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr);} 10677adad957SLisandro Dalcin ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %g \n",lu->myid,lu->id.RINFO(1));CHKERRQ(ierr); 10687adad957SLisandro Dalcin ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 1069f6c57405SHong Zhang if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, " RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr);} 10707adad957SLisandro Dalcin ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %g \n",lu->myid,lu->id.RINFO(2));CHKERRQ(ierr); 10717adad957SLisandro Dalcin ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 1072f6c57405SHong Zhang if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, " RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr);} 10737adad957SLisandro Dalcin ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %g \n",lu->myid,lu->id.RINFO(3));CHKERRQ(ierr); 10747adad957SLisandro Dalcin ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 1075f6c57405SHong Zhang 1076f6c57405SHong Zhang if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, " INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization): \n");CHKERRQ(ierr);} 10777adad957SLisandro Dalcin ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %d \n",lu->myid,lu->id.INFO(15));CHKERRQ(ierr); 10787adad957SLisandro Dalcin ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 1079f6c57405SHong Zhang 1080f6c57405SHong Zhang if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, " INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization): \n");CHKERRQ(ierr);} 10817adad957SLisandro Dalcin ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %d \n",lu->myid,lu->id.INFO(16));CHKERRQ(ierr); 10827adad957SLisandro Dalcin ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 1083f6c57405SHong Zhang 1084f6c57405SHong Zhang if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, " INFO(23) (num of pivots eliminated on this processor after factorization): \n");CHKERRQ(ierr);} 10857adad957SLisandro Dalcin ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %d \n",lu->myid,lu->id.INFO(23));CHKERRQ(ierr); 10867adad957SLisandro Dalcin ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 1087f6c57405SHong Zhang 1088f6c57405SHong Zhang if (!lu->myid){ /* information from the host */ 1089f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",lu->id.RINFOG(1));CHKERRQ(ierr); 1090f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",lu->id.RINFOG(2));CHKERRQ(ierr); 1091f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",lu->id.RINFOG(3));CHKERRQ(ierr); 1092f6c57405SHong Zhang 1093f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(3));CHKERRQ(ierr); 1094f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(4));CHKERRQ(ierr); 1095f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(5) (estimated maximum front size in the complete tree): %d \n",lu->id.INFOG(5));CHKERRQ(ierr); 1096f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(6) (number of nodes in the complete tree): %d \n",lu->id.INFOG(6));CHKERRQ(ierr); 1097f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(7) (ordering option effectively uese after analysis): %d \n",lu->id.INFOG(7));CHKERRQ(ierr); 1098f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",lu->id.INFOG(8));CHKERRQ(ierr); 1099f6c57405SHong 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); 1100f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(10) (total integer space store the matrix factors after factorization): %d \n",lu->id.INFOG(10));CHKERRQ(ierr); 1101f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(11) (order of largest frontal matrix after factorization): %d \n",lu->id.INFOG(11));CHKERRQ(ierr); 1102f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(12) (number of off-diagonal pivots): %d \n",lu->id.INFOG(12));CHKERRQ(ierr); 1103f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(13) (number of delayed pivots after factorization): %d \n",lu->id.INFOG(13));CHKERRQ(ierr); 1104f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(14) (number of memory compress after factorization): %d \n",lu->id.INFOG(14));CHKERRQ(ierr); 1105f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(15) (number of steps of iterative refinement after solution): %d \n",lu->id.INFOG(15));CHKERRQ(ierr); 1106f6c57405SHong 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); 1107f6c57405SHong 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); 1108f6c57405SHong 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); 1109f6c57405SHong 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); 1110f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(20) (estimated number of entries in the factors): %d \n",lu->id.INFOG(20));CHKERRQ(ierr); 1111f6c57405SHong 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); 1112f6c57405SHong 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); 1113f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",lu->id.INFOG(23));CHKERRQ(ierr); 1114f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",lu->id.INFOG(24));CHKERRQ(ierr); 1115f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",lu->id.INFOG(25));CHKERRQ(ierr); 1116f6c57405SHong Zhang } 1117f6c57405SHong Zhang PetscFunctionReturn(0); 1118f6c57405SHong Zhang } 1119f6c57405SHong Zhang 1120f6c57405SHong Zhang #undef __FUNCT__ 1121f6c57405SHong Zhang #define __FUNCT__ "MatView_MUMPS" 1122b24902e0SBarry Smith PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer) 1123b24902e0SBarry Smith { 1124f6c57405SHong Zhang PetscErrorCode ierr; 1125f6c57405SHong Zhang PetscTruth iascii; 1126f6c57405SHong Zhang PetscViewerFormat format; 1127cb828f0fSHong Zhang Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 1128f6c57405SHong Zhang 1129f6c57405SHong Zhang PetscFunctionBegin; 1130cb828f0fSHong Zhang /* check if matrix is mumps type */ 1131cb828f0fSHong Zhang if (A->ops->solve != MatSolve_MUMPS) PetscFunctionReturn(0); 1132cb828f0fSHong Zhang 1133f6c57405SHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 1134f6c57405SHong Zhang if (iascii) { 1135f6c57405SHong Zhang ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 1136cb828f0fSHong Zhang if (format == PETSC_VIEWER_ASCII_INFO || mumps->mumpsview){ 1137cb828f0fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr); 1138cb828f0fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," SYM (matrix type): %d \n",mumps->id.sym);CHKERRQ(ierr); 1139cb828f0fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," PAR (host participation): %d \n",mumps->id.par);CHKERRQ(ierr); 1140cb828f0fSHong Zhang if (mumps->mumpsview){ /* View all MUMPS parameters */ 1141f6c57405SHong Zhang ierr = MatFactorInfo_MUMPS(A,viewer);CHKERRQ(ierr); 1142f6c57405SHong Zhang } 1143f6c57405SHong Zhang } 1144cb828f0fSHong Zhang } 1145f6c57405SHong Zhang PetscFunctionReturn(0); 1146f6c57405SHong Zhang } 1147f6c57405SHong Zhang 114835bd34faSBarry Smith #undef __FUNCT__ 114935bd34faSBarry Smith #define __FUNCT__ "MatGetInfo_MUMPS" 115035bd34faSBarry Smith PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info) 115135bd34faSBarry Smith { 1152cb828f0fSHong Zhang Mat_MUMPS *mumps =(Mat_MUMPS*)A->spptr; 115335bd34faSBarry Smith 115435bd34faSBarry Smith PetscFunctionBegin; 115535bd34faSBarry Smith info->block_size = 1.0; 1156cb828f0fSHong Zhang info->nz_allocated = mumps->id.INFOG(20); 1157cb828f0fSHong Zhang info->nz_used = mumps->id.INFOG(20); 115835bd34faSBarry Smith info->nz_unneeded = 0.0; 115935bd34faSBarry Smith info->assemblies = 0.0; 116035bd34faSBarry Smith info->mallocs = 0.0; 116135bd34faSBarry Smith info->memory = 0.0; 116235bd34faSBarry Smith info->fill_ratio_given = 0; 116335bd34faSBarry Smith info->fill_ratio_needed = 0; 116435bd34faSBarry Smith info->factor_mallocs = 0; 116535bd34faSBarry Smith PetscFunctionReturn(0); 116635bd34faSBarry Smith } 116735bd34faSBarry Smith 116824b6179bSKris Buschelman /*MC 116941c8de11SBarry Smith MAT_SOLVER_MUMPS - A matrix type providing direct solvers (LU and Cholesky) for 117024b6179bSKris Buschelman distributed and sequential matrices via the external package MUMPS. 117124b6179bSKris Buschelman 117241c8de11SBarry Smith Works with MATAIJ and MATSBAIJ matrices 117324b6179bSKris Buschelman 117424b6179bSKris Buschelman Options Database Keys: 117541c8de11SBarry Smith + -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric 117624b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,...,4> - print level 117724b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide) 117824b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide) 117924b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T 118024b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements 118194b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view 118224b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide) 118324b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide) 118424b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide) 118524b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide) 118624b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold 118724b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement 118824b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold 118924b6179bSKris Buschelman 119024b6179bSKris Buschelman Level: beginner 119124b6179bSKris Buschelman 119241c8de11SBarry Smith .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage 119341c8de11SBarry Smith 119424b6179bSKris Buschelman M*/ 119524b6179bSKris Buschelman 11962877fffaSHong Zhang EXTERN_C_BEGIN 119735bd34faSBarry Smith #undef __FUNCT__ 119835bd34faSBarry Smith #define __FUNCT__ "MatFactorGetSolverPackage_mumps" 119935bd34faSBarry Smith PetscErrorCode MatFactorGetSolverPackage_mumps(Mat A,const MatSolverPackage *type) 120035bd34faSBarry Smith { 120135bd34faSBarry Smith PetscFunctionBegin; 120235bd34faSBarry Smith *type = MAT_SOLVER_MUMPS; 120335bd34faSBarry Smith PetscFunctionReturn(0); 120435bd34faSBarry Smith } 120535bd34faSBarry Smith EXTERN_C_END 120635bd34faSBarry Smith 120735bd34faSBarry Smith EXTERN_C_BEGIN 1208bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI AIJ matrices */ 12092877fffaSHong Zhang #undef __FUNCT__ 1210bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_aij_mumps" 1211bccb9932SShri Abhyankar PetscErrorCode MatGetFactor_aij_mumps(Mat A,MatFactorType ftype,Mat *F) 12122877fffaSHong Zhang { 12132877fffaSHong Zhang Mat B; 12142877fffaSHong Zhang PetscErrorCode ierr; 12152877fffaSHong Zhang Mat_MUMPS *mumps; 1216bccb9932SShri Abhyankar PetscTruth isSeqAIJ; 12172877fffaSHong Zhang 12182877fffaSHong Zhang PetscFunctionBegin; 12192877fffaSHong Zhang /* Create the factorization matrix */ 1220bccb9932SShri Abhyankar ierr = PetscTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr); 12212877fffaSHong Zhang ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 12222877fffaSHong Zhang ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 12232877fffaSHong Zhang ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 1224bccb9932SShri Abhyankar if (isSeqAIJ) { 12252877fffaSHong Zhang ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 1226bccb9932SShri Abhyankar } else { 1227bccb9932SShri Abhyankar ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 1228bccb9932SShri Abhyankar } 12292877fffaSHong Zhang 123016ebf90aSShri Abhyankar ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr); 12312877fffaSHong Zhang B->ops->view = MatView_MUMPS; 123235bd34faSBarry Smith B->ops->getinfo = MatGetInfo_MUMPS; 123335bd34faSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr); 123497969023SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMumpsSetIcntl_C","MatMumpsSetIcntl",MatMumpsSetIcntl);CHKERRQ(ierr); 1235450b117fSShri Abhyankar if (ftype == MAT_FACTOR_LU) { 1236450b117fSShri Abhyankar B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS; 1237d5f3da31SBarry Smith B->factortype = MAT_FACTOR_LU; 1238bccb9932SShri Abhyankar if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqaij; 1239bccb9932SShri Abhyankar else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpiaij; 1240dcd589f8SShri Abhyankar } else { 124167877ebaSShri Abhyankar B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS; 1242450b117fSShri Abhyankar B->factortype = MAT_FACTOR_CHOLESKY; 1243bccb9932SShri Abhyankar if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqsbaij; 1244bccb9932SShri Abhyankar else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpisbaij; 1245450b117fSShri Abhyankar } 12462877fffaSHong Zhang 12472877fffaSHong Zhang mumps->CleanUpMUMPS = PETSC_FALSE; 12482877fffaSHong Zhang mumps->isAIJ = PETSC_TRUE; 12492877fffaSHong Zhang mumps->scat_rhs = PETSC_NULL; 12502877fffaSHong Zhang mumps->scat_sol = PETSC_NULL; 12512877fffaSHong Zhang mumps->nSolve = 0; 12522877fffaSHong Zhang mumps->MatDestroy = B->ops->destroy; 12532877fffaSHong Zhang B->ops->destroy = MatDestroy_MUMPS; 12542877fffaSHong Zhang B->spptr = (void*)mumps; 12552877fffaSHong Zhang 12562877fffaSHong Zhang *F = B; 12572877fffaSHong Zhang PetscFunctionReturn(0); 12582877fffaSHong Zhang } 12592877fffaSHong Zhang EXTERN_C_END 12602877fffaSHong Zhang 1261bccb9932SShri Abhyankar 12622877fffaSHong Zhang EXTERN_C_BEGIN 1263bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI SBAIJ matrices */ 12642877fffaSHong Zhang #undef __FUNCT__ 1265bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_sbaij_mumps" 1266bccb9932SShri Abhyankar PetscErrorCode MatGetFactor_sbaij_mumps(Mat A,MatFactorType ftype,Mat *F) 12672877fffaSHong Zhang { 12682877fffaSHong Zhang Mat B; 12692877fffaSHong Zhang PetscErrorCode ierr; 12702877fffaSHong Zhang Mat_MUMPS *mumps; 1271bccb9932SShri Abhyankar PetscTruth isSeqSBAIJ; 12722877fffaSHong Zhang 12732877fffaSHong Zhang PetscFunctionBegin; 1274e7e72b3dSBarry Smith if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix"); 1275bccb9932SShri 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"); 1276bccb9932SShri Abhyankar ierr = PetscTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr); 12772877fffaSHong Zhang /* Create the factorization matrix */ 12782877fffaSHong Zhang ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 12792877fffaSHong Zhang ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 12802877fffaSHong Zhang ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 128116ebf90aSShri Abhyankar ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr); 1282bccb9932SShri Abhyankar if (isSeqSBAIJ) { 1283bccb9932SShri Abhyankar ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr); 128416ebf90aSShri Abhyankar mumps->ConvertToTriples = MatConvertToTriples_seqsbaij_seqsbaij; 1285dcd589f8SShri Abhyankar } else { 1286bccb9932SShri Abhyankar ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 1287bccb9932SShri Abhyankar mumps->ConvertToTriples = MatConvertToTriples_mpisbaij_mpisbaij; 1288bccb9932SShri Abhyankar } 1289bccb9932SShri Abhyankar 129067877ebaSShri Abhyankar B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS; 1291bccb9932SShri Abhyankar B->ops->view = MatView_MUMPS; 1292bccb9932SShri Abhyankar ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr); 1293bccb9932SShri Abhyankar ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMumpsSetIcntl_C","MatMumpsSetIcntl",MatMumpsSetIcntl);CHKERRQ(ierr); 1294f4762488SHong Zhang B->factortype = MAT_FACTOR_CHOLESKY; 1295a214ac2aSShri Abhyankar 1296a214ac2aSShri Abhyankar mumps->CleanUpMUMPS = PETSC_FALSE; 1297bccb9932SShri Abhyankar mumps->isAIJ = PETSC_FALSE; 12982877fffaSHong Zhang mumps->scat_rhs = PETSC_NULL; 12992877fffaSHong Zhang mumps->scat_sol = PETSC_NULL; 13002877fffaSHong Zhang mumps->nSolve = 0; 1301f3c0ef26SHong Zhang mumps->MatDestroy = B->ops->destroy; 1302f3c0ef26SHong Zhang B->ops->destroy = MatDestroy_MUMPS; 13032877fffaSHong Zhang B->spptr = (void*)mumps; 1304f3c0ef26SHong Zhang 13052877fffaSHong Zhang *F = B; 13062877fffaSHong Zhang PetscFunctionReturn(0); 13072877fffaSHong Zhang } 13082877fffaSHong Zhang EXTERN_C_END 130997969023SHong Zhang 1310450b117fSShri Abhyankar EXTERN_C_BEGIN 1311450b117fSShri Abhyankar #undef __FUNCT__ 1312bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_baij_mumps" 1313bccb9932SShri Abhyankar PetscErrorCode MatGetFactor_baij_mumps(Mat A,MatFactorType ftype,Mat *F) 131467877ebaSShri Abhyankar { 131567877ebaSShri Abhyankar Mat B; 131667877ebaSShri Abhyankar PetscErrorCode ierr; 131767877ebaSShri Abhyankar Mat_MUMPS *mumps; 1318bccb9932SShri Abhyankar PetscTruth isSeqBAIJ; 131967877ebaSShri Abhyankar 132067877ebaSShri Abhyankar PetscFunctionBegin; 132167877ebaSShri Abhyankar /* Create the factorization matrix */ 1322bccb9932SShri Abhyankar ierr = PetscTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);CHKERRQ(ierr); 132367877ebaSShri Abhyankar ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 132467877ebaSShri Abhyankar ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 132567877ebaSShri Abhyankar ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 1326bccb9932SShri Abhyankar if (isSeqBAIJ) { 132767877ebaSShri Abhyankar ierr = MatSeqBAIJSetPreallocation(B,A->rmap->bs,0,PETSC_NULL);CHKERRQ(ierr); 1328bccb9932SShri Abhyankar } else { 132967877ebaSShri Abhyankar ierr = MatMPIBAIJSetPreallocation(B,A->rmap->bs,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 1330bccb9932SShri Abhyankar } 1331450b117fSShri Abhyankar 133267877ebaSShri Abhyankar ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr); 1333450b117fSShri Abhyankar if (ftype == MAT_FACTOR_LU) { 1334450b117fSShri Abhyankar B->ops->lufactorsymbolic = MatLUFactorSymbolic_BAIJMUMPS; 1335450b117fSShri Abhyankar B->factortype = MAT_FACTOR_LU; 1336bccb9932SShri Abhyankar if (isSeqBAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqbaij_seqaij; 1337bccb9932SShri Abhyankar else mumps->ConvertToTriples = MatConvertToTriples_mpibaij_mpiaij; 1338450b117fSShri Abhyankar } 1339bccb9932SShri Abhyankar else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc BAIJ matrices with MUMPS Cholesky, use SBAIJ or AIJ matrix instead\n"); 1340bccb9932SShri Abhyankar 1341450b117fSShri Abhyankar B->ops->view = MatView_MUMPS; 1342450b117fSShri Abhyankar ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr); 1343450b117fSShri Abhyankar ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMumpsSetIcntl_C","MatMumpsSetIcntl",MatMumpsSetIcntl);CHKERRQ(ierr); 1344450b117fSShri Abhyankar 1345450b117fSShri Abhyankar mumps->CleanUpMUMPS = PETSC_FALSE; 1346450b117fSShri Abhyankar mumps->isAIJ = PETSC_TRUE; 1347450b117fSShri Abhyankar mumps->scat_rhs = PETSC_NULL; 1348450b117fSShri Abhyankar mumps->scat_sol = PETSC_NULL; 1349450b117fSShri Abhyankar mumps->nSolve = 0; 1350450b117fSShri Abhyankar mumps->MatDestroy = B->ops->destroy; 1351450b117fSShri Abhyankar B->ops->destroy = MatDestroy_MUMPS; 1352450b117fSShri Abhyankar B->spptr = (void*)mumps; 1353450b117fSShri Abhyankar 1354450b117fSShri Abhyankar *F = B; 1355450b117fSShri Abhyankar PetscFunctionReturn(0); 1356450b117fSShri Abhyankar } 1357450b117fSShri Abhyankar EXTERN_C_END 1358a214ac2aSShri Abhyankar 135961288eaaSHong Zhang /* -------------------------------------------------------------------------------------------*/ 136061288eaaSHong Zhang /*@ 136161288eaaSHong Zhang MatMumpsSetIcntl - Set MUMPS parameter ICNTL() 136261288eaaSHong Zhang 136361288eaaSHong Zhang Collective on Mat 136461288eaaSHong Zhang 136561288eaaSHong Zhang Input Parameters: 136661288eaaSHong Zhang + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 136761288eaaSHong Zhang . idx - index of MUMPS parameter array ICNTL() 136861288eaaSHong Zhang - icntl - value of MUMPS ICNTL(imumps) 136961288eaaSHong Zhang 137061288eaaSHong Zhang Options Database: 137161288eaaSHong Zhang . -mat_mumps_icntl_<idx> <icntl> 137261288eaaSHong Zhang 137361288eaaSHong Zhang Level: beginner 137461288eaaSHong Zhang 137561288eaaSHong Zhang References: MUMPS Users' Guide 137661288eaaSHong Zhang 137761288eaaSHong Zhang .seealso: MatGetFactor() 137861288eaaSHong Zhang @*/ 137997969023SHong Zhang #undef __FUNCT__ 138097969023SHong Zhang #define __FUNCT__ "MatMumpsSetIcntl" 138186e5884dSMatthew Knepley PetscErrorCode MatMumpsSetIcntl(Mat F,PetscInt idx,PetscInt icntl) 138297969023SHong Zhang { 1383dcd589f8SShri Abhyankar Mat_MUMPS *lu =(Mat_MUMPS*)(F)->spptr; 138497969023SHong Zhang 138597969023SHong Zhang PetscFunctionBegin; 138661288eaaSHong Zhang lu->id.ICNTL(idx) = icntl; 138797969023SHong Zhang PetscFunctionReturn(0); 138897969023SHong Zhang } 138997969023SHong Zhang 1390