11c2a3de1SBarry Smith 2397b6df1SKris Buschelman /* 3c2b5dc30SHong Zhang Provides an interface to the MUMPS sparse solver 4397b6df1SKris Buschelman */ 551d5961aSHong Zhang 6c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h> /*I "petscmat.h" I*/ 7c6db04a5SJed Brown #include <../src/mat/impls/sbaij/mpi/mpisbaij.h> 8397b6df1SKris Buschelman 9397b6df1SKris Buschelman EXTERN_C_BEGIN 10397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 112907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 122907cef9SHong Zhang #include <cmumps_c.h> 132907cef9SHong Zhang #else 14c6db04a5SJed Brown #include <zmumps_c.h> 152907cef9SHong Zhang #endif 162907cef9SHong Zhang #else 172907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 182907cef9SHong Zhang #include <smumps_c.h> 19397b6df1SKris Buschelman #else 20c6db04a5SJed Brown #include <dmumps_c.h> 21397b6df1SKris Buschelman #endif 222907cef9SHong Zhang #endif 23397b6df1SKris Buschelman EXTERN_C_END 24397b6df1SKris Buschelman #define JOB_INIT -1 253d472b54SHong Zhang #define JOB_FACTSYMBOLIC 1 263d472b54SHong Zhang #define JOB_FACTNUMERIC 2 273d472b54SHong Zhang #define JOB_SOLVE 3 28397b6df1SKris Buschelman #define JOB_END -2 293d472b54SHong Zhang 302907cef9SHong Zhang /* calls to MUMPS */ 312907cef9SHong Zhang #if defined(PETSC_USE_COMPLEX) 322907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 332907cef9SHong Zhang #define PetscMUMPS_c cmumps_c 342907cef9SHong Zhang #else 352907cef9SHong Zhang #define PetscMUMPS_c zmumps_c 362907cef9SHong Zhang #endif 372907cef9SHong Zhang #else 382907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 392907cef9SHong Zhang #define PetscMUMPS_c smumps_c 402907cef9SHong Zhang #else 412907cef9SHong Zhang #define PetscMUMPS_c dmumps_c 422907cef9SHong Zhang #endif 432907cef9SHong Zhang #endif 442907cef9SHong Zhang 453d472b54SHong Zhang 46397b6df1SKris Buschelman /* macros s.t. indices match MUMPS documentation */ 47397b6df1SKris Buschelman #define ICNTL(I) icntl[(I)-1] 48397b6df1SKris Buschelman #define CNTL(I) cntl[(I)-1] 49397b6df1SKris Buschelman #define INFOG(I) infog[(I)-1] 50a7aca84bSHong Zhang #define INFO(I) info[(I)-1] 51397b6df1SKris Buschelman #define RINFOG(I) rinfog[(I)-1] 52adc1d99fSHong Zhang #define RINFO(I) rinfo[(I)-1] 53397b6df1SKris Buschelman 54397b6df1SKris Buschelman typedef struct { 55397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 562907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 572907cef9SHong Zhang CMUMPS_STRUC_C id; 582907cef9SHong Zhang #else 59397b6df1SKris Buschelman ZMUMPS_STRUC_C id; 602907cef9SHong Zhang #endif 612907cef9SHong Zhang #else 622907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 632907cef9SHong Zhang SMUMPS_STRUC_C id; 64397b6df1SKris Buschelman #else 65397b6df1SKris Buschelman DMUMPS_STRUC_C id; 66397b6df1SKris Buschelman #endif 672907cef9SHong Zhang #endif 682907cef9SHong Zhang 69397b6df1SKris Buschelman MatStructure matstruc; 70c1490034SHong Zhang PetscMPIInt myid,size; 71a5e57a09SHong Zhang PetscInt *irn,*jcn,nz,sym; 72397b6df1SKris Buschelman PetscScalar *val; 73397b6df1SKris Buschelman MPI_Comm comm_mumps; 7464e6c443SBarry Smith PetscBool isAIJ,CleanUpMUMPS; 75a5e57a09SHong Zhang PetscInt ICNTL9_pre; /* check if ICNTL(9) is changed from previous MatSolve */ 76801fbe65SHong Zhang VecScatter scat_rhs, scat_sol; /* used by MatSolve() */ 77801fbe65SHong Zhang Vec b_seq,x_seq; 78*b34f08ffSHong Zhang PetscInt ninfo,*info; /* display INFO */ 792205254eSKarl Rupp 80bf0cc555SLisandro Dalcin PetscErrorCode (*Destroy)(Mat); 81bccb9932SShri Abhyankar PetscErrorCode (*ConvertToTriples)(Mat, int, MatReuse, int*, int**, int**, PetscScalar**); 82f0c56d0fSKris Buschelman } Mat_MUMPS; 83f0c56d0fSKris Buschelman 8409573ac7SBarry Smith extern PetscErrorCode MatDuplicate_MUMPS(Mat,MatDuplicateOption,Mat*); 85b24902e0SBarry Smith 86397b6df1SKris Buschelman /* 87d341cd04SHong Zhang MatConvertToTriples_A_B - convert Petsc matrix to triples: row[nz], col[nz], val[nz] 88d341cd04SHong Zhang 89397b6df1SKris Buschelman input: 9067877ebaSShri Abhyankar A - matrix in aij,baij or sbaij (bs=1) format 91397b6df1SKris Buschelman shift - 0: C style output triple; 1: Fortran style output triple. 92bccb9932SShri Abhyankar reuse - MAT_INITIAL_MATRIX: spaces are allocated and values are set for the triple 93bccb9932SShri Abhyankar MAT_REUSE_MATRIX: only the values in v array are updated 94397b6df1SKris Buschelman output: 95397b6df1SKris Buschelman nnz - dim of r, c, and v (number of local nonzero entries of A) 96397b6df1SKris Buschelman r, c, v - row and col index, matrix values (matrix triples) 97eb9baa12SBarry Smith 98eb9baa12SBarry Smith The returned values r, c, and sometimes v are obtained in a single PetscMalloc(). Then in MatDestroy_MUMPS() it is 99eb9baa12SBarry Smith freed with PetscFree((mumps->irn); This is not ideal code, the fact that v is ONLY sometimes part of mumps->irn means 100eb9baa12SBarry Smith that the PetscMalloc() cannot easily be replaced with a PetscMalloc3(). 101eb9baa12SBarry Smith 102397b6df1SKris Buschelman */ 10316ebf90aSShri Abhyankar 10416ebf90aSShri Abhyankar #undef __FUNCT__ 10516ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqaij_seqaij" 106bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 107b24902e0SBarry Smith { 108185f6596SHong Zhang const PetscInt *ai,*aj,*ajj,M=A->rmap->n; 10967877ebaSShri Abhyankar PetscInt nz,rnz,i,j; 110dfbe8321SBarry Smith PetscErrorCode ierr; 111c1490034SHong Zhang PetscInt *row,*col; 11216ebf90aSShri Abhyankar Mat_SeqAIJ *aa=(Mat_SeqAIJ*)A->data; 113397b6df1SKris Buschelman 114397b6df1SKris Buschelman PetscFunctionBegin; 11516ebf90aSShri Abhyankar *v=aa->a; 116bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 1172205254eSKarl Rupp nz = aa->nz; 1182205254eSKarl Rupp ai = aa->i; 1192205254eSKarl Rupp aj = aa->j; 12016ebf90aSShri Abhyankar *nnz = nz; 121785e854fSJed Brown ierr = PetscMalloc1(2*nz, &row);CHKERRQ(ierr); 122185f6596SHong Zhang col = row + nz; 123185f6596SHong Zhang 12416ebf90aSShri Abhyankar nz = 0; 12516ebf90aSShri Abhyankar for (i=0; i<M; i++) { 12616ebf90aSShri Abhyankar rnz = ai[i+1] - ai[i]; 12767877ebaSShri Abhyankar ajj = aj + ai[i]; 12867877ebaSShri Abhyankar for (j=0; j<rnz; j++) { 12967877ebaSShri Abhyankar row[nz] = i+shift; col[nz++] = ajj[j] + shift; 13016ebf90aSShri Abhyankar } 13116ebf90aSShri Abhyankar } 13216ebf90aSShri Abhyankar *r = row; *c = col; 13316ebf90aSShri Abhyankar } 13416ebf90aSShri Abhyankar PetscFunctionReturn(0); 13516ebf90aSShri Abhyankar } 136397b6df1SKris Buschelman 13716ebf90aSShri Abhyankar #undef __FUNCT__ 13867877ebaSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqbaij_seqaij" 139bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqbaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 14067877ebaSShri Abhyankar { 14167877ebaSShri Abhyankar Mat_SeqBAIJ *aa=(Mat_SeqBAIJ*)A->data; 14233d57670SJed Brown const PetscInt *ai,*aj,*ajj,bs2 = aa->bs2; 14333d57670SJed Brown PetscInt bs,M,nz,idx=0,rnz,i,j,k,m; 14467877ebaSShri Abhyankar PetscErrorCode ierr; 14567877ebaSShri Abhyankar PetscInt *row,*col; 14667877ebaSShri Abhyankar 14767877ebaSShri Abhyankar PetscFunctionBegin; 14833d57670SJed Brown ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 14933d57670SJed Brown M = A->rmap->N/bs; 150cf3759fdSShri Abhyankar *v = aa->a; 151bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 152cf3759fdSShri Abhyankar ai = aa->i; aj = aa->j; 15367877ebaSShri Abhyankar nz = bs2*aa->nz; 15467877ebaSShri Abhyankar *nnz = nz; 155785e854fSJed Brown ierr = PetscMalloc1(2*nz, &row);CHKERRQ(ierr); 156185f6596SHong Zhang col = row + nz; 157185f6596SHong Zhang 15867877ebaSShri Abhyankar for (i=0; i<M; i++) { 15967877ebaSShri Abhyankar ajj = aj + ai[i]; 16067877ebaSShri Abhyankar rnz = ai[i+1] - ai[i]; 16167877ebaSShri Abhyankar for (k=0; k<rnz; k++) { 16267877ebaSShri Abhyankar for (j=0; j<bs; j++) { 16367877ebaSShri Abhyankar for (m=0; m<bs; m++) { 16467877ebaSShri Abhyankar row[idx] = i*bs + m + shift; 165cf3759fdSShri Abhyankar col[idx++] = bs*(ajj[k]) + j + shift; 16667877ebaSShri Abhyankar } 16767877ebaSShri Abhyankar } 16867877ebaSShri Abhyankar } 16967877ebaSShri Abhyankar } 170cf3759fdSShri Abhyankar *r = row; *c = col; 17167877ebaSShri Abhyankar } 17267877ebaSShri Abhyankar PetscFunctionReturn(0); 17367877ebaSShri Abhyankar } 17467877ebaSShri Abhyankar 17567877ebaSShri Abhyankar #undef __FUNCT__ 17616ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqsbaij_seqsbaij" 177bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqsbaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 17816ebf90aSShri Abhyankar { 17967877ebaSShri Abhyankar const PetscInt *ai, *aj,*ajj,M=A->rmap->n; 18067877ebaSShri Abhyankar PetscInt nz,rnz,i,j; 18116ebf90aSShri Abhyankar PetscErrorCode ierr; 18216ebf90aSShri Abhyankar PetscInt *row,*col; 18316ebf90aSShri Abhyankar Mat_SeqSBAIJ *aa=(Mat_SeqSBAIJ*)A->data; 18416ebf90aSShri Abhyankar 18516ebf90aSShri Abhyankar PetscFunctionBegin; 186882afa5aSHong Zhang *v = aa->a; 187bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 1882205254eSKarl Rupp nz = aa->nz; 1892205254eSKarl Rupp ai = aa->i; 1902205254eSKarl Rupp aj = aa->j; 1912205254eSKarl Rupp *v = aa->a; 19216ebf90aSShri Abhyankar *nnz = nz; 193785e854fSJed Brown ierr = PetscMalloc1(2*nz, &row);CHKERRQ(ierr); 194185f6596SHong Zhang col = row + nz; 195185f6596SHong Zhang 19616ebf90aSShri Abhyankar nz = 0; 19716ebf90aSShri Abhyankar for (i=0; i<M; i++) { 19816ebf90aSShri Abhyankar rnz = ai[i+1] - ai[i]; 19967877ebaSShri Abhyankar ajj = aj + ai[i]; 20067877ebaSShri Abhyankar for (j=0; j<rnz; j++) { 20167877ebaSShri Abhyankar row[nz] = i+shift; col[nz++] = ajj[j] + shift; 20216ebf90aSShri Abhyankar } 20316ebf90aSShri Abhyankar } 20416ebf90aSShri Abhyankar *r = row; *c = col; 20516ebf90aSShri Abhyankar } 20616ebf90aSShri Abhyankar PetscFunctionReturn(0); 20716ebf90aSShri Abhyankar } 20816ebf90aSShri Abhyankar 20916ebf90aSShri Abhyankar #undef __FUNCT__ 21016ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqaij_seqsbaij" 211bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 21216ebf90aSShri Abhyankar { 21367877ebaSShri Abhyankar const PetscInt *ai,*aj,*ajj,*adiag,M=A->rmap->n; 21467877ebaSShri Abhyankar PetscInt nz,rnz,i,j; 21567877ebaSShri Abhyankar const PetscScalar *av,*v1; 21616ebf90aSShri Abhyankar PetscScalar *val; 21716ebf90aSShri Abhyankar PetscErrorCode ierr; 21816ebf90aSShri Abhyankar PetscInt *row,*col; 219829b1710SHong Zhang Mat_SeqAIJ *aa=(Mat_SeqAIJ*)A->data; 22016ebf90aSShri Abhyankar 22116ebf90aSShri Abhyankar PetscFunctionBegin; 22216ebf90aSShri Abhyankar ai =aa->i; aj=aa->j;av=aa->a; 22316ebf90aSShri Abhyankar adiag=aa->diag; 224bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 225829b1710SHong Zhang /* count nz in the uppper triangular part of A */ 226829b1710SHong Zhang nz = 0; 227829b1710SHong Zhang for (i=0; i<M; i++) nz += ai[i+1] - adiag[i]; 22816ebf90aSShri Abhyankar *nnz = nz; 229829b1710SHong Zhang 230185f6596SHong Zhang ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr); 231185f6596SHong Zhang col = row + nz; 232185f6596SHong Zhang val = (PetscScalar*)(col + nz); 233185f6596SHong Zhang 23416ebf90aSShri Abhyankar nz = 0; 23516ebf90aSShri Abhyankar for (i=0; i<M; i++) { 23616ebf90aSShri Abhyankar rnz = ai[i+1] - adiag[i]; 23767877ebaSShri Abhyankar ajj = aj + adiag[i]; 238cf3759fdSShri Abhyankar v1 = av + adiag[i]; 23967877ebaSShri Abhyankar for (j=0; j<rnz; j++) { 24067877ebaSShri Abhyankar row[nz] = i+shift; col[nz] = ajj[j] + shift; val[nz++] = v1[j]; 24116ebf90aSShri Abhyankar } 24216ebf90aSShri Abhyankar } 24316ebf90aSShri Abhyankar *r = row; *c = col; *v = val; 244397b6df1SKris Buschelman } else { 24516ebf90aSShri Abhyankar nz = 0; val = *v; 24616ebf90aSShri Abhyankar for (i=0; i <M; i++) { 24716ebf90aSShri Abhyankar rnz = ai[i+1] - adiag[i]; 24867877ebaSShri Abhyankar ajj = aj + adiag[i]; 24967877ebaSShri Abhyankar v1 = av + adiag[i]; 25067877ebaSShri Abhyankar for (j=0; j<rnz; j++) { 25167877ebaSShri Abhyankar val[nz++] = v1[j]; 25216ebf90aSShri Abhyankar } 25316ebf90aSShri Abhyankar } 25416ebf90aSShri Abhyankar } 25516ebf90aSShri Abhyankar PetscFunctionReturn(0); 25616ebf90aSShri Abhyankar } 25716ebf90aSShri Abhyankar 25816ebf90aSShri Abhyankar #undef __FUNCT__ 25916ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpisbaij_mpisbaij" 260bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpisbaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 26116ebf90aSShri Abhyankar { 26216ebf90aSShri Abhyankar const PetscInt *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj; 26316ebf90aSShri Abhyankar PetscErrorCode ierr; 26416ebf90aSShri Abhyankar PetscInt rstart,nz,i,j,jj,irow,countA,countB; 26516ebf90aSShri Abhyankar PetscInt *row,*col; 26616ebf90aSShri Abhyankar const PetscScalar *av, *bv,*v1,*v2; 26716ebf90aSShri Abhyankar PetscScalar *val; 268397b6df1SKris Buschelman Mat_MPISBAIJ *mat = (Mat_MPISBAIJ*)A->data; 269397b6df1SKris Buschelman Mat_SeqSBAIJ *aa = (Mat_SeqSBAIJ*)(mat->A)->data; 270397b6df1SKris Buschelman Mat_SeqBAIJ *bb = (Mat_SeqBAIJ*)(mat->B)->data; 27116ebf90aSShri Abhyankar 27216ebf90aSShri Abhyankar PetscFunctionBegin; 273d0f46423SBarry Smith ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart; 274397b6df1SKris Buschelman av=aa->a; bv=bb->a; 275397b6df1SKris Buschelman 2762205254eSKarl Rupp garray = mat->garray; 2772205254eSKarl Rupp 278bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 27916ebf90aSShri Abhyankar nz = aa->nz + bb->nz; 28016ebf90aSShri Abhyankar *nnz = nz; 281185f6596SHong Zhang ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr); 282185f6596SHong Zhang col = row + nz; 283185f6596SHong Zhang val = (PetscScalar*)(col + nz); 284185f6596SHong Zhang 285397b6df1SKris Buschelman *r = row; *c = col; *v = val; 286397b6df1SKris Buschelman } else { 287397b6df1SKris Buschelman row = *r; col = *c; val = *v; 288397b6df1SKris Buschelman } 289397b6df1SKris Buschelman 290028e57e8SHong Zhang jj = 0; irow = rstart; 291397b6df1SKris Buschelman for (i=0; i<m; i++) { 292397b6df1SKris Buschelman ajj = aj + ai[i]; /* ptr to the beginning of this row */ 293397b6df1SKris Buschelman countA = ai[i+1] - ai[i]; 294397b6df1SKris Buschelman countB = bi[i+1] - bi[i]; 295397b6df1SKris Buschelman bjj = bj + bi[i]; 29616ebf90aSShri Abhyankar v1 = av + ai[i]; 29716ebf90aSShri Abhyankar v2 = bv + bi[i]; 298397b6df1SKris Buschelman 299397b6df1SKris Buschelman /* A-part */ 300397b6df1SKris Buschelman for (j=0; j<countA; j++) { 301bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 302397b6df1SKris Buschelman row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift; 303397b6df1SKris Buschelman } 30416ebf90aSShri Abhyankar val[jj++] = v1[j]; 305397b6df1SKris Buschelman } 30616ebf90aSShri Abhyankar 30716ebf90aSShri Abhyankar /* B-part */ 30816ebf90aSShri Abhyankar for (j=0; j < countB; j++) { 309bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 310397b6df1SKris Buschelman row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift; 311397b6df1SKris Buschelman } 31216ebf90aSShri Abhyankar val[jj++] = v2[j]; 31316ebf90aSShri Abhyankar } 31416ebf90aSShri Abhyankar irow++; 31516ebf90aSShri Abhyankar } 31616ebf90aSShri Abhyankar PetscFunctionReturn(0); 31716ebf90aSShri Abhyankar } 31816ebf90aSShri Abhyankar 31916ebf90aSShri Abhyankar #undef __FUNCT__ 32016ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpiaij_mpiaij" 321bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpiaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 32216ebf90aSShri Abhyankar { 32316ebf90aSShri Abhyankar const PetscInt *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj; 32416ebf90aSShri Abhyankar PetscErrorCode ierr; 32516ebf90aSShri Abhyankar PetscInt rstart,nz,i,j,jj,irow,countA,countB; 32616ebf90aSShri Abhyankar PetscInt *row,*col; 32716ebf90aSShri Abhyankar const PetscScalar *av, *bv,*v1,*v2; 32816ebf90aSShri Abhyankar PetscScalar *val; 32916ebf90aSShri Abhyankar Mat_MPIAIJ *mat = (Mat_MPIAIJ*)A->data; 33016ebf90aSShri Abhyankar Mat_SeqAIJ *aa = (Mat_SeqAIJ*)(mat->A)->data; 33116ebf90aSShri Abhyankar Mat_SeqAIJ *bb = (Mat_SeqAIJ*)(mat->B)->data; 33216ebf90aSShri Abhyankar 33316ebf90aSShri Abhyankar PetscFunctionBegin; 33416ebf90aSShri Abhyankar ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart; 33516ebf90aSShri Abhyankar av=aa->a; bv=bb->a; 33616ebf90aSShri Abhyankar 3372205254eSKarl Rupp garray = mat->garray; 3382205254eSKarl Rupp 339bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 34016ebf90aSShri Abhyankar nz = aa->nz + bb->nz; 34116ebf90aSShri Abhyankar *nnz = nz; 342185f6596SHong Zhang ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr); 343185f6596SHong Zhang col = row + nz; 344185f6596SHong Zhang val = (PetscScalar*)(col + nz); 345185f6596SHong Zhang 34616ebf90aSShri Abhyankar *r = row; *c = col; *v = val; 34716ebf90aSShri Abhyankar } else { 34816ebf90aSShri Abhyankar row = *r; col = *c; val = *v; 34916ebf90aSShri Abhyankar } 35016ebf90aSShri Abhyankar 35116ebf90aSShri Abhyankar jj = 0; irow = rstart; 35216ebf90aSShri Abhyankar for (i=0; i<m; i++) { 35316ebf90aSShri Abhyankar ajj = aj + ai[i]; /* ptr to the beginning of this row */ 35416ebf90aSShri Abhyankar countA = ai[i+1] - ai[i]; 35516ebf90aSShri Abhyankar countB = bi[i+1] - bi[i]; 35616ebf90aSShri Abhyankar bjj = bj + bi[i]; 35716ebf90aSShri Abhyankar v1 = av + ai[i]; 35816ebf90aSShri Abhyankar v2 = bv + bi[i]; 35916ebf90aSShri Abhyankar 36016ebf90aSShri Abhyankar /* A-part */ 36116ebf90aSShri Abhyankar for (j=0; j<countA; j++) { 362bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 36316ebf90aSShri Abhyankar row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift; 36416ebf90aSShri Abhyankar } 36516ebf90aSShri Abhyankar val[jj++] = v1[j]; 36616ebf90aSShri Abhyankar } 36716ebf90aSShri Abhyankar 36816ebf90aSShri Abhyankar /* B-part */ 36916ebf90aSShri Abhyankar for (j=0; j < countB; j++) { 370bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 37116ebf90aSShri Abhyankar row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift; 37216ebf90aSShri Abhyankar } 37316ebf90aSShri Abhyankar val[jj++] = v2[j]; 37416ebf90aSShri Abhyankar } 37516ebf90aSShri Abhyankar irow++; 37616ebf90aSShri Abhyankar } 37716ebf90aSShri Abhyankar PetscFunctionReturn(0); 37816ebf90aSShri Abhyankar } 37916ebf90aSShri Abhyankar 38016ebf90aSShri Abhyankar #undef __FUNCT__ 38167877ebaSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpibaij_mpiaij" 382bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpibaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 38367877ebaSShri Abhyankar { 38467877ebaSShri Abhyankar Mat_MPIBAIJ *mat = (Mat_MPIBAIJ*)A->data; 38567877ebaSShri Abhyankar Mat_SeqBAIJ *aa = (Mat_SeqBAIJ*)(mat->A)->data; 38667877ebaSShri Abhyankar Mat_SeqBAIJ *bb = (Mat_SeqBAIJ*)(mat->B)->data; 38767877ebaSShri Abhyankar const PetscInt *ai = aa->i, *bi = bb->i, *aj = aa->j, *bj = bb->j,*ajj, *bjj; 388d985c460SShri Abhyankar const PetscInt *garray = mat->garray,mbs=mat->mbs,rstart=A->rmap->rstart; 38933d57670SJed Brown const PetscInt bs2=mat->bs2; 39067877ebaSShri Abhyankar PetscErrorCode ierr; 39133d57670SJed Brown PetscInt bs,nz,i,j,k,n,jj,irow,countA,countB,idx; 39267877ebaSShri Abhyankar PetscInt *row,*col; 39367877ebaSShri Abhyankar const PetscScalar *av=aa->a, *bv=bb->a,*v1,*v2; 39467877ebaSShri Abhyankar PetscScalar *val; 39567877ebaSShri Abhyankar 39667877ebaSShri Abhyankar PetscFunctionBegin; 39733d57670SJed Brown ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 398bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 39967877ebaSShri Abhyankar nz = bs2*(aa->nz + bb->nz); 40067877ebaSShri Abhyankar *nnz = nz; 401185f6596SHong Zhang ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr); 402185f6596SHong Zhang col = row + nz; 403185f6596SHong Zhang val = (PetscScalar*)(col + nz); 404185f6596SHong Zhang 40567877ebaSShri Abhyankar *r = row; *c = col; *v = val; 40667877ebaSShri Abhyankar } else { 40767877ebaSShri Abhyankar row = *r; col = *c; val = *v; 40867877ebaSShri Abhyankar } 40967877ebaSShri Abhyankar 410d985c460SShri Abhyankar jj = 0; irow = rstart; 41167877ebaSShri Abhyankar for (i=0; i<mbs; i++) { 41267877ebaSShri Abhyankar countA = ai[i+1] - ai[i]; 41367877ebaSShri Abhyankar countB = bi[i+1] - bi[i]; 41467877ebaSShri Abhyankar ajj = aj + ai[i]; 41567877ebaSShri Abhyankar bjj = bj + bi[i]; 41667877ebaSShri Abhyankar v1 = av + bs2*ai[i]; 41767877ebaSShri Abhyankar v2 = bv + bs2*bi[i]; 41867877ebaSShri Abhyankar 41967877ebaSShri Abhyankar idx = 0; 42067877ebaSShri Abhyankar /* A-part */ 42167877ebaSShri Abhyankar for (k=0; k<countA; k++) { 42267877ebaSShri Abhyankar for (j=0; j<bs; j++) { 42367877ebaSShri Abhyankar for (n=0; n<bs; n++) { 424bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 425d985c460SShri Abhyankar row[jj] = irow + n + shift; 426d985c460SShri Abhyankar col[jj] = rstart + bs*ajj[k] + j + shift; 42767877ebaSShri Abhyankar } 42867877ebaSShri Abhyankar val[jj++] = v1[idx++]; 42967877ebaSShri Abhyankar } 43067877ebaSShri Abhyankar } 43167877ebaSShri Abhyankar } 43267877ebaSShri Abhyankar 43367877ebaSShri Abhyankar idx = 0; 43467877ebaSShri Abhyankar /* B-part */ 43567877ebaSShri Abhyankar for (k=0; k<countB; k++) { 43667877ebaSShri Abhyankar for (j=0; j<bs; j++) { 43767877ebaSShri Abhyankar for (n=0; n<bs; n++) { 438bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 439d985c460SShri Abhyankar row[jj] = irow + n + shift; 440d985c460SShri Abhyankar col[jj] = bs*garray[bjj[k]] + j + shift; 44167877ebaSShri Abhyankar } 442d985c460SShri Abhyankar val[jj++] = v2[idx++]; 44367877ebaSShri Abhyankar } 44467877ebaSShri Abhyankar } 44567877ebaSShri Abhyankar } 446d985c460SShri Abhyankar irow += bs; 44767877ebaSShri Abhyankar } 44867877ebaSShri Abhyankar PetscFunctionReturn(0); 44967877ebaSShri Abhyankar } 45067877ebaSShri Abhyankar 45167877ebaSShri Abhyankar #undef __FUNCT__ 45216ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpiaij_mpisbaij" 453bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpiaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 45416ebf90aSShri Abhyankar { 45516ebf90aSShri Abhyankar const PetscInt *ai, *aj,*adiag, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj; 45616ebf90aSShri Abhyankar PetscErrorCode ierr; 457e0bace9bSHong Zhang PetscInt rstart,nz,nza,nzb,i,j,jj,irow,countA,countB; 45816ebf90aSShri Abhyankar PetscInt *row,*col; 45916ebf90aSShri Abhyankar const PetscScalar *av, *bv,*v1,*v2; 46016ebf90aSShri Abhyankar PetscScalar *val; 46116ebf90aSShri Abhyankar Mat_MPIAIJ *mat = (Mat_MPIAIJ*)A->data; 46216ebf90aSShri Abhyankar Mat_SeqAIJ *aa =(Mat_SeqAIJ*)(mat->A)->data; 46316ebf90aSShri Abhyankar Mat_SeqAIJ *bb =(Mat_SeqAIJ*)(mat->B)->data; 46416ebf90aSShri Abhyankar 46516ebf90aSShri Abhyankar PetscFunctionBegin; 46616ebf90aSShri Abhyankar ai=aa->i; aj=aa->j; adiag=aa->diag; 46716ebf90aSShri Abhyankar bi=bb->i; bj=bb->j; garray = mat->garray; 46816ebf90aSShri Abhyankar av=aa->a; bv=bb->a; 4692205254eSKarl Rupp 47016ebf90aSShri Abhyankar rstart = A->rmap->rstart; 47116ebf90aSShri Abhyankar 472bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 473e0bace9bSHong Zhang nza = 0; /* num of upper triangular entries in mat->A, including diagonals */ 474e0bace9bSHong Zhang nzb = 0; /* num of upper triangular entries in mat->B */ 47516ebf90aSShri Abhyankar for (i=0; i<m; i++) { 476e0bace9bSHong Zhang nza += (ai[i+1] - adiag[i]); 47716ebf90aSShri Abhyankar countB = bi[i+1] - bi[i]; 47816ebf90aSShri Abhyankar bjj = bj + bi[i]; 479e0bace9bSHong Zhang for (j=0; j<countB; j++) { 480e0bace9bSHong Zhang if (garray[bjj[j]] > rstart) nzb++; 481e0bace9bSHong Zhang } 482e0bace9bSHong Zhang } 48316ebf90aSShri Abhyankar 484e0bace9bSHong Zhang nz = nza + nzb; /* total nz of upper triangular part of mat */ 48516ebf90aSShri Abhyankar *nnz = nz; 486185f6596SHong Zhang ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr); 487185f6596SHong Zhang col = row + nz; 488185f6596SHong Zhang val = (PetscScalar*)(col + nz); 489185f6596SHong Zhang 49016ebf90aSShri Abhyankar *r = row; *c = col; *v = val; 49116ebf90aSShri Abhyankar } else { 49216ebf90aSShri Abhyankar row = *r; col = *c; val = *v; 49316ebf90aSShri Abhyankar } 49416ebf90aSShri Abhyankar 49516ebf90aSShri Abhyankar jj = 0; irow = rstart; 49616ebf90aSShri Abhyankar for (i=0; i<m; i++) { 49716ebf90aSShri Abhyankar ajj = aj + adiag[i]; /* ptr to the beginning of the diagonal of this row */ 49816ebf90aSShri Abhyankar v1 = av + adiag[i]; 49916ebf90aSShri Abhyankar countA = ai[i+1] - adiag[i]; 50016ebf90aSShri Abhyankar countB = bi[i+1] - bi[i]; 50116ebf90aSShri Abhyankar bjj = bj + bi[i]; 50216ebf90aSShri Abhyankar v2 = bv + bi[i]; 50316ebf90aSShri Abhyankar 50416ebf90aSShri Abhyankar /* A-part */ 50516ebf90aSShri Abhyankar for (j=0; j<countA; j++) { 506bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 50716ebf90aSShri Abhyankar row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift; 50816ebf90aSShri Abhyankar } 50916ebf90aSShri Abhyankar val[jj++] = v1[j]; 51016ebf90aSShri Abhyankar } 51116ebf90aSShri Abhyankar 51216ebf90aSShri Abhyankar /* B-part */ 51316ebf90aSShri Abhyankar for (j=0; j < countB; j++) { 51416ebf90aSShri Abhyankar if (garray[bjj[j]] > rstart) { 515bccb9932SShri Abhyankar if (reuse == MAT_INITIAL_MATRIX) { 51616ebf90aSShri Abhyankar row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift; 51716ebf90aSShri Abhyankar } 51816ebf90aSShri Abhyankar val[jj++] = v2[j]; 51916ebf90aSShri Abhyankar } 520397b6df1SKris Buschelman } 521397b6df1SKris Buschelman irow++; 522397b6df1SKris Buschelman } 523397b6df1SKris Buschelman PetscFunctionReturn(0); 524397b6df1SKris Buschelman } 525397b6df1SKris Buschelman 526397b6df1SKris Buschelman #undef __FUNCT__ 52720be8e61SHong Zhang #define __FUNCT__ "MatGetDiagonal_MUMPS" 52820be8e61SHong Zhang PetscErrorCode MatGetDiagonal_MUMPS(Mat A,Vec v) 52920be8e61SHong Zhang { 53020be8e61SHong Zhang PetscFunctionBegin; 53120be8e61SHong Zhang SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Mat type: MUMPS factor"); 53220be8e61SHong Zhang PetscFunctionReturn(0); 53320be8e61SHong Zhang } 53420be8e61SHong Zhang 53520be8e61SHong Zhang #undef __FUNCT__ 5363924e44cSKris Buschelman #define __FUNCT__ "MatDestroy_MUMPS" 537dfbe8321SBarry Smith PetscErrorCode MatDestroy_MUMPS(Mat A) 538dfbe8321SBarry Smith { 539a5e57a09SHong Zhang Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 540dfbe8321SBarry Smith PetscErrorCode ierr; 541b24902e0SBarry Smith 542397b6df1SKris Buschelman PetscFunctionBegin; 543a5e57a09SHong Zhang if (mumps->CleanUpMUMPS) { 544397b6df1SKris Buschelman /* Terminate instance, deallocate memories */ 545a5e57a09SHong Zhang ierr = PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);CHKERRQ(ierr); 546a5e57a09SHong Zhang ierr = VecScatterDestroy(&mumps->scat_rhs);CHKERRQ(ierr); 547a5e57a09SHong Zhang ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr); 548801fbe65SHong Zhang ierr = VecDestroy(&mumps->b_seq);CHKERRQ(ierr); 549a5e57a09SHong Zhang ierr = VecDestroy(&mumps->x_seq);CHKERRQ(ierr); 550a5e57a09SHong Zhang ierr = PetscFree(mumps->id.perm_in);CHKERRQ(ierr); 551a5e57a09SHong Zhang ierr = PetscFree(mumps->irn);CHKERRQ(ierr); 552*b34f08ffSHong Zhang ierr = PetscFree(mumps->info);CHKERRQ(ierr); 5532205254eSKarl Rupp 554a5e57a09SHong Zhang mumps->id.job = JOB_END; 555a5e57a09SHong Zhang PetscMUMPS_c(&mumps->id); 556a5e57a09SHong Zhang ierr = MPI_Comm_free(&(mumps->comm_mumps));CHKERRQ(ierr); 557397b6df1SKris Buschelman } 558a5e57a09SHong Zhang if (mumps->Destroy) { 559a5e57a09SHong Zhang ierr = (mumps->Destroy)(A);CHKERRQ(ierr); 560bf0cc555SLisandro Dalcin } 561bf0cc555SLisandro Dalcin ierr = PetscFree(A->spptr);CHKERRQ(ierr); 562bf0cc555SLisandro Dalcin 56397969023SHong Zhang /* clear composed functions */ 564bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverPackage_C",NULL);CHKERRQ(ierr); 565bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetIcntl_C",NULL);CHKERRQ(ierr); 566bc6112feSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetIcntl_C",NULL);CHKERRQ(ierr); 567bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetCntl_C",NULL);CHKERRQ(ierr); 568bc6112feSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetCntl_C",NULL);CHKERRQ(ierr); 569bc6112feSHong Zhang 570ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfo_C",NULL);CHKERRQ(ierr); 571ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfog_C",NULL);CHKERRQ(ierr); 572ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfo_C",NULL);CHKERRQ(ierr); 573ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfog_C",NULL);CHKERRQ(ierr); 574397b6df1SKris Buschelman PetscFunctionReturn(0); 575397b6df1SKris Buschelman } 576397b6df1SKris Buschelman 577397b6df1SKris Buschelman #undef __FUNCT__ 578f6c57405SHong Zhang #define __FUNCT__ "MatSolve_MUMPS" 579b24902e0SBarry Smith PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x) 580b24902e0SBarry Smith { 581a5e57a09SHong Zhang Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 582d54de34fSKris Buschelman PetscScalar *array; 58367877ebaSShri Abhyankar Vec b_seq; 584329ec9b3SHong Zhang IS is_iden,is_petsc; 585dfbe8321SBarry Smith PetscErrorCode ierr; 586329ec9b3SHong Zhang PetscInt i; 587883f2eb9SBarry Smith static PetscBool cite1 = PETSC_FALSE,cite2 = PETSC_FALSE; 588397b6df1SKris Buschelman 589397b6df1SKris Buschelman PetscFunctionBegin; 590883f2eb9SBarry Smith ierr = PetscCitationsRegister("@article{MUMPS01,\n author = {P.~R. Amestoy and I.~S. Duff and J.-Y. L'Excellent and J. Koster},\n title = {A fully asynchronous multifrontal solver using distributed dynamic scheduling},\n journal = {SIAM Journal on Matrix Analysis and Applications},\n volume = {23},\n number = {1},\n pages = {15--41},\n year = {2001}\n}\n",&cite1);CHKERRQ(ierr); 591883f2eb9SBarry Smith ierr = PetscCitationsRegister("@article{MUMPS02,\n author = {P.~R. Amestoy and A. Guermouche and J.-Y. L'Excellent and S. Pralet},\n title = {Hybrid scheduling for the parallel solution of linear systems},\n journal = {Parallel Computing},\n volume = {32},\n number = {2},\n pages = {136--156},\n year = {2006}\n}\n",&cite2);CHKERRQ(ierr); 592a5e57a09SHong Zhang mumps->id.nrhs = 1; 593a5e57a09SHong Zhang b_seq = mumps->b_seq; 594a5e57a09SHong Zhang if (mumps->size > 1) { 595329ec9b3SHong Zhang /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */ 596a5e57a09SHong Zhang ierr = VecScatterBegin(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 597a5e57a09SHong Zhang ierr = VecScatterEnd(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 598a5e57a09SHong Zhang if (!mumps->myid) {ierr = VecGetArray(b_seq,&array);CHKERRQ(ierr);} 599397b6df1SKris Buschelman } else { /* size == 1 */ 600397b6df1SKris Buschelman ierr = VecCopy(b,x);CHKERRQ(ierr); 601397b6df1SKris Buschelman ierr = VecGetArray(x,&array);CHKERRQ(ierr); 602397b6df1SKris Buschelman } 603a5e57a09SHong Zhang if (!mumps->myid) { /* define rhs on the host */ 604a5e57a09SHong Zhang mumps->id.nrhs = 1; 605397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 6062907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 607a5e57a09SHong Zhang mumps->id.rhs = (mumps_complex*)array; 6082907cef9SHong Zhang #else 609a5e57a09SHong Zhang mumps->id.rhs = (mumps_double_complex*)array; 6102907cef9SHong Zhang #endif 611397b6df1SKris Buschelman #else 612a5e57a09SHong Zhang mumps->id.rhs = array; 613397b6df1SKris Buschelman #endif 614397b6df1SKris Buschelman } 615397b6df1SKris Buschelman 616397b6df1SKris Buschelman /* solve phase */ 617329ec9b3SHong Zhang /*-------------*/ 618a5e57a09SHong Zhang mumps->id.job = JOB_SOLVE; 619a5e57a09SHong Zhang PetscMUMPS_c(&mumps->id); 620a5e57a09SHong Zhang if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",mumps->id.INFOG(1)); 621397b6df1SKris Buschelman 622a5e57a09SHong Zhang if (mumps->size > 1) { /* convert mumps distributed solution to petsc mpi x */ 623a5e57a09SHong Zhang if (mumps->scat_sol && mumps->ICNTL9_pre != mumps->id.ICNTL(9)) { 624a5e57a09SHong Zhang /* when id.ICNTL(9) changes, the contents of lsol_loc may change (not its size, lsol_loc), recreates scat_sol */ 625a5e57a09SHong Zhang ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr); 626397b6df1SKris Buschelman } 627a5e57a09SHong Zhang if (!mumps->scat_sol) { /* create scatter scat_sol */ 628a5e57a09SHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,mumps->id.lsol_loc,0,1,&is_iden);CHKERRQ(ierr); /* from */ 629a5e57a09SHong Zhang for (i=0; i<mumps->id.lsol_loc; i++) { 630a5e57a09SHong Zhang mumps->id.isol_loc[i] -= 1; /* change Fortran style to C style */ 631a5e57a09SHong Zhang } 632a5e57a09SHong Zhang ierr = ISCreateGeneral(PETSC_COMM_SELF,mumps->id.lsol_loc,mumps->id.isol_loc,PETSC_COPY_VALUES,&is_petsc);CHKERRQ(ierr); /* to */ 633a5e57a09SHong Zhang ierr = VecScatterCreate(mumps->x_seq,is_iden,x,is_petsc,&mumps->scat_sol);CHKERRQ(ierr); 6346bf464f9SBarry Smith ierr = ISDestroy(&is_iden);CHKERRQ(ierr); 6356bf464f9SBarry Smith ierr = ISDestroy(&is_petsc);CHKERRQ(ierr); 6362205254eSKarl Rupp 637a5e57a09SHong Zhang mumps->ICNTL9_pre = mumps->id.ICNTL(9); /* save current value of id.ICNTL(9) */ 638397b6df1SKris Buschelman } 639a5e57a09SHong Zhang 640a5e57a09SHong Zhang ierr = VecScatterBegin(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 641a5e57a09SHong Zhang ierr = VecScatterEnd(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 642329ec9b3SHong Zhang } 643397b6df1SKris Buschelman PetscFunctionReturn(0); 644397b6df1SKris Buschelman } 645397b6df1SKris Buschelman 64651d5961aSHong Zhang #undef __FUNCT__ 64751d5961aSHong Zhang #define __FUNCT__ "MatSolveTranspose_MUMPS" 64851d5961aSHong Zhang PetscErrorCode MatSolveTranspose_MUMPS(Mat A,Vec b,Vec x) 64951d5961aSHong Zhang { 650a5e57a09SHong Zhang Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 65151d5961aSHong Zhang PetscErrorCode ierr; 65251d5961aSHong Zhang 65351d5961aSHong Zhang PetscFunctionBegin; 654a5e57a09SHong Zhang mumps->id.ICNTL(9) = 0; 6550ad0caddSJed Brown ierr = MatSolve_MUMPS(A,b,x);CHKERRQ(ierr); 656a5e57a09SHong Zhang mumps->id.ICNTL(9) = 1; 65751d5961aSHong Zhang PetscFunctionReturn(0); 65851d5961aSHong Zhang } 65951d5961aSHong Zhang 660e0b74bf9SHong Zhang #undef __FUNCT__ 661e0b74bf9SHong Zhang #define __FUNCT__ "MatMatSolve_MUMPS" 662e0b74bf9SHong Zhang PetscErrorCode MatMatSolve_MUMPS(Mat A,Mat B,Mat X) 663e0b74bf9SHong Zhang { 664bda8bf91SBarry Smith PetscErrorCode ierr; 665bda8bf91SBarry Smith PetscBool flg; 6664e34a73bSHong Zhang Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 667334c5f61SHong Zhang PetscInt i,nrhs,M; 6682cd7d884SHong Zhang PetscScalar *array,*bray; 669bda8bf91SBarry Smith 670e0b74bf9SHong Zhang PetscFunctionBegin; 6710298fd71SBarry Smith ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 672801fbe65SHong Zhang if (!flg) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix"); 6730298fd71SBarry Smith ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 674801fbe65SHong Zhang if (!flg) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix"); 675801fbe65SHong Zhang if (B->rmap->n != X->rmap->n) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONG,"Matrix B and X must have same row distribution"); 6764e34a73bSHong Zhang 6772cd7d884SHong Zhang ierr = MatGetSize(B,&M,&nrhs);CHKERRQ(ierr); 678334c5f61SHong Zhang mumps->id.nrhs = nrhs; 679334c5f61SHong Zhang mumps->id.lrhs = M; 6804e34a73bSHong Zhang 6812cd7d884SHong Zhang if (mumps->size == 1) { 6822cd7d884SHong Zhang /* copy B to X */ 6832cd7d884SHong Zhang ierr = MatDenseGetArray(B,&bray);CHKERRQ(ierr); 6842cd7d884SHong Zhang ierr = MatDenseGetArray(X,&array);CHKERRQ(ierr); 6852cd7d884SHong Zhang for (i=0; i<M*nrhs; i++) array[i] = bray[i]; 6862cd7d884SHong Zhang ierr = MatDenseRestoreArray(B,&bray);CHKERRQ(ierr); 6872cd7d884SHong Zhang 6882cd7d884SHong Zhang #if defined(PETSC_USE_COMPLEX) 6892cd7d884SHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 6902cd7d884SHong Zhang mumps->id.rhs = (mumps_complex*)array; 6912cd7d884SHong Zhang #else 6922cd7d884SHong Zhang mumps->id.rhs = (mumps_double_complex*)array; 6932cd7d884SHong Zhang #endif 6942cd7d884SHong Zhang #else 6952cd7d884SHong Zhang mumps->id.rhs = array; 6962cd7d884SHong Zhang #endif 697801fbe65SHong Zhang 6982cd7d884SHong Zhang /* solve phase */ 6992cd7d884SHong Zhang /*-------------*/ 7002cd7d884SHong Zhang mumps->id.job = JOB_SOLVE; 7012cd7d884SHong Zhang PetscMUMPS_c(&mumps->id); 7022cd7d884SHong Zhang if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",mumps->id.INFOG(1)); 7032cd7d884SHong Zhang ierr = MatDenseRestoreArray(X,&array);CHKERRQ(ierr); 704334c5f61SHong Zhang } else { /*--------- parallel case --------*/ 70571aed81dSHong Zhang PetscInt lsol_loc,nlsol_loc,*isol_loc,*idx,*iidx,*idxx,*isol_loc_save; 70671aed81dSHong Zhang PetscScalar *sol_loc,*sol_loc_save; 707801fbe65SHong Zhang IS is_to,is_from; 708334c5f61SHong Zhang PetscInt k,proc,j,m; 709801fbe65SHong Zhang const PetscInt *rstart; 710334c5f61SHong Zhang Vec v_mpi,b_seq,x_seq; 711334c5f61SHong Zhang VecScatter scat_rhs,scat_sol; 712801fbe65SHong Zhang 713801fbe65SHong Zhang /* create x_seq to hold local solution */ 71471aed81dSHong Zhang isol_loc_save = mumps->id.isol_loc; /* save it for MatSovle() */ 71571aed81dSHong Zhang sol_loc_save = mumps->id.sol_loc; 716801fbe65SHong Zhang 71771aed81dSHong Zhang lsol_loc = mumps->id.INFO(23); 71871aed81dSHong Zhang nlsol_loc = nrhs*lsol_loc; /* length of sol_loc */ 71971aed81dSHong Zhang ierr = PetscMalloc2(nlsol_loc,&sol_loc,nlsol_loc,&isol_loc);CHKERRQ(ierr); 720801fbe65SHong Zhang #if defined(PETSC_USE_COMPLEX) 721801fbe65SHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 722801fbe65SHong Zhang mumps->id.sol_loc = (mumps_complex*)sol_loc; 723801fbe65SHong Zhang #else 724801fbe65SHong Zhang mumps->id.sol_loc = (mumps_double_complex*)sol_loc; 725801fbe65SHong Zhang #endif 726801fbe65SHong Zhang #else 727801fbe65SHong Zhang mumps->id.sol_loc = sol_loc; 728801fbe65SHong Zhang #endif 729801fbe65SHong Zhang mumps->id.isol_loc = isol_loc; 730801fbe65SHong Zhang 73171aed81dSHong Zhang ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,nlsol_loc,sol_loc,&x_seq);CHKERRQ(ierr); 7322cd7d884SHong Zhang 73374f0fcc7SHong Zhang /* copy rhs matrix B into vector v_mpi */ 734334c5f61SHong Zhang ierr = MatGetLocalSize(B,&m,NULL);CHKERRQ(ierr); 735801fbe65SHong Zhang ierr = MatDenseGetArray(B,&bray);CHKERRQ(ierr); 73674f0fcc7SHong Zhang ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)B),1,nrhs*m,nrhs*M,(const PetscScalar*)bray,&v_mpi);CHKERRQ(ierr); 737801fbe65SHong Zhang ierr = MatDenseRestoreArray(B,&bray);CHKERRQ(ierr); 738801fbe65SHong Zhang 739334c5f61SHong Zhang /* scatter v_mpi to b_seq because MUMPS only supports centralized rhs */ 74074f0fcc7SHong Zhang /* idx: maps from k-th index of v_mpi to (i,j)-th global entry of B; 741801fbe65SHong Zhang iidx: inverse of idx, will be used by scattering xx_seq -> X */ 742801fbe65SHong Zhang ierr = PetscMalloc2(nrhs*M,&idx,nrhs*M,&iidx);CHKERRQ(ierr); 743801fbe65SHong Zhang ierr = MatGetOwnershipRanges(B,&rstart);CHKERRQ(ierr); 744801fbe65SHong Zhang k = 0; 745801fbe65SHong Zhang for (proc=0; proc<mumps->size; proc++){ 746801fbe65SHong Zhang for (j=0; j<nrhs; j++){ 747801fbe65SHong Zhang for (i=rstart[proc]; i<rstart[proc+1]; i++){ 748801fbe65SHong Zhang iidx[j*M + i] = k; 749801fbe65SHong Zhang idx[k++] = j*M + i; 750801fbe65SHong Zhang } 751801fbe65SHong Zhang } 7522cd7d884SHong Zhang } 7532cd7d884SHong Zhang 754801fbe65SHong Zhang if (!mumps->myid) { 755334c5f61SHong Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,nrhs*M,&b_seq);CHKERRQ(ierr); 756801fbe65SHong Zhang ierr = ISCreateGeneral(PETSC_COMM_SELF,nrhs*M,idx,PETSC_COPY_VALUES,&is_to);CHKERRQ(ierr); 757801fbe65SHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,nrhs*M,0,1,&is_from);CHKERRQ(ierr); 758801fbe65SHong Zhang } else { 759334c5f61SHong Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,0,&b_seq);CHKERRQ(ierr); 760801fbe65SHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_to);CHKERRQ(ierr); 761801fbe65SHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_from);CHKERRQ(ierr); 762801fbe65SHong Zhang } 763334c5f61SHong Zhang ierr = VecScatterCreate(v_mpi,is_from,b_seq,is_to,&scat_rhs);CHKERRQ(ierr); 764334c5f61SHong Zhang ierr = VecScatterBegin(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 765801fbe65SHong Zhang ierr = ISDestroy(&is_to);CHKERRQ(ierr); 766801fbe65SHong Zhang ierr = ISDestroy(&is_from);CHKERRQ(ierr); 767334c5f61SHong Zhang ierr = VecScatterEnd(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 768801fbe65SHong Zhang 769801fbe65SHong Zhang if (!mumps->myid) { /* define rhs on the host */ 770334c5f61SHong Zhang ierr = VecGetArray(b_seq,&bray);CHKERRQ(ierr); 771801fbe65SHong Zhang #if defined(PETSC_USE_COMPLEX) 772801fbe65SHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 773801fbe65SHong Zhang mumps->id.rhs = (mumps_complex*)bray; 774801fbe65SHong Zhang #else 775801fbe65SHong Zhang mumps->id.rhs = (mumps_double_complex*)bray; 776801fbe65SHong Zhang #endif 777801fbe65SHong Zhang #else 778801fbe65SHong Zhang mumps->id.rhs = bray; 779801fbe65SHong Zhang #endif 780334c5f61SHong Zhang ierr = VecRestoreArray(b_seq,&bray);CHKERRQ(ierr); 781801fbe65SHong Zhang } 782801fbe65SHong Zhang 783801fbe65SHong Zhang /* solve phase */ 784801fbe65SHong Zhang /*-------------*/ 785801fbe65SHong Zhang mumps->id.job = JOB_SOLVE; 786801fbe65SHong Zhang PetscMUMPS_c(&mumps->id); 787801fbe65SHong Zhang if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",mumps->id.INFOG(1)); 788801fbe65SHong Zhang 789334c5f61SHong Zhang /* scatter mumps distributed solution to petsc vector v_mpi, which shares local arrays with solution matrix X */ 79074f0fcc7SHong Zhang ierr = MatDenseGetArray(X,&array);CHKERRQ(ierr); 79174f0fcc7SHong Zhang ierr = VecPlaceArray(v_mpi,array);CHKERRQ(ierr); 792801fbe65SHong Zhang 793334c5f61SHong Zhang /* create scatter scat_sol */ 79471aed81dSHong Zhang ierr = PetscMalloc1(nlsol_loc,&idxx);CHKERRQ(ierr); 79571aed81dSHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,nlsol_loc,0,1,&is_from);CHKERRQ(ierr); 79671aed81dSHong Zhang for (i=0; i<lsol_loc; i++) { 797334c5f61SHong Zhang isol_loc[i] -= 1; /* change Fortran style to C style */ 798334c5f61SHong Zhang idxx[i] = iidx[isol_loc[i]]; 799801fbe65SHong Zhang for (j=1; j<nrhs; j++){ 800334c5f61SHong Zhang idxx[j*lsol_loc+i] = iidx[isol_loc[i]+j*M]; 801801fbe65SHong Zhang } 802801fbe65SHong Zhang } 80371aed81dSHong Zhang ierr = ISCreateGeneral(PETSC_COMM_SELF,nlsol_loc,idxx,PETSC_COPY_VALUES,&is_to);CHKERRQ(ierr); 804334c5f61SHong Zhang ierr = VecScatterCreate(x_seq,is_from,v_mpi,is_to,&scat_sol);CHKERRQ(ierr); 805334c5f61SHong Zhang ierr = VecScatterBegin(scat_sol,x_seq,v_mpi,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 806801fbe65SHong Zhang ierr = ISDestroy(&is_from);CHKERRQ(ierr); 807801fbe65SHong Zhang ierr = ISDestroy(&is_to);CHKERRQ(ierr); 808334c5f61SHong Zhang ierr = VecScatterEnd(scat_sol,x_seq,v_mpi,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 809801fbe65SHong Zhang ierr = MatDenseRestoreArray(X,&array);CHKERRQ(ierr); 81071aed81dSHong Zhang 81171aed81dSHong Zhang /* free spaces */ 81271aed81dSHong Zhang mumps->id.sol_loc = sol_loc_save; 81371aed81dSHong Zhang mumps->id.isol_loc = isol_loc_save; 81471aed81dSHong Zhang 81571aed81dSHong Zhang ierr = PetscFree2(sol_loc,isol_loc);CHKERRQ(ierr); 816801fbe65SHong Zhang ierr = PetscFree2(idx,iidx);CHKERRQ(ierr); 817801fbe65SHong Zhang ierr = PetscFree(idxx);CHKERRQ(ierr); 81871aed81dSHong Zhang ierr = VecDestroy(&x_seq);CHKERRQ(ierr); 81974f0fcc7SHong Zhang ierr = VecDestroy(&v_mpi);CHKERRQ(ierr); 820334c5f61SHong Zhang ierr = VecDestroy(&b_seq);CHKERRQ(ierr); 821334c5f61SHong Zhang ierr = VecScatterDestroy(&scat_rhs);CHKERRQ(ierr); 822334c5f61SHong Zhang ierr = VecScatterDestroy(&scat_sol);CHKERRQ(ierr); 823801fbe65SHong Zhang } 824e0b74bf9SHong Zhang PetscFunctionReturn(0); 825e0b74bf9SHong Zhang } 826e0b74bf9SHong Zhang 827ace3df97SHong Zhang #if !defined(PETSC_USE_COMPLEX) 828a58c3f20SHong Zhang /* 829a58c3f20SHong Zhang input: 830a58c3f20SHong Zhang F: numeric factor 831a58c3f20SHong Zhang output: 832a58c3f20SHong Zhang nneg: total number of negative pivots 833a58c3f20SHong Zhang nzero: 0 834a58c3f20SHong Zhang npos: (global dimension of F) - nneg 835a58c3f20SHong Zhang */ 836a58c3f20SHong Zhang 837a58c3f20SHong Zhang #undef __FUNCT__ 838a58c3f20SHong Zhang #define __FUNCT__ "MatGetInertia_SBAIJMUMPS" 839dfbe8321SBarry Smith PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos) 840a58c3f20SHong Zhang { 841a5e57a09SHong Zhang Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 842dfbe8321SBarry Smith PetscErrorCode ierr; 843c1490034SHong Zhang PetscMPIInt size; 844a58c3f20SHong Zhang 845a58c3f20SHong Zhang PetscFunctionBegin; 846ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)F),&size);CHKERRQ(ierr); 847bcb30aebSHong 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 */ 848a5e57a09SHong Zhang if (size > 1 && mumps->id.ICNTL(13) != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"ICNTL(13)=%d. -mat_mumps_icntl_13 must be set as 1 for correct global matrix inertia\n",mumps->id.INFOG(13)); 849ed85ac9fSHong Zhang 850710ac8efSHong Zhang if (nneg) *nneg = mumps->id.INFOG(12); 851ed85ac9fSHong Zhang if (nzero || npos) { 852ed85ac9fSHong Zhang if (mumps->id.ICNTL(24) != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"-mat_mumps_icntl_24 must be set as 1 for null pivot row detection"); 853710ac8efSHong Zhang if (nzero) *nzero = mumps->id.INFOG(28); 854710ac8efSHong Zhang if (npos) *npos = F->rmap->N - (mumps->id.INFOG(12) + mumps->id.INFOG(28)); 855a58c3f20SHong Zhang } 856a58c3f20SHong Zhang PetscFunctionReturn(0); 857a58c3f20SHong Zhang } 858ace3df97SHong Zhang #endif /* !defined(PETSC_USE_COMPLEX) */ 859a58c3f20SHong Zhang 860397b6df1SKris Buschelman #undef __FUNCT__ 861f6c57405SHong Zhang #define __FUNCT__ "MatFactorNumeric_MUMPS" 8620481f469SBarry Smith PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info) 863af281ebdSHong Zhang { 864a5e57a09SHong Zhang Mat_MUMPS *mumps =(Mat_MUMPS*)(F)->spptr; 8656849ba73SBarry Smith PetscErrorCode ierr; 866e09efc27SHong Zhang Mat F_diag; 867ace3abfcSBarry Smith PetscBool isMPIAIJ; 868397b6df1SKris Buschelman 869397b6df1SKris Buschelman PetscFunctionBegin; 870a5e57a09SHong Zhang ierr = (*mumps->ConvertToTriples)(A, 1, MAT_REUSE_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr); 871397b6df1SKris Buschelman 872397b6df1SKris Buschelman /* numerical factorization phase */ 873329ec9b3SHong Zhang /*-------------------------------*/ 874a5e57a09SHong Zhang mumps->id.job = JOB_FACTNUMERIC; 8754e34a73bSHong Zhang if (!mumps->id.ICNTL(18)) { /* A is centralized */ 876a5e57a09SHong Zhang if (!mumps->myid) { 877397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 8782907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 879a5e57a09SHong Zhang mumps->id.a = (mumps_complex*)mumps->val; 8802907cef9SHong Zhang #else 881a5e57a09SHong Zhang mumps->id.a = (mumps_double_complex*)mumps->val; 8822907cef9SHong Zhang #endif 883397b6df1SKris Buschelman #else 884a5e57a09SHong Zhang mumps->id.a = mumps->val; 885397b6df1SKris Buschelman #endif 886397b6df1SKris Buschelman } 887397b6df1SKris Buschelman } else { 888397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 8892907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 890a5e57a09SHong Zhang mumps->id.a_loc = (mumps_complex*)mumps->val; 8912907cef9SHong Zhang #else 892a5e57a09SHong Zhang mumps->id.a_loc = (mumps_double_complex*)mumps->val; 8932907cef9SHong Zhang #endif 894397b6df1SKris Buschelman #else 895a5e57a09SHong Zhang mumps->id.a_loc = mumps->val; 896397b6df1SKris Buschelman #endif 897397b6df1SKris Buschelman } 898a5e57a09SHong Zhang PetscMUMPS_c(&mumps->id); 899a5e57a09SHong Zhang if (mumps->id.INFOG(1) < 0) { 900151787a6SHong Zhang if (mumps->id.INFO(1) == -13) { 901151787a6SHong Zhang if (mumps->id.INFO(2) < 0) { 902151787a6SHong Zhang SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: Cannot allocate required memory %d megabytes\n",-mumps->id.INFO(2)); 903151787a6SHong Zhang } else { 904151787a6SHong Zhang SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: Cannot allocate required memory %d bytes\n",mumps->id.INFO(2)); 905151787a6SHong Zhang } 906151787a6SHong Zhang } else SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: INFO(1)=%d, INFO(2)=%d\n",mumps->id.INFO(1),mumps->id.INFO(2)); 907397b6df1SKris Buschelman } 908a5e57a09SHong Zhang if (!mumps->myid && mumps->id.ICNTL(16) > 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB," mumps->id.ICNTL(16):=%d\n",mumps->id.INFOG(16)); 909397b6df1SKris Buschelman 910dcd589f8SShri Abhyankar (F)->assembled = PETSC_TRUE; 911a5e57a09SHong Zhang mumps->matstruc = SAME_NONZERO_PATTERN; 912a5e57a09SHong Zhang mumps->CleanUpMUMPS = PETSC_TRUE; 91367877ebaSShri Abhyankar 914a5e57a09SHong Zhang if (mumps->size > 1) { 91567877ebaSShri Abhyankar PetscInt lsol_loc; 91667877ebaSShri Abhyankar PetscScalar *sol_loc; 9172205254eSKarl Rupp 918c2093ab7SHong Zhang ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr); 919c2093ab7SHong Zhang if (isMPIAIJ) F_diag = ((Mat_MPIAIJ*)(F)->data)->A; 920c2093ab7SHong Zhang else F_diag = ((Mat_MPISBAIJ*)(F)->data)->A; 921c2093ab7SHong Zhang F_diag->assembled = PETSC_TRUE; 922c2093ab7SHong Zhang 923c2093ab7SHong Zhang /* distributed solution; Create x_seq=sol_loc for repeated use */ 924c2093ab7SHong Zhang if (mumps->x_seq) { 925c2093ab7SHong Zhang ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr); 926c2093ab7SHong Zhang ierr = PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);CHKERRQ(ierr); 927c2093ab7SHong Zhang ierr = VecDestroy(&mumps->x_seq);CHKERRQ(ierr); 928c2093ab7SHong Zhang } 929a5e57a09SHong Zhang lsol_loc = mumps->id.INFO(23); /* length of sol_loc */ 930dcca6d9dSJed Brown ierr = PetscMalloc2(lsol_loc,&sol_loc,lsol_loc,&mumps->id.isol_loc);CHKERRQ(ierr); 931a5e57a09SHong Zhang mumps->id.lsol_loc = lsol_loc; 93267877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 9332907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 934a5e57a09SHong Zhang mumps->id.sol_loc = (mumps_complex*)sol_loc; 9352907cef9SHong Zhang #else 936a5e57a09SHong Zhang mumps->id.sol_loc = (mumps_double_complex*)sol_loc; 9372907cef9SHong Zhang #endif 93867877ebaSShri Abhyankar #else 939a5e57a09SHong Zhang mumps->id.sol_loc = sol_loc; 94067877ebaSShri Abhyankar #endif 941a5e57a09SHong Zhang ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,lsol_loc,sol_loc,&mumps->x_seq);CHKERRQ(ierr); 94267877ebaSShri Abhyankar } 943397b6df1SKris Buschelman PetscFunctionReturn(0); 944397b6df1SKris Buschelman } 945397b6df1SKris Buschelman 9469a2535b5SHong Zhang /* Sets MUMPS options from the options database */ 947dcd589f8SShri Abhyankar #undef __FUNCT__ 9489a2535b5SHong Zhang #define __FUNCT__ "PetscSetMUMPSFromOptions" 9499a2535b5SHong Zhang PetscErrorCode PetscSetMUMPSFromOptions(Mat F, Mat A) 950dcd589f8SShri Abhyankar { 9519a2535b5SHong Zhang Mat_MUMPS *mumps = (Mat_MUMPS*)F->spptr; 952dcd589f8SShri Abhyankar PetscErrorCode ierr; 953*b34f08ffSHong Zhang PetscInt icntl,info[40],i,ninfo=40; 954ace3abfcSBarry Smith PetscBool flg; 955dcd589f8SShri Abhyankar 956dcd589f8SShri Abhyankar PetscFunctionBegin; 957ce94432eSBarry Smith ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MUMPS Options","Mat");CHKERRQ(ierr); 9589a2535b5SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_1","ICNTL(1): output stream for error messages","None",mumps->id.ICNTL(1),&icntl,&flg);CHKERRQ(ierr); 9599a2535b5SHong Zhang if (flg) mumps->id.ICNTL(1) = icntl; 9609a2535b5SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_2","ICNTL(2): output stream for diagnostic printing, statistics, and warning","None",mumps->id.ICNTL(2),&icntl,&flg);CHKERRQ(ierr); 9619a2535b5SHong Zhang if (flg) mumps->id.ICNTL(2) = icntl; 9629a2535b5SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_3","ICNTL(3): output stream for global information, collected on the host","None",mumps->id.ICNTL(3),&icntl,&flg);CHKERRQ(ierr); 9639a2535b5SHong Zhang if (flg) mumps->id.ICNTL(3) = icntl; 964dcd589f8SShri Abhyankar 9659a2535b5SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",mumps->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr); 9669a2535b5SHong Zhang if (flg) mumps->id.ICNTL(4) = icntl; 9679a2535b5SHong Zhang if (mumps->id.ICNTL(4) || PetscLogPrintInfo) mumps->id.ICNTL(3) = 6; /* resume MUMPS default id.ICNTL(3) = 6 */ 9689a2535b5SHong Zhang 969d341cd04SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_6","ICNTL(6): permutes to a zero-free diagonal and/or scale the matrix (0 to 7)","None",mumps->id.ICNTL(6),&icntl,&flg);CHKERRQ(ierr); 9709a2535b5SHong Zhang if (flg) mumps->id.ICNTL(6) = icntl; 9719a2535b5SHong Zhang 972d341cd04SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): computes a symmetric permutation in sequential analysis (0 to 7). 3=Scotch, 4=PORD, 5=Metis","None",mumps->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr); 973dcd589f8SShri Abhyankar if (flg) { 9742205254eSKarl Rupp if (icntl== 1 && mumps->size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"pivot order be set by the user in PERM_IN -- not supported by the PETSc/MUMPS interface\n"); 9752205254eSKarl Rupp else mumps->id.ICNTL(7) = icntl; 976dcd589f8SShri Abhyankar } 977e0b74bf9SHong Zhang 9780298fd71SBarry Smith ierr = PetscOptionsInt("-mat_mumps_icntl_8","ICNTL(8): scaling strategy (-2 to 8 or 77)","None",mumps->id.ICNTL(8),&mumps->id.ICNTL(8),NULL);CHKERRQ(ierr); 979d341cd04SHong Zhang /* ierr = PetscOptionsInt("-mat_mumps_icntl_9","ICNTL(9): computes the solution using A or A^T","None",mumps->id.ICNTL(9),&mumps->id.ICNTL(9),NULL);CHKERRQ(ierr); handled by MatSolveTranspose_MUMPS() */ 9800298fd71SBarry Smith ierr = PetscOptionsInt("-mat_mumps_icntl_10","ICNTL(10): max num of refinements","None",mumps->id.ICNTL(10),&mumps->id.ICNTL(10),NULL);CHKERRQ(ierr); 981d341cd04SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_11","ICNTL(11): statistics related to an error analysis (via -ksp_view)","None",mumps->id.ICNTL(11),&mumps->id.ICNTL(11),NULL);CHKERRQ(ierr); 982d341cd04SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): an ordering strategy for symmetric matrices (0 to 3)","None",mumps->id.ICNTL(12),&mumps->id.ICNTL(12),NULL);CHKERRQ(ierr); 983d341cd04SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): parallelism of the root node (enable ScaLAPACK) and its splitting","None",mumps->id.ICNTL(13),&mumps->id.ICNTL(13),NULL);CHKERRQ(ierr); 984d341cd04SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_14","ICNTL(14): percentage increase in the estimated working space","None",mumps->id.ICNTL(14),&mumps->id.ICNTL(14),NULL);CHKERRQ(ierr); 985d341cd04SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_19","ICNTL(19): computes the Schur complement","None",mumps->id.ICNTL(19),&mumps->id.ICNTL(19),NULL);CHKERRQ(ierr); 9864e34a73bSHong Zhang /* ierr = PetscOptionsInt("-mat_mumps_icntl_20","ICNTL(20): the format (dense or sparse) of the right-hand sides","None",mumps->id.ICNTL(20),&mumps->id.ICNTL(20),NULL);CHKERRQ(ierr); -- sparse rhs is not supported in PETSc API */ 987d341cd04SHong Zhang /* ierr = PetscOptionsInt("-mat_mumps_icntl_21","ICNTL(21): the distribution (centralized or distributed) of the solution vectors","None",mumps->id.ICNTL(21),&mumps->id.ICNTL(21),NULL);CHKERRQ(ierr); we only use distributed solution vector */ 9889a2535b5SHong Zhang 989d341cd04SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_22","ICNTL(22): in-core/out-of-core factorization and solve (0 or 1)","None",mumps->id.ICNTL(22),&mumps->id.ICNTL(22),NULL);CHKERRQ(ierr); 9900298fd71SBarry Smith ierr = PetscOptionsInt("-mat_mumps_icntl_23","ICNTL(23): max size of the working memory (MB) that can allocate per processor","None",mumps->id.ICNTL(23),&mumps->id.ICNTL(23),NULL);CHKERRQ(ierr); 9910298fd71SBarry Smith ierr = PetscOptionsInt("-mat_mumps_icntl_24","ICNTL(24): detection of null pivot rows (0 or 1)","None",mumps->id.ICNTL(24),&mumps->id.ICNTL(24),NULL);CHKERRQ(ierr); 9929a2535b5SHong Zhang if (mumps->id.ICNTL(24)) { 9939a2535b5SHong Zhang mumps->id.ICNTL(13) = 1; /* turn-off ScaLAPACK to help with the correct detection of null pivots */ 994d7ebd59bSHong Zhang } 995d7ebd59bSHong Zhang 996d341cd04SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_25","ICNTL(25): compute a solution of a deficient matrix and a null space basis","None",mumps->id.ICNTL(25),&mumps->id.ICNTL(25),NULL);CHKERRQ(ierr); 997d341cd04SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_26","ICNTL(26): drives the solution phase if a Schur complement matrix","None",mumps->id.ICNTL(26),&mumps->id.ICNTL(26),NULL);CHKERRQ(ierr); 9982cd7d884SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_27","ICNTL(27): the blocking size for multiple right-hand sides","None",mumps->id.ICNTL(27),&mumps->id.ICNTL(27),NULL);CHKERRQ(ierr); 9990298fd71SBarry Smith ierr = PetscOptionsInt("-mat_mumps_icntl_28","ICNTL(28): use 1 for sequential analysis and ictnl(7) ordering, or 2 for parallel analysis and ictnl(29) ordering","None",mumps->id.ICNTL(28),&mumps->id.ICNTL(28),NULL);CHKERRQ(ierr); 1000d341cd04SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_29","ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis","None",mumps->id.ICNTL(29),&mumps->id.ICNTL(29),NULL);CHKERRQ(ierr); 10010298fd71SBarry Smith ierr = PetscOptionsInt("-mat_mumps_icntl_30","ICNTL(30): compute user-specified set of entries in inv(A)","None",mumps->id.ICNTL(30),&mumps->id.ICNTL(30),NULL);CHKERRQ(ierr); 1002d341cd04SHong Zhang ierr = PetscOptionsInt("-mat_mumps_icntl_31","ICNTL(31): indicates which factors may be discarded during factorization","None",mumps->id.ICNTL(31),&mumps->id.ICNTL(31),NULL);CHKERRQ(ierr); 10034e34a73bSHong Zhang /* ierr = PetscOptionsInt("-mat_mumps_icntl_32","ICNTL(32): performs the forward elemination of the right-hand sides during factorization","None",mumps->id.ICNTL(32),&mumps->id.ICNTL(32),NULL);CHKERRQ(ierr); -- not supported by PETSc API */ 10040298fd71SBarry Smith ierr = PetscOptionsInt("-mat_mumps_icntl_33","ICNTL(33): compute determinant","None",mumps->id.ICNTL(33),&mumps->id.ICNTL(33),NULL);CHKERRQ(ierr); 1005dcd589f8SShri Abhyankar 10060298fd71SBarry Smith ierr = PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",mumps->id.CNTL(1),&mumps->id.CNTL(1),NULL);CHKERRQ(ierr); 10070298fd71SBarry Smith ierr = PetscOptionsReal("-mat_mumps_cntl_2","CNTL(2): stopping criterion of refinement","None",mumps->id.CNTL(2),&mumps->id.CNTL(2),NULL);CHKERRQ(ierr); 10080298fd71SBarry Smith ierr = PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",mumps->id.CNTL(3),&mumps->id.CNTL(3),NULL);CHKERRQ(ierr); 10090298fd71SBarry Smith ierr = PetscOptionsReal("-mat_mumps_cntl_4","CNTL(4): value for static pivoting","None",mumps->id.CNTL(4),&mumps->id.CNTL(4),NULL);CHKERRQ(ierr); 10100298fd71SBarry Smith ierr = PetscOptionsReal("-mat_mumps_cntl_5","CNTL(5): fixation for null pivots","None",mumps->id.CNTL(5),&mumps->id.CNTL(5),NULL);CHKERRQ(ierr); 1011e5bb22a1SHong Zhang 10120298fd71SBarry Smith ierr = PetscOptionsString("-mat_mumps_ooc_tmpdir", "out of core directory", "None", mumps->id.ooc_tmpdir, mumps->id.ooc_tmpdir, 256, NULL); 1013*b34f08ffSHong Zhang 1014*b34f08ffSHong Zhang ierr = PetscOptionsGetIntArray(NULL,"-mat_mumps_view_info",info,&ninfo,NULL);CHKERRQ(ierr); /* why does not show with '-help'??? */ 1015*b34f08ffSHong Zhang if (ninfo) { 1016*b34f08ffSHong Zhang if (ninfo > 40) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"number of INFO %d must <= 40\n",ninfo); 1017*b34f08ffSHong Zhang ierr = PetscMalloc1(ninfo,&mumps->info);CHKERRQ(ierr); 1018*b34f08ffSHong Zhang mumps->ninfo = ninfo; 1019*b34f08ffSHong Zhang for (i=0; i<ninfo; i++) { 1020*b34f08ffSHong Zhang if (info[i] < 0 || info[i]>40) { 1021*b34f08ffSHong Zhang SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index of INFO %d must between 1 and 40\n",ninfo); 1022*b34f08ffSHong Zhang } else { 1023*b34f08ffSHong Zhang mumps->info[i] = info[i]; 1024*b34f08ffSHong Zhang } 1025*b34f08ffSHong Zhang } 1026*b34f08ffSHong Zhang } 1027*b34f08ffSHong Zhang 1028dcd589f8SShri Abhyankar PetscOptionsEnd(); 1029dcd589f8SShri Abhyankar PetscFunctionReturn(0); 1030dcd589f8SShri Abhyankar } 1031dcd589f8SShri Abhyankar 1032dcd589f8SShri Abhyankar #undef __FUNCT__ 1033dcd589f8SShri Abhyankar #define __FUNCT__ "PetscInitializeMUMPS" 1034f697e70eSHong Zhang PetscErrorCode PetscInitializeMUMPS(Mat A,Mat_MUMPS *mumps) 1035dcd589f8SShri Abhyankar { 1036dcd589f8SShri Abhyankar PetscErrorCode ierr; 1037dcd589f8SShri Abhyankar 1038dcd589f8SShri Abhyankar PetscFunctionBegin; 1039ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)A), &mumps->myid); 1040ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&mumps->size);CHKERRQ(ierr); 1041ce94432eSBarry Smith ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)A),&(mumps->comm_mumps));CHKERRQ(ierr); 10422205254eSKarl Rupp 1043f697e70eSHong Zhang mumps->id.comm_fortran = MPI_Comm_c2f(mumps->comm_mumps); 1044f697e70eSHong Zhang 1045f697e70eSHong Zhang mumps->id.job = JOB_INIT; 1046f697e70eSHong Zhang mumps->id.par = 1; /* host participates factorizaton and solve */ 1047f697e70eSHong Zhang mumps->id.sym = mumps->sym; 10482907cef9SHong Zhang PetscMUMPS_c(&mumps->id); 1049f697e70eSHong Zhang 1050f697e70eSHong Zhang mumps->CleanUpMUMPS = PETSC_FALSE; 10510298fd71SBarry Smith mumps->scat_rhs = NULL; 10520298fd71SBarry Smith mumps->scat_sol = NULL; 10539a2535b5SHong Zhang 105470544d5fSHong Zhang /* set PETSc-MUMPS default options - override MUMPS default */ 10559a2535b5SHong Zhang mumps->id.ICNTL(3) = 0; 10569a2535b5SHong Zhang mumps->id.ICNTL(4) = 0; 10579a2535b5SHong Zhang if (mumps->size == 1) { 10589a2535b5SHong Zhang mumps->id.ICNTL(18) = 0; /* centralized assembled matrix input */ 10599a2535b5SHong Zhang } else { 10609a2535b5SHong Zhang mumps->id.ICNTL(18) = 3; /* distributed assembled matrix input */ 10614e34a73bSHong Zhang mumps->id.ICNTL(20) = 0; /* rhs is in dense format */ 106270544d5fSHong Zhang mumps->id.ICNTL(21) = 1; /* distributed solution */ 10639a2535b5SHong Zhang } 1064dcd589f8SShri Abhyankar PetscFunctionReturn(0); 1065dcd589f8SShri Abhyankar } 1066dcd589f8SShri Abhyankar 1067a5e57a09SHong Zhang /* Note Petsc r(=c) permutation is used when mumps->id.ICNTL(7)==1 with centralized assembled matrix input; otherwise r and c are ignored */ 1068397b6df1SKris Buschelman #undef __FUNCT__ 1069f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS" 10700481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 1071b24902e0SBarry Smith { 1072a5e57a09SHong Zhang Mat_MUMPS *mumps = (Mat_MUMPS*)F->spptr; 1073dcd589f8SShri Abhyankar PetscErrorCode ierr; 107467877ebaSShri Abhyankar Vec b; 107567877ebaSShri Abhyankar IS is_iden; 107667877ebaSShri Abhyankar const PetscInt M = A->rmap->N; 1077397b6df1SKris Buschelman 1078397b6df1SKris Buschelman PetscFunctionBegin; 1079a5e57a09SHong Zhang mumps->matstruc = DIFFERENT_NONZERO_PATTERN; 1080dcd589f8SShri Abhyankar 10819a2535b5SHong Zhang /* Set MUMPS options from the options database */ 10829a2535b5SHong Zhang ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr); 1083dcd589f8SShri Abhyankar 1084a5e57a09SHong Zhang ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr); 1085dcd589f8SShri Abhyankar 108667877ebaSShri Abhyankar /* analysis phase */ 108767877ebaSShri Abhyankar /*----------------*/ 1088a5e57a09SHong Zhang mumps->id.job = JOB_FACTSYMBOLIC; 1089a5e57a09SHong Zhang mumps->id.n = M; 1090a5e57a09SHong Zhang switch (mumps->id.ICNTL(18)) { 109167877ebaSShri Abhyankar case 0: /* centralized assembled matrix input */ 1092a5e57a09SHong Zhang if (!mumps->myid) { 1093a5e57a09SHong Zhang mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn; 1094a5e57a09SHong Zhang if (mumps->id.ICNTL(6)>1) { 109567877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 10962907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 1097a5e57a09SHong Zhang mumps->id.a = (mumps_complex*)mumps->val; 10982907cef9SHong Zhang #else 1099a5e57a09SHong Zhang mumps->id.a = (mumps_double_complex*)mumps->val; 11002907cef9SHong Zhang #endif 110167877ebaSShri Abhyankar #else 1102a5e57a09SHong Zhang mumps->id.a = mumps->val; 110367877ebaSShri Abhyankar #endif 110467877ebaSShri Abhyankar } 1105a5e57a09SHong Zhang if (mumps->id.ICNTL(7) == 1) { /* use user-provide matrix ordering - assuming r = c ordering */ 11065248a706SHong Zhang /* 11075248a706SHong Zhang PetscBool flag; 11085248a706SHong Zhang ierr = ISEqual(r,c,&flag);CHKERRQ(ierr); 11095248a706SHong Zhang if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"row_perm != col_perm"); 11105248a706SHong Zhang ierr = ISView(r,PETSC_VIEWER_STDOUT_SELF); 11115248a706SHong Zhang */ 1112a5e57a09SHong Zhang if (!mumps->myid) { 1113e0b74bf9SHong Zhang const PetscInt *idx; 1114e0b74bf9SHong Zhang PetscInt i,*perm_in; 11152205254eSKarl Rupp 1116785e854fSJed Brown ierr = PetscMalloc1(M,&perm_in);CHKERRQ(ierr); 1117e0b74bf9SHong Zhang ierr = ISGetIndices(r,&idx);CHKERRQ(ierr); 11182205254eSKarl Rupp 1119a5e57a09SHong Zhang mumps->id.perm_in = perm_in; 1120e0b74bf9SHong Zhang for (i=0; i<M; i++) perm_in[i] = idx[i]+1; /* perm_in[]: start from 1, not 0! */ 1121e0b74bf9SHong Zhang ierr = ISRestoreIndices(r,&idx);CHKERRQ(ierr); 1122e0b74bf9SHong Zhang } 1123e0b74bf9SHong Zhang } 112467877ebaSShri Abhyankar } 112567877ebaSShri Abhyankar break; 112667877ebaSShri Abhyankar case 3: /* distributed assembled matrix input (size>1) */ 1127a5e57a09SHong Zhang mumps->id.nz_loc = mumps->nz; 1128a5e57a09SHong Zhang mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn; 1129a5e57a09SHong Zhang if (mumps->id.ICNTL(6)>1) { 113067877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 11312907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 1132a5e57a09SHong Zhang mumps->id.a_loc = (mumps_complex*)mumps->val; 11332907cef9SHong Zhang #else 1134a5e57a09SHong Zhang mumps->id.a_loc = (mumps_double_complex*)mumps->val; 11352907cef9SHong Zhang #endif 113667877ebaSShri Abhyankar #else 1137a5e57a09SHong Zhang mumps->id.a_loc = mumps->val; 113867877ebaSShri Abhyankar #endif 113967877ebaSShri Abhyankar } 114067877ebaSShri Abhyankar /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */ 1141a5e57a09SHong Zhang if (!mumps->myid) { 11422cd7d884SHong Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,A->rmap->N,&mumps->b_seq);CHKERRQ(ierr); 11432cd7d884SHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,A->rmap->N,0,1,&is_iden);CHKERRQ(ierr); 114467877ebaSShri Abhyankar } else { 1145a5e57a09SHong Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr); 114667877ebaSShri Abhyankar ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr); 114767877ebaSShri Abhyankar } 11482a7a6963SBarry Smith ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr); 1149a5e57a09SHong Zhang ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr); 11506bf464f9SBarry Smith ierr = ISDestroy(&is_iden);CHKERRQ(ierr); 11516bf464f9SBarry Smith ierr = VecDestroy(&b);CHKERRQ(ierr); 115267877ebaSShri Abhyankar break; 115367877ebaSShri Abhyankar } 1154a5e57a09SHong Zhang PetscMUMPS_c(&mumps->id); 1155a5e57a09SHong Zhang if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",mumps->id.INFOG(1)); 115667877ebaSShri Abhyankar 1157719d5645SBarry Smith F->ops->lufactornumeric = MatFactorNumeric_MUMPS; 1158dcd589f8SShri Abhyankar F->ops->solve = MatSolve_MUMPS; 115951d5961aSHong Zhang F->ops->solvetranspose = MatSolveTranspose_MUMPS; 11604e34a73bSHong Zhang F->ops->matsolve = MatMatSolve_MUMPS; 1161b24902e0SBarry Smith PetscFunctionReturn(0); 1162b24902e0SBarry Smith } 1163b24902e0SBarry Smith 1164450b117fSShri Abhyankar /* Note the Petsc r and c permutations are ignored */ 1165450b117fSShri Abhyankar #undef __FUNCT__ 1166450b117fSShri Abhyankar #define __FUNCT__ "MatLUFactorSymbolic_BAIJMUMPS" 1167450b117fSShri Abhyankar PetscErrorCode MatLUFactorSymbolic_BAIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 1168450b117fSShri Abhyankar { 1169a5e57a09SHong Zhang Mat_MUMPS *mumps = (Mat_MUMPS*)F->spptr; 1170dcd589f8SShri Abhyankar PetscErrorCode ierr; 117167877ebaSShri Abhyankar Vec b; 117267877ebaSShri Abhyankar IS is_iden; 117367877ebaSShri Abhyankar const PetscInt M = A->rmap->N; 1174450b117fSShri Abhyankar 1175450b117fSShri Abhyankar PetscFunctionBegin; 1176a5e57a09SHong Zhang mumps->matstruc = DIFFERENT_NONZERO_PATTERN; 1177dcd589f8SShri Abhyankar 11789a2535b5SHong Zhang /* Set MUMPS options from the options database */ 11799a2535b5SHong Zhang ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr); 1180dcd589f8SShri Abhyankar 1181a5e57a09SHong Zhang ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr); 118267877ebaSShri Abhyankar 118367877ebaSShri Abhyankar /* analysis phase */ 118467877ebaSShri Abhyankar /*----------------*/ 1185a5e57a09SHong Zhang mumps->id.job = JOB_FACTSYMBOLIC; 1186a5e57a09SHong Zhang mumps->id.n = M; 1187a5e57a09SHong Zhang switch (mumps->id.ICNTL(18)) { 118867877ebaSShri Abhyankar case 0: /* centralized assembled matrix input */ 1189a5e57a09SHong Zhang if (!mumps->myid) { 1190a5e57a09SHong Zhang mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn; 1191a5e57a09SHong Zhang if (mumps->id.ICNTL(6)>1) { 119267877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 11932907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 1194a5e57a09SHong Zhang mumps->id.a = (mumps_complex*)mumps->val; 11952907cef9SHong Zhang #else 1196a5e57a09SHong Zhang mumps->id.a = (mumps_double_complex*)mumps->val; 11972907cef9SHong Zhang #endif 119867877ebaSShri Abhyankar #else 1199a5e57a09SHong Zhang mumps->id.a = mumps->val; 120067877ebaSShri Abhyankar #endif 120167877ebaSShri Abhyankar } 120267877ebaSShri Abhyankar } 120367877ebaSShri Abhyankar break; 120467877ebaSShri Abhyankar case 3: /* distributed assembled matrix input (size>1) */ 1205a5e57a09SHong Zhang mumps->id.nz_loc = mumps->nz; 1206a5e57a09SHong Zhang mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn; 1207a5e57a09SHong Zhang if (mumps->id.ICNTL(6)>1) { 120867877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 12092907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 1210a5e57a09SHong Zhang mumps->id.a_loc = (mumps_complex*)mumps->val; 12112907cef9SHong Zhang #else 1212a5e57a09SHong Zhang mumps->id.a_loc = (mumps_double_complex*)mumps->val; 12132907cef9SHong Zhang #endif 121467877ebaSShri Abhyankar #else 1215a5e57a09SHong Zhang mumps->id.a_loc = mumps->val; 121667877ebaSShri Abhyankar #endif 121767877ebaSShri Abhyankar } 121867877ebaSShri Abhyankar /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */ 1219a5e57a09SHong Zhang if (!mumps->myid) { 1220a5e57a09SHong Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr); 122167877ebaSShri Abhyankar ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr); 122267877ebaSShri Abhyankar } else { 1223a5e57a09SHong Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr); 122467877ebaSShri Abhyankar ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr); 122567877ebaSShri Abhyankar } 12262a7a6963SBarry Smith ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr); 1227a5e57a09SHong Zhang ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr); 12286bf464f9SBarry Smith ierr = ISDestroy(&is_iden);CHKERRQ(ierr); 12296bf464f9SBarry Smith ierr = VecDestroy(&b);CHKERRQ(ierr); 123067877ebaSShri Abhyankar break; 123167877ebaSShri Abhyankar } 1232a5e57a09SHong Zhang PetscMUMPS_c(&mumps->id); 1233a5e57a09SHong Zhang if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",mumps->id.INFOG(1)); 123467877ebaSShri Abhyankar 1235450b117fSShri Abhyankar F->ops->lufactornumeric = MatFactorNumeric_MUMPS; 1236dcd589f8SShri Abhyankar F->ops->solve = MatSolve_MUMPS; 123751d5961aSHong Zhang F->ops->solvetranspose = MatSolveTranspose_MUMPS; 1238450b117fSShri Abhyankar PetscFunctionReturn(0); 1239450b117fSShri Abhyankar } 1240b24902e0SBarry Smith 1241141f4205SHong Zhang /* Note the Petsc r permutation and factor info are ignored */ 1242397b6df1SKris Buschelman #undef __FUNCT__ 124367877ebaSShri Abhyankar #define __FUNCT__ "MatCholeskyFactorSymbolic_MUMPS" 124467877ebaSShri Abhyankar PetscErrorCode MatCholeskyFactorSymbolic_MUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info) 1245b24902e0SBarry Smith { 1246a5e57a09SHong Zhang Mat_MUMPS *mumps = (Mat_MUMPS*)F->spptr; 1247dcd589f8SShri Abhyankar PetscErrorCode ierr; 124867877ebaSShri Abhyankar Vec b; 124967877ebaSShri Abhyankar IS is_iden; 125067877ebaSShri Abhyankar const PetscInt M = A->rmap->N; 1251397b6df1SKris Buschelman 1252397b6df1SKris Buschelman PetscFunctionBegin; 1253a5e57a09SHong Zhang mumps->matstruc = DIFFERENT_NONZERO_PATTERN; 1254dcd589f8SShri Abhyankar 12559a2535b5SHong Zhang /* Set MUMPS options from the options database */ 12569a2535b5SHong Zhang ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr); 1257dcd589f8SShri Abhyankar 1258a5e57a09SHong Zhang ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr); 1259dcd589f8SShri Abhyankar 126067877ebaSShri Abhyankar /* analysis phase */ 126167877ebaSShri Abhyankar /*----------------*/ 1262a5e57a09SHong Zhang mumps->id.job = JOB_FACTSYMBOLIC; 1263a5e57a09SHong Zhang mumps->id.n = M; 1264a5e57a09SHong Zhang switch (mumps->id.ICNTL(18)) { 126567877ebaSShri Abhyankar case 0: /* centralized assembled matrix input */ 1266a5e57a09SHong Zhang if (!mumps->myid) { 1267a5e57a09SHong Zhang mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn; 1268a5e57a09SHong Zhang if (mumps->id.ICNTL(6)>1) { 126967877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 12702907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 1271a5e57a09SHong Zhang mumps->id.a = (mumps_complex*)mumps->val; 12722907cef9SHong Zhang #else 1273a5e57a09SHong Zhang mumps->id.a = (mumps_double_complex*)mumps->val; 12742907cef9SHong Zhang #endif 127567877ebaSShri Abhyankar #else 1276a5e57a09SHong Zhang mumps->id.a = mumps->val; 127767877ebaSShri Abhyankar #endif 127867877ebaSShri Abhyankar } 127967877ebaSShri Abhyankar } 128067877ebaSShri Abhyankar break; 128167877ebaSShri Abhyankar case 3: /* distributed assembled matrix input (size>1) */ 1282a5e57a09SHong Zhang mumps->id.nz_loc = mumps->nz; 1283a5e57a09SHong Zhang mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn; 1284a5e57a09SHong Zhang if (mumps->id.ICNTL(6)>1) { 128567877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX) 12862907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE) 1287a5e57a09SHong Zhang mumps->id.a_loc = (mumps_complex*)mumps->val; 12882907cef9SHong Zhang #else 1289a5e57a09SHong Zhang mumps->id.a_loc = (mumps_double_complex*)mumps->val; 12902907cef9SHong Zhang #endif 129167877ebaSShri Abhyankar #else 1292a5e57a09SHong Zhang mumps->id.a_loc = mumps->val; 129367877ebaSShri Abhyankar #endif 129467877ebaSShri Abhyankar } 129567877ebaSShri Abhyankar /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */ 1296a5e57a09SHong Zhang if (!mumps->myid) { 1297a5e57a09SHong Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr); 129867877ebaSShri Abhyankar ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr); 129967877ebaSShri Abhyankar } else { 1300a5e57a09SHong Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr); 130167877ebaSShri Abhyankar ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr); 130267877ebaSShri Abhyankar } 13032a7a6963SBarry Smith ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr); 1304a5e57a09SHong Zhang ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr); 13056bf464f9SBarry Smith ierr = ISDestroy(&is_iden);CHKERRQ(ierr); 13066bf464f9SBarry Smith ierr = VecDestroy(&b);CHKERRQ(ierr); 130767877ebaSShri Abhyankar break; 130867877ebaSShri Abhyankar } 1309a5e57a09SHong Zhang PetscMUMPS_c(&mumps->id); 1310a5e57a09SHong Zhang if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",mumps->id.INFOG(1)); 131167877ebaSShri Abhyankar 13122792810eSHong Zhang F->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS; 1313dcd589f8SShri Abhyankar F->ops->solve = MatSolve_MUMPS; 131451d5961aSHong Zhang F->ops->solvetranspose = MatSolve_MUMPS; 13154e34a73bSHong Zhang F->ops->matsolve = MatMatSolve_MUMPS; 13164e34a73bSHong Zhang #if defined(PETSC_USE_COMPLEX) 13170298fd71SBarry Smith F->ops->getinertia = NULL; 13184e34a73bSHong Zhang #else 13194e34a73bSHong Zhang F->ops->getinertia = MatGetInertia_SBAIJMUMPS; 1320db4efbfdSBarry Smith #endif 1321b24902e0SBarry Smith PetscFunctionReturn(0); 1322b24902e0SBarry Smith } 1323b24902e0SBarry Smith 1324397b6df1SKris Buschelman #undef __FUNCT__ 132564e6c443SBarry Smith #define __FUNCT__ "MatView_MUMPS" 132664e6c443SBarry Smith PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer) 132774ed9c26SBarry Smith { 1328f6c57405SHong Zhang PetscErrorCode ierr; 132964e6c443SBarry Smith PetscBool iascii; 133064e6c443SBarry Smith PetscViewerFormat format; 1331a5e57a09SHong Zhang Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 1332f6c57405SHong Zhang 1333f6c57405SHong Zhang PetscFunctionBegin; 133464e6c443SBarry Smith /* check if matrix is mumps type */ 133564e6c443SBarry Smith if (A->ops->solve != MatSolve_MUMPS) PetscFunctionReturn(0); 133664e6c443SBarry Smith 1337251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 133864e6c443SBarry Smith if (iascii) { 133964e6c443SBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 134064e6c443SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO) { 134164e6c443SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr); 1342a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," SYM (matrix type): %d \n",mumps->id.sym);CHKERRQ(ierr); 1343a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," PAR (host participation): %d \n",mumps->id.par);CHKERRQ(ierr); 1344a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(1) (output for error): %d \n",mumps->id.ICNTL(1));CHKERRQ(ierr); 1345a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(2) (output of diagnostic msg): %d \n",mumps->id.ICNTL(2));CHKERRQ(ierr); 1346a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(3) (output for global info): %d \n",mumps->id.ICNTL(3));CHKERRQ(ierr); 1347a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(4) (level of printing): %d \n",mumps->id.ICNTL(4));CHKERRQ(ierr); 1348a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(5) (input mat struct): %d \n",mumps->id.ICNTL(5));CHKERRQ(ierr); 1349a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(6) (matrix prescaling): %d \n",mumps->id.ICNTL(6));CHKERRQ(ierr); 1350a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(7) (sequentia matrix ordering):%d \n",mumps->id.ICNTL(7));CHKERRQ(ierr); 1351a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(8) (scalling strategy): %d \n",mumps->id.ICNTL(8));CHKERRQ(ierr); 1352a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(10) (max num of refinements): %d \n",mumps->id.ICNTL(10));CHKERRQ(ierr); 1353a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(11) (error analysis): %d \n",mumps->id.ICNTL(11));CHKERRQ(ierr); 1354a5e57a09SHong Zhang if (mumps->id.ICNTL(11)>0) { 1355a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(4) (inf norm of input mat): %g\n",mumps->id.RINFOG(4));CHKERRQ(ierr); 1356a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(5) (inf norm of solution): %g\n",mumps->id.RINFOG(5));CHKERRQ(ierr); 1357a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(6) (inf norm of residual): %g\n",mumps->id.RINFOG(6));CHKERRQ(ierr); 1358a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",mumps->id.RINFOG(7),mumps->id.RINFOG(8));CHKERRQ(ierr); 1359a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(9) (error estimate): %g \n",mumps->id.RINFOG(9));CHKERRQ(ierr); 1360a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",mumps->id.RINFOG(10),mumps->id.RINFOG(11));CHKERRQ(ierr); 1361f6c57405SHong Zhang } 1362a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(12) (efficiency control): %d \n",mumps->id.ICNTL(12));CHKERRQ(ierr); 1363a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(13) (efficiency control): %d \n",mumps->id.ICNTL(13));CHKERRQ(ierr); 1364a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(14) (percentage of estimated workspace increase): %d \n",mumps->id.ICNTL(14));CHKERRQ(ierr); 1365f6c57405SHong Zhang /* ICNTL(15-17) not used */ 1366a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(18) (input mat struct): %d \n",mumps->id.ICNTL(18));CHKERRQ(ierr); 1367a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(19) (Shur complement info): %d \n",mumps->id.ICNTL(19));CHKERRQ(ierr); 1368a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(20) (rhs sparse pattern): %d \n",mumps->id.ICNTL(20));CHKERRQ(ierr); 1369a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(21) (somumpstion struct): %d \n",mumps->id.ICNTL(21));CHKERRQ(ierr); 1370a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(22) (in-core/out-of-core facility): %d \n",mumps->id.ICNTL(22));CHKERRQ(ierr); 1371a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(23) (max size of memory can be allocated locally):%d \n",mumps->id.ICNTL(23));CHKERRQ(ierr); 1372c0165424SHong Zhang 1373a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(24) (detection of null pivot rows): %d \n",mumps->id.ICNTL(24));CHKERRQ(ierr); 1374a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(25) (computation of a null space basis): %d \n",mumps->id.ICNTL(25));CHKERRQ(ierr); 1375a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(26) (Schur options for rhs or solution): %d \n",mumps->id.ICNTL(26));CHKERRQ(ierr); 1376a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(27) (experimental parameter): %d \n",mumps->id.ICNTL(27));CHKERRQ(ierr); 1377a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(28) (use parallel or sequential ordering): %d \n",mumps->id.ICNTL(28));CHKERRQ(ierr); 1378a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(29) (parallel ordering): %d \n",mumps->id.ICNTL(29));CHKERRQ(ierr); 137942179a6aSHong Zhang 1380a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(30) (user-specified set of entries in inv(A)): %d \n",mumps->id.ICNTL(30));CHKERRQ(ierr); 1381a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(31) (factors is discarded in the solve phase): %d \n",mumps->id.ICNTL(31));CHKERRQ(ierr); 1382a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(33) (compute determinant): %d \n",mumps->id.ICNTL(33));CHKERRQ(ierr); 1383f6c57405SHong Zhang 1384a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(1) (relative pivoting threshold): %g \n",mumps->id.CNTL(1));CHKERRQ(ierr); 1385a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(2) (stopping criterion of refinement): %g \n",mumps->id.CNTL(2));CHKERRQ(ierr); 1386a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(3) (absomumpste pivoting threshold): %g \n",mumps->id.CNTL(3));CHKERRQ(ierr); 1387a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(4) (vamumpse of static pivoting): %g \n",mumps->id.CNTL(4));CHKERRQ(ierr); 1388a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(5) (fixation for null pivots): %g \n",mumps->id.CNTL(5));CHKERRQ(ierr); 1389f6c57405SHong Zhang 1390f6c57405SHong Zhang /* infomation local to each processor */ 139134ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer, " RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr); 13927b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 1393a5e57a09SHong Zhang ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %g \n",mumps->myid,mumps->id.RINFO(1));CHKERRQ(ierr); 139434ed7027SBarry Smith ierr = PetscViewerFlush(viewer); 139534ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer, " RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr); 1396a5e57a09SHong Zhang ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %g \n",mumps->myid,mumps->id.RINFO(2));CHKERRQ(ierr); 139734ed7027SBarry Smith ierr = PetscViewerFlush(viewer); 139834ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer, " RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr); 1399a5e57a09SHong Zhang ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %g \n",mumps->myid,mumps->id.RINFO(3));CHKERRQ(ierr); 140034ed7027SBarry Smith ierr = PetscViewerFlush(viewer); 1401f6c57405SHong Zhang 140234ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer, " INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization): \n");CHKERRQ(ierr); 1403a5e57a09SHong Zhang ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %d \n",mumps->myid,mumps->id.INFO(15));CHKERRQ(ierr); 140434ed7027SBarry Smith ierr = PetscViewerFlush(viewer); 1405f6c57405SHong Zhang 140634ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer, " INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization): \n");CHKERRQ(ierr); 1407a5e57a09SHong Zhang ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %d \n",mumps->myid,mumps->id.INFO(16));CHKERRQ(ierr); 140834ed7027SBarry Smith ierr = PetscViewerFlush(viewer); 1409f6c57405SHong Zhang 141034ed7027SBarry Smith ierr = PetscViewerASCIIPrintf(viewer, " INFO(23) (num of pivots eliminated on this processor after factorization): \n");CHKERRQ(ierr); 1411a5e57a09SHong Zhang ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %d \n",mumps->myid,mumps->id.INFO(23));CHKERRQ(ierr); 141234ed7027SBarry Smith ierr = PetscViewerFlush(viewer); 1413*b34f08ffSHong Zhang 1414*b34f08ffSHong Zhang if (mumps->ninfo && mumps->ninfo <= 40){ 1415*b34f08ffSHong Zhang PetscInt i; 1416*b34f08ffSHong Zhang for (i=0; i<mumps->ninfo; i++){ 1417*b34f08ffSHong Zhang ierr = PetscViewerASCIIPrintf(viewer, " INFO(%d): \n",mumps->info[i]);CHKERRQ(ierr); 1418*b34f08ffSHong Zhang ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %d \n",mumps->myid,mumps->id.INFO(mumps->info[i]));CHKERRQ(ierr); 1419*b34f08ffSHong Zhang ierr = PetscViewerFlush(viewer); 1420*b34f08ffSHong Zhang } 1421*b34f08ffSHong Zhang } 1422*b34f08ffSHong Zhang 1423*b34f08ffSHong Zhang 14247b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr); 1425f6c57405SHong Zhang 1426a5e57a09SHong Zhang if (!mumps->myid) { /* information from the host */ 1427a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",mumps->id.RINFOG(1));CHKERRQ(ierr); 1428a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",mumps->id.RINFOG(2));CHKERRQ(ierr); 1429a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",mumps->id.RINFOG(3));CHKERRQ(ierr); 1430a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," (RINFOG(12) RINFOG(13))*2^INFOG(34) (determinant): (%g,%g)*(2^%d)\n",mumps->id.RINFOG(12),mumps->id.RINFOG(13),mumps->id.INFOG(34));CHKERRQ(ierr); 1431f6c57405SHong Zhang 1432a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(3));CHKERRQ(ierr); 1433a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(4));CHKERRQ(ierr); 1434a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(5) (estimated maximum front size in the complete tree): %d \n",mumps->id.INFOG(5));CHKERRQ(ierr); 1435a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(6) (number of nodes in the complete tree): %d \n",mumps->id.INFOG(6));CHKERRQ(ierr); 1436a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(7) (ordering option effectively use after analysis): %d \n",mumps->id.INFOG(7));CHKERRQ(ierr); 1437a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",mumps->id.INFOG(8));CHKERRQ(ierr); 1438a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(9) (total real/complex workspace to store the matrix factors after factorization): %d \n",mumps->id.INFOG(9));CHKERRQ(ierr); 1439a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(10) (total integer space store the matrix factors after factorization): %d \n",mumps->id.INFOG(10));CHKERRQ(ierr); 1440a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(11) (order of largest frontal matrix after factorization): %d \n",mumps->id.INFOG(11));CHKERRQ(ierr); 1441a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(12) (number of off-diagonal pivots): %d \n",mumps->id.INFOG(12));CHKERRQ(ierr); 1442a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(13) (number of delayed pivots after factorization): %d \n",mumps->id.INFOG(13));CHKERRQ(ierr); 1443a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(14) (number of memory compress after factorization): %d \n",mumps->id.INFOG(14));CHKERRQ(ierr); 1444a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(15) (number of steps of iterative refinement after solution): %d \n",mumps->id.INFOG(15));CHKERRQ(ierr); 1445a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(16) (estimated size (in MB) of all MUMPS internal data for factorization after analysis: value on the most memory consuming processor): %d \n",mumps->id.INFOG(16));CHKERRQ(ierr); 1446a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(17) (estimated size of all MUMPS internal data for factorization after analysis: sum over all processors): %d \n",mumps->id.INFOG(17));CHKERRQ(ierr); 1447a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(18) (size of all MUMPS internal data allocated during factorization: value on the most memory consuming processor): %d \n",mumps->id.INFOG(18));CHKERRQ(ierr); 1448a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(19) (size of all MUMPS internal data allocated during factorization: sum over all processors): %d \n",mumps->id.INFOG(19));CHKERRQ(ierr); 1449a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(20) (estimated number of entries in the factors): %d \n",mumps->id.INFOG(20));CHKERRQ(ierr); 1450a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(21) (size in MB of memory effectively used during factorization - value on the most memory consuming processor): %d \n",mumps->id.INFOG(21));CHKERRQ(ierr); 1451a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(22) (size in MB of memory effectively used during factorization - sum over all processors): %d \n",mumps->id.INFOG(22));CHKERRQ(ierr); 1452a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",mumps->id.INFOG(23));CHKERRQ(ierr); 1453a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",mumps->id.INFOG(24));CHKERRQ(ierr); 1454a5e57a09SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",mumps->id.INFOG(25));CHKERRQ(ierr); 145540d435e3SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(28) (after factorization: number of null pivots encountered): %d\n",mumps->id.INFOG(28));CHKERRQ(ierr); 145640d435e3SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(29) (after factorization: effective number of entries in the factors (sum over all processors)): %d\n",mumps->id.INFOG(29));CHKERRQ(ierr); 145740d435e3SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(30, 31) (after solution: size in Mbytes of memory used during solution phase): %d, %d\n",mumps->id.INFOG(30),mumps->id.INFOG(31));CHKERRQ(ierr); 145840d435e3SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(32) (after analysis: type of analysis done): %d\n",mumps->id.INFOG(32));CHKERRQ(ierr); 145940d435e3SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(33) (value used for ICNTL(8)): %d\n",mumps->id.INFOG(33));CHKERRQ(ierr); 146040d435e3SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(34) (exponent of the determinant if determinant is requested): %d\n",mumps->id.INFOG(34));CHKERRQ(ierr); 1461f6c57405SHong Zhang } 1462f6c57405SHong Zhang } 1463cb828f0fSHong Zhang } 1464f6c57405SHong Zhang PetscFunctionReturn(0); 1465f6c57405SHong Zhang } 1466f6c57405SHong Zhang 146735bd34faSBarry Smith #undef __FUNCT__ 146835bd34faSBarry Smith #define __FUNCT__ "MatGetInfo_MUMPS" 146935bd34faSBarry Smith PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info) 147035bd34faSBarry Smith { 1471cb828f0fSHong Zhang Mat_MUMPS *mumps =(Mat_MUMPS*)A->spptr; 147235bd34faSBarry Smith 147335bd34faSBarry Smith PetscFunctionBegin; 147435bd34faSBarry Smith info->block_size = 1.0; 1475cb828f0fSHong Zhang info->nz_allocated = mumps->id.INFOG(20); 1476cb828f0fSHong Zhang info->nz_used = mumps->id.INFOG(20); 147735bd34faSBarry Smith info->nz_unneeded = 0.0; 147835bd34faSBarry Smith info->assemblies = 0.0; 147935bd34faSBarry Smith info->mallocs = 0.0; 148035bd34faSBarry Smith info->memory = 0.0; 148135bd34faSBarry Smith info->fill_ratio_given = 0; 148235bd34faSBarry Smith info->fill_ratio_needed = 0; 148335bd34faSBarry Smith info->factor_mallocs = 0; 148435bd34faSBarry Smith PetscFunctionReturn(0); 148535bd34faSBarry Smith } 148635bd34faSBarry Smith 14875ccb76cbSHong Zhang /* -------------------------------------------------------------------------------------------*/ 14885ccb76cbSHong Zhang #undef __FUNCT__ 14895ccb76cbSHong Zhang #define __FUNCT__ "MatMumpsSetIcntl_MUMPS" 14905ccb76cbSHong Zhang PetscErrorCode MatMumpsSetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt ival) 14915ccb76cbSHong Zhang { 1492a5e57a09SHong Zhang Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 14935ccb76cbSHong Zhang 14945ccb76cbSHong Zhang PetscFunctionBegin; 1495a5e57a09SHong Zhang mumps->id.ICNTL(icntl) = ival; 14965ccb76cbSHong Zhang PetscFunctionReturn(0); 14975ccb76cbSHong Zhang } 14985ccb76cbSHong Zhang 14995ccb76cbSHong Zhang #undef __FUNCT__ 1500bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetIcntl_MUMPS" 1501bc6112feSHong Zhang PetscErrorCode MatMumpsGetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt *ival) 1502bc6112feSHong Zhang { 1503bc6112feSHong Zhang Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 1504bc6112feSHong Zhang 1505bc6112feSHong Zhang PetscFunctionBegin; 1506bc6112feSHong Zhang *ival = mumps->id.ICNTL(icntl); 1507bc6112feSHong Zhang PetscFunctionReturn(0); 1508bc6112feSHong Zhang } 1509bc6112feSHong Zhang 1510bc6112feSHong Zhang #undef __FUNCT__ 15115ccb76cbSHong Zhang #define __FUNCT__ "MatMumpsSetIcntl" 15125ccb76cbSHong Zhang /*@ 15135ccb76cbSHong Zhang MatMumpsSetIcntl - Set MUMPS parameter ICNTL() 15145ccb76cbSHong Zhang 15155ccb76cbSHong Zhang Logically Collective on Mat 15165ccb76cbSHong Zhang 15175ccb76cbSHong Zhang Input Parameters: 15185ccb76cbSHong Zhang + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 15195ccb76cbSHong Zhang . icntl - index of MUMPS parameter array ICNTL() 15205ccb76cbSHong Zhang - ival - value of MUMPS ICNTL(icntl) 15215ccb76cbSHong Zhang 15225ccb76cbSHong Zhang Options Database: 15235ccb76cbSHong Zhang . -mat_mumps_icntl_<icntl> <ival> 15245ccb76cbSHong Zhang 15255ccb76cbSHong Zhang Level: beginner 15265ccb76cbSHong Zhang 15275ccb76cbSHong Zhang References: MUMPS Users' Guide 15285ccb76cbSHong Zhang 15295ccb76cbSHong Zhang .seealso: MatGetFactor() 15305ccb76cbSHong Zhang @*/ 15315ccb76cbSHong Zhang PetscErrorCode MatMumpsSetIcntl(Mat F,PetscInt icntl,PetscInt ival) 15325ccb76cbSHong Zhang { 15335ccb76cbSHong Zhang PetscErrorCode ierr; 15345ccb76cbSHong Zhang 15355ccb76cbSHong Zhang PetscFunctionBegin; 15365ccb76cbSHong Zhang PetscValidLogicalCollectiveInt(F,icntl,2); 15375ccb76cbSHong Zhang PetscValidLogicalCollectiveInt(F,ival,3); 15385ccb76cbSHong Zhang ierr = PetscTryMethod(F,"MatMumpsSetIcntl_C",(Mat,PetscInt,PetscInt),(F,icntl,ival));CHKERRQ(ierr); 15395ccb76cbSHong Zhang PetscFunctionReturn(0); 15405ccb76cbSHong Zhang } 15415ccb76cbSHong Zhang 1542bc6112feSHong Zhang #undef __FUNCT__ 1543bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetIcntl" 1544a21f80fcSHong Zhang /*@ 1545a21f80fcSHong Zhang MatMumpsGetIcntl - Get MUMPS parameter ICNTL() 1546a21f80fcSHong Zhang 1547a21f80fcSHong Zhang Logically Collective on Mat 1548a21f80fcSHong Zhang 1549a21f80fcSHong Zhang Input Parameters: 1550a21f80fcSHong Zhang + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 1551a21f80fcSHong Zhang - icntl - index of MUMPS parameter array ICNTL() 1552a21f80fcSHong Zhang 1553a21f80fcSHong Zhang Output Parameter: 1554a21f80fcSHong Zhang . ival - value of MUMPS ICNTL(icntl) 1555a21f80fcSHong Zhang 1556a21f80fcSHong Zhang Level: beginner 1557a21f80fcSHong Zhang 1558a21f80fcSHong Zhang References: MUMPS Users' Guide 1559a21f80fcSHong Zhang 1560a21f80fcSHong Zhang .seealso: MatGetFactor() 1561a21f80fcSHong Zhang @*/ 1562bc6112feSHong Zhang PetscErrorCode MatMumpsGetIcntl(Mat F,PetscInt icntl,PetscInt *ival) 1563bc6112feSHong Zhang { 1564bc6112feSHong Zhang PetscErrorCode ierr; 1565bc6112feSHong Zhang 1566bc6112feSHong Zhang PetscFunctionBegin; 1567bc6112feSHong Zhang PetscValidLogicalCollectiveInt(F,icntl,2); 1568bc6112feSHong Zhang PetscValidIntPointer(ival,3); 1569bc6112feSHong Zhang ierr = PetscTryMethod(F,"MatMumpsGetIcntl_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr); 1570bc6112feSHong Zhang PetscFunctionReturn(0); 1571bc6112feSHong Zhang } 1572bc6112feSHong Zhang 15738928b65cSHong Zhang /* -------------------------------------------------------------------------------------------*/ 15748928b65cSHong Zhang #undef __FUNCT__ 15758928b65cSHong Zhang #define __FUNCT__ "MatMumpsSetCntl_MUMPS" 15768928b65cSHong Zhang PetscErrorCode MatMumpsSetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal val) 15778928b65cSHong Zhang { 15788928b65cSHong Zhang Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 15798928b65cSHong Zhang 15808928b65cSHong Zhang PetscFunctionBegin; 15818928b65cSHong Zhang mumps->id.CNTL(icntl) = val; 15828928b65cSHong Zhang PetscFunctionReturn(0); 15838928b65cSHong Zhang } 15848928b65cSHong Zhang 15858928b65cSHong Zhang #undef __FUNCT__ 1586bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetCntl_MUMPS" 1587bc6112feSHong Zhang PetscErrorCode MatMumpsGetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal *val) 1588bc6112feSHong Zhang { 1589bc6112feSHong Zhang Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 1590bc6112feSHong Zhang 1591bc6112feSHong Zhang PetscFunctionBegin; 1592bc6112feSHong Zhang *val = mumps->id.CNTL(icntl); 1593bc6112feSHong Zhang PetscFunctionReturn(0); 1594bc6112feSHong Zhang } 1595bc6112feSHong Zhang 1596bc6112feSHong Zhang #undef __FUNCT__ 15978928b65cSHong Zhang #define __FUNCT__ "MatMumpsSetCntl" 15988928b65cSHong Zhang /*@ 15998928b65cSHong Zhang MatMumpsSetCntl - Set MUMPS parameter CNTL() 16008928b65cSHong Zhang 16018928b65cSHong Zhang Logically Collective on Mat 16028928b65cSHong Zhang 16038928b65cSHong Zhang Input Parameters: 16048928b65cSHong Zhang + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 16058928b65cSHong Zhang . icntl - index of MUMPS parameter array CNTL() 16068928b65cSHong Zhang - val - value of MUMPS CNTL(icntl) 16078928b65cSHong Zhang 16088928b65cSHong Zhang Options Database: 16098928b65cSHong Zhang . -mat_mumps_cntl_<icntl> <val> 16108928b65cSHong Zhang 16118928b65cSHong Zhang Level: beginner 16128928b65cSHong Zhang 16138928b65cSHong Zhang References: MUMPS Users' Guide 16148928b65cSHong Zhang 16158928b65cSHong Zhang .seealso: MatGetFactor() 16168928b65cSHong Zhang @*/ 16178928b65cSHong Zhang PetscErrorCode MatMumpsSetCntl(Mat F,PetscInt icntl,PetscReal val) 16188928b65cSHong Zhang { 16198928b65cSHong Zhang PetscErrorCode ierr; 16208928b65cSHong Zhang 16218928b65cSHong Zhang PetscFunctionBegin; 16228928b65cSHong Zhang PetscValidLogicalCollectiveInt(F,icntl,2); 1623bc6112feSHong Zhang PetscValidLogicalCollectiveReal(F,val,3); 16248928b65cSHong Zhang ierr = PetscTryMethod(F,"MatMumpsSetCntl_C",(Mat,PetscInt,PetscReal),(F,icntl,val));CHKERRQ(ierr); 16258928b65cSHong Zhang PetscFunctionReturn(0); 16268928b65cSHong Zhang } 16278928b65cSHong Zhang 1628bc6112feSHong Zhang #undef __FUNCT__ 1629bc6112feSHong Zhang #define __FUNCT__ "MatMumpsGetCntl" 1630a21f80fcSHong Zhang /*@ 1631a21f80fcSHong Zhang MatMumpsGetCntl - Get MUMPS parameter CNTL() 1632a21f80fcSHong Zhang 1633a21f80fcSHong Zhang Logically Collective on Mat 1634a21f80fcSHong Zhang 1635a21f80fcSHong Zhang Input Parameters: 1636a21f80fcSHong Zhang + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 1637a21f80fcSHong Zhang - icntl - index of MUMPS parameter array CNTL() 1638a21f80fcSHong Zhang 1639a21f80fcSHong Zhang Output Parameter: 1640a21f80fcSHong Zhang . val - value of MUMPS CNTL(icntl) 1641a21f80fcSHong Zhang 1642a21f80fcSHong Zhang Level: beginner 1643a21f80fcSHong Zhang 1644a21f80fcSHong Zhang References: MUMPS Users' Guide 1645a21f80fcSHong Zhang 1646a21f80fcSHong Zhang .seealso: MatGetFactor() 1647a21f80fcSHong Zhang @*/ 1648bc6112feSHong Zhang PetscErrorCode MatMumpsGetCntl(Mat F,PetscInt icntl,PetscReal *val) 1649bc6112feSHong Zhang { 1650bc6112feSHong Zhang PetscErrorCode ierr; 1651bc6112feSHong Zhang 1652bc6112feSHong Zhang PetscFunctionBegin; 1653bc6112feSHong Zhang PetscValidLogicalCollectiveInt(F,icntl,2); 1654bc6112feSHong Zhang PetscValidRealPointer(val,3); 1655bc6112feSHong Zhang ierr = PetscTryMethod(F,"MatMumpsGetCntl_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr); 1656bc6112feSHong Zhang PetscFunctionReturn(0); 1657bc6112feSHong Zhang } 1658bc6112feSHong Zhang 1659bc6112feSHong Zhang #undef __FUNCT__ 1660ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfo_MUMPS" 1661ca810319SHong Zhang PetscErrorCode MatMumpsGetInfo_MUMPS(Mat F,PetscInt icntl,PetscInt *info) 1662bc6112feSHong Zhang { 1663bc6112feSHong Zhang Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 1664bc6112feSHong Zhang 1665bc6112feSHong Zhang PetscFunctionBegin; 1666bc6112feSHong Zhang *info = mumps->id.INFO(icntl); 1667bc6112feSHong Zhang PetscFunctionReturn(0); 1668bc6112feSHong Zhang } 1669bc6112feSHong Zhang 1670bc6112feSHong Zhang #undef __FUNCT__ 1671ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfog_MUMPS" 1672ca810319SHong Zhang PetscErrorCode MatMumpsGetInfog_MUMPS(Mat F,PetscInt icntl,PetscInt *infog) 1673bc6112feSHong Zhang { 1674bc6112feSHong Zhang Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 1675bc6112feSHong Zhang 1676bc6112feSHong Zhang PetscFunctionBegin; 1677bc6112feSHong Zhang *infog = mumps->id.INFOG(icntl); 1678bc6112feSHong Zhang PetscFunctionReturn(0); 1679bc6112feSHong Zhang } 1680bc6112feSHong Zhang 1681bc6112feSHong Zhang #undef __FUNCT__ 1682ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfo_MUMPS" 1683ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfo_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfo) 1684bc6112feSHong Zhang { 1685bc6112feSHong Zhang Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 1686bc6112feSHong Zhang 1687bc6112feSHong Zhang PetscFunctionBegin; 1688bc6112feSHong Zhang *rinfo = mumps->id.RINFO(icntl); 1689bc6112feSHong Zhang PetscFunctionReturn(0); 1690bc6112feSHong Zhang } 1691bc6112feSHong Zhang 1692bc6112feSHong Zhang #undef __FUNCT__ 1693ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfog_MUMPS" 1694ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfog_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfog) 1695bc6112feSHong Zhang { 1696bc6112feSHong Zhang Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 1697bc6112feSHong Zhang 1698bc6112feSHong Zhang PetscFunctionBegin; 1699bc6112feSHong Zhang *rinfog = mumps->id.RINFOG(icntl); 1700bc6112feSHong Zhang PetscFunctionReturn(0); 1701bc6112feSHong Zhang } 1702bc6112feSHong Zhang 1703bc6112feSHong Zhang #undef __FUNCT__ 1704ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfo" 1705a21f80fcSHong Zhang /*@ 1706a21f80fcSHong Zhang MatMumpsGetInfo - Get MUMPS parameter INFO() 1707a21f80fcSHong Zhang 1708a21f80fcSHong Zhang Logically Collective on Mat 1709a21f80fcSHong Zhang 1710a21f80fcSHong Zhang Input Parameters: 1711a21f80fcSHong Zhang + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 1712a21f80fcSHong Zhang - icntl - index of MUMPS parameter array INFO() 1713a21f80fcSHong Zhang 1714a21f80fcSHong Zhang Output Parameter: 1715a21f80fcSHong Zhang . ival - value of MUMPS INFO(icntl) 1716a21f80fcSHong Zhang 1717a21f80fcSHong Zhang Level: beginner 1718a21f80fcSHong Zhang 1719a21f80fcSHong Zhang References: MUMPS Users' Guide 1720a21f80fcSHong Zhang 1721a21f80fcSHong Zhang .seealso: MatGetFactor() 1722a21f80fcSHong Zhang @*/ 1723ca810319SHong Zhang PetscErrorCode MatMumpsGetInfo(Mat F,PetscInt icntl,PetscInt *ival) 1724bc6112feSHong Zhang { 1725bc6112feSHong Zhang PetscErrorCode ierr; 1726bc6112feSHong Zhang 1727bc6112feSHong Zhang PetscFunctionBegin; 1728ca810319SHong Zhang PetscValidIntPointer(ival,3); 1729ca810319SHong Zhang ierr = PetscTryMethod(F,"MatMumpsGetInfo_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr); 1730bc6112feSHong Zhang PetscFunctionReturn(0); 1731bc6112feSHong Zhang } 1732bc6112feSHong Zhang 1733bc6112feSHong Zhang #undef __FUNCT__ 1734ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetInfog" 1735a21f80fcSHong Zhang /*@ 1736a21f80fcSHong Zhang MatMumpsGetInfog - Get MUMPS parameter INFOG() 1737a21f80fcSHong Zhang 1738a21f80fcSHong Zhang Logically Collective on Mat 1739a21f80fcSHong Zhang 1740a21f80fcSHong Zhang Input Parameters: 1741a21f80fcSHong Zhang + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 1742a21f80fcSHong Zhang - icntl - index of MUMPS parameter array INFOG() 1743a21f80fcSHong Zhang 1744a21f80fcSHong Zhang Output Parameter: 1745a21f80fcSHong Zhang . ival - value of MUMPS INFOG(icntl) 1746a21f80fcSHong Zhang 1747a21f80fcSHong Zhang Level: beginner 1748a21f80fcSHong Zhang 1749a21f80fcSHong Zhang References: MUMPS Users' Guide 1750a21f80fcSHong Zhang 1751a21f80fcSHong Zhang .seealso: MatGetFactor() 1752a21f80fcSHong Zhang @*/ 1753ca810319SHong Zhang PetscErrorCode MatMumpsGetInfog(Mat F,PetscInt icntl,PetscInt *ival) 1754bc6112feSHong Zhang { 1755bc6112feSHong Zhang PetscErrorCode ierr; 1756bc6112feSHong Zhang 1757bc6112feSHong Zhang PetscFunctionBegin; 1758ca810319SHong Zhang PetscValidIntPointer(ival,3); 1759ca810319SHong Zhang ierr = PetscTryMethod(F,"MatMumpsGetInfog_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr); 1760bc6112feSHong Zhang PetscFunctionReturn(0); 1761bc6112feSHong Zhang } 1762bc6112feSHong Zhang 1763bc6112feSHong Zhang #undef __FUNCT__ 1764ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfo" 1765a21f80fcSHong Zhang /*@ 1766a21f80fcSHong Zhang MatMumpsGetRinfo - Get MUMPS parameter RINFO() 1767a21f80fcSHong Zhang 1768a21f80fcSHong Zhang Logically Collective on Mat 1769a21f80fcSHong Zhang 1770a21f80fcSHong Zhang Input Parameters: 1771a21f80fcSHong Zhang + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 1772a21f80fcSHong Zhang - icntl - index of MUMPS parameter array RINFO() 1773a21f80fcSHong Zhang 1774a21f80fcSHong Zhang Output Parameter: 1775a21f80fcSHong Zhang . val - value of MUMPS RINFO(icntl) 1776a21f80fcSHong Zhang 1777a21f80fcSHong Zhang Level: beginner 1778a21f80fcSHong Zhang 1779a21f80fcSHong Zhang References: MUMPS Users' Guide 1780a21f80fcSHong Zhang 1781a21f80fcSHong Zhang .seealso: MatGetFactor() 1782a21f80fcSHong Zhang @*/ 1783ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfo(Mat F,PetscInt icntl,PetscReal *val) 1784bc6112feSHong Zhang { 1785bc6112feSHong Zhang PetscErrorCode ierr; 1786bc6112feSHong Zhang 1787bc6112feSHong Zhang PetscFunctionBegin; 1788bc6112feSHong Zhang PetscValidRealPointer(val,3); 1789ca810319SHong Zhang ierr = PetscTryMethod(F,"MatMumpsGetRinfo_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr); 1790bc6112feSHong Zhang PetscFunctionReturn(0); 1791bc6112feSHong Zhang } 1792bc6112feSHong Zhang 1793bc6112feSHong Zhang #undef __FUNCT__ 1794ca810319SHong Zhang #define __FUNCT__ "MatMumpsGetRinfog" 1795a21f80fcSHong Zhang /*@ 1796a21f80fcSHong Zhang MatMumpsGetRinfog - Get MUMPS parameter RINFOG() 1797a21f80fcSHong Zhang 1798a21f80fcSHong Zhang Logically Collective on Mat 1799a21f80fcSHong Zhang 1800a21f80fcSHong Zhang Input Parameters: 1801a21f80fcSHong Zhang + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 1802a21f80fcSHong Zhang - icntl - index of MUMPS parameter array RINFOG() 1803a21f80fcSHong Zhang 1804a21f80fcSHong Zhang Output Parameter: 1805a21f80fcSHong Zhang . val - value of MUMPS RINFOG(icntl) 1806a21f80fcSHong Zhang 1807a21f80fcSHong Zhang Level: beginner 1808a21f80fcSHong Zhang 1809a21f80fcSHong Zhang References: MUMPS Users' Guide 1810a21f80fcSHong Zhang 1811a21f80fcSHong Zhang .seealso: MatGetFactor() 1812a21f80fcSHong Zhang @*/ 1813ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfog(Mat F,PetscInt icntl,PetscReal *val) 1814bc6112feSHong Zhang { 1815bc6112feSHong Zhang PetscErrorCode ierr; 1816bc6112feSHong Zhang 1817bc6112feSHong Zhang PetscFunctionBegin; 1818bc6112feSHong Zhang PetscValidRealPointer(val,3); 1819ca810319SHong Zhang ierr = PetscTryMethod(F,"MatMumpsGetRinfog_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr); 1820bc6112feSHong Zhang PetscFunctionReturn(0); 1821bc6112feSHong Zhang } 1822bc6112feSHong Zhang 182324b6179bSKris Buschelman /*MC 18242692d6eeSBarry Smith MATSOLVERMUMPS - A matrix type providing direct solvers (LU and Cholesky) for 182524b6179bSKris Buschelman distributed and sequential matrices via the external package MUMPS. 182624b6179bSKris Buschelman 182741c8de11SBarry Smith Works with MATAIJ and MATSBAIJ matrices 182824b6179bSKris Buschelman 182924b6179bSKris Buschelman Options Database Keys: 18304e34a73bSHong Zhang + -mat_mumps_icntl_1 <6>: ICNTL(1): output stream for error messages (None) 18314e34a73bSHong Zhang . -mat_mumps_icntl_2 <0>: ICNTL(2): output stream for diagnostic printing, statistics, and warning (None) 18324e34a73bSHong Zhang . -mat_mumps_icntl_3 <0>: ICNTL(3): output stream for global information, collected on the host (None) 18334e34a73bSHong Zhang . -mat_mumps_icntl_4 <0>: ICNTL(4): level of printing (0 to 4) (None) 18344e34a73bSHong Zhang . -mat_mumps_icntl_6 <7>: ICNTL(6): permutes to a zero-free diagonal and/or scale the matrix (0 to 7) (None) 18354e34a73bSHong Zhang . -mat_mumps_icntl_7 <7>: ICNTL(7): computes a symmetric permutation in sequential analysis (0 to 7). 3=Scotch, 4=PORD, 5=Metis (None) 18364e34a73bSHong Zhang . -mat_mumps_icntl_8 <77>: ICNTL(8): scaling strategy (-2 to 8 or 77) (None) 18374e34a73bSHong Zhang . -mat_mumps_icntl_10 <0>: ICNTL(10): max num of refinements (None) 18384e34a73bSHong Zhang . -mat_mumps_icntl_11 <0>: ICNTL(11): statistics related to an error analysis (via -ksp_view) (None) 18394e34a73bSHong Zhang . -mat_mumps_icntl_12 <1>: ICNTL(12): an ordering strategy for symmetric matrices (0 to 3) (None) 18404e34a73bSHong Zhang . -mat_mumps_icntl_13 <0>: ICNTL(13): parallelism of the root node (enable ScaLAPACK) and its splitting (None) 18414e34a73bSHong Zhang . -mat_mumps_icntl_14 <20>: ICNTL(14): percentage increase in the estimated working space (None) 18424e34a73bSHong Zhang . -mat_mumps_icntl_19 <0>: ICNTL(19): computes the Schur complement (None) 18434e34a73bSHong Zhang . -mat_mumps_icntl_22 <0>: ICNTL(22): in-core/out-of-core factorization and solve (0 or 1) (None) 18444e34a73bSHong Zhang . -mat_mumps_icntl_23 <0>: ICNTL(23): max size of the working memory (MB) that can allocate per processor (None) 18454e34a73bSHong Zhang . -mat_mumps_icntl_24 <0>: ICNTL(24): detection of null pivot rows (0 or 1) (None) 18464e34a73bSHong Zhang . -mat_mumps_icntl_25 <0>: ICNTL(25): compute a solution of a deficient matrix and a null space basis (None) 18474e34a73bSHong Zhang . -mat_mumps_icntl_26 <0>: ICNTL(26): drives the solution phase if a Schur complement matrix (None) 18484e34a73bSHong Zhang . -mat_mumps_icntl_28 <1>: ICNTL(28): use 1 for sequential analysis and ictnl(7) ordering, or 2 for parallel analysis and ictnl(29) ordering (None) 18494e34a73bSHong Zhang . -mat_mumps_icntl_29 <0>: ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis (None) 18504e34a73bSHong Zhang . -mat_mumps_icntl_30 <0>: ICNTL(30): compute user-specified set of entries in inv(A) (None) 18514e34a73bSHong Zhang . -mat_mumps_icntl_31 <0>: ICNTL(31): indicates which factors may be discarded during factorization (None) 18524e34a73bSHong Zhang . -mat_mumps_icntl_33 <0>: ICNTL(33): compute determinant (None) 18534e34a73bSHong Zhang . -mat_mumps_cntl_1 <0.01>: CNTL(1): relative pivoting threshold (None) 18544e34a73bSHong Zhang . -mat_mumps_cntl_2 <1.49012e-08>: CNTL(2): stopping criterion of refinement (None) 18554e34a73bSHong Zhang . -mat_mumps_cntl_3 <0>: CNTL(3): absolute pivoting threshold (None) 18564e34a73bSHong Zhang . -mat_mumps_cntl_4 <-1>: CNTL(4): value for static pivoting (None) 18574e34a73bSHong Zhang - -mat_mumps_cntl_5 <0>: CNTL(5): fixation for null pivots (None) 185824b6179bSKris Buschelman 185924b6179bSKris Buschelman Level: beginner 186024b6179bSKris Buschelman 186141c8de11SBarry Smith .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage 186241c8de11SBarry Smith 186324b6179bSKris Buschelman M*/ 186424b6179bSKris Buschelman 186535bd34faSBarry Smith #undef __FUNCT__ 186635bd34faSBarry Smith #define __FUNCT__ "MatFactorGetSolverPackage_mumps" 1867f7a08781SBarry Smith static PetscErrorCode MatFactorGetSolverPackage_mumps(Mat A,const MatSolverPackage *type) 186835bd34faSBarry Smith { 186935bd34faSBarry Smith PetscFunctionBegin; 18702692d6eeSBarry Smith *type = MATSOLVERMUMPS; 187135bd34faSBarry Smith PetscFunctionReturn(0); 187235bd34faSBarry Smith } 187335bd34faSBarry Smith 1874bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI AIJ matrices */ 18752877fffaSHong Zhang #undef __FUNCT__ 1876bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_aij_mumps" 18778cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_aij_mumps(Mat A,MatFactorType ftype,Mat *F) 18782877fffaSHong Zhang { 18792877fffaSHong Zhang Mat B; 18802877fffaSHong Zhang PetscErrorCode ierr; 18812877fffaSHong Zhang Mat_MUMPS *mumps; 1882ace3abfcSBarry Smith PetscBool isSeqAIJ; 18832877fffaSHong Zhang 18842877fffaSHong Zhang PetscFunctionBegin; 18852877fffaSHong Zhang /* Create the factorization matrix */ 1886251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr); 1887ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 18882877fffaSHong Zhang ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 18892877fffaSHong Zhang ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 1890bccb9932SShri Abhyankar if (isSeqAIJ) { 18910298fd71SBarry Smith ierr = MatSeqAIJSetPreallocation(B,0,NULL);CHKERRQ(ierr); 1892bccb9932SShri Abhyankar } else { 18930298fd71SBarry Smith ierr = MatMPIAIJSetPreallocation(B,0,NULL,0,NULL);CHKERRQ(ierr); 1894bccb9932SShri Abhyankar } 18952877fffaSHong Zhang 1896b00a9115SJed Brown ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr); 18972205254eSKarl Rupp 18982877fffaSHong Zhang B->ops->view = MatView_MUMPS; 189935bd34faSBarry Smith B->ops->getinfo = MatGetInfo_MUMPS; 190020be8e61SHong Zhang B->ops->getdiagonal = MatGetDiagonal_MUMPS; 19012205254eSKarl Rupp 1902bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr); 1903bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr); 1904bc6112feSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr); 1905bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr); 1906bc6112feSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr); 1907bc6112feSHong Zhang 1908ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr); 1909ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr); 1910ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr); 1911ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr); 1912450b117fSShri Abhyankar if (ftype == MAT_FACTOR_LU) { 1913450b117fSShri Abhyankar B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS; 1914d5f3da31SBarry Smith B->factortype = MAT_FACTOR_LU; 1915bccb9932SShri Abhyankar if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqaij; 1916bccb9932SShri Abhyankar else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpiaij; 1917746480a1SHong Zhang mumps->sym = 0; 1918dcd589f8SShri Abhyankar } else { 191967877ebaSShri Abhyankar B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS; 1920450b117fSShri Abhyankar B->factortype = MAT_FACTOR_CHOLESKY; 1921bccb9932SShri Abhyankar if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqsbaij; 1922bccb9932SShri Abhyankar else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpisbaij; 19236fdc2a6dSBarry Smith if (A->spd_set && A->spd) mumps->sym = 1; 19246fdc2a6dSBarry Smith else mumps->sym = 2; 1925450b117fSShri Abhyankar } 19262877fffaSHong Zhang 19272877fffaSHong Zhang mumps->isAIJ = PETSC_TRUE; 1928bf0cc555SLisandro Dalcin mumps->Destroy = B->ops->destroy; 19292877fffaSHong Zhang B->ops->destroy = MatDestroy_MUMPS; 19302877fffaSHong Zhang B->spptr = (void*)mumps; 19312205254eSKarl Rupp 1932f697e70eSHong Zhang ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr); 1933746480a1SHong Zhang 19342877fffaSHong Zhang *F = B; 19352877fffaSHong Zhang PetscFunctionReturn(0); 19362877fffaSHong Zhang } 19372877fffaSHong Zhang 1938bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI SBAIJ matrices */ 19392877fffaSHong Zhang #undef __FUNCT__ 1940bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_sbaij_mumps" 19418cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_sbaij_mumps(Mat A,MatFactorType ftype,Mat *F) 19422877fffaSHong Zhang { 19432877fffaSHong Zhang Mat B; 19442877fffaSHong Zhang PetscErrorCode ierr; 19452877fffaSHong Zhang Mat_MUMPS *mumps; 1946ace3abfcSBarry Smith PetscBool isSeqSBAIJ; 19472877fffaSHong Zhang 19482877fffaSHong Zhang PetscFunctionBegin; 1949ce94432eSBarry Smith if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix"); 1950ce94432eSBarry Smith if (A->rmap->bs > 1) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with block size > 1 with MUMPS Cholesky, use AIJ matrix instead"); 1951251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr); 19522877fffaSHong Zhang /* Create the factorization matrix */ 1953ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 19542877fffaSHong Zhang ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 19552877fffaSHong Zhang ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 1956b00a9115SJed Brown ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr); 1957bccb9932SShri Abhyankar if (isSeqSBAIJ) { 19580298fd71SBarry Smith ierr = MatSeqSBAIJSetPreallocation(B,1,0,NULL);CHKERRQ(ierr); 19592205254eSKarl Rupp 196016ebf90aSShri Abhyankar mumps->ConvertToTriples = MatConvertToTriples_seqsbaij_seqsbaij; 1961dcd589f8SShri Abhyankar } else { 19620298fd71SBarry Smith ierr = MatMPISBAIJSetPreallocation(B,1,0,NULL,0,NULL);CHKERRQ(ierr); 19632205254eSKarl Rupp 1964bccb9932SShri Abhyankar mumps->ConvertToTriples = MatConvertToTriples_mpisbaij_mpisbaij; 1965bccb9932SShri Abhyankar } 1966bccb9932SShri Abhyankar 196767877ebaSShri Abhyankar B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS; 1968bccb9932SShri Abhyankar B->ops->view = MatView_MUMPS; 196920be8e61SHong Zhang B->ops->getdiagonal = MatGetDiagonal_MUMPS; 19702205254eSKarl Rupp 1971bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr); 1972b13644aeSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr); 1973b13644aeSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr); 1974b13644aeSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr); 1975b13644aeSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr); 1976bc6112feSHong Zhang 1977ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr); 1978ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr); 1979ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr); 1980ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr); 19812205254eSKarl Rupp 1982f4762488SHong Zhang B->factortype = MAT_FACTOR_CHOLESKY; 19836fdc2a6dSBarry Smith if (A->spd_set && A->spd) mumps->sym = 1; 19846fdc2a6dSBarry Smith else mumps->sym = 2; 1985a214ac2aSShri Abhyankar 1986bccb9932SShri Abhyankar mumps->isAIJ = PETSC_FALSE; 1987bf0cc555SLisandro Dalcin mumps->Destroy = B->ops->destroy; 1988f3c0ef26SHong Zhang B->ops->destroy = MatDestroy_MUMPS; 19892877fffaSHong Zhang B->spptr = (void*)mumps; 19902205254eSKarl Rupp 1991f697e70eSHong Zhang ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr); 1992746480a1SHong Zhang 19932877fffaSHong Zhang *F = B; 19942877fffaSHong Zhang PetscFunctionReturn(0); 19952877fffaSHong Zhang } 199697969023SHong Zhang 1997450b117fSShri Abhyankar #undef __FUNCT__ 1998bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_baij_mumps" 19998cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_baij_mumps(Mat A,MatFactorType ftype,Mat *F) 200067877ebaSShri Abhyankar { 200167877ebaSShri Abhyankar Mat B; 200267877ebaSShri Abhyankar PetscErrorCode ierr; 200367877ebaSShri Abhyankar Mat_MUMPS *mumps; 2004ace3abfcSBarry Smith PetscBool isSeqBAIJ; 200567877ebaSShri Abhyankar 200667877ebaSShri Abhyankar PetscFunctionBegin; 200767877ebaSShri Abhyankar /* Create the factorization matrix */ 2008251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);CHKERRQ(ierr); 2009ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 201067877ebaSShri Abhyankar ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 201167877ebaSShri Abhyankar ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 2012bccb9932SShri Abhyankar if (isSeqBAIJ) { 20130298fd71SBarry Smith ierr = MatSeqBAIJSetPreallocation(B,A->rmap->bs,0,NULL);CHKERRQ(ierr); 2014bccb9932SShri Abhyankar } else { 20150298fd71SBarry Smith ierr = MatMPIBAIJSetPreallocation(B,A->rmap->bs,0,NULL,0,NULL);CHKERRQ(ierr); 2016bccb9932SShri Abhyankar } 2017450b117fSShri Abhyankar 2018b00a9115SJed Brown ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr); 2019450b117fSShri Abhyankar if (ftype == MAT_FACTOR_LU) { 2020450b117fSShri Abhyankar B->ops->lufactorsymbolic = MatLUFactorSymbolic_BAIJMUMPS; 2021450b117fSShri Abhyankar B->factortype = MAT_FACTOR_LU; 2022bccb9932SShri Abhyankar if (isSeqBAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqbaij_seqaij; 2023bccb9932SShri Abhyankar else mumps->ConvertToTriples = MatConvertToTriples_mpibaij_mpiaij; 2024746480a1SHong Zhang mumps->sym = 0; 2025f23aa3ddSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc BAIJ matrices with MUMPS Cholesky, use SBAIJ or AIJ matrix instead\n"); 2026bccb9932SShri Abhyankar 2027450b117fSShri Abhyankar B->ops->view = MatView_MUMPS; 202820be8e61SHong Zhang B->ops->getdiagonal = MatGetDiagonal_MUMPS; 20292205254eSKarl Rupp 2030bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr); 2031bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr); 2032bc6112feSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr); 2033bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr); 2034bc6112feSHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr); 2035bc6112feSHong Zhang 2036ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr); 2037ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr); 2038ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr); 2039ca810319SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr); 2040450b117fSShri Abhyankar 2041450b117fSShri Abhyankar mumps->isAIJ = PETSC_TRUE; 2042bf0cc555SLisandro Dalcin mumps->Destroy = B->ops->destroy; 2043450b117fSShri Abhyankar B->ops->destroy = MatDestroy_MUMPS; 2044450b117fSShri Abhyankar B->spptr = (void*)mumps; 20452205254eSKarl Rupp 2046f697e70eSHong Zhang ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr); 2047746480a1SHong Zhang 2048450b117fSShri Abhyankar *F = B; 2049450b117fSShri Abhyankar PetscFunctionReturn(0); 2050450b117fSShri Abhyankar } 205142c9c57cSBarry Smith 205242c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor_aij_mumps(Mat,MatFactorType,Mat*); 205342c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor_baij_mumps(Mat,MatFactorType,Mat*); 205442c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor_sbaij_mumps(Mat,MatFactorType,Mat*); 205542c9c57cSBarry Smith 205642c9c57cSBarry Smith #undef __FUNCT__ 205742c9c57cSBarry Smith #define __FUNCT__ "MatSolverPackageRegister_MUMPS" 205829b38603SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_MUMPS(void) 205942c9c57cSBarry Smith { 206042c9c57cSBarry Smith PetscErrorCode ierr; 206142c9c57cSBarry Smith 206242c9c57cSBarry Smith PetscFunctionBegin; 206342c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIAIJ, MAT_FACTOR_LU,MatGetFactor_aij_mumps);CHKERRQ(ierr); 206442c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIAIJ, MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);CHKERRQ(ierr); 206542c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIBAIJ, MAT_FACTOR_LU,MatGetFactor_baij_mumps);CHKERRQ(ierr); 206642c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIBAIJ, MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);CHKERRQ(ierr); 206742c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPISBAIJ, MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);CHKERRQ(ierr); 206842c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQAIJ, MAT_FACTOR_LU,MatGetFactor_aij_mumps);CHKERRQ(ierr); 206942c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQAIJ, MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);CHKERRQ(ierr); 207042c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQBAIJ, MAT_FACTOR_LU,MatGetFactor_baij_mumps);CHKERRQ(ierr); 207142c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQBAIJ, MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);CHKERRQ(ierr); 207242c9c57cSBarry Smith ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQSBAIJ, MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);CHKERRQ(ierr); 207342c9c57cSBarry Smith PetscFunctionReturn(0); 207442c9c57cSBarry Smith } 207542c9c57cSBarry Smith 2076