xref: /petsc/src/mat/impls/aij/mpi/mumps/mumps.c (revision 2205254efee3a00a594e5e2a3a70f74dcb40bc03)
11c2a3de1SBarry Smith 
2397b6df1SKris Buschelman /*
3c2b5dc30SHong Zhang     Provides an interface to the MUMPS sparse solver
4397b6df1SKris Buschelman */
551d5961aSHong Zhang 
6c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h> /*I  "petscmat.h"  I*/
7c6db04a5SJed Brown #include <../src/mat/impls/sbaij/mpi/mpisbaij.h>
8397b6df1SKris Buschelman 
9397b6df1SKris Buschelman EXTERN_C_BEGIN
10397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
112907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
122907cef9SHong Zhang #include <cmumps_c.h>
132907cef9SHong Zhang #else
14c6db04a5SJed Brown #include <zmumps_c.h>
152907cef9SHong Zhang #endif
162907cef9SHong Zhang #else
172907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
182907cef9SHong Zhang #include <smumps_c.h>
19397b6df1SKris Buschelman #else
20c6db04a5SJed Brown #include <dmumps_c.h>
21397b6df1SKris Buschelman #endif
222907cef9SHong Zhang #endif
23397b6df1SKris Buschelman EXTERN_C_END
24397b6df1SKris Buschelman #define JOB_INIT -1
253d472b54SHong Zhang #define JOB_FACTSYMBOLIC 1
263d472b54SHong Zhang #define JOB_FACTNUMERIC 2
273d472b54SHong Zhang #define JOB_SOLVE 3
28397b6df1SKris Buschelman #define JOB_END -2
293d472b54SHong Zhang 
302907cef9SHong Zhang /* calls to MUMPS */
312907cef9SHong Zhang #if defined(PETSC_USE_COMPLEX)
322907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
332907cef9SHong Zhang #define PetscMUMPS_c cmumps_c
342907cef9SHong Zhang #else
352907cef9SHong Zhang #define PetscMUMPS_c zmumps_c
362907cef9SHong Zhang #endif
372907cef9SHong Zhang #else
382907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
392907cef9SHong Zhang #define PetscMUMPS_c smumps_c
402907cef9SHong Zhang #else
412907cef9SHong Zhang #define PetscMUMPS_c dmumps_c
422907cef9SHong Zhang #endif
432907cef9SHong Zhang #endif
442907cef9SHong Zhang 
453d472b54SHong Zhang 
46397b6df1SKris Buschelman /* macros s.t. indices match MUMPS documentation */
47397b6df1SKris Buschelman #define ICNTL(I) icntl[(I)-1]
48397b6df1SKris Buschelman #define CNTL(I) cntl[(I)-1]
49397b6df1SKris Buschelman #define INFOG(I) infog[(I)-1]
50a7aca84bSHong Zhang #define INFO(I) info[(I)-1]
51397b6df1SKris Buschelman #define RINFOG(I) rinfog[(I)-1]
52adc1d99fSHong Zhang #define RINFO(I) rinfo[(I)-1]
53397b6df1SKris Buschelman 
54397b6df1SKris Buschelman typedef struct {
55397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
562907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
572907cef9SHong Zhang   CMUMPS_STRUC_C id;
582907cef9SHong Zhang #else
59397b6df1SKris Buschelman   ZMUMPS_STRUC_C id;
602907cef9SHong Zhang #endif
612907cef9SHong Zhang #else
622907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
632907cef9SHong Zhang   SMUMPS_STRUC_C id;
64397b6df1SKris Buschelman #else
65397b6df1SKris Buschelman   DMUMPS_STRUC_C id;
66397b6df1SKris Buschelman #endif
672907cef9SHong Zhang #endif
682907cef9SHong Zhang 
69397b6df1SKris Buschelman   MatStructure matstruc;
70c1490034SHong Zhang   PetscMPIInt  myid,size;
71a5e57a09SHong Zhang   PetscInt     *irn,*jcn,nz,sym;
72397b6df1SKris Buschelman   PetscScalar  *val;
73397b6df1SKris Buschelman   MPI_Comm     comm_mumps;
74329ec9b3SHong Zhang   VecScatter   scat_rhs, scat_sol;
7564e6c443SBarry Smith   PetscBool    isAIJ,CleanUpMUMPS;
76329ec9b3SHong Zhang   Vec          b_seq,x_seq;
77a5e57a09SHong Zhang   PetscInt     ICNTL9_pre;   /* check if ICNTL(9) is changed from previous MatSolve */
78*2205254eSKarl Rupp 
79bf0cc555SLisandro Dalcin   PetscErrorCode (*Destroy)(Mat);
80bccb9932SShri Abhyankar   PetscErrorCode (*ConvertToTriples)(Mat, int, MatReuse, int*, int**, int**, PetscScalar**);
81f0c56d0fSKris Buschelman } Mat_MUMPS;
82f0c56d0fSKris Buschelman 
8309573ac7SBarry Smith extern PetscErrorCode MatDuplicate_MUMPS(Mat,MatDuplicateOption,Mat*);
84b24902e0SBarry Smith 
8567877ebaSShri Abhyankar 
8667877ebaSShri Abhyankar /* MatConvertToTriples_A_B */
8767877ebaSShri Abhyankar /*convert Petsc matrix to triples: row[nz], col[nz], val[nz] */
88397b6df1SKris Buschelman /*
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)
97397b6df1SKris Buschelman  */
9816ebf90aSShri Abhyankar 
9916ebf90aSShri Abhyankar #undef __FUNCT__
10016ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqaij_seqaij"
101bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
102b24902e0SBarry Smith {
103185f6596SHong Zhang   const PetscInt *ai,*aj,*ajj,M=A->rmap->n;
10467877ebaSShri Abhyankar   PetscInt       nz,rnz,i,j;
105dfbe8321SBarry Smith   PetscErrorCode ierr;
106c1490034SHong Zhang   PetscInt       *row,*col;
10716ebf90aSShri Abhyankar   Mat_SeqAIJ     *aa=(Mat_SeqAIJ*)A->data;
108397b6df1SKris Buschelman 
109397b6df1SKris Buschelman   PetscFunctionBegin;
11016ebf90aSShri Abhyankar   *v=aa->a;
111bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
112*2205254eSKarl Rupp     nz   = aa->nz;
113*2205254eSKarl Rupp     ai   = aa->i;
114*2205254eSKarl Rupp     aj   = aa->j;
11516ebf90aSShri Abhyankar     *nnz = nz;
116185f6596SHong Zhang     ierr = PetscMalloc(2*nz*sizeof(PetscInt), &row);CHKERRQ(ierr);
117185f6596SHong Zhang     col  = row + nz;
118185f6596SHong Zhang 
11916ebf90aSShri Abhyankar     nz = 0;
12016ebf90aSShri Abhyankar     for (i=0; i<M; i++) {
12116ebf90aSShri Abhyankar       rnz = ai[i+1] - ai[i];
12267877ebaSShri Abhyankar       ajj = aj + ai[i];
12367877ebaSShri Abhyankar       for (j=0; j<rnz; j++) {
12467877ebaSShri Abhyankar         row[nz] = i+shift; col[nz++] = ajj[j] + shift;
12516ebf90aSShri Abhyankar       }
12616ebf90aSShri Abhyankar     }
12716ebf90aSShri Abhyankar     *r = row; *c = col;
12816ebf90aSShri Abhyankar   }
12916ebf90aSShri Abhyankar   PetscFunctionReturn(0);
13016ebf90aSShri Abhyankar }
131397b6df1SKris Buschelman 
13216ebf90aSShri Abhyankar #undef __FUNCT__
13367877ebaSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqbaij_seqaij"
134bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqbaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
13567877ebaSShri Abhyankar {
13667877ebaSShri Abhyankar   Mat_SeqBAIJ    *aa=(Mat_SeqBAIJ*)A->data;
13767877ebaSShri Abhyankar   const PetscInt *ai,*aj,*ajj,bs=A->rmap->bs,bs2=aa->bs2,M=A->rmap->N/bs;
1380ad0caddSJed Brown   PetscInt       nz,idx=0,rnz,i,j,k,m;
13967877ebaSShri Abhyankar   PetscErrorCode ierr;
14067877ebaSShri Abhyankar   PetscInt       *row,*col;
14167877ebaSShri Abhyankar 
14267877ebaSShri Abhyankar   PetscFunctionBegin;
143cf3759fdSShri Abhyankar   *v = aa->a;
144bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
145cf3759fdSShri Abhyankar     ai   = aa->i; aj = aa->j;
14667877ebaSShri Abhyankar     nz   = bs2*aa->nz;
14767877ebaSShri Abhyankar     *nnz = nz;
148185f6596SHong Zhang     ierr = PetscMalloc(2*nz*sizeof(PetscInt), &row);CHKERRQ(ierr);
149185f6596SHong Zhang     col  = row + nz;
150185f6596SHong Zhang 
15167877ebaSShri Abhyankar     for (i=0; i<M; i++) {
15267877ebaSShri Abhyankar       ajj = aj + ai[i];
15367877ebaSShri Abhyankar       rnz = ai[i+1] - ai[i];
15467877ebaSShri Abhyankar       for (k=0; k<rnz; k++) {
15567877ebaSShri Abhyankar         for (j=0; j<bs; j++) {
15667877ebaSShri Abhyankar           for (m=0; m<bs; m++) {
15767877ebaSShri Abhyankar             row[idx]   = i*bs + m + shift;
158cf3759fdSShri Abhyankar             col[idx++] = bs*(ajj[k]) + j + shift;
15967877ebaSShri Abhyankar           }
16067877ebaSShri Abhyankar         }
16167877ebaSShri Abhyankar       }
16267877ebaSShri Abhyankar     }
163cf3759fdSShri Abhyankar     *r = row; *c = col;
16467877ebaSShri Abhyankar   }
16567877ebaSShri Abhyankar   PetscFunctionReturn(0);
16667877ebaSShri Abhyankar }
16767877ebaSShri Abhyankar 
16867877ebaSShri Abhyankar #undef __FUNCT__
16916ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqsbaij_seqsbaij"
170bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqsbaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
17116ebf90aSShri Abhyankar {
17267877ebaSShri Abhyankar   const PetscInt *ai, *aj,*ajj,M=A->rmap->n;
17367877ebaSShri Abhyankar   PetscInt       nz,rnz,i,j;
17416ebf90aSShri Abhyankar   PetscErrorCode ierr;
17516ebf90aSShri Abhyankar   PetscInt       *row,*col;
17616ebf90aSShri Abhyankar   Mat_SeqSBAIJ   *aa=(Mat_SeqSBAIJ*)A->data;
17716ebf90aSShri Abhyankar 
17816ebf90aSShri Abhyankar   PetscFunctionBegin;
179882afa5aSHong Zhang   *v = aa->a;
180bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
181*2205254eSKarl Rupp     nz   = aa->nz;
182*2205254eSKarl Rupp     ai   = aa->i;
183*2205254eSKarl Rupp     aj   = aa->j;
184*2205254eSKarl Rupp     *v   = aa->a;
18516ebf90aSShri Abhyankar     *nnz = nz;
186185f6596SHong Zhang     ierr = PetscMalloc(2*nz*sizeof(PetscInt), &row);CHKERRQ(ierr);
187185f6596SHong Zhang     col  = row + nz;
188185f6596SHong Zhang 
18916ebf90aSShri Abhyankar     nz = 0;
19016ebf90aSShri Abhyankar     for (i=0; i<M; i++) {
19116ebf90aSShri Abhyankar       rnz = ai[i+1] - ai[i];
19267877ebaSShri Abhyankar       ajj = aj + ai[i];
19367877ebaSShri Abhyankar       for (j=0; j<rnz; j++) {
19467877ebaSShri Abhyankar         row[nz] = i+shift; col[nz++] = ajj[j] + shift;
19516ebf90aSShri Abhyankar       }
19616ebf90aSShri Abhyankar     }
19716ebf90aSShri Abhyankar     *r = row; *c = col;
19816ebf90aSShri Abhyankar   }
19916ebf90aSShri Abhyankar   PetscFunctionReturn(0);
20016ebf90aSShri Abhyankar }
20116ebf90aSShri Abhyankar 
20216ebf90aSShri Abhyankar #undef __FUNCT__
20316ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_seqaij_seqsbaij"
204bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
20516ebf90aSShri Abhyankar {
20667877ebaSShri Abhyankar   const PetscInt    *ai,*aj,*ajj,*adiag,M=A->rmap->n;
20767877ebaSShri Abhyankar   PetscInt          nz,rnz,i,j;
20867877ebaSShri Abhyankar   const PetscScalar *av,*v1;
20916ebf90aSShri Abhyankar   PetscScalar       *val;
21016ebf90aSShri Abhyankar   PetscErrorCode    ierr;
21116ebf90aSShri Abhyankar   PetscInt          *row,*col;
21216ebf90aSShri Abhyankar   Mat_SeqSBAIJ      *aa=(Mat_SeqSBAIJ*)A->data;
21316ebf90aSShri Abhyankar 
21416ebf90aSShri Abhyankar   PetscFunctionBegin;
21516ebf90aSShri Abhyankar   ai   =aa->i; aj=aa->j;av=aa->a;
21616ebf90aSShri Abhyankar   adiag=aa->diag;
217bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
21816ebf90aSShri Abhyankar     nz   = M + (aa->nz-M)/2;
21916ebf90aSShri Abhyankar     *nnz = nz;
220185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
221185f6596SHong Zhang     col  = row + nz;
222185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
223185f6596SHong Zhang 
22416ebf90aSShri Abhyankar     nz = 0;
22516ebf90aSShri Abhyankar     for (i=0; i<M; i++) {
22616ebf90aSShri Abhyankar       rnz = ai[i+1] - adiag[i];
22767877ebaSShri Abhyankar       ajj = aj + adiag[i];
228cf3759fdSShri Abhyankar       v1  = av + adiag[i];
22967877ebaSShri Abhyankar       for (j=0; j<rnz; j++) {
23067877ebaSShri Abhyankar         row[nz] = i+shift; col[nz] = ajj[j] + shift; val[nz++] = v1[j];
23116ebf90aSShri Abhyankar       }
23216ebf90aSShri Abhyankar     }
23316ebf90aSShri Abhyankar     *r = row; *c = col; *v = val;
234397b6df1SKris Buschelman   } else {
23516ebf90aSShri Abhyankar     nz = 0; val = *v;
23616ebf90aSShri Abhyankar     for (i=0; i <M; i++) {
23716ebf90aSShri Abhyankar       rnz = ai[i+1] - adiag[i];
23867877ebaSShri Abhyankar       ajj = aj + adiag[i];
23967877ebaSShri Abhyankar       v1  = av + adiag[i];
24067877ebaSShri Abhyankar       for (j=0; j<rnz; j++) {
24167877ebaSShri Abhyankar         val[nz++] = v1[j];
24216ebf90aSShri Abhyankar       }
24316ebf90aSShri Abhyankar     }
24416ebf90aSShri Abhyankar   }
24516ebf90aSShri Abhyankar   PetscFunctionReturn(0);
24616ebf90aSShri Abhyankar }
24716ebf90aSShri Abhyankar 
24816ebf90aSShri Abhyankar #undef __FUNCT__
24916ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpisbaij_mpisbaij"
250bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpisbaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
25116ebf90aSShri Abhyankar {
25216ebf90aSShri Abhyankar   const PetscInt    *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
25316ebf90aSShri Abhyankar   PetscErrorCode    ierr;
25416ebf90aSShri Abhyankar   PetscInt          rstart,nz,i,j,jj,irow,countA,countB;
25516ebf90aSShri Abhyankar   PetscInt          *row,*col;
25616ebf90aSShri Abhyankar   const PetscScalar *av, *bv,*v1,*v2;
25716ebf90aSShri Abhyankar   PetscScalar       *val;
258397b6df1SKris Buschelman   Mat_MPISBAIJ      *mat = (Mat_MPISBAIJ*)A->data;
259397b6df1SKris Buschelman   Mat_SeqSBAIJ      *aa  = (Mat_SeqSBAIJ*)(mat->A)->data;
260397b6df1SKris Buschelman   Mat_SeqBAIJ       *bb  = (Mat_SeqBAIJ*)(mat->B)->data;
26116ebf90aSShri Abhyankar 
26216ebf90aSShri Abhyankar   PetscFunctionBegin;
263d0f46423SBarry Smith   ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart;
264397b6df1SKris Buschelman   av=aa->a; bv=bb->a;
265397b6df1SKris Buschelman 
266*2205254eSKarl Rupp   garray = mat->garray;
267*2205254eSKarl Rupp 
268bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
26916ebf90aSShri Abhyankar     nz   = aa->nz + bb->nz;
27016ebf90aSShri Abhyankar     *nnz = nz;
271185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
272185f6596SHong Zhang     col  = row + nz;
273185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
274185f6596SHong Zhang 
275397b6df1SKris Buschelman     *r = row; *c = col; *v = val;
276397b6df1SKris Buschelman   } else {
277397b6df1SKris Buschelman     row = *r; col = *c; val = *v;
278397b6df1SKris Buschelman   }
279397b6df1SKris Buschelman 
280028e57e8SHong Zhang   jj = 0; irow = rstart;
281397b6df1SKris Buschelman   for (i=0; i<m; i++) {
282397b6df1SKris Buschelman     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
283397b6df1SKris Buschelman     countA = ai[i+1] - ai[i];
284397b6df1SKris Buschelman     countB = bi[i+1] - bi[i];
285397b6df1SKris Buschelman     bjj    = bj + bi[i];
28616ebf90aSShri Abhyankar     v1     = av + ai[i];
28716ebf90aSShri Abhyankar     v2     = bv + bi[i];
288397b6df1SKris Buschelman 
289397b6df1SKris Buschelman     /* A-part */
290397b6df1SKris Buschelman     for (j=0; j<countA; j++) {
291bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
292397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
293397b6df1SKris Buschelman       }
29416ebf90aSShri Abhyankar       val[jj++] = v1[j];
295397b6df1SKris Buschelman     }
29616ebf90aSShri Abhyankar 
29716ebf90aSShri Abhyankar     /* B-part */
29816ebf90aSShri Abhyankar     for (j=0; j < countB; j++) {
299bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
300397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
301397b6df1SKris Buschelman       }
30216ebf90aSShri Abhyankar       val[jj++] = v2[j];
30316ebf90aSShri Abhyankar     }
30416ebf90aSShri Abhyankar     irow++;
30516ebf90aSShri Abhyankar   }
30616ebf90aSShri Abhyankar   PetscFunctionReturn(0);
30716ebf90aSShri Abhyankar }
30816ebf90aSShri Abhyankar 
30916ebf90aSShri Abhyankar #undef __FUNCT__
31016ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpiaij_mpiaij"
311bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpiaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
31216ebf90aSShri Abhyankar {
31316ebf90aSShri Abhyankar   const PetscInt    *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
31416ebf90aSShri Abhyankar   PetscErrorCode    ierr;
31516ebf90aSShri Abhyankar   PetscInt          rstart,nz,i,j,jj,irow,countA,countB;
31616ebf90aSShri Abhyankar   PetscInt          *row,*col;
31716ebf90aSShri Abhyankar   const PetscScalar *av, *bv,*v1,*v2;
31816ebf90aSShri Abhyankar   PetscScalar       *val;
31916ebf90aSShri Abhyankar   Mat_MPIAIJ        *mat = (Mat_MPIAIJ*)A->data;
32016ebf90aSShri Abhyankar   Mat_SeqAIJ        *aa  = (Mat_SeqAIJ*)(mat->A)->data;
32116ebf90aSShri Abhyankar   Mat_SeqAIJ        *bb  = (Mat_SeqAIJ*)(mat->B)->data;
32216ebf90aSShri Abhyankar 
32316ebf90aSShri Abhyankar   PetscFunctionBegin;
32416ebf90aSShri Abhyankar   ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart;
32516ebf90aSShri Abhyankar   av=aa->a; bv=bb->a;
32616ebf90aSShri Abhyankar 
327*2205254eSKarl Rupp   garray = mat->garray;
328*2205254eSKarl Rupp 
329bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
33016ebf90aSShri Abhyankar     nz   = aa->nz + bb->nz;
33116ebf90aSShri Abhyankar     *nnz = nz;
332185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
333185f6596SHong Zhang     col  = row + nz;
334185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
335185f6596SHong Zhang 
33616ebf90aSShri Abhyankar     *r = row; *c = col; *v = val;
33716ebf90aSShri Abhyankar   } else {
33816ebf90aSShri Abhyankar     row = *r; col = *c; val = *v;
33916ebf90aSShri Abhyankar   }
34016ebf90aSShri Abhyankar 
34116ebf90aSShri Abhyankar   jj = 0; irow = rstart;
34216ebf90aSShri Abhyankar   for (i=0; i<m; i++) {
34316ebf90aSShri Abhyankar     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
34416ebf90aSShri Abhyankar     countA = ai[i+1] - ai[i];
34516ebf90aSShri Abhyankar     countB = bi[i+1] - bi[i];
34616ebf90aSShri Abhyankar     bjj    = bj + bi[i];
34716ebf90aSShri Abhyankar     v1     = av + ai[i];
34816ebf90aSShri Abhyankar     v2     = bv + bi[i];
34916ebf90aSShri Abhyankar 
35016ebf90aSShri Abhyankar     /* A-part */
35116ebf90aSShri Abhyankar     for (j=0; j<countA; j++) {
352bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
35316ebf90aSShri Abhyankar         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
35416ebf90aSShri Abhyankar       }
35516ebf90aSShri Abhyankar       val[jj++] = v1[j];
35616ebf90aSShri Abhyankar     }
35716ebf90aSShri Abhyankar 
35816ebf90aSShri Abhyankar     /* B-part */
35916ebf90aSShri Abhyankar     for (j=0; j < countB; j++) {
360bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
36116ebf90aSShri Abhyankar         row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
36216ebf90aSShri Abhyankar       }
36316ebf90aSShri Abhyankar       val[jj++] = v2[j];
36416ebf90aSShri Abhyankar     }
36516ebf90aSShri Abhyankar     irow++;
36616ebf90aSShri Abhyankar   }
36716ebf90aSShri Abhyankar   PetscFunctionReturn(0);
36816ebf90aSShri Abhyankar }
36916ebf90aSShri Abhyankar 
37016ebf90aSShri Abhyankar #undef __FUNCT__
37167877ebaSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpibaij_mpiaij"
372bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpibaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
37367877ebaSShri Abhyankar {
37467877ebaSShri Abhyankar   Mat_MPIBAIJ       *mat    = (Mat_MPIBAIJ*)A->data;
37567877ebaSShri Abhyankar   Mat_SeqBAIJ       *aa     = (Mat_SeqBAIJ*)(mat->A)->data;
37667877ebaSShri Abhyankar   Mat_SeqBAIJ       *bb     = (Mat_SeqBAIJ*)(mat->B)->data;
37767877ebaSShri Abhyankar   const PetscInt    *ai     = aa->i, *bi = bb->i, *aj = aa->j, *bj = bb->j,*ajj, *bjj;
378d985c460SShri Abhyankar   const PetscInt    *garray = mat->garray,mbs=mat->mbs,rstart=A->rmap->rstart;
37967877ebaSShri Abhyankar   const PetscInt    bs      = A->rmap->bs,bs2=mat->bs2;
38067877ebaSShri Abhyankar   PetscErrorCode    ierr;
38167877ebaSShri Abhyankar   PetscInt          nz,i,j,k,n,jj,irow,countA,countB,idx;
38267877ebaSShri Abhyankar   PetscInt          *row,*col;
38367877ebaSShri Abhyankar   const PetscScalar *av=aa->a, *bv=bb->a,*v1,*v2;
38467877ebaSShri Abhyankar   PetscScalar       *val;
38567877ebaSShri Abhyankar 
38667877ebaSShri Abhyankar   PetscFunctionBegin;
387bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
38867877ebaSShri Abhyankar     nz   = bs2*(aa->nz + bb->nz);
38967877ebaSShri Abhyankar     *nnz = nz;
390185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
391185f6596SHong Zhang     col  = row + nz;
392185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
393185f6596SHong Zhang 
39467877ebaSShri Abhyankar     *r = row; *c = col; *v = val;
39567877ebaSShri Abhyankar   } else {
39667877ebaSShri Abhyankar     row = *r; col = *c; val = *v;
39767877ebaSShri Abhyankar   }
39867877ebaSShri Abhyankar 
399d985c460SShri Abhyankar   jj = 0; irow = rstart;
40067877ebaSShri Abhyankar   for (i=0; i<mbs; i++) {
40167877ebaSShri Abhyankar     countA = ai[i+1] - ai[i];
40267877ebaSShri Abhyankar     countB = bi[i+1] - bi[i];
40367877ebaSShri Abhyankar     ajj    = aj + ai[i];
40467877ebaSShri Abhyankar     bjj    = bj + bi[i];
40567877ebaSShri Abhyankar     v1     = av + bs2*ai[i];
40667877ebaSShri Abhyankar     v2     = bv + bs2*bi[i];
40767877ebaSShri Abhyankar 
40867877ebaSShri Abhyankar     idx = 0;
40967877ebaSShri Abhyankar     /* A-part */
41067877ebaSShri Abhyankar     for (k=0; k<countA; k++) {
41167877ebaSShri Abhyankar       for (j=0; j<bs; j++) {
41267877ebaSShri Abhyankar         for (n=0; n<bs; n++) {
413bccb9932SShri Abhyankar           if (reuse == MAT_INITIAL_MATRIX) {
414d985c460SShri Abhyankar             row[jj] = irow + n + shift;
415d985c460SShri Abhyankar             col[jj] = rstart + bs*ajj[k] + j + shift;
41667877ebaSShri Abhyankar           }
41767877ebaSShri Abhyankar           val[jj++] = v1[idx++];
41867877ebaSShri Abhyankar         }
41967877ebaSShri Abhyankar       }
42067877ebaSShri Abhyankar     }
42167877ebaSShri Abhyankar 
42267877ebaSShri Abhyankar     idx = 0;
42367877ebaSShri Abhyankar     /* B-part */
42467877ebaSShri Abhyankar     for (k=0; k<countB; k++) {
42567877ebaSShri Abhyankar       for (j=0; j<bs; j++) {
42667877ebaSShri Abhyankar         for (n=0; n<bs; n++) {
427bccb9932SShri Abhyankar           if (reuse == MAT_INITIAL_MATRIX) {
428d985c460SShri Abhyankar             row[jj] = irow + n + shift;
429d985c460SShri Abhyankar             col[jj] = bs*garray[bjj[k]] + j + shift;
43067877ebaSShri Abhyankar           }
431d985c460SShri Abhyankar           val[jj++] = v2[idx++];
43267877ebaSShri Abhyankar         }
43367877ebaSShri Abhyankar       }
43467877ebaSShri Abhyankar     }
435d985c460SShri Abhyankar     irow += bs;
43667877ebaSShri Abhyankar   }
43767877ebaSShri Abhyankar   PetscFunctionReturn(0);
43867877ebaSShri Abhyankar }
43967877ebaSShri Abhyankar 
44067877ebaSShri Abhyankar #undef __FUNCT__
44116ebf90aSShri Abhyankar #define __FUNCT__ "MatConvertToTriples_mpiaij_mpisbaij"
442bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpiaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
44316ebf90aSShri Abhyankar {
44416ebf90aSShri Abhyankar   const PetscInt    *ai, *aj,*adiag, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
44516ebf90aSShri Abhyankar   PetscErrorCode    ierr;
446e0bace9bSHong Zhang   PetscInt          rstart,nz,nza,nzb,i,j,jj,irow,countA,countB;
44716ebf90aSShri Abhyankar   PetscInt          *row,*col;
44816ebf90aSShri Abhyankar   const PetscScalar *av, *bv,*v1,*v2;
44916ebf90aSShri Abhyankar   PetscScalar       *val;
45016ebf90aSShri Abhyankar   Mat_MPIAIJ        *mat =  (Mat_MPIAIJ*)A->data;
45116ebf90aSShri Abhyankar   Mat_SeqAIJ        *aa  =(Mat_SeqAIJ*)(mat->A)->data;
45216ebf90aSShri Abhyankar   Mat_SeqAIJ        *bb  =(Mat_SeqAIJ*)(mat->B)->data;
45316ebf90aSShri Abhyankar 
45416ebf90aSShri Abhyankar   PetscFunctionBegin;
45516ebf90aSShri Abhyankar   ai=aa->i; aj=aa->j; adiag=aa->diag;
45616ebf90aSShri Abhyankar   bi=bb->i; bj=bb->j; garray = mat->garray;
45716ebf90aSShri Abhyankar   av=aa->a; bv=bb->a;
458*2205254eSKarl Rupp 
45916ebf90aSShri Abhyankar   rstart = A->rmap->rstart;
46016ebf90aSShri Abhyankar 
461bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
462e0bace9bSHong Zhang     nza = 0;    /* num of upper triangular entries in mat->A, including diagonals */
463e0bace9bSHong Zhang     nzb = 0;    /* num of upper triangular entries in mat->B */
46416ebf90aSShri Abhyankar     for (i=0; i<m; i++) {
465e0bace9bSHong Zhang       nza   += (ai[i+1] - adiag[i]);
46616ebf90aSShri Abhyankar       countB = bi[i+1] - bi[i];
46716ebf90aSShri Abhyankar       bjj    = bj + bi[i];
468e0bace9bSHong Zhang       for (j=0; j<countB; j++) {
469e0bace9bSHong Zhang         if (garray[bjj[j]] > rstart) nzb++;
470e0bace9bSHong Zhang       }
471e0bace9bSHong Zhang     }
47216ebf90aSShri Abhyankar 
473e0bace9bSHong Zhang     nz = nza + nzb; /* total nz of upper triangular part of mat */
47416ebf90aSShri Abhyankar     *nnz = nz;
475185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
476185f6596SHong Zhang     col  = row + nz;
477185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
478185f6596SHong Zhang 
47916ebf90aSShri Abhyankar     *r = row; *c = col; *v = val;
48016ebf90aSShri Abhyankar   } else {
48116ebf90aSShri Abhyankar     row = *r; col = *c; val = *v;
48216ebf90aSShri Abhyankar   }
48316ebf90aSShri Abhyankar 
48416ebf90aSShri Abhyankar   jj = 0; irow = rstart;
48516ebf90aSShri Abhyankar   for (i=0; i<m; i++) {
48616ebf90aSShri Abhyankar     ajj    = aj + adiag[i];                 /* ptr to the beginning of the diagonal of this row */
48716ebf90aSShri Abhyankar     v1     = av + adiag[i];
48816ebf90aSShri Abhyankar     countA = ai[i+1] - adiag[i];
48916ebf90aSShri Abhyankar     countB = bi[i+1] - bi[i];
49016ebf90aSShri Abhyankar     bjj    = bj + bi[i];
49116ebf90aSShri Abhyankar     v2     = bv + bi[i];
49216ebf90aSShri Abhyankar 
49316ebf90aSShri Abhyankar     /* A-part */
49416ebf90aSShri Abhyankar     for (j=0; j<countA; j++) {
495bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
49616ebf90aSShri Abhyankar         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
49716ebf90aSShri Abhyankar       }
49816ebf90aSShri Abhyankar       val[jj++] = v1[j];
49916ebf90aSShri Abhyankar     }
50016ebf90aSShri Abhyankar 
50116ebf90aSShri Abhyankar     /* B-part */
50216ebf90aSShri Abhyankar     for (j=0; j < countB; j++) {
50316ebf90aSShri Abhyankar       if (garray[bjj[j]] > rstart) {
504bccb9932SShri Abhyankar         if (reuse == MAT_INITIAL_MATRIX) {
50516ebf90aSShri Abhyankar           row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
50616ebf90aSShri Abhyankar         }
50716ebf90aSShri Abhyankar         val[jj++] = v2[j];
50816ebf90aSShri Abhyankar       }
509397b6df1SKris Buschelman     }
510397b6df1SKris Buschelman     irow++;
511397b6df1SKris Buschelman   }
512397b6df1SKris Buschelman   PetscFunctionReturn(0);
513397b6df1SKris Buschelman }
514397b6df1SKris Buschelman 
515397b6df1SKris Buschelman #undef __FUNCT__
5163924e44cSKris Buschelman #define __FUNCT__ "MatDestroy_MUMPS"
517dfbe8321SBarry Smith PetscErrorCode MatDestroy_MUMPS(Mat A)
518dfbe8321SBarry Smith {
519a5e57a09SHong Zhang   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->spptr;
520dfbe8321SBarry Smith   PetscErrorCode ierr;
521b24902e0SBarry Smith 
522397b6df1SKris Buschelman   PetscFunctionBegin;
523a5e57a09SHong Zhang   if (mumps->CleanUpMUMPS) {
524397b6df1SKris Buschelman     /* Terminate instance, deallocate memories */
525a5e57a09SHong Zhang     ierr = PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);CHKERRQ(ierr);
526a5e57a09SHong Zhang     ierr = VecScatterDestroy(&mumps->scat_rhs);CHKERRQ(ierr);
527a5e57a09SHong Zhang     ierr = VecDestroy(&mumps->b_seq);CHKERRQ(ierr);
528a5e57a09SHong Zhang     ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr);
529a5e57a09SHong Zhang     ierr = VecDestroy(&mumps->x_seq);CHKERRQ(ierr);
530a5e57a09SHong Zhang     ierr = PetscFree(mumps->id.perm_in);CHKERRQ(ierr);
531a5e57a09SHong Zhang     ierr = PetscFree(mumps->irn);CHKERRQ(ierr);
532*2205254eSKarl Rupp 
533a5e57a09SHong Zhang     mumps->id.job = JOB_END;
534a5e57a09SHong Zhang     PetscMUMPS_c(&mumps->id);
535a5e57a09SHong Zhang     ierr = MPI_Comm_free(&(mumps->comm_mumps));CHKERRQ(ierr);
536397b6df1SKris Buschelman   }
537a5e57a09SHong Zhang   if (mumps->Destroy) {
538a5e57a09SHong Zhang     ierr = (mumps->Destroy)(A);CHKERRQ(ierr);
539bf0cc555SLisandro Dalcin   }
540bf0cc555SLisandro Dalcin   ierr = PetscFree(A->spptr);CHKERRQ(ierr);
541bf0cc555SLisandro Dalcin 
54297969023SHong Zhang   /* clear composed functions */
54397969023SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatFactorGetSolverPackage_C","",PETSC_NULL);CHKERRQ(ierr);
544f250808bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMumpsSetIcntl_C","",PETSC_NULL);CHKERRQ(ierr);
545397b6df1SKris Buschelman   PetscFunctionReturn(0);
546397b6df1SKris Buschelman }
547397b6df1SKris Buschelman 
548397b6df1SKris Buschelman #undef __FUNCT__
549f6c57405SHong Zhang #define __FUNCT__ "MatSolve_MUMPS"
550b24902e0SBarry Smith PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x)
551b24902e0SBarry Smith {
552a5e57a09SHong Zhang   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->spptr;
553d54de34fSKris Buschelman   PetscScalar    *array;
55467877ebaSShri Abhyankar   Vec            b_seq;
555329ec9b3SHong Zhang   IS             is_iden,is_petsc;
556dfbe8321SBarry Smith   PetscErrorCode ierr;
557329ec9b3SHong Zhang   PetscInt       i;
558397b6df1SKris Buschelman 
559397b6df1SKris Buschelman   PetscFunctionBegin;
560a5e57a09SHong Zhang   mumps->id.nrhs = 1;
561a5e57a09SHong Zhang   b_seq          = mumps->b_seq;
562a5e57a09SHong Zhang   if (mumps->size > 1) {
563329ec9b3SHong Zhang     /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */
564a5e57a09SHong Zhang     ierr = VecScatterBegin(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
565a5e57a09SHong Zhang     ierr = VecScatterEnd(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
566a5e57a09SHong Zhang     if (!mumps->myid) {ierr = VecGetArray(b_seq,&array);CHKERRQ(ierr);}
567397b6df1SKris Buschelman   } else {  /* size == 1 */
568397b6df1SKris Buschelman     ierr = VecCopy(b,x);CHKERRQ(ierr);
569397b6df1SKris Buschelman     ierr = VecGetArray(x,&array);CHKERRQ(ierr);
570397b6df1SKris Buschelman   }
571a5e57a09SHong Zhang   if (!mumps->myid) { /* define rhs on the host */
572a5e57a09SHong Zhang     mumps->id.nrhs = 1;
573397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
5742907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
575a5e57a09SHong Zhang     mumps->id.rhs = (mumps_complex*)array;
5762907cef9SHong Zhang #else
577a5e57a09SHong Zhang     mumps->id.rhs = (mumps_double_complex*)array;
5782907cef9SHong Zhang #endif
579397b6df1SKris Buschelman #else
580a5e57a09SHong Zhang     mumps->id.rhs = array;
581397b6df1SKris Buschelman #endif
582397b6df1SKris Buschelman   }
583397b6df1SKris Buschelman 
584397b6df1SKris Buschelman   /* solve phase */
585329ec9b3SHong Zhang   /*-------------*/
586a5e57a09SHong Zhang   mumps->id.job = JOB_SOLVE;
587a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
588a5e57a09SHong 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));
589397b6df1SKris Buschelman 
590a5e57a09SHong Zhang   if (mumps->size > 1) { /* convert mumps distributed solution to petsc mpi x */
591a5e57a09SHong Zhang     if (mumps->scat_sol && mumps->ICNTL9_pre != mumps->id.ICNTL(9)) {
592a5e57a09SHong Zhang       /* when id.ICNTL(9) changes, the contents of lsol_loc may change (not its size, lsol_loc), recreates scat_sol */
593a5e57a09SHong Zhang       ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr);
594397b6df1SKris Buschelman     }
595a5e57a09SHong Zhang     if (!mumps->scat_sol) { /* create scatter scat_sol */
596a5e57a09SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,mumps->id.lsol_loc,0,1,&is_iden);CHKERRQ(ierr); /* from */
597a5e57a09SHong Zhang       for (i=0; i<mumps->id.lsol_loc; i++) {
598a5e57a09SHong Zhang         mumps->id.isol_loc[i] -= 1; /* change Fortran style to C style */
599a5e57a09SHong Zhang       }
600a5e57a09SHong Zhang       ierr = ISCreateGeneral(PETSC_COMM_SELF,mumps->id.lsol_loc,mumps->id.isol_loc,PETSC_COPY_VALUES,&is_petsc);CHKERRQ(ierr);  /* to */
601a5e57a09SHong Zhang       ierr = VecScatterCreate(mumps->x_seq,is_iden,x,is_petsc,&mumps->scat_sol);CHKERRQ(ierr);
6026bf464f9SBarry Smith       ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
6036bf464f9SBarry Smith       ierr = ISDestroy(&is_petsc);CHKERRQ(ierr);
604*2205254eSKarl Rupp 
605a5e57a09SHong Zhang       mumps->ICNTL9_pre = mumps->id.ICNTL(9); /* save current value of id.ICNTL(9) */
606397b6df1SKris Buschelman     }
607a5e57a09SHong Zhang 
608a5e57a09SHong Zhang     ierr = VecScatterBegin(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
609a5e57a09SHong Zhang     ierr = VecScatterEnd(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
610329ec9b3SHong Zhang   }
611397b6df1SKris Buschelman   PetscFunctionReturn(0);
612397b6df1SKris Buschelman }
613397b6df1SKris Buschelman 
61451d5961aSHong Zhang #undef __FUNCT__
61551d5961aSHong Zhang #define __FUNCT__ "MatSolveTranspose_MUMPS"
61651d5961aSHong Zhang PetscErrorCode MatSolveTranspose_MUMPS(Mat A,Vec b,Vec x)
61751d5961aSHong Zhang {
618a5e57a09SHong Zhang   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->spptr;
61951d5961aSHong Zhang   PetscErrorCode ierr;
62051d5961aSHong Zhang 
62151d5961aSHong Zhang   PetscFunctionBegin;
622a5e57a09SHong Zhang   mumps->id.ICNTL(9) = 0;
623*2205254eSKarl Rupp 
6240ad0caddSJed Brown   ierr = MatSolve_MUMPS(A,b,x);CHKERRQ(ierr);
625*2205254eSKarl Rupp 
626a5e57a09SHong Zhang   mumps->id.ICNTL(9) = 1;
62751d5961aSHong Zhang   PetscFunctionReturn(0);
62851d5961aSHong Zhang }
62951d5961aSHong Zhang 
630e0b74bf9SHong Zhang #undef __FUNCT__
631e0b74bf9SHong Zhang #define __FUNCT__ "MatMatSolve_MUMPS"
632e0b74bf9SHong Zhang PetscErrorCode MatMatSolve_MUMPS(Mat A,Mat B,Mat X)
633e0b74bf9SHong Zhang {
634bda8bf91SBarry Smith   PetscErrorCode ierr;
635bda8bf91SBarry Smith   PetscBool      flg;
636bda8bf91SBarry Smith 
637e0b74bf9SHong Zhang   PetscFunctionBegin;
638251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,PETSC_NULL);CHKERRQ(ierr);
639bda8bf91SBarry Smith   if (!flg) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
640251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,PETSC_NULL);CHKERRQ(ierr);
641*2205254eSKarl Rupp   if (!flg) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");
642*2205254eSKarl Rupp   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatSolve_MUMPS() is not implemented yet");
643e0b74bf9SHong Zhang   PetscFunctionReturn(0);
644e0b74bf9SHong Zhang }
645e0b74bf9SHong Zhang 
646ace3df97SHong Zhang #if !defined(PETSC_USE_COMPLEX)
647a58c3f20SHong Zhang /*
648a58c3f20SHong Zhang   input:
649a58c3f20SHong Zhang    F:        numeric factor
650a58c3f20SHong Zhang   output:
651a58c3f20SHong Zhang    nneg:     total number of negative pivots
652a58c3f20SHong Zhang    nzero:    0
653a58c3f20SHong Zhang    npos:     (global dimension of F) - nneg
654a58c3f20SHong Zhang */
655a58c3f20SHong Zhang 
656a58c3f20SHong Zhang #undef __FUNCT__
657a58c3f20SHong Zhang #define __FUNCT__ "MatGetInertia_SBAIJMUMPS"
658dfbe8321SBarry Smith PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
659a58c3f20SHong Zhang {
660a5e57a09SHong Zhang   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->spptr;
661dfbe8321SBarry Smith   PetscErrorCode ierr;
662c1490034SHong Zhang   PetscMPIInt    size;
663a58c3f20SHong Zhang 
664a58c3f20SHong Zhang   PetscFunctionBegin;
6657adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)F)->comm,&size);CHKERRQ(ierr);
666bcb30aebSHong 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 */
667a5e57a09SHong 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));
668a58c3f20SHong Zhang   if (nneg) {
669a5e57a09SHong Zhang     if (!mumps->myid) {
670a5e57a09SHong Zhang       *nneg = mumps->id.INFOG(12);
671a58c3f20SHong Zhang     }
672a5e57a09SHong Zhang     ierr = MPI_Bcast(nneg,1,MPI_INT,0,mumps->comm_mumps);CHKERRQ(ierr);
673a58c3f20SHong Zhang   }
674a58c3f20SHong Zhang   if (nzero) *nzero = 0;
675d0f46423SBarry Smith   if (npos)  *npos  = F->rmap->N - (*nneg);
676a58c3f20SHong Zhang   PetscFunctionReturn(0);
677a58c3f20SHong Zhang }
678ace3df97SHong Zhang #endif /* !defined(PETSC_USE_COMPLEX) */
679a58c3f20SHong Zhang 
680397b6df1SKris Buschelman #undef __FUNCT__
681f6c57405SHong Zhang #define __FUNCT__ "MatFactorNumeric_MUMPS"
6820481f469SBarry Smith PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info)
683af281ebdSHong Zhang {
684a5e57a09SHong Zhang   Mat_MUMPS      *mumps =(Mat_MUMPS*)(F)->spptr;
6856849ba73SBarry Smith   PetscErrorCode ierr;
686e09efc27SHong Zhang   Mat            F_diag;
687ace3abfcSBarry Smith   PetscBool      isMPIAIJ;
688397b6df1SKris Buschelman 
689397b6df1SKris Buschelman   PetscFunctionBegin;
690a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_REUSE_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
691397b6df1SKris Buschelman 
692397b6df1SKris Buschelman   /* numerical factorization phase */
693329ec9b3SHong Zhang   /*-------------------------------*/
694a5e57a09SHong Zhang   mumps->id.job = JOB_FACTNUMERIC;
695a5e57a09SHong Zhang   if (!mumps->id.ICNTL(18)) {
696a5e57a09SHong Zhang     if (!mumps->myid) {
697397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
6982907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
699a5e57a09SHong Zhang       mumps->id.a = (mumps_complex*)mumps->val;
7002907cef9SHong Zhang #else
701a5e57a09SHong Zhang       mumps->id.a = (mumps_double_complex*)mumps->val;
7022907cef9SHong Zhang #endif
703397b6df1SKris Buschelman #else
704a5e57a09SHong Zhang       mumps->id.a = mumps->val;
705397b6df1SKris Buschelman #endif
706397b6df1SKris Buschelman     }
707397b6df1SKris Buschelman   } else {
708397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
7092907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
710a5e57a09SHong Zhang     mumps->id.a_loc = (mumps_complex*)mumps->val;
7112907cef9SHong Zhang #else
712a5e57a09SHong Zhang     mumps->id.a_loc = (mumps_double_complex*)mumps->val;
7132907cef9SHong Zhang #endif
714397b6df1SKris Buschelman #else
715a5e57a09SHong Zhang     mumps->id.a_loc = mumps->val;
716397b6df1SKris Buschelman #endif
717397b6df1SKris Buschelman   }
718a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
719a5e57a09SHong Zhang   if (mumps->id.INFOG(1) < 0) {
720a5e57a09SHong Zhang     if (mumps->id.INFO(1) == -13) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: Cannot allocate required memory %d megabytes\n",mumps->id.INFO(2));
721a5e57a09SHong 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));
722397b6df1SKris Buschelman   }
723a5e57a09SHong 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));
724397b6df1SKris Buschelman 
725a5e57a09SHong Zhang   if (mumps->size > 1) {
726251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr);
727*2205254eSKarl Rupp     if (isMPIAIJ) F_diag = ((Mat_MPIAIJ*)(F)->data)->A;
728*2205254eSKarl Rupp     else F_diag = ((Mat_MPISBAIJ*)(F)->data)->A;
729e09efc27SHong Zhang     F_diag->assembled = PETSC_TRUE;
730a5e57a09SHong Zhang     if (mumps->scat_sol) {
731a5e57a09SHong Zhang       ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr);
732a5e57a09SHong Zhang       ierr = PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);CHKERRQ(ierr);
733a5e57a09SHong Zhang       ierr = VecDestroy(&mumps->x_seq);CHKERRQ(ierr);
734329ec9b3SHong Zhang     }
7358ada1bb4SHong Zhang   }
736dcd589f8SShri Abhyankar   (F)->assembled      = PETSC_TRUE;
737a5e57a09SHong Zhang   mumps->matstruc     = SAME_NONZERO_PATTERN;
738a5e57a09SHong Zhang   mumps->CleanUpMUMPS = PETSC_TRUE;
73967877ebaSShri Abhyankar 
740a5e57a09SHong Zhang   if (mumps->size > 1) {
74167877ebaSShri Abhyankar     /* distributed solution */
742a5e57a09SHong Zhang     if (!mumps->scat_sol) {
74367877ebaSShri Abhyankar       /* Create x_seq=sol_loc for repeated use */
74467877ebaSShri Abhyankar       PetscInt    lsol_loc;
74567877ebaSShri Abhyankar       PetscScalar *sol_loc;
746*2205254eSKarl Rupp 
747a5e57a09SHong Zhang       lsol_loc = mumps->id.INFO(23); /* length of sol_loc */
748*2205254eSKarl Rupp 
749a5e57a09SHong Zhang       ierr = PetscMalloc2(lsol_loc,PetscScalar,&sol_loc,lsol_loc,PetscInt,&mumps->id.isol_loc);CHKERRQ(ierr);
750*2205254eSKarl Rupp 
751a5e57a09SHong Zhang       mumps->id.lsol_loc = lsol_loc;
75267877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
7532907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
754a5e57a09SHong Zhang       mumps->id.sol_loc = (mumps_complex*)sol_loc;
7552907cef9SHong Zhang #else
756a5e57a09SHong Zhang       mumps->id.sol_loc = (mumps_double_complex*)sol_loc;
7572907cef9SHong Zhang #endif
75867877ebaSShri Abhyankar #else
759a5e57a09SHong Zhang       mumps->id.sol_loc = sol_loc;
76067877ebaSShri Abhyankar #endif
761a5e57a09SHong Zhang       ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,lsol_loc,sol_loc,&mumps->x_seq);CHKERRQ(ierr);
76267877ebaSShri Abhyankar     }
76367877ebaSShri Abhyankar   }
764397b6df1SKris Buschelman   PetscFunctionReturn(0);
765397b6df1SKris Buschelman }
766397b6df1SKris Buschelman 
7679a2535b5SHong Zhang /* Sets MUMPS options from the options database */
768dcd589f8SShri Abhyankar #undef __FUNCT__
7699a2535b5SHong Zhang #define __FUNCT__ "PetscSetMUMPSFromOptions"
7709a2535b5SHong Zhang PetscErrorCode PetscSetMUMPSFromOptions(Mat F, Mat A)
771dcd589f8SShri Abhyankar {
7729a2535b5SHong Zhang   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->spptr;
773dcd589f8SShri Abhyankar   PetscErrorCode ierr;
774dcd589f8SShri Abhyankar   PetscInt       icntl;
775ace3abfcSBarry Smith   PetscBool      flg;
776dcd589f8SShri Abhyankar 
777dcd589f8SShri Abhyankar   PetscFunctionBegin;
778dcd589f8SShri Abhyankar   ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"MUMPS Options","Mat");CHKERRQ(ierr);
7799a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_1","ICNTL(1): output stream for error messages","None",mumps->id.ICNTL(1),&icntl,&flg);CHKERRQ(ierr);
7809a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(1) = icntl;
7819a2535b5SHong 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);
7829a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(2) = icntl;
7839a2535b5SHong 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);
7849a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(3) = icntl;
785dcd589f8SShri Abhyankar 
7869a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",mumps->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr);
7879a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(4) = icntl;
7889a2535b5SHong Zhang   if (mumps->id.ICNTL(4) || PetscLogPrintInfo) mumps->id.ICNTL(3) = 6; /* resume MUMPS default id.ICNTL(3) = 6 */
7899a2535b5SHong Zhang 
7909a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_6","ICNTL(6): permuting and/or scaling the matrix (0 to 7)","None",mumps->id.ICNTL(6),&icntl,&flg);CHKERRQ(ierr);
7919a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(6) = icntl;
7929a2535b5SHong Zhang 
7939a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): matrix ordering (0 to 7). 3=Scotch, 4=PORD, 5=Metis","None",mumps->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr);
794dcd589f8SShri Abhyankar   if (flg) {
795*2205254eSKarl 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");
796*2205254eSKarl Rupp     else mumps->id.ICNTL(7) = icntl;
797dcd589f8SShri Abhyankar   }
798e0b74bf9SHong Zhang 
79970544d5fSHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_8","ICNTL(8): scaling strategy (-2 to 8 or 77)","None",mumps->id.ICNTL(8),&mumps->id.ICNTL(8),PETSC_NULL);CHKERRQ(ierr);
8009a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_10","ICNTL(10): max num of refinements","None",mumps->id.ICNTL(10),&mumps->id.ICNTL(10),PETSC_NULL);CHKERRQ(ierr);
8019a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_11","ICNTL(11): statistics related to the linear system solved (via -ksp_view)","None",mumps->id.ICNTL(11),&mumps->id.ICNTL(11),PETSC_NULL);CHKERRQ(ierr);
80270544d5fSHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): efficiency control: defines the ordering strategy with scaling constraints (0 to 3)","None",mumps->id.ICNTL(12),&mumps->id.ICNTL(12),PETSC_NULL);CHKERRQ(ierr);
8039a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): efficiency control: with or without ScaLAPACK","None",mumps->id.ICNTL(13),&mumps->id.ICNTL(13),PETSC_NULL);CHKERRQ(ierr);
8049a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_14","ICNTL(14): percentage of estimated workspace increase","None",mumps->id.ICNTL(14),&mumps->id.ICNTL(14),PETSC_NULL);CHKERRQ(ierr);
8059a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_19","ICNTL(19): Schur complement","None",mumps->id.ICNTL(19),&mumps->id.ICNTL(19),PETSC_NULL);CHKERRQ(ierr);
8069a2535b5SHong Zhang 
8079a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_22","ICNTL(22): in-core/out-of-core facility (0 or 1)","None",mumps->id.ICNTL(22),&mumps->id.ICNTL(22),PETSC_NULL);CHKERRQ(ierr);
8089a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_23","ICNTL(23): max size of the working memory (MB) that can allocate per processor","None",mumps->id.ICNTL(23),&mumps->id.ICNTL(23),PETSC_NULL);CHKERRQ(ierr);
8099a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_24","ICNTL(24): detection of null pivot rows (0 or 1)","None",mumps->id.ICNTL(24),&mumps->id.ICNTL(24),PETSC_NULL);CHKERRQ(ierr);
8109a2535b5SHong Zhang   if (mumps->id.ICNTL(24)) {
8119a2535b5SHong Zhang     mumps->id.ICNTL(13) = 1; /* turn-off ScaLAPACK to help with the correct detection of null pivots */
812d7ebd59bSHong Zhang   }
813d7ebd59bSHong Zhang 
8149a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_25","ICNTL(25): computation of a null space basis","None",mumps->id.ICNTL(25),&mumps->id.ICNTL(25),PETSC_NULL);CHKERRQ(ierr);
8159a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_26","ICNTL(26): Schur options for right-hand side or solution vector","None",mumps->id.ICNTL(26),&mumps->id.ICNTL(26),PETSC_NULL);CHKERRQ(ierr);
8169a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_27","ICNTL(27): experimental parameter","None",mumps->id.ICNTL(27),&mumps->id.ICNTL(27),PETSC_NULL);CHKERRQ(ierr);
8179a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_28","ICNTL(28): use 1 for sequential analysis and ictnl(7) ordering, or 2 for parallel analysis and ictnl(29) ordering","None",mumps->id.ICNTL(28),&mumps->id.ICNTL(28),PETSC_NULL);CHKERRQ(ierr);
8189a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_29","ICNTL(29): parallel ordering 1 = ptscotch 2 = parmetis","None",mumps->id.ICNTL(29),&mumps->id.ICNTL(29),PETSC_NULL);CHKERRQ(ierr);
8199a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_30","ICNTL(30): compute user-specified set of entries in inv(A)","None",mumps->id.ICNTL(30),&mumps->id.ICNTL(30),PETSC_NULL);CHKERRQ(ierr);
8209a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_31","ICNTL(31): factors can be discarded in the solve phase","None",mumps->id.ICNTL(31),&mumps->id.ICNTL(31),PETSC_NULL);CHKERRQ(ierr);
8219a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_33","ICNTL(33): compute determinant","None",mumps->id.ICNTL(33),&mumps->id.ICNTL(33),PETSC_NULL);CHKERRQ(ierr);
822dcd589f8SShri Abhyankar 
8239a2535b5SHong Zhang   ierr = PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",mumps->id.CNTL(1),&mumps->id.CNTL(1),PETSC_NULL);CHKERRQ(ierr);
8249a2535b5SHong Zhang   ierr = PetscOptionsReal("-mat_mumps_cntl_2","CNTL(2): stopping criterion of refinement","None",mumps->id.CNTL(2),&mumps->id.CNTL(2),PETSC_NULL);CHKERRQ(ierr);
8259a2535b5SHong Zhang   ierr = PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",mumps->id.CNTL(3),&mumps->id.CNTL(3),PETSC_NULL);CHKERRQ(ierr);
8269a2535b5SHong Zhang   ierr = PetscOptionsReal("-mat_mumps_cntl_4","CNTL(4): value for static pivoting","None",mumps->id.CNTL(4),&mumps->id.CNTL(4),PETSC_NULL);CHKERRQ(ierr);
8279a2535b5SHong Zhang   ierr = PetscOptionsReal("-mat_mumps_cntl_5","CNTL(5): fixation for null pivots","None",mumps->id.CNTL(5),&mumps->id.CNTL(5),PETSC_NULL);CHKERRQ(ierr);
828e5bb22a1SHong Zhang 
829e5bb22a1SHong Zhang   ierr = PetscOptionsString("-mat_mumps_ooc_tmpdir", "out of core directory", "None", mumps->id.ooc_tmpdir, mumps->id.ooc_tmpdir, 256, PETSC_NULL);
830dcd589f8SShri Abhyankar   PetscOptionsEnd();
831dcd589f8SShri Abhyankar   PetscFunctionReturn(0);
832dcd589f8SShri Abhyankar }
833dcd589f8SShri Abhyankar 
834dcd589f8SShri Abhyankar #undef __FUNCT__
835dcd589f8SShri Abhyankar #define __FUNCT__ "PetscInitializeMUMPS"
836f697e70eSHong Zhang PetscErrorCode PetscInitializeMUMPS(Mat A,Mat_MUMPS *mumps)
837dcd589f8SShri Abhyankar {
838dcd589f8SShri Abhyankar   PetscErrorCode ierr;
839dcd589f8SShri Abhyankar 
840dcd589f8SShri Abhyankar   PetscFunctionBegin;
841f697e70eSHong Zhang   ierr = MPI_Comm_rank(((PetscObject)A)->comm, &mumps->myid);
842f697e70eSHong Zhang   ierr = MPI_Comm_size(((PetscObject)A)->comm,&mumps->size);CHKERRQ(ierr);
843f697e70eSHong Zhang   ierr = MPI_Comm_dup(((PetscObject)A)->comm,&(mumps->comm_mumps));CHKERRQ(ierr);
844*2205254eSKarl Rupp 
845f697e70eSHong Zhang   mumps->id.comm_fortran = MPI_Comm_c2f(mumps->comm_mumps);
846f697e70eSHong Zhang 
847f697e70eSHong Zhang   mumps->id.job = JOB_INIT;
848f697e70eSHong Zhang   mumps->id.par = 1;  /* host participates factorizaton and solve */
849f697e70eSHong Zhang   mumps->id.sym = mumps->sym;
8502907cef9SHong Zhang   PetscMUMPS_c(&mumps->id);
851f697e70eSHong Zhang 
852f697e70eSHong Zhang   mumps->CleanUpMUMPS = PETSC_FALSE;
853f697e70eSHong Zhang   mumps->scat_rhs     = PETSC_NULL;
854f697e70eSHong Zhang   mumps->scat_sol     = PETSC_NULL;
8559a2535b5SHong Zhang 
85670544d5fSHong Zhang   /* set PETSc-MUMPS default options - override MUMPS default */
8579a2535b5SHong Zhang   mumps->id.ICNTL(3) = 0;
8589a2535b5SHong Zhang   mumps->id.ICNTL(4) = 0;
8599a2535b5SHong Zhang   if (mumps->size == 1) {
8609a2535b5SHong Zhang     mumps->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
8619a2535b5SHong Zhang   } else {
8629a2535b5SHong Zhang     mumps->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
86370544d5fSHong Zhang     mumps->id.ICNTL(21) = 1;   /* distributed solution */
8649a2535b5SHong Zhang   }
865dcd589f8SShri Abhyankar   PetscFunctionReturn(0);
866dcd589f8SShri Abhyankar }
867dcd589f8SShri Abhyankar 
868a5e57a09SHong 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 */
869397b6df1SKris Buschelman #undef __FUNCT__
870f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS"
8710481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
872b24902e0SBarry Smith {
873a5e57a09SHong Zhang   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->spptr;
874dcd589f8SShri Abhyankar   PetscErrorCode ierr;
87567877ebaSShri Abhyankar   Vec            b;
87667877ebaSShri Abhyankar   IS             is_iden;
87767877ebaSShri Abhyankar   const PetscInt M = A->rmap->N;
878397b6df1SKris Buschelman 
879397b6df1SKris Buschelman   PetscFunctionBegin;
880a5e57a09SHong Zhang   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
881dcd589f8SShri Abhyankar 
8829a2535b5SHong Zhang   /* Set MUMPS options from the options database */
8839a2535b5SHong Zhang   ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr);
884dcd589f8SShri Abhyankar 
885a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
886dcd589f8SShri Abhyankar 
88767877ebaSShri Abhyankar   /* analysis phase */
88867877ebaSShri Abhyankar   /*----------------*/
889a5e57a09SHong Zhang   mumps->id.job = JOB_FACTSYMBOLIC;
890a5e57a09SHong Zhang   mumps->id.n   = M;
891a5e57a09SHong Zhang   switch (mumps->id.ICNTL(18)) {
89267877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
893a5e57a09SHong Zhang     if (!mumps->myid) {
894a5e57a09SHong Zhang       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
895a5e57a09SHong Zhang       if (mumps->id.ICNTL(6)>1) {
89667877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
8972907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
898a5e57a09SHong Zhang         mumps->id.a = (mumps_complex*)mumps->val;
8992907cef9SHong Zhang #else
900a5e57a09SHong Zhang         mumps->id.a = (mumps_double_complex*)mumps->val;
9012907cef9SHong Zhang #endif
90267877ebaSShri Abhyankar #else
903a5e57a09SHong Zhang         mumps->id.a = mumps->val;
90467877ebaSShri Abhyankar #endif
90567877ebaSShri Abhyankar       }
906a5e57a09SHong Zhang       if (mumps->id.ICNTL(7) == 1) { /* use user-provide matrix ordering - assuming r = c ordering */
9075248a706SHong Zhang         /*
9085248a706SHong Zhang         PetscBool      flag;
9095248a706SHong Zhang         ierr = ISEqual(r,c,&flag);CHKERRQ(ierr);
9105248a706SHong Zhang         if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"row_perm != col_perm");
9115248a706SHong Zhang         ierr = ISView(r,PETSC_VIEWER_STDOUT_SELF);
9125248a706SHong Zhang          */
913a5e57a09SHong Zhang         if (!mumps->myid) {
914e0b74bf9SHong Zhang           const PetscInt *idx;
915e0b74bf9SHong Zhang           PetscInt i,*perm_in;
916*2205254eSKarl Rupp 
917e0b74bf9SHong Zhang           ierr = PetscMalloc(M*sizeof(PetscInt),&perm_in);CHKERRQ(ierr);
918e0b74bf9SHong Zhang           ierr = ISGetIndices(r,&idx);CHKERRQ(ierr);
919*2205254eSKarl Rupp 
920a5e57a09SHong Zhang           mumps->id.perm_in = perm_in;
921e0b74bf9SHong Zhang           for (i=0; i<M; i++) perm_in[i] = idx[i]+1; /* perm_in[]: start from 1, not 0! */
922e0b74bf9SHong Zhang           ierr = ISRestoreIndices(r,&idx);CHKERRQ(ierr);
923e0b74bf9SHong Zhang         }
924e0b74bf9SHong Zhang       }
92567877ebaSShri Abhyankar     }
92667877ebaSShri Abhyankar     break;
92767877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
928a5e57a09SHong Zhang     mumps->id.nz_loc = mumps->nz;
929a5e57a09SHong Zhang     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
930a5e57a09SHong Zhang     if (mumps->id.ICNTL(6)>1) {
93167877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
9322907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
933a5e57a09SHong Zhang       mumps->id.a_loc = (mumps_complex*)mumps->val;
9342907cef9SHong Zhang #else
935a5e57a09SHong Zhang       mumps->id.a_loc = (mumps_double_complex*)mumps->val;
9362907cef9SHong Zhang #endif
93767877ebaSShri Abhyankar #else
938a5e57a09SHong Zhang       mumps->id.a_loc = mumps->val;
93967877ebaSShri Abhyankar #endif
94067877ebaSShri Abhyankar     }
94167877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
942a5e57a09SHong Zhang     if (!mumps->myid) {
943a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr);
94467877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
94567877ebaSShri Abhyankar     } else {
946a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr);
94767877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
94867877ebaSShri Abhyankar     }
94967877ebaSShri Abhyankar     ierr = VecCreate(((PetscObject)A)->comm,&b);CHKERRQ(ierr);
95067877ebaSShri Abhyankar     ierr = VecSetSizes(b,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
95167877ebaSShri Abhyankar     ierr = VecSetFromOptions(b);CHKERRQ(ierr);
95267877ebaSShri Abhyankar 
953a5e57a09SHong Zhang     ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr);
9546bf464f9SBarry Smith     ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
9556bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
95667877ebaSShri Abhyankar     break;
95767877ebaSShri Abhyankar   }
958a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
959a5e57a09SHong 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));
96067877ebaSShri Abhyankar 
961719d5645SBarry Smith   F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
962dcd589f8SShri Abhyankar   F->ops->solve           = MatSolve_MUMPS;
96351d5961aSHong Zhang   F->ops->solvetranspose  = MatSolveTranspose_MUMPS;
96417f96c7aSHong Zhang   F->ops->matsolve        = 0;  /* use MatMatSolve_Basic() until mumps supports distributed rhs */
965b24902e0SBarry Smith   PetscFunctionReturn(0);
966b24902e0SBarry Smith }
967b24902e0SBarry Smith 
968450b117fSShri Abhyankar /* Note the Petsc r and c permutations are ignored */
969450b117fSShri Abhyankar #undef __FUNCT__
970450b117fSShri Abhyankar #define __FUNCT__ "MatLUFactorSymbolic_BAIJMUMPS"
971450b117fSShri Abhyankar PetscErrorCode MatLUFactorSymbolic_BAIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
972450b117fSShri Abhyankar {
973a5e57a09SHong Zhang   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->spptr;
974dcd589f8SShri Abhyankar   PetscErrorCode ierr;
97567877ebaSShri Abhyankar   Vec            b;
97667877ebaSShri Abhyankar   IS             is_iden;
97767877ebaSShri Abhyankar   const PetscInt M = A->rmap->N;
978450b117fSShri Abhyankar 
979450b117fSShri Abhyankar   PetscFunctionBegin;
980a5e57a09SHong Zhang   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
981dcd589f8SShri Abhyankar 
9829a2535b5SHong Zhang   /* Set MUMPS options from the options database */
9839a2535b5SHong Zhang   ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr);
984dcd589f8SShri Abhyankar 
985a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
98667877ebaSShri Abhyankar 
98767877ebaSShri Abhyankar   /* analysis phase */
98867877ebaSShri Abhyankar   /*----------------*/
989a5e57a09SHong Zhang   mumps->id.job = JOB_FACTSYMBOLIC;
990a5e57a09SHong Zhang   mumps->id.n   = M;
991a5e57a09SHong Zhang   switch (mumps->id.ICNTL(18)) {
99267877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
993a5e57a09SHong Zhang     if (!mumps->myid) {
994a5e57a09SHong Zhang       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
995a5e57a09SHong Zhang       if (mumps->id.ICNTL(6)>1) {
99667877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
9972907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
998a5e57a09SHong Zhang         mumps->id.a = (mumps_complex*)mumps->val;
9992907cef9SHong Zhang #else
1000a5e57a09SHong Zhang         mumps->id.a = (mumps_double_complex*)mumps->val;
10012907cef9SHong Zhang #endif
100267877ebaSShri Abhyankar #else
1003a5e57a09SHong Zhang         mumps->id.a = mumps->val;
100467877ebaSShri Abhyankar #endif
100567877ebaSShri Abhyankar       }
100667877ebaSShri Abhyankar     }
100767877ebaSShri Abhyankar     break;
100867877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
1009a5e57a09SHong Zhang     mumps->id.nz_loc = mumps->nz;
1010a5e57a09SHong Zhang     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1011a5e57a09SHong Zhang     if (mumps->id.ICNTL(6)>1) {
101267877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
10132907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
1014a5e57a09SHong Zhang       mumps->id.a_loc = (mumps_complex*)mumps->val;
10152907cef9SHong Zhang #else
1016a5e57a09SHong Zhang       mumps->id.a_loc = (mumps_double_complex*)mumps->val;
10172907cef9SHong Zhang #endif
101867877ebaSShri Abhyankar #else
1019a5e57a09SHong Zhang       mumps->id.a_loc = mumps->val;
102067877ebaSShri Abhyankar #endif
102167877ebaSShri Abhyankar     }
102267877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1023a5e57a09SHong Zhang     if (!mumps->myid) {
1024a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr);
102567877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
102667877ebaSShri Abhyankar     } else {
1027a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr);
102867877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
102967877ebaSShri Abhyankar     }
103067877ebaSShri Abhyankar     ierr = VecCreate(((PetscObject)A)->comm,&b);CHKERRQ(ierr);
103167877ebaSShri Abhyankar     ierr = VecSetSizes(b,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
103267877ebaSShri Abhyankar     ierr = VecSetFromOptions(b);CHKERRQ(ierr);
103367877ebaSShri Abhyankar 
1034a5e57a09SHong Zhang     ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr);
10356bf464f9SBarry Smith     ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
10366bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
103767877ebaSShri Abhyankar     break;
103867877ebaSShri Abhyankar   }
1039a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
1040a5e57a09SHong 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));
104167877ebaSShri Abhyankar 
1042450b117fSShri Abhyankar   F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
1043dcd589f8SShri Abhyankar   F->ops->solve           = MatSolve_MUMPS;
104451d5961aSHong Zhang   F->ops->solvetranspose  = MatSolveTranspose_MUMPS;
1045450b117fSShri Abhyankar   PetscFunctionReturn(0);
1046450b117fSShri Abhyankar }
1047b24902e0SBarry Smith 
1048141f4205SHong Zhang /* Note the Petsc r permutation and factor info are ignored */
1049397b6df1SKris Buschelman #undef __FUNCT__
105067877ebaSShri Abhyankar #define __FUNCT__ "MatCholeskyFactorSymbolic_MUMPS"
105167877ebaSShri Abhyankar PetscErrorCode MatCholeskyFactorSymbolic_MUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info)
1052b24902e0SBarry Smith {
1053a5e57a09SHong Zhang   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->spptr;
1054dcd589f8SShri Abhyankar   PetscErrorCode ierr;
105567877ebaSShri Abhyankar   Vec            b;
105667877ebaSShri Abhyankar   IS             is_iden;
105767877ebaSShri Abhyankar   const PetscInt M = A->rmap->N;
1058397b6df1SKris Buschelman 
1059397b6df1SKris Buschelman   PetscFunctionBegin;
1060a5e57a09SHong Zhang   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
1061dcd589f8SShri Abhyankar 
10629a2535b5SHong Zhang   /* Set MUMPS options from the options database */
10639a2535b5SHong Zhang   ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr);
1064dcd589f8SShri Abhyankar 
1065a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
1066dcd589f8SShri Abhyankar 
106767877ebaSShri Abhyankar   /* analysis phase */
106867877ebaSShri Abhyankar   /*----------------*/
1069a5e57a09SHong Zhang   mumps->id.job = JOB_FACTSYMBOLIC;
1070a5e57a09SHong Zhang   mumps->id.n   = M;
1071a5e57a09SHong Zhang   switch (mumps->id.ICNTL(18)) {
107267877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
1073a5e57a09SHong Zhang     if (!mumps->myid) {
1074a5e57a09SHong Zhang       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1075a5e57a09SHong Zhang       if (mumps->id.ICNTL(6)>1) {
107667877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
10772907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
1078a5e57a09SHong Zhang         mumps->id.a = (mumps_complex*)mumps->val;
10792907cef9SHong Zhang #else
1080a5e57a09SHong Zhang         mumps->id.a = (mumps_double_complex*)mumps->val;
10812907cef9SHong Zhang #endif
108267877ebaSShri Abhyankar #else
1083a5e57a09SHong Zhang         mumps->id.a = mumps->val;
108467877ebaSShri Abhyankar #endif
108567877ebaSShri Abhyankar       }
108667877ebaSShri Abhyankar     }
108767877ebaSShri Abhyankar     break;
108867877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
1089a5e57a09SHong Zhang     mumps->id.nz_loc = mumps->nz;
1090a5e57a09SHong Zhang     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1091a5e57a09SHong Zhang     if (mumps->id.ICNTL(6)>1) {
109267877ebaSShri Abhyankar #if defined(PETSC_USE_COMPLEX)
10932907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
1094a5e57a09SHong Zhang       mumps->id.a_loc = (mumps_complex*)mumps->val;
10952907cef9SHong Zhang #else
1096a5e57a09SHong Zhang       mumps->id.a_loc = (mumps_double_complex*)mumps->val;
10972907cef9SHong Zhang #endif
109867877ebaSShri Abhyankar #else
1099a5e57a09SHong Zhang       mumps->id.a_loc = mumps->val;
110067877ebaSShri Abhyankar #endif
110167877ebaSShri Abhyankar     }
110267877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1103a5e57a09SHong Zhang     if (!mumps->myid) {
1104a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr);
110567877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
110667877ebaSShri Abhyankar     } else {
1107a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr);
110867877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
110967877ebaSShri Abhyankar     }
111067877ebaSShri Abhyankar     ierr = VecCreate(((PetscObject)A)->comm,&b);CHKERRQ(ierr);
111167877ebaSShri Abhyankar     ierr = VecSetSizes(b,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
111267877ebaSShri Abhyankar     ierr = VecSetFromOptions(b);CHKERRQ(ierr);
111367877ebaSShri Abhyankar 
1114a5e57a09SHong Zhang     ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr);
11156bf464f9SBarry Smith     ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
11166bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
111767877ebaSShri Abhyankar     break;
111867877ebaSShri Abhyankar   }
1119a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
1120a5e57a09SHong 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));
112167877ebaSShri Abhyankar 
11222792810eSHong Zhang   F->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS;
1123dcd589f8SShri Abhyankar   F->ops->solve                 = MatSolve_MUMPS;
112451d5961aSHong Zhang   F->ops->solvetranspose        = MatSolve_MUMPS;
112530c107b7SHong Zhang   F->ops->matsolve              = 0; /* use MatMatSolve_Basic() until mumps supports distributed rhs */
1126db4efbfdSBarry Smith #if !defined(PETSC_USE_COMPLEX)
112705aa0992SJose Roman   F->ops->getinertia = MatGetInertia_SBAIJMUMPS;
112805aa0992SJose Roman #else
112905aa0992SJose Roman   F->ops->getinertia = PETSC_NULL;
1130db4efbfdSBarry Smith #endif
1131b24902e0SBarry Smith   PetscFunctionReturn(0);
1132b24902e0SBarry Smith }
1133b24902e0SBarry Smith 
1134397b6df1SKris Buschelman #undef __FUNCT__
113564e6c443SBarry Smith #define __FUNCT__ "MatView_MUMPS"
113664e6c443SBarry Smith PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer)
113774ed9c26SBarry Smith {
1138f6c57405SHong Zhang   PetscErrorCode    ierr;
113964e6c443SBarry Smith   PetscBool         iascii;
114064e6c443SBarry Smith   PetscViewerFormat format;
1141a5e57a09SHong Zhang   Mat_MUMPS         *mumps=(Mat_MUMPS*)A->spptr;
1142f6c57405SHong Zhang 
1143f6c57405SHong Zhang   PetscFunctionBegin;
114464e6c443SBarry Smith   /* check if matrix is mumps type */
114564e6c443SBarry Smith   if (A->ops->solve != MatSolve_MUMPS) PetscFunctionReturn(0);
114664e6c443SBarry Smith 
1147251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
114864e6c443SBarry Smith   if (iascii) {
114964e6c443SBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
115064e6c443SBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO) {
115164e6c443SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr);
1152a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                   %d \n",mumps->id.sym);CHKERRQ(ierr);
1153a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  PAR (host participation):            %d \n",mumps->id.par);CHKERRQ(ierr);
1154a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(1) (output for error):         %d \n",mumps->id.ICNTL(1));CHKERRQ(ierr);
1155a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(2) (output of diagnostic msg): %d \n",mumps->id.ICNTL(2));CHKERRQ(ierr);
1156a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(3) (output for global info):   %d \n",mumps->id.ICNTL(3));CHKERRQ(ierr);
1157a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):        %d \n",mumps->id.ICNTL(4));CHKERRQ(ierr);
1158a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):         %d \n",mumps->id.ICNTL(5));CHKERRQ(ierr);
1159a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):        %d \n",mumps->id.ICNTL(6));CHKERRQ(ierr);
1160a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (sequentia matrix ordering):%d \n",mumps->id.ICNTL(7));CHKERRQ(ierr);
1161a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(8) (scalling strategy):        %d \n",mumps->id.ICNTL(8));CHKERRQ(ierr);
1162a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements):  %d \n",mumps->id.ICNTL(10));CHKERRQ(ierr);
1163a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):          %d \n",mumps->id.ICNTL(11));CHKERRQ(ierr);
1164a5e57a09SHong Zhang       if (mumps->id.ICNTL(11)>0) {
1165a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(4) (inf norm of input mat):        %g\n",mumps->id.RINFOG(4));CHKERRQ(ierr);
1166a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(5) (inf norm of solution):         %g\n",mumps->id.RINFOG(5));CHKERRQ(ierr);
1167a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(6) (inf norm of residual):         %g\n",mumps->id.RINFOG(6));CHKERRQ(ierr);
1168a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",mumps->id.RINFOG(7),mumps->id.RINFOG(8));CHKERRQ(ierr);
1169a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(9) (error estimate):               %g \n",mumps->id.RINFOG(9));CHKERRQ(ierr);
1170a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",mumps->id.RINFOG(10),mumps->id.RINFOG(11));CHKERRQ(ierr);
1171f6c57405SHong Zhang       }
1172a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):                         %d \n",mumps->id.ICNTL(12));CHKERRQ(ierr);
1173a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):                         %d \n",mumps->id.ICNTL(13));CHKERRQ(ierr);
1174a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (percentage of estimated workspace increase): %d \n",mumps->id.ICNTL(14));CHKERRQ(ierr);
1175f6c57405SHong Zhang       /* ICNTL(15-17) not used */
1176a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):                           %d \n",mumps->id.ICNTL(18));CHKERRQ(ierr);
1177a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(19) (Shur complement info):                       %d \n",mumps->id.ICNTL(19));CHKERRQ(ierr);
1178a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(20) (rhs sparse pattern):                         %d \n",mumps->id.ICNTL(20));CHKERRQ(ierr);
1179a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(21) (somumpstion struct):                            %d \n",mumps->id.ICNTL(21));CHKERRQ(ierr);
1180a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(22) (in-core/out-of-core facility):               %d \n",mumps->id.ICNTL(22));CHKERRQ(ierr);
1181a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(23) (max size of memory can be allocated locally):%d \n",mumps->id.ICNTL(23));CHKERRQ(ierr);
1182c0165424SHong Zhang 
1183a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(24) (detection of null pivot rows):               %d \n",mumps->id.ICNTL(24));CHKERRQ(ierr);
1184a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(25) (computation of a null space basis):          %d \n",mumps->id.ICNTL(25));CHKERRQ(ierr);
1185a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(26) (Schur options for rhs or solution):          %d \n",mumps->id.ICNTL(26));CHKERRQ(ierr);
1186a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(27) (experimental parameter):                     %d \n",mumps->id.ICNTL(27));CHKERRQ(ierr);
1187a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(28) (use parallel or sequential ordering):        %d \n",mumps->id.ICNTL(28));CHKERRQ(ierr);
1188a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(29) (parallel ordering):                          %d \n",mumps->id.ICNTL(29));CHKERRQ(ierr);
118942179a6aSHong Zhang 
1190a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(30) (user-specified set of entries in inv(A)):    %d \n",mumps->id.ICNTL(30));CHKERRQ(ierr);
1191a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(31) (factors is discarded in the solve phase):    %d \n",mumps->id.ICNTL(31));CHKERRQ(ierr);
1192a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(33) (compute determinant):                        %d \n",mumps->id.ICNTL(33));CHKERRQ(ierr);
1193f6c57405SHong Zhang 
1194a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",mumps->id.CNTL(1));CHKERRQ(ierr);
1195a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",mumps->id.CNTL(2));CHKERRQ(ierr);
1196a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absomumpste pivoting threshold):      %g \n",mumps->id.CNTL(3));CHKERRQ(ierr);
1197a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(4) (vamumpse of static pivoting):         %g \n",mumps->id.CNTL(4));CHKERRQ(ierr);
1198a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(5) (fixation for null pivots):         %g \n",mumps->id.CNTL(5));CHKERRQ(ierr);
1199f6c57405SHong Zhang 
1200f6c57405SHong Zhang       /* infomation local to each processor */
120134ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr);
12027b23a99aSBarry Smith       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
1203a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %g \n",mumps->myid,mumps->id.RINFO(1));CHKERRQ(ierr);
120434ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
120534ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr);
1206a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",mumps->myid,mumps->id.RINFO(2));CHKERRQ(ierr);
120734ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
120834ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr);
1209a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",mumps->myid,mumps->id.RINFO(3));CHKERRQ(ierr);
121034ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
1211f6c57405SHong Zhang 
121234ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization): \n");CHKERRQ(ierr);
1213a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] %d \n",mumps->myid,mumps->id.INFO(15));CHKERRQ(ierr);
121434ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
1215f6c57405SHong Zhang 
121634ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization): \n");CHKERRQ(ierr);
1217a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",mumps->myid,mumps->id.INFO(16));CHKERRQ(ierr);
121834ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
1219f6c57405SHong Zhang 
122034ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(23) (num of pivots eliminated on this processor after factorization): \n");CHKERRQ(ierr);
1221a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",mumps->myid,mumps->id.INFO(23));CHKERRQ(ierr);
122234ed7027SBarry Smith       ierr = PetscViewerFlush(viewer);
12237b23a99aSBarry Smith       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
1224f6c57405SHong Zhang 
1225a5e57a09SHong Zhang       if (!mumps->myid) { /* information from the host */
1226a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",mumps->id.RINFOG(1));CHKERRQ(ierr);
1227a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",mumps->id.RINFOG(2));CHKERRQ(ierr);
1228a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",mumps->id.RINFOG(3));CHKERRQ(ierr);
1229a5e57a09SHong 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);
1230f6c57405SHong Zhang 
1231a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(3));CHKERRQ(ierr);
1232a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(4));CHKERRQ(ierr);
1233a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(5) (estimated maximum front size in the complete tree): %d \n",mumps->id.INFOG(5));CHKERRQ(ierr);
1234a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(6) (number of nodes in the complete tree): %d \n",mumps->id.INFOG(6));CHKERRQ(ierr);
1235a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(7) (ordering option effectively use after analysis): %d \n",mumps->id.INFOG(7));CHKERRQ(ierr);
1236a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",mumps->id.INFOG(8));CHKERRQ(ierr);
1237a5e57a09SHong 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);
1238a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(10) (total integer space store the matrix factors after factorization): %d \n",mumps->id.INFOG(10));CHKERRQ(ierr);
1239a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(11) (order of largest frontal matrix after factorization): %d \n",mumps->id.INFOG(11));CHKERRQ(ierr);
1240a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(12) (number of off-diagonal pivots): %d \n",mumps->id.INFOG(12));CHKERRQ(ierr);
1241a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(13) (number of delayed pivots after factorization): %d \n",mumps->id.INFOG(13));CHKERRQ(ierr);
1242a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(14) (number of memory compress after factorization): %d \n",mumps->id.INFOG(14));CHKERRQ(ierr);
1243a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(15) (number of steps of iterative refinement after solution): %d \n",mumps->id.INFOG(15));CHKERRQ(ierr);
1244a5e57a09SHong 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);
1245a5e57a09SHong 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);
1246a5e57a09SHong 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);
1247a5e57a09SHong 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);
1248a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(20) (estimated number of entries in the factors): %d \n",mumps->id.INFOG(20));CHKERRQ(ierr);
1249a5e57a09SHong 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);
1250a5e57a09SHong 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);
1251a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",mumps->id.INFOG(23));CHKERRQ(ierr);
1252a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",mumps->id.INFOG(24));CHKERRQ(ierr);
1253a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",mumps->id.INFOG(25));CHKERRQ(ierr);
1254f6c57405SHong Zhang       }
1255f6c57405SHong Zhang     }
1256cb828f0fSHong Zhang   }
1257f6c57405SHong Zhang   PetscFunctionReturn(0);
1258f6c57405SHong Zhang }
1259f6c57405SHong Zhang 
126035bd34faSBarry Smith #undef __FUNCT__
126135bd34faSBarry Smith #define __FUNCT__ "MatGetInfo_MUMPS"
126235bd34faSBarry Smith PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info)
126335bd34faSBarry Smith {
1264cb828f0fSHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)A->spptr;
126535bd34faSBarry Smith 
126635bd34faSBarry Smith   PetscFunctionBegin;
126735bd34faSBarry Smith   info->block_size        = 1.0;
1268cb828f0fSHong Zhang   info->nz_allocated      = mumps->id.INFOG(20);
1269cb828f0fSHong Zhang   info->nz_used           = mumps->id.INFOG(20);
127035bd34faSBarry Smith   info->nz_unneeded       = 0.0;
127135bd34faSBarry Smith   info->assemblies        = 0.0;
127235bd34faSBarry Smith   info->mallocs           = 0.0;
127335bd34faSBarry Smith   info->memory            = 0.0;
127435bd34faSBarry Smith   info->fill_ratio_given  = 0;
127535bd34faSBarry Smith   info->fill_ratio_needed = 0;
127635bd34faSBarry Smith   info->factor_mallocs    = 0;
127735bd34faSBarry Smith   PetscFunctionReturn(0);
127835bd34faSBarry Smith }
127935bd34faSBarry Smith 
12805ccb76cbSHong Zhang /* -------------------------------------------------------------------------------------------*/
12815ccb76cbSHong Zhang #undef __FUNCT__
12825ccb76cbSHong Zhang #define __FUNCT__ "MatMumpsSetIcntl_MUMPS"
12835ccb76cbSHong Zhang PetscErrorCode MatMumpsSetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt ival)
12845ccb76cbSHong Zhang {
1285a5e57a09SHong Zhang   Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr;
12865ccb76cbSHong Zhang 
12875ccb76cbSHong Zhang   PetscFunctionBegin;
1288a5e57a09SHong Zhang   mumps->id.ICNTL(icntl) = ival;
12895ccb76cbSHong Zhang   PetscFunctionReturn(0);
12905ccb76cbSHong Zhang }
12915ccb76cbSHong Zhang 
12925ccb76cbSHong Zhang #undef __FUNCT__
12935ccb76cbSHong Zhang #define __FUNCT__ "MatMumpsSetIcntl"
12945ccb76cbSHong Zhang /*@
12955ccb76cbSHong Zhang   MatMumpsSetIcntl - Set MUMPS parameter ICNTL()
12965ccb76cbSHong Zhang 
12975ccb76cbSHong Zhang    Logically Collective on Mat
12985ccb76cbSHong Zhang 
12995ccb76cbSHong Zhang    Input Parameters:
13005ccb76cbSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
13015ccb76cbSHong Zhang .  icntl - index of MUMPS parameter array ICNTL()
13025ccb76cbSHong Zhang -  ival - value of MUMPS ICNTL(icntl)
13035ccb76cbSHong Zhang 
13045ccb76cbSHong Zhang   Options Database:
13055ccb76cbSHong Zhang .   -mat_mumps_icntl_<icntl> <ival>
13065ccb76cbSHong Zhang 
13075ccb76cbSHong Zhang    Level: beginner
13085ccb76cbSHong Zhang 
13095ccb76cbSHong Zhang    References: MUMPS Users' Guide
13105ccb76cbSHong Zhang 
13115ccb76cbSHong Zhang .seealso: MatGetFactor()
13125ccb76cbSHong Zhang @*/
13135ccb76cbSHong Zhang PetscErrorCode MatMumpsSetIcntl(Mat F,PetscInt icntl,PetscInt ival)
13145ccb76cbSHong Zhang {
13155ccb76cbSHong Zhang   PetscErrorCode ierr;
13165ccb76cbSHong Zhang 
13175ccb76cbSHong Zhang   PetscFunctionBegin;
13185ccb76cbSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
13195ccb76cbSHong Zhang   PetscValidLogicalCollectiveInt(F,ival,3);
13205ccb76cbSHong Zhang   ierr = PetscTryMethod(F,"MatMumpsSetIcntl_C",(Mat,PetscInt,PetscInt),(F,icntl,ival));CHKERRQ(ierr);
13215ccb76cbSHong Zhang   PetscFunctionReturn(0);
13225ccb76cbSHong Zhang }
13235ccb76cbSHong Zhang 
132424b6179bSKris Buschelman /*MC
13252692d6eeSBarry Smith   MATSOLVERMUMPS -  A matrix type providing direct solvers (LU and Cholesky) for
132624b6179bSKris Buschelman   distributed and sequential matrices via the external package MUMPS.
132724b6179bSKris Buschelman 
132841c8de11SBarry Smith   Works with MATAIJ and MATSBAIJ matrices
132924b6179bSKris Buschelman 
133024b6179bSKris Buschelman   Options Database Keys:
1331fb8376fbSHong Zhang + -mat_mumps_icntl_4 <0,...,4> - print level
133224b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
133364e6c443SBarry Smith . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guidec)
133424b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
133524b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
133694b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view
133724b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
133824b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
133924b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
134024b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
134124b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
134224b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
134324b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
134424b6179bSKris Buschelman 
134524b6179bSKris Buschelman   Level: beginner
134624b6179bSKris Buschelman 
134741c8de11SBarry Smith .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage
134841c8de11SBarry Smith 
134924b6179bSKris Buschelman M*/
135024b6179bSKris Buschelman 
13512877fffaSHong Zhang EXTERN_C_BEGIN
135235bd34faSBarry Smith #undef __FUNCT__
135335bd34faSBarry Smith #define __FUNCT__ "MatFactorGetSolverPackage_mumps"
135435bd34faSBarry Smith PetscErrorCode MatFactorGetSolverPackage_mumps(Mat A,const MatSolverPackage *type)
135535bd34faSBarry Smith {
135635bd34faSBarry Smith   PetscFunctionBegin;
13572692d6eeSBarry Smith   *type = MATSOLVERMUMPS;
135835bd34faSBarry Smith   PetscFunctionReturn(0);
135935bd34faSBarry Smith }
136035bd34faSBarry Smith EXTERN_C_END
136135bd34faSBarry Smith 
136235bd34faSBarry Smith EXTERN_C_BEGIN
1363bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI AIJ matrices */
13642877fffaSHong Zhang #undef __FUNCT__
1365bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_aij_mumps"
1366bccb9932SShri Abhyankar PetscErrorCode MatGetFactor_aij_mumps(Mat A,MatFactorType ftype,Mat *F)
13672877fffaSHong Zhang {
13682877fffaSHong Zhang   Mat            B;
13692877fffaSHong Zhang   PetscErrorCode ierr;
13702877fffaSHong Zhang   Mat_MUMPS      *mumps;
1371ace3abfcSBarry Smith   PetscBool      isSeqAIJ;
13722877fffaSHong Zhang 
13732877fffaSHong Zhang   PetscFunctionBegin;
13742877fffaSHong Zhang   /* Create the factorization matrix */
1375251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr);
13762877fffaSHong Zhang   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
13772877fffaSHong Zhang   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
13782877fffaSHong Zhang   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
1379bccb9932SShri Abhyankar   if (isSeqAIJ) {
13802877fffaSHong Zhang     ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
1381bccb9932SShri Abhyankar   } else {
1382bccb9932SShri Abhyankar     ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
1383bccb9932SShri Abhyankar   }
13842877fffaSHong Zhang 
138516ebf90aSShri Abhyankar   ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr);
1386*2205254eSKarl Rupp 
13872877fffaSHong Zhang   B->ops->view    = MatView_MUMPS;
138835bd34faSBarry Smith   B->ops->getinfo = MatGetInfo_MUMPS;
1389*2205254eSKarl Rupp 
139035bd34faSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
13915ccb76cbSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMumpsSetIcntl_C","MatMumpsSetIcntl_MUMPS",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
1392450b117fSShri Abhyankar   if (ftype == MAT_FACTOR_LU) {
1393450b117fSShri Abhyankar     B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
1394d5f3da31SBarry Smith     B->factortype            = MAT_FACTOR_LU;
1395bccb9932SShri Abhyankar     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqaij;
1396bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpiaij;
1397746480a1SHong Zhang     mumps->sym = 0;
1398dcd589f8SShri Abhyankar   } else {
139967877ebaSShri Abhyankar     B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
1400450b117fSShri Abhyankar     B->factortype                  = MAT_FACTOR_CHOLESKY;
1401bccb9932SShri Abhyankar     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqsbaij;
1402bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpisbaij;
14036fdc2a6dSBarry Smith     if (A->spd_set && A->spd) mumps->sym = 1;
14046fdc2a6dSBarry Smith     else                      mumps->sym = 2;
1405450b117fSShri Abhyankar   }
14062877fffaSHong Zhang 
14072877fffaSHong Zhang   mumps->isAIJ    = PETSC_TRUE;
1408bf0cc555SLisandro Dalcin   mumps->Destroy  = B->ops->destroy;
14092877fffaSHong Zhang   B->ops->destroy = MatDestroy_MUMPS;
14102877fffaSHong Zhang   B->spptr        = (void*)mumps;
1411*2205254eSKarl Rupp 
1412f697e70eSHong Zhang   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
1413746480a1SHong Zhang 
14142877fffaSHong Zhang   *F = B;
14152877fffaSHong Zhang   PetscFunctionReturn(0);
14162877fffaSHong Zhang }
14172877fffaSHong Zhang EXTERN_C_END
14182877fffaSHong Zhang 
1419bccb9932SShri Abhyankar 
14202877fffaSHong Zhang EXTERN_C_BEGIN
1421bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI SBAIJ matrices */
14222877fffaSHong Zhang #undef __FUNCT__
1423bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_sbaij_mumps"
1424bccb9932SShri Abhyankar PetscErrorCode MatGetFactor_sbaij_mumps(Mat A,MatFactorType ftype,Mat *F)
14252877fffaSHong Zhang {
14262877fffaSHong Zhang   Mat            B;
14272877fffaSHong Zhang   PetscErrorCode ierr;
14282877fffaSHong Zhang   Mat_MUMPS      *mumps;
1429ace3abfcSBarry Smith   PetscBool      isSeqSBAIJ;
14302877fffaSHong Zhang 
14312877fffaSHong Zhang   PetscFunctionBegin;
1432e7e72b3dSBarry Smith   if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix");
1433bccb9932SShri Abhyankar   if (A->rmap->bs > 1) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with block size > 1 with MUMPS Cholesky, use AIJ matrix instead");
1434251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr);
14352877fffaSHong Zhang   /* Create the factorization matrix */
14362877fffaSHong Zhang   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
14372877fffaSHong Zhang   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
14382877fffaSHong Zhang   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
143916ebf90aSShri Abhyankar   ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr);
1440bccb9932SShri Abhyankar   if (isSeqSBAIJ) {
1441bccb9932SShri Abhyankar     ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
1442*2205254eSKarl Rupp 
144316ebf90aSShri Abhyankar     mumps->ConvertToTriples = MatConvertToTriples_seqsbaij_seqsbaij;
1444dcd589f8SShri Abhyankar   } else {
1445bccb9932SShri Abhyankar     ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
1446*2205254eSKarl Rupp 
1447bccb9932SShri Abhyankar     mumps->ConvertToTriples = MatConvertToTriples_mpisbaij_mpisbaij;
1448bccb9932SShri Abhyankar   }
1449bccb9932SShri Abhyankar 
145067877ebaSShri Abhyankar   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
1451bccb9932SShri Abhyankar   B->ops->view                   = MatView_MUMPS;
1452*2205254eSKarl Rupp 
1453bccb9932SShri Abhyankar   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
1454f250808bSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMumpsSetIcntl_C","MatMumpsSetIcntl",MatMumpsSetIcntl);CHKERRQ(ierr);
1455*2205254eSKarl Rupp 
1456f4762488SHong Zhang   B->factortype = MAT_FACTOR_CHOLESKY;
14576fdc2a6dSBarry Smith   if (A->spd_set && A->spd) mumps->sym = 1;
14586fdc2a6dSBarry Smith   else                      mumps->sym = 2;
1459a214ac2aSShri Abhyankar 
1460bccb9932SShri Abhyankar   mumps->isAIJ    = PETSC_FALSE;
1461bf0cc555SLisandro Dalcin   mumps->Destroy  = B->ops->destroy;
1462f3c0ef26SHong Zhang   B->ops->destroy = MatDestroy_MUMPS;
14632877fffaSHong Zhang   B->spptr        = (void*)mumps;
1464*2205254eSKarl Rupp 
1465f697e70eSHong Zhang   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
1466746480a1SHong Zhang 
14672877fffaSHong Zhang   *F = B;
14682877fffaSHong Zhang   PetscFunctionReturn(0);
14692877fffaSHong Zhang }
14702877fffaSHong Zhang EXTERN_C_END
147197969023SHong Zhang 
1472450b117fSShri Abhyankar EXTERN_C_BEGIN
1473450b117fSShri Abhyankar #undef __FUNCT__
1474bccb9932SShri Abhyankar #define __FUNCT__ "MatGetFactor_baij_mumps"
1475bccb9932SShri Abhyankar PetscErrorCode MatGetFactor_baij_mumps(Mat A,MatFactorType ftype,Mat *F)
147667877ebaSShri Abhyankar {
147767877ebaSShri Abhyankar   Mat            B;
147867877ebaSShri Abhyankar   PetscErrorCode ierr;
147967877ebaSShri Abhyankar   Mat_MUMPS      *mumps;
1480ace3abfcSBarry Smith   PetscBool      isSeqBAIJ;
148167877ebaSShri Abhyankar 
148267877ebaSShri Abhyankar   PetscFunctionBegin;
148367877ebaSShri Abhyankar   /* Create the factorization matrix */
1484251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);CHKERRQ(ierr);
148567877ebaSShri Abhyankar   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
148667877ebaSShri Abhyankar   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
148767877ebaSShri Abhyankar   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
1488bccb9932SShri Abhyankar   if (isSeqBAIJ) {
148967877ebaSShri Abhyankar     ierr = MatSeqBAIJSetPreallocation(B,A->rmap->bs,0,PETSC_NULL);CHKERRQ(ierr);
1490bccb9932SShri Abhyankar   } else {
149167877ebaSShri Abhyankar     ierr = MatMPIBAIJSetPreallocation(B,A->rmap->bs,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
1492bccb9932SShri Abhyankar   }
1493450b117fSShri Abhyankar 
149467877ebaSShri Abhyankar   ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr);
1495450b117fSShri Abhyankar   if (ftype == MAT_FACTOR_LU) {
1496450b117fSShri Abhyankar     B->ops->lufactorsymbolic = MatLUFactorSymbolic_BAIJMUMPS;
1497450b117fSShri Abhyankar     B->factortype            = MAT_FACTOR_LU;
1498bccb9932SShri Abhyankar     if (isSeqBAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqbaij_seqaij;
1499bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpibaij_mpiaij;
1500746480a1SHong Zhang     mumps->sym = 0;
1501f23aa3ddSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc BAIJ matrices with MUMPS Cholesky, use SBAIJ or AIJ matrix instead\n");
1502bccb9932SShri Abhyankar 
1503450b117fSShri Abhyankar   B->ops->view = MatView_MUMPS;
1504*2205254eSKarl Rupp 
1505450b117fSShri Abhyankar   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
15065ccb76cbSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMumpsSetIcntl_C","MatMumpsSetIcntl_MUMPS",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
1507450b117fSShri Abhyankar 
1508450b117fSShri Abhyankar   mumps->isAIJ    = PETSC_TRUE;
1509bf0cc555SLisandro Dalcin   mumps->Destroy  = B->ops->destroy;
1510450b117fSShri Abhyankar   B->ops->destroy = MatDestroy_MUMPS;
1511450b117fSShri Abhyankar   B->spptr        = (void*)mumps;
1512*2205254eSKarl Rupp 
1513f697e70eSHong Zhang   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
1514746480a1SHong Zhang 
1515450b117fSShri Abhyankar   *F = B;
1516450b117fSShri Abhyankar   PetscFunctionReturn(0);
1517450b117fSShri Abhyankar }
1518450b117fSShri Abhyankar EXTERN_C_END
1519a214ac2aSShri Abhyankar 
1520